mappedPatchBase.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2015-2023 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "mappedPatchBase.H"
31 #include "ListListOps.H"
32 #include "meshSearchMeshObject.H"
34 #include "meshTools.H"
35 #include "OFstream.H"
36 #include "Random.H"
37 #include "treeDataFace.H"
38 #include "treeDataPoint.H"
39 #include "indexedOctree.H"
40 #include "polyMesh.H"
41 #include "polyPatch.H"
42 #include "Time.H"
43 #include "mapDistribute.H"
44 #include "SubField.H"
45 #include "triangle.H"
46 #include "syncTools.H"
47 #include "treeDataCell.H"
48 #include "DynamicField.H"
49 #include "faceAreaWeightAMI.H"
50 #include "OTstream.H"
51 
52 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
53 
54 namespace Foam
55 {
56  defineTypeNameAndDebug(mappedPatchBase, 0);
57 }
58 
59 
60 const Foam::Enum
61 <
63 >
65 ({
66  { sampleMode::NEARESTCELL, "nearestCell" },
67  { sampleMode::NEARESTPATCHFACE, "nearestPatchFace" },
68  { sampleMode::NEARESTPATCHFACEAMI, "nearestPatchFaceAMI" },
69  { sampleMode::NEARESTPATCHPOINT, "nearestPatchPoint" },
70  { sampleMode::NEARESTFACE, "nearestFace" },
71  { sampleMode::NEARESTONLYCELL, "nearestOnlyCell" },
72 });
73 
74 
75 const Foam::Enum
76 <
78 >
80 ({
81  { offsetMode::UNIFORM, "uniform" },
82  { offsetMode::NONUNIFORM, "nonuniform" },
83  { offsetMode::NORMAL, "normal" },
84 });
85 
86 
87 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
88 
90 (
91  const dictionary& dict
92 )
93 {
94  if (dict.found("sampleDatabase"))
95  {
96  if (dict.get<bool>("sampleDatabase"))
97  {
99  (
100  dict.getOrDefault<fileName>
101  (
102  "sampleDatabasePath",
104  )
105  );
106  }
107  }
108  else if (dict.found("sampleDatabasePath"))
109  {
110  return autoPtr<fileName>::New(dict.get<fileName>("sampleDatabasePath"));
111  }
112 
113  return nullptr;
114 }
115 
116 
118 {
119  if (sameWorld())
120  {
121  return true;
122  }
123 
124  const Time& runTime = patch_.boundaryMesh().mesh().time();
125  return const_cast<multiWorldConnections&>
126  (
128  ).addConnectionByName(sampleWorld_);
129 }
130 
131 
133 {
134  if (sameWorld())
135  {
136  return UPstream::commWorld();
137  }
138 
139  const Time& runTime = patch_.boundaryMesh().mesh().time();
140  return
142 }
143 
144 
146 (
147  const polyPatch& pp
148 ) const
149 {
150  const polyMesh& mesh = pp.boundaryMesh().mesh();
151 
152  // Force construction of min-tet decomp
153  (void)mesh.tetBasePtIs();
154 
155  // Initialise to face-centre
156  auto tfacePoints = tmp<pointField>::New(patch_.size());
157  auto& facePoints = tfacePoints.ref();
158 
159  forAll(pp, facei)
160  {
161  facePoints[facei] = facePoint
162  (
163  mesh,
164  pp.start()+facei,
166  ).point();
167  }
168 
169  return tfacePoints;
170 }
171 
172 
174 (
175  const label mySampleWorld, // Wanted world
176  const pointField& facePoints,
177  pointField& samples, // All samples
178  labelList& patchFaceWorlds, // Per sample: wanted world
179  labelList& patchFaceProcs, // Per sample: originating processor
180  labelList& patchFaces, // Per sample: originating patchFace index
181  pointField& patchFc // Per sample: originating centre
182 ) const
183 {
184  DebugInFunction << nl;
185 
186  const label myComm = getCommunicator(); // Get or create
187  const label myRank = Pstream::myProcNo(myComm);
188  const label nProcs = Pstream::nProcs(myComm);
189 
190  const label oldWarnComm = UPstream::commWarn(myComm);
191 
192  if (debug & 2)
193  {
194  Perr<< "patch: " << patch_.name()
195  << "[rank=" << myRank << " procs=" << nProcs
196  << " comm=" << myComm << "] collect samples" << endl;
197  }
198 
199  // Collect all sample points and the faces they come from.
200  {
201  List<pointField> globalFc(nProcs);
202  globalFc[myRank] = facePoints;
203 
204  Pstream::allGatherList(globalFc, Pstream::msgType(), myComm);
205 
206  // Rework into straight list
207  patchFc = ListListOps::combine<pointField>
208  (
209  globalFc,
210  accessOp<pointField>()
211  );
212  }
213 
214  {
215  List<pointField> globalSamples(nProcs);
216  globalSamples[myRank] = samplePoints(facePoints);
217  Pstream::allGatherList(globalSamples, Pstream::msgType(), myComm);
218 
219  // Rework into straight list
220  samples = ListListOps::combine<pointField>
221  (
222  globalSamples,
223  accessOp<pointField>()
224  );
225  }
226 
227  {
228  labelListList globalFaces(nProcs);
229  globalFaces[myRank] = identity(patch_.size());
230  // Distribute to all processors
231  Pstream::allGatherList(globalFaces, Pstream::msgType(), myComm);
232 
233  patchFaces = ListListOps::combine<labelList>
234  (
235  globalFaces,
236  accessOp<labelList>()
237  );
238  }
239 
240  {
241  labelList procToWorldIndex
242  (
243  UPstream::allGatherValues<label>(mySampleWorld, myComm)
244  );
245  labelList nPerProc
246  (
247  UPstream::allGatherValues<label>(patch_.size(), myComm)
248  );
249 
250  patchFaceWorlds.resize(patchFaces.size());
251  patchFaceProcs.resize(patchFaces.size());
252 
253  label sampleI = 0;
254  forAll(nPerProc, proci)
255  {
256  for (label i = 0; i < nPerProc[proci]; i++)
257  {
258  patchFaceWorlds[sampleI] = procToWorldIndex[proci];
259  patchFaceProcs[sampleI] = proci;
260  sampleI++;
261  }
262  }
263  }
265  // Restore communicator settings
266  UPstream::commWarn(oldWarnComm);
267 }
268 
269 
271 (
272  const sampleMode mode,
273 
274  const label mySampleWorld, // local world to sample == my own world
275  const word& sampleRegion, // local region to sample
276  const word& samplePatch, // local patch to sample
277 
278  const pointField& samples,
279  List<nearInfoWorld>& nearest
280 ) const
281 {
282  DebugInFunction << nl;
283 
284  const label myComm = getCommunicator(); // Get or create
285 
286  // Find the local cell containing the samples
287  const label myRank = Pstream::myProcNo(myComm);
288 
289  // Lookup the correct region
290  const polyMesh& mesh = lookupMesh(sampleRegion);
291 
292  // All the info for nearest. Construct to miss
293  nearest.setSize(samples.size());
294  nearInfoWorld miss;
295  {
296  miss.first().second() = Tuple2<scalar, label>(Foam::sqr(GREAT), -1);
297  miss.second() = -1; // set world to be ignored
298  }
299  nearest = miss;
300 
301  switch (mode)
302  {
303  case NEARESTCELL:
304  {
305  if (samplePatch.size() && samplePatch != "none")
306  {
308  << "No need to supply a patch name when in "
309  << sampleModeNames_[mode] << " mode." << exit(FatalError);
310  }
311 
312  //- Note: face-diagonal decomposition
313  const indexedOctree<Foam::treeDataCell>& tree = mesh.cellTree();
314 
315  forAll(samples, sampleI)
316  {
317  const point& sample = samples[sampleI];
318  nearInfoWorld& near = nearest[sampleI];
319 
320  label celli = tree.findInside(sample);
321 
322  if (celli == -1)
323  {
324  near.first().second().first() = Foam::sqr(GREAT);
325  near.first().second().second() = myRank;
326  near.second() = mySampleWorld;
327  }
328  else
329  {
330  const point& cc = mesh.cellCentres()[celli];
331 
332  near.first().first() = pointIndexHit
333  (
334  true,
335  cc,
336  celli
337  );
338  near.first().second().first() = sample.distSqr(cc);
339  near.first().second().second() = myRank;
340  near.second() = mySampleWorld;
341  }
342  }
343  break;
344  }
345 
346  case NEARESTONLYCELL:
347  {
348  if (samplePatch.size() && samplePatch != "none")
349  {
351  << "No need to supply a patch name when in "
352  << sampleModeNames_[mode] << " mode." << exit(FatalError);
353  }
354 
355  //- Note: face-diagonal decomposition
356  const indexedOctree<Foam::treeDataCell>& tree = mesh.cellTree();
357 
358  forAll(samples, sampleI)
359  {
360  const point& sample = samples[sampleI];
361  nearInfoWorld& near = nearest[sampleI];
362 
363  near.first().first() = tree.findNearest(sample, sqr(GREAT));
364 
365  near.first().second().first() =
366  near.first().first().hitPoint().distSqr(sample);
367 
368  near.first().second().second() = myRank;
369  near.second() = mySampleWorld;
370  }
371  break;
372  }
373 
374  case NEARESTPATCHFACE:
375  {
376  Random rndGen(123456);
377 
378  const polyPatch& pp = lookupPatch(sampleRegion, samplePatch);
379 
380  if (pp.empty())
381  {
382  forAll(samples, sampleI)
383  {
384  nearInfoWorld& near = nearest[sampleI];
385  near.first().second().first() = Foam::sqr(GREAT);
386  near.first().second().second() = myRank;
387  near.second() = mySampleWorld;
388  }
389  }
390  else
391  {
392  treeBoundBox patchBb
393  (
394  treeBoundBox(pp.points(), pp.meshPoints()).extend
395  (
396  rndGen,
397  1e-4,
398  ROOTVSMALL
399  )
400  );
401 
402  indexedOctree<treeDataFace> boundaryTree
403  (
404  treeDataFace(mesh, pp.range()), // Patch faces
405 
406  patchBb, // overall search domain
407  8, // maxLevel
408  10, // leafsize
409  3.0 // duplicity
410  );
411  const auto& treeData = boundaryTree.shapes();
412 
413  forAll(samples, sampleI)
414  {
415  const point& sample = samples[sampleI];
416 
417  nearInfoWorld& near = nearest[sampleI];
418  pointIndexHit& nearInfo = near.first().first();
419  nearInfo = boundaryTree.findNearest
420  (
421  sample,
422  patchBb.magSqr()
423  );
424 
425  if (!nearInfo.hit())
426  {
427  near.first().second().first() = Foam::sqr(GREAT);
428  near.first().second().second() = myRank;
429  near.second() = mySampleWorld;
430  }
431  else
432  {
433  const point& fc = treeData.centre(nearInfo.index());
434 
435  nearInfo.setPoint(fc);
436  near.first().second().first() = sample.distSqr(fc);
437  near.first().second().second() = myRank;
438  near.second() = mySampleWorld;
439  }
440  }
441  }
442  break;
443  }
444 
445  case NEARESTPATCHPOINT:
446  {
447  Random rndGen(123456);
448 
449  const polyPatch& pp = lookupPatch(sampleRegion, samplePatch);
450 
451  if (pp.empty())
452  {
453  forAll(samples, sampleI)
454  {
455  nearInfoWorld& near = nearest[sampleI];
456  near.first().second().first() = Foam::sqr(GREAT);
457  near.first().second().second() = myRank;
458  near.second() = mySampleWorld;
459  }
460  }
461  else
462  {
463  // patch (local) points
464  treeBoundBox patchBb
465  (
466  treeBoundBox(pp.points(), pp.meshPoints()).extend
467  (
468  rndGen,
469  1e-4,
470  ROOTVSMALL
471  )
472  );
473 
474  indexedOctree<treeDataPoint> boundaryTree
475  (
476  // Patch points
477  treeDataPoint(mesh.points(), pp.meshPoints()),
478 
479  patchBb, // overall search domain
480  8, // maxLevel
481  10, // leafsize
482  3.0 // duplicity
483  );
484 
485  forAll(samples, sampleI)
486  {
487  const point& sample = samples[sampleI];
488 
489  nearInfoWorld& near = nearest[sampleI];
490  pointIndexHit& nearInfo = near.first().first();
491  nearInfo = boundaryTree.findNearest
492  (
493  sample,
494  patchBb.magSqr()
495  );
496 
497  if (!nearInfo.hit())
498  {
499  near.first().second().first() = Foam::sqr(GREAT);
500  near.first().second().second() = myRank;
501  near.second() = mySampleWorld;
502  }
503  else
504  {
505  const point& pt = nearInfo.point();
506 
507  near.first().second().first() = sample.distSqr(pt);
508  near.first().second().second() = myRank;
509  near.second() = mySampleWorld;
510  }
511  }
512  }
513  break;
514  }
515 
516  case NEARESTFACE:
517  {
518  if (samplePatch.size() && samplePatch != "none")
519  {
521  << "No need to supply a patch name when in "
522  << sampleModeNames_[mode] << " mode." << exit(FatalError);
523  }
524 
525  //- Note: face-diagonal decomposition
526  const meshSearchMeshObject& meshSearchEngine =
528 
529  forAll(samples, sampleI)
530  {
531  const point& sample = samples[sampleI];
532  nearInfoWorld& near = nearest[sampleI];
533 
534  label facei = meshSearchEngine.findNearestFace(sample);
535 
536  if (facei == -1)
537  {
538  near.first().second().first() = Foam::sqr(GREAT);
539  near.first().second().second() = myRank;
540  near.second() = mySampleWorld;
541  }
542  else
543  {
544  const point& fc = mesh.faceCentres()[facei];
545 
546  near.first().first() = pointIndexHit(true, fc, facei);
547  near.first().second().first() = sample.distSqr(fc);
548  near.first().second().second() = myRank;
549  near.second() = mySampleWorld;
550  }
551  }
552  break;
553  }
554 
555  case NEARESTPATCHFACEAMI:
556  {
557  // nothing to do here
558  return;
559  }
560 
561  default:
562  {
564  << "problem." << abort(FatalError);
565  }
566  }
567 }
568 
569 
571 (
572  const sampleMode mode,
573  const label myWorld,
574  const pointField& samples,
575  const labelList& wantedWorlds,
576  const labelList& origProcs,
577 
578  labelList& sampleProcs,
579  labelList& sampleIndices,
580  pointField& sampleLocations
581 ) const
582 {
583  DebugInFunction << nl;
584 
585  // Find the processor/cell containing the samples. Does not account
586  // for samples being found in two processors.
587 
588  const label myComm = getCommunicator(); // Get or create
589  const label myRank = Pstream::myProcNo(myComm);
590  const label nProcs = Pstream::nProcs(myComm);
591 
592  const label oldWarnComm = UPstream::commWarn(myComm);
593 
594  wordList samplePatches(nProcs);
595  {
596  samplePatches[myRank] = samplePatch_;
597  Pstream::allGatherList(samplePatches, Pstream::msgType(), myComm);
598  }
599  wordList sampleRegions(nProcs);
600  {
601  sampleRegions[myRank] = sampleRegion_;
602  Pstream::allGatherList(sampleRegions, Pstream::msgType(), myComm);
603  }
604 
605 
606  // Find all the info for nearest
607  List<nearInfoWorld> nearest(samples.size());
608  forAll(nearest, samplei)
609  {
610  nearest[samplei].first() = nearInfo
611  (
612  pointIndexHit(),
613  Tuple2<scalar, label>(Foam::sqr(GREAT), -1)
614  );
615  nearest[samplei].second() = wantedWorlds[samplei];
616  }
617 
618 
619  // Extract samples to search for locally
620  {
621  DynamicList<label> localMap(samples.size());
622  forAll(wantedWorlds, samplei)
623  {
624  if (wantedWorlds[samplei] == myWorld)
625  {
626  localMap.append(samplei);
627  }
628  }
629 
630  if (localMap.size())
631  {
632  pointField localSamples(samples, localMap);
633  labelList localOrigProcs(origProcs, localMap);
634 
635  //Assume single patch to sample for now
636  const word localOrigPatch(samplePatches[localOrigProcs[0]]);
637  const word localOrigRegion(sampleRegions[localOrigProcs[0]]);
638  List<nearInfoWorld> localNearest(localSamples.size());
639 
640  if (debug)
641  {
642  Pout<< "Searching locally for " << localSamples.size()
643  << " samples on region:" << localOrigRegion
644  << " on patch:" << localOrigPatch << endl;
645  }
646  findLocalSamples
647  (
648  mode,
649  myWorld,
650  localOrigRegion,
651  localOrigPatch,
652  localSamples,
653  localNearest
654  );
655  UIndirectList<nearInfoWorld>(nearest, localMap) = localNearest;
656  }
657  }
658 
659 
660  // Find nearest - globally consistent
662  (
663  nearest,
666  myComm
667  );
668 
669  //if (debug)
670  //{
671  // Pout<< "** After combining:" << endl;
672  // forAll(nearest, samplei)
673  // {
674  // Pout<< " sample:" << samples[samplei]
675  // << " origating from proc:" << origProcs[samplei]
676  // << " to be found on world:" << wantedWorlds[samplei] << nl
677  // << " found on world:" << nearest[samplei].second() << nl
678  // << " found on proc:"
679  // << nearest[samplei].first().second().second() << nl
680  // << " found on patchfacei:"
681  // << nearest[samplei].first().first().index() << nl
682  // << " found at location:"
683  // << nearest[samplei].first().first().point() << nl;
684  // }
685  // Pout<< endl;
686  //}
687 
688  // Convert back into proc+local index
689  sampleProcs.setSize(samples.size());
690  sampleIndices.setSize(samples.size());
691  sampleLocations.setSize(samples.size());
692 
693  forAll(nearest, sampleI)
694  {
695  const nearInfo& ni = nearest[sampleI].first();
696 
697  if (!ni.first().hit())
698  {
699  sampleProcs[sampleI] = -1;
700  sampleIndices[sampleI] = -1;
701  sampleLocations[sampleI] = vector::max;
702  }
703  else
704  {
705  sampleProcs[sampleI] = ni.second().second();
706  sampleIndices[sampleI] = ni.first().index();
707  sampleLocations[sampleI] = ni.first().point();
708  }
709  }
710 
711  // Return communicator settings
712  UPstream::commWarn(oldWarnComm);
713 }
714 
717 {
719 
720  static bool hasWarned = false;
721  if (mapPtr_)
722  {
724  << "Mapping already calculated" << exit(FatalError);
725  }
726 
727  DebugInFunction << nl;
728 
729 
730  // Early construction of tetBasePtIs since does parallel
731  // communication which might conflict with inter-world communicator
732  // (for multi-world operation)
733  (void)patch_.boundaryMesh().mesh().tetBasePtIs();
734 
735  const label myComm = getCommunicator(); // Get or create
736 
738  //if (sampleDatabase() && !sameWorld())
739  //{
740  // const word syncName("syncObjects");
741  // const polyMesh& mesh = patch_.boundaryMesh().mesh();
742  // const Time& runTime = mesh.time();
743  // const functionObjectList& fos = runTime.functionObjects();
744  //
745  // forAll(fos, i)
746  // {
747  // Pout<< "** FO:" << fos[i].name() << " tpye:" << fos[i].type()
748  // << endl;
749  // }
750  //
751  // if (fos.findObjectID(syncName) == -1)
752  // {
753  // //FatalErrorInFunction
754  // WarningInFunction
755  // << "Did not detect functionObject " << syncName
756  // << ". This is used to synchronise data inbetween worlds"
757  // << endl;
758  // }
759  //}
760 
761 
762  // Get points on face (since cannot use face-centres - might be off
763  // face-diagonal decomposed tets.
764  tmp<pointField> patchPoints(facePoints(patch_));
765 
766  // Get offsetted points
767  const pointField offsettedPoints(samplePoints(patchPoints()));
768 
769  // Do a sanity check - am I sampling my own patch?
770  // This only makes sense for a non-zero offset.
771  bool sampleMyself =
772  (
773  mode_ == NEARESTPATCHFACE
774  && sampleWorld() == UPstream::myWorld()
775  && sampleRegion() == patch_.boundaryMesh().mesh().name()
776  && samplePatch() == patch_.name()
777  );
778 
779  if (sampleMyself)
780  {
781  // Check offset
782  vectorField d(offsettedPoints-patchPoints());
783  bool coincident = (gAverage(mag(d)) <= ROOTVSMALL);
784 
785  if (sampleMyself && coincident)
786  {
788  << "Invalid offset " << d << endl
789  << "Offset is the vector added to the patch face centres to"
790  << " find the patch face supplying the data." << endl
791  << "Setting it to " << d
792  << " on the same patch, on the same region"
793  << " will find the faces themselves which does not make sense"
794  << " for anything but testing." << endl
795  << "patch:" << patch_.name() << endl
796  << "sampleRegion:" << sampleRegion() << endl
797  << "mode:" << sampleModeNames_[mode_] << endl
798  << "samplePatch:" << samplePatch() << endl
799  << "offsetMode:" << offsetModeNames_[offsetMode_] << endl;
800  }
801  }
802 
803  // Get local world and world-to-sample in index form
804  const label myWorld = Pstream::myWorldID();
805  const label mySampleWorld = Pstream::allWorlds().find(sampleWorld_);
806 
807 
808  // Get global list of all samples and the processor and face they come from.
809  pointField samples; // coordinates
810  labelList patchFaceWorlds; // world to sample
811  labelList patchFaceProcs; // originating processor
812  labelList patchFaces; // originating face on processor patch
813  pointField patchFc; // originating face centre
814  collectSamples
815  (
816  mySampleWorld, // world I want to sample
817  patchPoints,
818 
819  samples,
820  patchFaceWorlds,
821  patchFaceProcs,
822  patchFaces,
823  patchFc
824  );
825 
826  //if (debug)
827  //{
828  // forAll(samples, samplei)
829  // {
830  // Pout<< " sample:" << samples[samplei]
831  // << " origating from proc:" << patchFaceProcs[samplei]
832  // << " face:" << patchFaces[samplei]
833  // << " to be found on world:" << patchFaceWorlds[samplei] << nl;
834  // }
835  //}
836 
837  // Find processor and cell/face samples are in and actual location.
838  labelList sampleProcs;
839  labelList sampleIndices;
840  pointField sampleLocations;
841  findSamples
842  (
843  mode_,
844  myWorld, // my world (in index form)
845  samples,
846  patchFaceWorlds,
847  patchFaceProcs,
848 
849  sampleProcs,
850  sampleIndices,
851  sampleLocations
852  );
853 
854  // Check for samples that were not found. This will only happen for
855  // NEARESTCELL since finds cell containing a location
856  if (mode_ == NEARESTCELL)
857  {
858  label nNotFound = 0;
859  forAll(sampleProcs, sampleI)
860  {
861  if (sampleProcs[sampleI] == -1)
862  {
863  nNotFound++;
864  }
865  }
866  reduce(nNotFound, sumOp<label>(), Pstream::msgType(), myComm);
867 
868  if (nNotFound > 0)
869  {
870  if (!hasWarned)
871  {
873  << "Did not find " << nNotFound
874  << " out of " << sampleProcs.size() << " total samples."
875  << " Sampling these on owner cell centre instead." << endl
876  << "On patch " << patch_.name()
877  << " on region " << sampleRegion()
878  << " in mode " << sampleModeNames_[mode_] << endl
879  << "with offset mode " << offsetModeNames_[offsetMode_]
880  << ". Suppressing further warnings from " << type() << endl;
881 
882  hasWarned = true;
883  }
884 
885  // Collect the samples that cannot be found
886  DynamicList<label> subMap;
887  forAll(sampleProcs, sampleI)
888  {
889  if (sampleProcs[sampleI] == -1)
890  {
891  subMap.append(sampleI);
892  }
893  }
894 
895  // And re-search for pure nearest (should not fail)
896  labelList subSampleProcs;
897  labelList subSampleIndices;
898  pointField subSampleLocations;
899  findSamples
900  (
901  NEARESTONLYCELL,
902  myWorld, // my world (in index form)
903 
904  pointField(samples, subMap),
905  UIndirectList<label>(patchFaceWorlds, subMap)(),
906  UIndirectList<label>(patchFaceProcs, subMap)(),
907 
908  subSampleProcs,
909  subSampleIndices,
910  subSampleLocations
911  );
912 
913  // Insert
914  labelUIndList(sampleProcs, subMap) = subSampleProcs;
915  labelUIndList(sampleIndices, subMap) = subSampleIndices;
916  UIndirectList<point>(sampleLocations, subMap) = subSampleLocations;
917  }
918  }
919 
920  // Now we have all the data we need:
921  // - where sample originates from (so destination when mapping):
922  // patchFaces, patchFaceProcs.
923  // - cell/face sample is in (so source when mapping)
924  // sampleIndices, sampleProcs.
925 
926  if (Pstream::master(myComm))
927  {
928  forAll(samples, i)
929  {
930  if (sampleProcs[i] == -1)
931  {
932  FatalErrorInFunction << "did not find sample "
933  << patchFc[i] << " on patch " << patch_.name()
934  << " on region "
935  << patch_.boundaryMesh().mesh().name()
936  << " on processor " << patchFaceProcs[i]
937  << exit(FatalError);
938  }
939  }
940  }
941 
942 
943  if (debug && Pstream::master(myComm))
944  {
945  //forAll(samples, i)
946  //{
947  // Pout<< i << " need data in region "
948  // << patch_.boundaryMesh().mesh().name()
949  // << " for proc:" << patchFaceProcs[i]
950  // << " face:" << patchFaces[i]
951  // << " at:" << patchFc[i] << endl
952  // << "Found data in region " << sampleRegion()
953  // << " at proc:" << sampleProcs[i]
954  // << " face:" << sampleIndices[i]
955  // << " at:" << sampleLocations[i]
956  // << nl << endl;
957  //}
958 
959  OFstream str
960  (
961  patch_.boundaryMesh().mesh().time().path()
962  / patch_.name()
963  + "_mapped.obj"
964  );
965  Pout<< "Dumping mapping as lines from patch faceCentres to"
966  << " sampled cell/faceCentres/points to file " << str.name()
967  << endl;
968 
969  label vertI = 0;
970 
971  forAll(patchFc, i)
972  {
973  meshTools::writeOBJ(str, patchFc[i]);
974  vertI++;
975  meshTools::writeOBJ(str, sampleLocations[i]);
976  vertI++;
977  str << "l " << vertI-1 << ' ' << vertI << nl;
978  }
979  }
980 
981  // Determine schedule.
982  mapPtr_.reset(new mapDistribute(sampleProcs, patchFaceProcs, myComm));
983 
984  // Rework the schedule from indices into samples to cell data to send,
985  // face data to receive.
986 
987  labelListList& subMap = mapPtr_().subMap();
988  labelListList& constructMap = mapPtr_().constructMap();
989 
990  forAll(subMap, proci)
991  {
992  subMap[proci] = labelUIndList(sampleIndices, subMap[proci]);
993  constructMap[proci] = labelUIndList(patchFaces, constructMap[proci]);
994 
995  if (debug)
996  {
997  Pout<< "To proc:" << proci << " sending values of cells/faces:"
998  << subMap[proci] << endl;
999  Pout<< "From proc:" << proci
1000  << " receiving values of patch faces:"
1001  << constructMap[proci] << endl;
1002  }
1003  }
1004 
1005 
1006  // Redo constructSize
1007  mapPtr_().constructSize() = patch_.size();
1008 
1009  if (debug)
1010  {
1011  // Check that all elements get a value.
1012  bitSet used(patch_.size());
1013  forAll(constructMap, proci)
1014  {
1015  const labelList& map = constructMap[proci];
1016 
1017  forAll(map, i)
1018  {
1019  label facei = map[i];
1020 
1021  if (used.test(facei))
1022  {
1024  << "On patch " << patch_.name()
1025  << " patchface " << facei
1026  << " is assigned to more than once."
1027  << abort(FatalError);
1028  }
1029  else
1030  {
1031  used.set(facei);
1032  }
1033  }
1034  }
1035  forAll(used, facei)
1036  {
1037  if (!used.test(facei))
1038  {
1040  << "On patch " << patch_.name()
1041  << " patchface " << facei
1042  << " is never assigned to."
1043  << abort(FatalError);
1044  }
1045  }
1046  }
1047 
1048  updateMeshTime().setUpToDate();
1049  if (sameWorld())
1050  {
1051  updateSampleMeshTime().setUpToDate();
1052  }
1053 }
1054 
1057 const
1058 {
1059  const word surfType(surfDict_.getOrDefault<word>("type", "none"));
1060 
1061  if (!surfPtr_ && surfType != "none")
1062  {
1063  word surfName(surfDict_.getOrDefault("name", patch_.name()));
1064 
1065  const polyMesh& mesh = patch_.boundaryMesh().mesh();
1066 
1067  surfPtr_ =
1069  (
1070  surfType,
1071  IOobject
1072  (
1073  surfName,
1074  mesh.time().constant(),
1075  "triSurface",
1076  mesh,
1079  ),
1080  surfDict_
1081  );
1082  }
1083 
1084  return surfPtr_;
1085 }
1086 
1088 void Foam::mappedPatchBase::calcAMI() const
1089 {
1090  if (AMIPtr_->upToDate())
1091  {
1093  << "AMI already up-to-date"
1094  << endl;
1095 
1096  return;
1097  }
1098 
1099  DebugInFunction << nl;
1100 
1101  const label myComm = getCommunicator(); // Get or create
1102 
1103  // Pre-calculate surface (if any)
1104  const auto& surf = surfPtr();
1105 
1106  // Check if running locally
1107  if (sampleWorld_.empty() || sameWorld())
1108  {
1109  const polyPatch& nbr = samplePolyPatch();
1110 
1111  // Transform neighbour patch to local system
1112  const pointField nbrPoints(samplePoints(nbr.localPoints()));
1113 
1114  const primitivePatch nbrPatch0
1115  (
1117  (
1118  nbr.localFaces(),
1119  nbr.size()
1120  ),
1121  nbrPoints
1122  );
1123 
1124 
1125  if (debug)
1126  {
1127  OFstream os(patch_.name() + "_neighbourPatch-org.obj");
1128  meshTools::writeOBJ(os, samplePolyPatch().localFaces(), nbrPoints);
1129 
1130  OFstream osN(patch_.name() + "_neighbourPatch-trans.obj");
1131  meshTools::writeOBJ(osN, nbrPatch0, nbrPoints);
1132 
1133  OFstream osO(patch_.name() + "_ownerPatch.obj");
1134  meshTools::writeOBJ(osO, patch_.localFaces(), patch_.localPoints());
1135  }
1136 
1137  // Construct/apply AMI interpolation to determine addressing and
1138  // weights.
1139 
1140  // Change to use inter-world communicator
1141  const label oldWarnComm = UPstream::commWarn(myComm);
1142  const label oldWorldComm = UPstream::commWorld(myComm);
1143 
1144  AMIPtr_->calculate(patch_, nbrPatch0, surf);
1145 
1146  // Restore communicator settings
1147  UPstream::commWarn(oldWarnComm);
1148  UPstream::commWorld(oldWorldComm);
1149  }
1150  else
1151  {
1152  faceList dummyFaces;
1153  pointField dummyPoints;
1154  const primitivePatch dummyPatch
1155  (
1156  SubList<face>(dummyFaces),
1157  dummyPoints
1158  );
1159 
1160  // Change to use inter-world communicator
1161  AMIPtr_->comm(myComm);
1162 
1163  const label oldWarnComm = UPstream::commWarn(myComm);
1164  const label oldWorldComm = UPstream::commWorld(myComm);
1165 
1166  if (masterWorld())
1167  {
1168  // Construct/apply AMI interpolation to determine addressing
1169  // and weights. Have patch_ for src faces, 0 faces for the
1170  // target side
1171  AMIPtr_->calculate(patch_, dummyPatch, surf);
1172  }
1173  else
1174  {
1175  // Construct/apply AMI interpolation to determine addressing
1176  // and weights. Have 0 faces for src side, patch_ for the tgt
1177  // side
1178  AMIPtr_->calculate(dummyPatch, patch_, surf);
1179  }
1180  // Now the AMI addressing/weights will be from src side (on masterWorld
1181  // processors) to tgt side (on other processors)
1182 
1183  // Restore communicator settings
1184  UPstream::commWarn(oldWarnComm);
1185  UPstream::commWorld(oldWorldComm);
1186  }
1187 }
1188 
1189 
1191 (
1192  const objectRegistry& obr,
1193  const wordList& names,
1194  const label index
1195 )
1196 {
1197  const objectRegistry& sub = obr.subRegistry(names[index], true);
1198  if (index == names.size()-1)
1199  {
1200  return sub;
1201  }
1202  else
1203  {
1204  return subRegistry(sub, names, index+1);
1205  }
1206 }
1207 
1208 
1209 // * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * * * * * //
1212 :
1213  patch_(pp),
1214  sampleWorld_(),
1215  sampleRegion_(patch_.boundaryMesh().mesh().name()),
1216  mode_(NEARESTPATCHFACE),
1217  samplePatch_(),
1218  coupleGroup_(),
1219  sampleDatabasePtr_(),
1220  offsetMode_(UNIFORM),
1221  offset_(Zero),
1222  offsets_(pp.size(), offset_),
1223  distance_(0),
1224  communicator_(-1), // Demand-driven (cached value)
1225  sameRegion_(true),
1226  mapPtr_(nullptr),
1227  AMIReverse_(false),
1228  AMIPtr_(new faceAreaWeightAMI(true, AMIReverse_)),
1229  surfPtr_(nullptr),
1230  surfDict_(fileName("surface"))
1231 {
1232  // NOTE: same region, no sample-world. Thus no world-world communication
1233 }
1234 
1235 
1237 (
1238  const polyPatch& pp,
1239  const word& sampleRegion,
1240  const sampleMode mode,
1241  const word& samplePatch,
1242  const vectorField& offsets
1243 )
1244 :
1246  (
1247  pp,
1248  sampleRegion,
1249  mode,
1250  samplePatch,
1251  scalar(0)
1252  )
1253 {
1255 }
1256 
1257 
1259 (
1260  const polyPatch& pp,
1261  const word& sampleRegion,
1262  const sampleMode mode,
1263  const word& samplePatch,
1264  const vector& uniformOffset
1265 )
1266 :
1268  (
1269  pp,
1270  sampleRegion,
1271  mode,
1272  samplePatch,
1273  scalar(0)
1274  )
1275 {
1276  mappedPatchBase::setOffset(uniformOffset);
1277 }
1278 
1279 
1281 (
1282  const polyPatch& pp,
1283  const word& sampleRegion,
1284  const sampleMode mode,
1285  const word& samplePatch,
1286  const scalar normalDistance
1287 )
1288 :
1289  patch_(pp),
1290  sampleWorld_(),
1291  sampleRegion_(sampleRegion),
1292  mode_(mode),
1293  samplePatch_(samplePatch),
1294  coupleGroup_(),
1295  sampleDatabasePtr_(),
1296  offsetMode_(NORMAL),
1297  offset_(Zero),
1298  offsets_(0),
1299  distance_(normalDistance),
1300  communicator_(-1), // Demand-driven (cached value)
1301  sameRegion_
1302  (
1303  sampleWorld_.empty()
1304  && sampleRegion_ == patch_.boundaryMesh().mesh().name()
1305  ),
1306  mapPtr_(nullptr),
1307  AMIReverse_(false),
1308  AMIPtr_(new faceAreaWeightAMI(true, AMIReverse_)),
1309  surfPtr_(nullptr),
1310  surfDict_(fileName("surface"))
1311 {
1313 }
1314 
1315 
1317 (
1318  const polyPatch& pp,
1319  const dictionary& dict
1320 )
1321 :
1322  patch_(pp),
1323  sampleWorld_(dict.getOrDefault<word>("sampleWorld", word::null)),
1324  sampleRegion_(dict.getOrDefault<word>("sampleRegion", word::null)),
1325  mode_(sampleModeNames_.get("sampleMode", dict)),
1326  samplePatch_(dict.getOrDefault<word>("samplePatch", word::null)),
1327  coupleGroup_(dict),
1328  sampleDatabasePtr_(readDatabase(dict)),
1329  offsetMode_(UNIFORM),
1330  offset_(Zero),
1331  offsets_(0),
1332  distance_(0),
1333  communicator_(-1), // Demand-driven (cached value)
1334  sameRegion_
1335  (
1336  sampleWorld_.empty()
1337  && sampleRegion_ == patch_.boundaryMesh().mesh().name()
1338  ),
1339  mapPtr_(nullptr),
1340  AMIReverse_
1341  (
1342  dict.getOrDefault
1343  (
1344  "reverseTarget", // AMIInterpolation uses this keyword
1345  dict.getOrDefault
1346  (
1347  "flipNormals",
1348  false
1349  )
1350  )
1351  ),
1352  AMIPtr_
1353  (
1355  (
1356  dict.getOrDefault("AMIMethod", faceAreaWeightAMI::typeName),
1357  dict,
1358  AMIReverse_
1359  )
1360  ),
1361  surfPtr_(nullptr),
1362  surfDict_(dict.subOrEmptyDict("surface"))
1363 {
1365 
1366  if (!coupleGroup_.good())
1367  {
1368  if (sampleWorld_.empty() && sampleRegion_.empty())
1369  {
1370  // If no coupleGroup and no sampleRegion assume local region
1372  sameRegion_ = true;
1373  }
1374  }
1375 
1376  if (offsetModeNames_.readIfPresent("offsetMode", dict, offsetMode_))
1377  {
1378  switch (offsetMode_)
1379  {
1380  case UNIFORM:
1381  {
1382  dict.readEntry("offset", offset_);
1383  }
1384  break;
1385 
1386  case NONUNIFORM:
1387  {
1388  offsets_ = pointField("offsets", dict, patch_.size());
1389  }
1390  break;
1391 
1392  case NORMAL:
1393  {
1394  dict.readEntry("distance", distance_);
1395  }
1396  break;
1397  }
1398  }
1399  else if (dict.readIfPresent("offset", offset_))
1400  {
1401  offsetMode_ = UNIFORM;
1402  }
1403  else if (dict.found("offsets"))
1404  {
1406  offsets_ = pointField("offsets", dict, patch_.size());
1407  }
1408  else if (mode_ != NEARESTPATCHFACE && mode_ != NEARESTPATCHFACEAMI)
1409  {
1411  << "Please supply the offsetMode as one of "
1412  << offsetModeNames_
1413  << exit(FatalIOError);
1414  }
1415 }
1416 
1417 
1419 (
1420  const polyPatch& pp,
1421  const sampleMode mode,
1422  const dictionary& dict
1423 )
1424 :
1425  patch_(pp),
1426  sampleWorld_(dict.getOrDefault<word>("sampleWorld", word::null)),
1427  sampleRegion_(dict.getOrDefault<word>("sampleRegion", word::null)),
1428  mode_(mode),
1429  samplePatch_(dict.getOrDefault<word>("samplePatch", word::null)),
1430  coupleGroup_(dict), //dict.getOrDefault<word>("coupleGroup", word::null)),
1431  sampleDatabasePtr_(readDatabase(dict)),
1432  offsetMode_(UNIFORM),
1433  offset_(Zero),
1434  offsets_(0),
1435  distance_(0),
1436  communicator_(-1), // Demand-driven (cached value)
1437  sameRegion_
1438  (
1439  sampleWorld_.empty()
1440  && sampleRegion_ == patch_.boundaryMesh().mesh().name()
1441  ),
1442  mapPtr_(nullptr),
1443  AMIReverse_
1444  (
1445  dict.getOrDefault
1446  (
1447  "reverseTarget", // AMIInterpolation uses this keyword
1448  dict.getOrDefault
1449  (
1450  "flipNormals",
1451  false
1452  )
1453  )
1454  ),
1455  AMIPtr_
1456  (
1458  (
1459  dict.getOrDefault("AMIMethod", faceAreaWeightAMI::typeName),
1460  dict,
1461  AMIReverse_
1462  )
1463  ),
1464  surfPtr_(nullptr),
1465  surfDict_(dict.subOrEmptyDict("surface"))
1466 {
1468 
1470  {
1472  << "Construct from sampleMode and dictionary only applicable for "
1473  << " collocated patches in modes "
1476  << exit(FatalIOError);
1477  }
1478 
1479  if (!coupleGroup_.good())
1480  {
1481  if (sampleWorld_.empty() && sampleRegion_.empty())
1482  {
1483  // If no coupleGroup and no sampleRegion assume local region
1485  sameRegion_ = true;
1486  }
1487  }
1488 }
1489 
1490 
1492 (
1493  const polyPatch& pp,
1494  const mappedPatchBase& mpb
1495 )
1496 :
1497  patch_(pp),
1498  sampleWorld_(mpb.sampleWorld_),
1499  sampleRegion_(mpb.sampleRegion_),
1500  mode_(mpb.mode_),
1501  samplePatch_(mpb.samplePatch_),
1502  coupleGroup_(mpb.coupleGroup_),
1503  sampleDatabasePtr_
1504  (
1505  mpb.sampleDatabasePtr_
1506  ? new fileName(mpb.sampleDatabasePtr_())
1507  : nullptr
1508  ),
1509  offsetMode_(mpb.offsetMode_),
1510  offset_(mpb.offset_),
1511  offsets_(mpb.offsets_),
1512  distance_(mpb.distance_),
1513  communicator_(mpb.communicator_),
1514  sameRegion_(mpb.sameRegion_),
1515  mapPtr_(nullptr),
1516  AMIReverse_(mpb.AMIReverse_),
1517  AMIPtr_(mpb.AMIPtr_->clone()),
1518  surfPtr_(nullptr),
1519  surfDict_(mpb.surfDict_)
1520 {}
1521 
1522 
1524 (
1525  const polyPatch& pp,
1526  const mappedPatchBase& mpb,
1527  const labelUList& mapAddressing
1528 )
1529 :
1530  patch_(pp),
1531  sampleWorld_(mpb.sampleWorld_),
1532  sampleRegion_(mpb.sampleRegion_),
1533  mode_(mpb.mode_),
1534  samplePatch_(mpb.samplePatch_),
1535  coupleGroup_(mpb.coupleGroup_),
1536  sampleDatabasePtr_
1537  (
1538  mpb.sampleDatabasePtr_
1539  ? new fileName(mpb.sampleDatabasePtr_())
1540  : nullptr
1541  ),
1542  offsetMode_(mpb.offsetMode_),
1543  offset_(mpb.offset_),
1544  offsets_
1545  (
1546  offsetMode_ == NONUNIFORM
1547  ? vectorField(mpb.offsets_, mapAddressing)
1548  : vectorField()
1549  ),
1550  distance_(mpb.distance_),
1551  communicator_(mpb.communicator_),
1552  sameRegion_(mpb.sameRegion_),
1553  mapPtr_(nullptr),
1554  AMIReverse_(mpb.AMIReverse_),
1555  AMIPtr_(mpb.AMIPtr_->clone()),
1556  surfPtr_(nullptr),
1557  surfDict_(mpb.surfDict_)
1558 {}
1559 
1560 
1561 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
1564 {
1565  clearOut();
1566 }
1567 
1570 {
1571  mapPtr_.reset(nullptr);
1572  surfPtr_.reset(nullptr);
1573  AMIPtr_->upToDate(false);
1574 }
1575 
1576 
1577 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
1579 void Foam::mappedPatchBase::setOffset(const scalar normalDist)
1580 {
1581  clearOut();
1582  offsetMode_ = offsetMode::NORMAL;
1583  offset_ = Zero;
1584  offsets_.clear();
1585  distance_ = normalDist;
1586 }
1587 
1589 void Foam::mappedPatchBase::setOffset(const vector& uniformOffset)
1590 {
1591  clearOut();
1592  offsetMode_ = offsetMode::UNIFORM;
1593  offset_ = uniformOffset;
1594  offsets_.clear();
1595  distance_ = Zero;
1596 }
1597 
1599 void Foam::mappedPatchBase::setOffset(const vectorField& offsets)
1600 {
1601  clearOut();
1602  offsetMode_ = offsetMode::NONUNIFORM;
1603  offset_ = Zero;
1604  offsets_ = offsets;
1605  distance_ = Zero;
1606 }
1607 
1608 
1610 (
1611  const word& sampleRegion
1612 ) const
1613 {
1614  const polyMesh& thisMesh = patch_.boundaryMesh().mesh();
1615  return
1616  (
1617  sampleRegion.empty() || sampleRegion == thisMesh.name()
1618  ? thisMesh
1619  : thisMesh.time().lookupObject<polyMesh>(sampleRegion)
1620  );
1621 }
1622 
1623 
1625 (
1626  const word& sampleRegion,
1627  const word& samplePatch
1628 ) const
1629 {
1630  const polyMesh& nbrMesh = lookupMesh(sampleRegion);
1631 
1632  const label patchi = nbrMesh.boundaryMesh().findPatchID(samplePatch);
1633 
1634  if (patchi == -1)
1635  {
1637  << "Cannot find patch " << samplePatch
1638  << " in region " << sampleRegion_ << endl
1639  << exit(FatalError);
1640  }
1641  return nbrMesh.boundaryMesh()[patchi];
1642 }
1643 
1646 {
1647  if (UPstream::myWorld() != sampleWorld_)
1648  {
1650  << "sampleWorld : " << sampleWorld_
1651  << " is not the current world : " << UPstream::myWorld()
1652  << exit(FatalError);
1653  }
1654  return lookupMesh(sampleRegion());
1655 }
1656 
1659 {
1660  const polyMesh& nbrMesh = sampleMesh();
1661 
1662  const label patchi = nbrMesh.boundaryMesh().findPatchID(samplePatch());
1663 
1664  if (patchi == -1)
1665  {
1667  << "Cannot find patch " << samplePatch()
1668  << " in region " << sampleRegion_ << endl
1669  << "Valid patches are " << nbrMesh.boundaryMesh().names()
1670  << exit(FatalError);
1671  }
1672 
1673  return nbrMesh.boundaryMesh()[patchi];
1674 }
1675 
1676 
1678 (
1679  const pointField& fc
1680 ) const
1681 {
1682  auto tfld = tmp<pointField>::New(fc);
1683  auto& fld = tfld.ref();
1684 
1685  switch (offsetMode_)
1686  {
1687  case UNIFORM:
1688  {
1689  fld += offset_;
1690  break;
1691  }
1692  case NONUNIFORM:
1693  {
1694  fld += offsets_;
1695  break;
1696  }
1697  case NORMAL:
1698  {
1699  // Get outwards pointing normal
1700  vectorField n(patch_.faceAreas());
1701  n /= mag(n);
1702 
1703  fld += distance_*n;
1704  break;
1705  }
1706  }
1707 
1708  return tfld;
1709 }
1710 
1713 {
1714  return samplePoints(facePoints(patch_));
1715 }
1716 
1717 
1719 (
1720  const polyMesh& mesh,
1721  const label facei,
1722  const polyMesh::cellDecomposition decompMode
1723 )
1724 {
1725  const point& fc = mesh.faceCentres()[facei];
1726 
1727  switch (decompMode)
1728  {
1729  case polyMesh::FACE_PLANES:
1731  {
1732  // For both decompositions the face centre is guaranteed to be
1733  // on the face
1734  return pointIndexHit(true, fc, facei);
1735  }
1736  break;
1737 
1739  case polyMesh::CELL_TETS:
1740  {
1741  // Find the intersection of a ray from face centre to cell centre
1742  // Find intersection of (face-centre-decomposition) centre to
1743  // cell-centre with face-diagonal-decomposition triangles.
1744 
1745  const pointField& p = mesh.points();
1746  const face& f = mesh.faces()[facei];
1747 
1748  if (f.size() <= 3)
1749  {
1750  // Return centre of triangle.
1751  return pointIndexHit(true, fc, 0);
1752  }
1753 
1754  const label celli = mesh.faceOwner()[facei];
1755  const point& cc = mesh.cellCentres()[celli];
1756  vector d = fc-cc;
1757 
1758  const label fp0 = mesh.tetBasePtIs()[facei];
1759  const point& basePoint = p[f[fp0]];
1760 
1761  label fp = f.fcIndex(fp0);
1762  for (label i = 2; i < f.size(); i++)
1763  {
1764  const point& thisPoint = p[f[fp]];
1765  label nextFp = f.fcIndex(fp);
1766  const point& nextPoint = p[f[nextFp]];
1767 
1768  const triPointRef tri(basePoint, thisPoint, nextPoint);
1769  pointHit hitInfo = tri.intersection
1770  (
1771  cc,
1772  d,
1774  );
1775 
1776  if (hitInfo.hit() && hitInfo.distance() > 0)
1777  {
1778  return pointIndexHit(true, hitInfo.point(), i-2);
1779  }
1780 
1781  fp = nextFp;
1782  }
1783 
1784  // Fall-back
1785  return pointIndexHit(false, fc, -1);
1786  }
1787  break;
1788 
1789  default:
1790  {
1792  << "problem" << abort(FatalError);
1793  return pointIndexHit();
1794  }
1795  }
1796 }
1797 
1798 
1800 (
1801  const objectRegistry& obr,
1802  const fileName& path
1803 )
1804 {
1805  // Lookup (and create if non-existing) a registry using
1806  // '/' separated path. Like 'mkdir -p'
1807 
1808  fileName cleanedPath(path);
1809  cleanedPath.clean(); // Remove unneeded ".."
1810  const wordList names(cleanedPath.components());
1811 
1812  if (names.empty())
1813  {
1814  return obr;
1815  }
1816  else
1817  {
1818  return subRegistry(obr, names, 0);
1819  }
1820 }
1821 
1822 
1824 (
1825  const fileName& root,
1826  const label proci
1827 )
1828 {
1829  const word processorName("processor"+Foam::name(proci));
1830  return root/"send"/processorName;
1831 }
1832 
1834 Foam::fileName Foam::mappedPatchBase::sendPath(const label proci) const
1835 {
1836  return sendPath(sampleDatabasePath(), proci);
1837 }
1838 
1839 
1841 (
1842  const fileName& root,
1843  const label proci
1844 )
1845 {
1846  const word processorName("processor"+Foam::name(proci));
1847  return root/"receive"/processorName;
1848 }
1849 
1851 Foam::fileName Foam::mappedPatchBase::receivePath(const label proci) const
1852 {
1853  return receivePath(sampleDatabasePath(), proci);
1854 }
1855 
1856 
1858 (
1859  const objectRegistry& obr,
1860  dictionary& dict
1861 )
1862 {
1863  forAllIters(obr, iter)
1864  {
1865  regIOobject* objPtr = iter.val();
1866  const regIOobject& obj = *objPtr;
1867 
1868  if (isA<objectRegistry>(obj))
1869  {
1870  dictionary& d = dict.subDictOrAdd(obj.name());
1871  writeDict(dynamic_cast<const objectRegistry&>(obj), d);
1872  }
1873  else if
1874  (
1875  writeIOField<scalar>(obj, dict)
1876  || writeIOField<vector>(obj, dict)
1877  || writeIOField<sphericalTensor>(obj, dict)
1878  || writeIOField<symmTensor>(obj, dict)
1879  || writeIOField<tensor>(obj, dict)
1880  )
1881  {
1882  // IOField specialisation
1883  }
1884  else
1885  {
1886  // Not tested. No way of retrieving data anyway ...
1887  OTstream os;
1888  obj.writeData(os);
1889 
1890  primitiveEntry* pePtr = new primitiveEntry(obj.name(), os.tokens());
1891  dict.add(pePtr);
1892  }
1893  }
1894 }
1895 
1898 {
1899  // Reverse of writeDict - reads fields from dictionary into objectRegistry
1900  for (const entry& e : d)
1901  {
1902  if (e.isDict())
1903  {
1904  // Add sub registry
1905  objectRegistry& sub = const_cast<objectRegistry&>
1906  (
1907  obr.subRegistry(e.keyword(), true)
1908  );
1909 
1910  readDict(e.dict(), sub);
1911  }
1912  else
1913  {
1914  ITstream& is = e.stream();
1915  token tok(is);
1916 
1917  if
1918  (
1919  constructIOField<scalar>(e.keyword(), tok, is, obr)
1920  || constructIOField<vector>(e.keyword(), tok, is, obr)
1921  || constructIOField<sphericalTensor>(e.keyword(), tok, is, obr)
1922  || constructIOField<symmTensor>(e.keyword(), tok, is, obr)
1923  || constructIOField<tensor>(e.keyword(), tok, is, obr)
1924  )
1925  {
1926  // Done storing field
1927  }
1928  else
1929  {
1930  FatalErrorInFunction << "Unsupported type " << e.keyword()
1931  << exit(FatalError);
1932  }
1933  }
1934  }
1935 }
1936 
1939 {
1940  os.writeEntry("sampleMode", sampleModeNames_[mode_]);
1941  os.writeEntryIfDifferent<word>("sampleWorld", word::null, sampleWorld_);
1942  os.writeEntryIfDifferent<word>("sampleRegion", word::null, sampleRegion_);
1943  os.writeEntryIfDifferent<word>("samplePatch", word::null, samplePatch_);
1944 
1945  if (sampleDatabasePtr_)
1946  {
1947  os.writeEntry("sampleDatabase", Switch::name(true));
1948  // Write database path if differing
1950  (
1951  "sampleDatabasePath",
1953  sampleDatabasePtr_()
1954  );
1955  }
1956  coupleGroup_.write(os);
1957 
1958  if
1959  (
1960  offsetMode_ == UNIFORM
1961  && offset_ == vector::zero
1962  && (mode_ == NEARESTPATCHFACE || mode_ == NEARESTPATCHFACEAMI)
1963  )
1964  {
1965  // Collocated mode. No need to write offset data
1966  }
1967  else
1968  {
1969  os.writeEntry("offsetMode", offsetModeNames_[offsetMode_]);
1970 
1971  switch (offsetMode_)
1972  {
1973  case UNIFORM:
1974  {
1975  os.writeEntry("offset", offset_);
1976  break;
1977  }
1978  case NONUNIFORM:
1979  {
1980  offsets_.writeEntry("offsets", os);
1981  break;
1982  }
1983  case NORMAL:
1984  {
1985  os.writeEntry("distance", distance_);
1986  break;
1987  }
1988  }
1989  }
1990 
1991  if (mode_ == NEARESTPATCHFACEAMI)
1992  {
1993  if (AMIPtr_)
1994  {
1995  // Use AMI to write itself. Problem: outputs:
1996  // - restartUncoveredSourceFace
1997  // - reverseTarget (instead of flipNormals)
1998  AMIPtr_->write(os);
1999  }
2000  if (!surfDict_.empty())
2001  {
2002  surfDict_.writeEntry(surfDict_.dictName(), os);
2003  }
2004  }
2005 }
2006 
2007 
2008 // ************************************************************************* //
label getWorldCommunicator() const
Get the communicator for the world-world connection.
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
label findPatchID(const word &patchName, const bool allowNotFound=true) const
Find patch index given a name, return -1 if not found.
static const word & myWorld()
My world.
Definition: UPstream.H:1165
use face normal + distance
prefixOSstream Perr
OSstream wrapped stderr (std::cerr) with parallel prefix.
const List< face_type > & localFaces() const
Return patch faces addressing into local point list.
bool set(const Key &key)
Same as insert (no value to overwrite)
Definition: HashSet.H:240
dictionary dict
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
const Type & lookupObject(const word &name, const bool recursive=false) const
Lookup and return const reference to the object of the given Type. Fatal if not found or the wrong ty...
const polyMesh & lookupMesh(const word &region) const
Lookup mesh.
A triangle primitive used to calculate face normals and swept volumes. Uses referred points...
Definition: triangle.H:217
nearest face on selected patch
vector offset_
Offset vector (uniform)
List< word > names(const UPtrList< T > &list, const UnaryMatchPredicate &matcher)
List of names generated by calling name() for each list item and filtered for matches.
A class for handling file names.
Definition: fileName.H:72
virtual const fileName & name() const
The name of the stream.
Definition: IOstream.C:33
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
static const Enum< offsetMode > offsetModeNames_
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:68
word sampleRegion_
Region to sample.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
const Field< point_type > & localPoints() const
Return pointField of points in patch.
labelHashSet used(const bitSet &select)
Convert a bitset to a labelHashSet of the indices used.
Definition: HashOps.C:26
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
tmp< pointField > samplePoints() const
Get the sample points.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
const polyPatch & patch_
Patch to sample.
A 2-tuple for storing two objects of dissimilar types. The container is similar in purpose to std::pa...
Definition: stringOps.H:54
const polyMesh & sampleMesh() const
Get the region mesh.
const labelIOList & tetBasePtIs() const
Return the tetBasePtIs.
Definition: polyMesh.C:899
const polyPatch & lookupPatch(const word &sampleRegion, const word &samplePatch) const
Lookup patch.
Centralized handling of multi-world MPI connections.
bool sameRegion_
Same region.
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:195
bool readIfPresent(const word &key, const dictionary &dict, EnumType &val, const bool warnOnly=false) const
Find an entry if present, and assign to T val.
Definition: EnumI.H:111
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: BitOps.H:56
bool test(const Key &key) const
Same as contains() - return true if key exists in the set.
Definition: HashSet.H:218
static const meshSearchMeshObject & New(const polyMesh &mesh, Args &&... args)
Get existing or create a new MeshObject. Registered with typeName.
Definition: MeshObject.C:53
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Output to file stream, using an OSstream.
Definition: OFstream.H:49
const coupleGroupIdentifier coupleGroup_
PatchGroup (if in sampleMode NEARESTPATCH*)
virtual ~mappedPatchBase()
Destructor.
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
bool empty() const noexcept
True if List is empty (ie, size() is zero)
Definition: UList.H:666
A token holds an item read from Istream.
Definition: token.H:65
static const fileName null
An empty fileName.
Definition: fileName.H:111
UIndirectList< label > labelUIndList
UIndirectList of labels.
Definition: IndirectList.H:65
Random rndGen
Definition: createFields.H:23
static fileName sendPath(const fileName &root, const label proci)
Helper: return path to store data to be sent to processor i.
T & first()
Access first element of the list, position [0].
Definition: UList.H:853
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface...
Definition: boundaryMesh.H:58
static const multiWorldConnections & New(const Time &runTime)
Access mesh object.
scalarField samples(nIntervals, Zero)
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
PointIndexHit< point > pointIndexHit
A PointIndexHit with a 3D point.
Definition: pointIndexHit.H:58
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
static const objectRegistry & subRegistry(const objectRegistry &obr, const wordList &names, const label index)
Lookup (sub)objectRegistry by following names of sub registries. Creates non-existing intermediate on...
void setOffset(const scalar normalDist)
Change to normal offset with given distance.
A simple output token stream that can be used to build token lists. Always UNCOMPRESSED.
Definition: OTstream.H:52
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:321
static int & msgType() noexcept
Message tag of standard messages.
Definition: UPstream.H:1229
This class describes the interaction of an object (often a face) and a point. It carries the info of ...
Definition: pointIndexHit.H:44
entry * add(entry *entryPtr, bool mergeEntry=false)
Add a new entry.
Definition: dictionary.C:625
void writeOBJ(Ostream &os, const point &pt)
Write obj representation of a point.
Definition: meshTools.C:196
Ignore writing from objectRegistry::writeObject()
static int myProcNo(const label communicator=worldComm)
Rank of this process in the communicator (starting from masterNo()). Can be negative if the process i...
Definition: UPstream.H:1074
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:360
const objectRegistry & subRegistry(const word &name, const bool forceCreate=false, const bool recursive=false) const
Lookup and return a const sub-objectRegistry.
virtual const fileName & name() const override
Get the name of the output serial stream. (eg, the name of the Fstream file name) ...
Definition: OSstream.H:128
List< labelList > labelListList
List of labelList.
Definition: labelList.H:38
cellDecomposition
Enumeration defining the decomposition of the cell for.
Definition: polyMesh.H:103
vectorField offsets_
Offset vector (nonuniform)
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:69
Macros for easy insertion into run-time selection tables.
const labelList & meshPoints() const
Return labelList of mesh points in patch.
scalar distance_
Offset distance (normal)
dictionary & subDictOrAdd(const word &keyword, enum keyType::option matchOpt=keyType::REGEX)
Find and return a sub-dictionary for manipulation.
Definition: dictionary.C:481
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1078
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
label fcIndex(const label i) const noexcept
The forward circular index. The next index in the list which returns to the first at the end of the l...
Definition: UListI.H:97
Face area weighted Arbitrary Mesh Interface (AMI) method.
A keyword and a list of tokens comprise a primitiveEntry. A primitiveEntry can be read...
Container (non-empty) with different values.
Definition: ListPolicy.H:107
nearest patch face + AMI interpolation
static void allGatherList(List< T > &values, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Gather data, but keep individual values separate. Uses linear/tree communication. ...
bool hit() const noexcept
Is there a hit.
Definition: pointHit.H:145
static label commWarn(const label communicator) noexcept
Alter communicator debugging setting. Warns for use of any communicator differing from specified...
Definition: UPstream.H:449
static void writeDict(const objectRegistry &obr, dictionary &dict)
Convert objectRegistry contents into dictionary.
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
A list of faces which address into the list of points.
A List obtained as a section of another List.
Definition: SubList.H:50
static const char * name(const bool b) noexcept
A string representation of bool as "false" / "true".
Definition: Switch.C:141
static label nProcs(const label communicator=worldComm)
Number of ranks in parallel run (for given communicator). It is 1 for serial run. ...
Definition: UPstream.H:1065
word sampleWorld_
World to sample.
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
void setSize(const label n)
Alias for resize()
Definition: List.H:316
static fileName receivePath(const fileName &root, const label proci)
Helper: return path to store data to be received from processor i.
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
bool good() const noexcept
The patchGroup has a non-empty name.
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundaryMesh reference.
Definition: polyPatch.C:309
const polyMesh & mesh() const noexcept
Return the mesh reference.
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i), works like std::iota() but returning a...
Definition: labelLists.C:44
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition: polyMesh.H:608
A class for handling words, derived from Foam::string.
Definition: word.H:63
const Time & time() const noexcept
Return time registry.
#define DebugInFunction
Report an information message using Foam::Info.
wordList names() const
Return a list of patch names.
Tree tree(triangles.begin(), triangles.end())
void findSamples(const sampleMode mode, const label myWorldIndex, const pointField &, const labelList &wantedWorlds, const labelList &origProcs, labelList &sampleProcs, labelList &sampleIndices, pointField &sampleLocations) const
Find (global) cells/faces containing samples.
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1116
#define forAllIters(container, iter)
Iterate across all elements in the container object.
Definition: stdFoam.H:336
static const word null
An empty word.
Definition: word.H:84
void calcMapping() const
Calculate mapping.
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition: tmp.H:206
Determines a mapping between patch face centres and mesh cell or face centres and processors they&#39;re ...
virtual bool write(const token &tok)=0
Write token to stream or otherwise handle it.
offsetMode offsetMode_
How to obtain samples.
Describes the interaction of a object and a (templated) point. It carries the info of a successful hi...
Definition: pointHit.H:43
const Field< point_type > & points() const noexcept
Return reference to global points.
void append(const T &val)
Copy append an element to the end of this list.
Definition: DynamicList.H:584
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1103
const vectorField & cellCentres() const
errorManip< error > abort(error &err)
Definition: errorManip.H:139
Ostream & writeEntryIfDifferent(const word &key, const T &value1, const T &value2)
Write a keyword/value entry only when the two values differ.
Definition: Ostream.H:336
label find(const T &val) const
Find index of the first occurrence of the value.
Definition: UList.C:173
tmp< pointField > facePoints(const polyPatch &) const
Get the points from face-centre-decomposition face centres and project them onto the face-diagonal-de...
const autoPtr< Foam::searchableSurface > & surfPtr() const
Return a pointer to the AMI projection surface.
static const Enum< sampleMode > sampleModeNames_
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
label getCommByName(const word &otherWorld) const
Get communicator for myWorld to other world connection by NAME.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
const word & constant() const noexcept
Return constant name.
Definition: TimePathsI.H:112
wordList components(const char delim='/') const
Return path components as wordList.
Definition: fileName.C:532
int debug
Static debugging option.
label facePoint(const int facei, const block &block, const label i, const label j)
OBJstream os(runTime.globalPath()/outputName)
fileName path(UMean.rootPath()/UMean.caseName()/"graphs"/UMean.instance())
virtual void write(Ostream &os) const
Write as a dictionary.
defineTypeNameAndDebug(combustionModel, 0)
static bool clean(std::string &str)
Cleanup filename string, possibly applies other transformations such as changing the path separator e...
Definition: fileName.C:192
static pointIndexHit facePoint(const polyMesh &, const label facei, const polyMesh::cellDecomposition)
Get a point on the face given a face decomposition method:
static autoPtr< fileName > readDatabase(const dictionary &dict)
Read optional database name from dictionary.
labelList f(nPoints)
offsetMode
How to project face centres.
static label commWorld() noexcept
Communicator for all ranks (respecting any local worlds)
Definition: UPstream.H:429
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< ' ';}gmvFile<< nl;for(const word &name :lagrangianScalarNames){ IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
const vectorField & faceCentres() const
const indexedOctree< treeDataCell > & cellTree() const
Return the cell search tree.
Definition: polyMesh.C:939
sampleMode mode() const noexcept
What to sample.
Class containing processor-to-processor mapping information.
void calcAMI() const
Calculate AMI interpolator.
vector point
Point is a vector.
Definition: point.H:37
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition: bitSet.H:59
#define WarningInFunction
Report a warning using Foam::Warning.
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition: error.H:64
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:627
Type gAverage(const FieldField< Field, Type > &f)
A List with indirect addressing. Like IndirectList but does not store addressing. ...
Definition: faMatrix.H:56
const T2 & second() const noexcept
Access the second element.
Definition: Tuple2.H:142
sampleMode
Mesh items to sample.
scalar distance() const noexcept
Return distance to hit.
Definition: pointHit.H:169
static bool master(const label communicator=worldComm)
True if process corresponds to the master rank in the communicator.
Definition: UPstream.H:1082
Interpolation class dealing with transfer of data between two primitive patches with an arbitrary mes...
static autoPtr< searchableSurface > New(const word &surfaceType, const IOobject &io, const dictionary &dict)
Return a reference to the selected searchableSurface.
void reduce(const List< UPstream::commsStruct > &comms, T &value, const BinaryOp &bop, const int tag, const label comm)
Reduce inplace (cf. MPI Allreduce) using specified communication schedule.
static void readDict(const dictionary &d, objectRegistry &obr)
(recursively) construct and register IOFields from dictionary
const point_type & point() const noexcept
Return the point, no checks.
Definition: pointHit.H:161
regIOobject is an abstract class derived from IOobject to handle automatic object registration with t...
Definition: regIOobject.H:66
label n
Container (non-empty) with identical values.
Definition: ListPolicy.H:106
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
void collectSamples(const label mySampleWorld, const pointField &facePoints, pointField &samples, labelList &patchFaceWorlds, labelList &patchFaceProcs, labelList &patchFaces, pointField &patchFc) const
Collect single list of samples and originating processor+face + wanted world.
List< label > labelList
A List of labels.
Definition: List.H:62
volScalarField & p
A class for managing temporary objects.
Definition: HashPtrTable.H:50
pointHit intersection(const point &p, const vector &q, const intersection::algorithm alg, const scalar tol=0.0) const
Fast intersection with a ray.
Definition: triangleI.H:630
Registry of regIOobjects.
static autoPtr< T > New(Args &&... args)
Construct autoPtr with forwarding arguments.
Definition: autoPtr.H:178
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:69
virtual bool writeData(Ostream &) const =0
Pure virtual writeData function.
const T1 & first() const noexcept
Access the first element.
Definition: Tuple2.H:132
static const wordList & allWorlds() noexcept
All worlds.
Definition: UPstream.H:1141
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
mode_t mode(const fileName &name, const bool followLink=true)
Return the file mode, normally following symbolic links.
Definition: POSIX.C:773
const polyPatch & samplePolyPatch() const
Get the patch on the region.
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
bool addWorldConnection()
Add a world-world connection.
mappedPatchBase(const polyPatch &)
Construct from patch.
void findLocalSamples(const sampleMode mode, const label sampleWorld, const word &sampleRegion, const word &samplePatch, const pointField &samplePoints, List< nearInfoWorld > &nearest) const
Find (local) cells/faces containing samples.
An input stream of tokens.
Definition: ITstream.H:52
const vectorField & offsets() const noexcept
Offset vectors (from patch faces to destination mesh objects)
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
Namespace for OpenFOAM.
A keyword and a list of tokens is an &#39;entry&#39;.
Definition: entry.H:63
static label myWorldID()
My worldID.
Definition: UPstream.H:1157
static void listCombineReduce(List< T > &values, const CombineOp &cop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
After completion all processors have the same data.
const sampleMode mode_
What to sample.
IOerror FatalIOError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL IO ERROR&#39; header text and ...
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127