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-2022 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.lookupOrDefault<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::worldComm;
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(Pstream::warnComm);
191  Pstream::warnComm = myComm;
192 
193  if (debug & 2)
194  {
195  Perr<< "patch: " << patch_.name()
196  << "[rank=" << myRank << " procs=" << nProcs
197  << " comm=" << myComm << "] collect samples" << endl;
198  }
199 
200  // Collect all sample points and the faces they come from.
201  {
202  List<pointField> globalFc(nProcs);
203  globalFc[myRank] = facePoints;
204 
205  Pstream::allGatherList(globalFc, Pstream::msgType(), myComm);
206 
207  // Rework into straight list
208  patchFc = ListListOps::combine<pointField>
209  (
210  globalFc,
211  accessOp<pointField>()
212  );
213  }
214 
215  {
216  List<pointField> globalSamples(nProcs);
217  globalSamples[myRank] = samplePoints(facePoints);
218  Pstream::allGatherList(globalSamples, Pstream::msgType(), myComm);
219 
220  // Rework into straight list
221  samples = ListListOps::combine<pointField>
222  (
223  globalSamples,
224  accessOp<pointField>()
225  );
226  }
227 
228  {
229  labelListList globalFaces(nProcs);
230  globalFaces[myRank] = identity(patch_.size());
231  // Distribute to all processors
232  Pstream::allGatherList(globalFaces, Pstream::msgType(), myComm);
233 
234  patchFaces = ListListOps::combine<labelList>
235  (
236  globalFaces,
237  accessOp<labelList>()
238  );
239  }
240 
241  {
242  labelList procToWorldIndex
243  (
244  UPstream::listGatherValues<label>(mySampleWorld, myComm)
245  );
246  labelList nPerProc
247  (
248  UPstream::listGatherValues<label>(patch_.size(), myComm)
249  );
250 
251  Pstream::broadcasts(myComm, procToWorldIndex, nPerProc);
252 
253  patchFaceWorlds.setSize(patchFaces.size());
254  patchFaceProcs.setSize(patchFaces.size());
255 
256  label sampleI = 0;
257  forAll(nPerProc, proci)
258  {
259  for (label i = 0; i < nPerProc[proci]; i++)
260  {
261  patchFaceWorlds[sampleI] = procToWorldIndex[proci];
262  patchFaceProcs[sampleI] = proci;
263  sampleI++;
264  }
265  }
266  }
267 
268  Pstream::warnComm = oldWarnComm;
269 }
270 
271 
273 (
274  const sampleMode mode,
275 
276  const label mySampleWorld, // local world to sample == my own world
277  const word& sampleRegion, // local region to sample
278  const word& samplePatch, // local patch to sample
279 
280  const pointField& samples,
281  List<nearInfoWorld>& nearest
282 ) const
283 {
284  DebugInFunction << nl;
285 
286  const label myComm = getCommunicator(); // Get or create
287 
288  // Find the local cell containing the samples
289  const label myRank = Pstream::myProcNo(myComm);
290 
291  // Lookup the correct region
292  const polyMesh& mesh = lookupMesh(sampleRegion);
293 
294  // All the info for nearest. Construct to miss
295  nearest.setSize(samples.size());
296  nearInfoWorld miss;
297  {
298  miss.first().second() = Tuple2<scalar, label>(Foam::sqr(GREAT), -1);
299  miss.second() = -1; // set world to be ignored
300  }
301  nearest = miss;
302 
303  switch (mode)
304  {
305  case NEARESTCELL:
306  {
307  if (samplePatch.size() && samplePatch != "none")
308  {
310  << "No need to supply a patch name when in "
311  << sampleModeNames_[mode] << " mode." << exit(FatalError);
312  }
313 
314  //- Note: face-diagonal decomposition
315  const indexedOctree<Foam::treeDataCell>& tree = mesh.cellTree();
316 
317  forAll(samples, sampleI)
318  {
319  const point& sample = samples[sampleI];
320  nearInfoWorld& near = nearest[sampleI];
321 
322  label celli = tree.findInside(sample);
323 
324  if (celli == -1)
325  {
326  near.first().second().first() = Foam::sqr(GREAT);
327  near.first().second().second() = myRank;
328  near.second() = mySampleWorld;
329  }
330  else
331  {
332  const point& cc = mesh.cellCentres()[celli];
333 
334  near.first().first() = pointIndexHit
335  (
336  true,
337  cc,
338  celli
339  );
340  near.first().second().first() = sample.distSqr(cc);
341  near.first().second().second() = myRank;
342  near.second() = mySampleWorld;
343  }
344  }
345  break;
346  }
347 
348  case NEARESTONLYCELL:
349  {
350  if (samplePatch.size() && samplePatch != "none")
351  {
353  << "No need to supply a patch name when in "
354  << sampleModeNames_[mode] << " mode." << exit(FatalError);
355  }
356 
357  //- Note: face-diagonal decomposition
358  const indexedOctree<Foam::treeDataCell>& tree = mesh.cellTree();
359 
360  forAll(samples, sampleI)
361  {
362  const point& sample = samples[sampleI];
363  nearInfoWorld& near = nearest[sampleI];
364 
365  near.first().first() = tree.findNearest(sample, sqr(GREAT));
366 
367  near.first().second().first() =
368  near.first().first().hitPoint().distSqr(sample);
369 
370  near.first().second().second() = myRank;
371  near.second() = mySampleWorld;
372  }
373  break;
374  }
375 
376  case NEARESTPATCHFACE:
377  {
378  Random rndGen(123456);
379 
380  const polyPatch& pp = lookupPatch(sampleRegion, samplePatch);
381 
382  if (pp.empty())
383  {
384  forAll(samples, sampleI)
385  {
386  nearInfoWorld& near = nearest[sampleI];
387  near.first().second().first() = Foam::sqr(GREAT);
388  near.first().second().second() = myRank;
389  near.second() = mySampleWorld;
390  }
391  }
392  else
393  {
394  treeBoundBox patchBb
395  (
396  treeBoundBox(pp.points(), pp.meshPoints()).extend
397  (
398  rndGen,
399  1e-4,
400  ROOTVSMALL
401  )
402  );
403 
404  indexedOctree<treeDataFace> boundaryTree
405  (
406  treeDataFace(mesh, pp.range()), // Patch faces
407 
408  patchBb, // overall search domain
409  8, // maxLevel
410  10, // leafsize
411  3.0 // duplicity
412  );
413  const auto& treeData = boundaryTree.shapes();
414 
415  forAll(samples, sampleI)
416  {
417  const point& sample = samples[sampleI];
418 
419  nearInfoWorld& near = nearest[sampleI];
420  pointIndexHit& nearInfo = near.first().first();
421  nearInfo = boundaryTree.findNearest
422  (
423  sample,
424  patchBb.magSqr()
425  );
426 
427  if (!nearInfo.hit())
428  {
429  near.first().second().first() = Foam::sqr(GREAT);
430  near.first().second().second() = myRank;
431  near.second() = mySampleWorld;
432  }
433  else
434  {
435  const point& fc = treeData.centre(nearInfo.index());
436 
437  nearInfo.setPoint(fc);
438  near.first().second().first() = sample.distSqr(fc);
439  near.first().second().second() = myRank;
440  near.second() = mySampleWorld;
441  }
442  }
443  }
444  break;
445  }
446 
447  case NEARESTPATCHPOINT:
448  {
449  Random rndGen(123456);
450 
451  const polyPatch& pp = lookupPatch(sampleRegion, samplePatch);
452 
453  if (pp.empty())
454  {
455  forAll(samples, sampleI)
456  {
457  nearInfoWorld& near = nearest[sampleI];
458  near.first().second().first() = Foam::sqr(GREAT);
459  near.first().second().second() = myRank;
460  near.second() = mySampleWorld;
461  }
462  }
463  else
464  {
465  // patch (local) points
466  treeBoundBox patchBb
467  (
468  treeBoundBox(pp.points(), pp.meshPoints()).extend
469  (
470  rndGen,
471  1e-4,
472  ROOTVSMALL
473  )
474  );
475 
476  indexedOctree<treeDataPoint> boundaryTree
477  (
478  // Patch points
479  treeDataPoint(mesh.points(), pp.meshPoints()),
480 
481  patchBb, // overall search domain
482  8, // maxLevel
483  10, // leafsize
484  3.0 // duplicity
485  );
486 
487  forAll(samples, sampleI)
488  {
489  const point& sample = samples[sampleI];
490 
491  nearInfoWorld& near = nearest[sampleI];
492  pointIndexHit& nearInfo = near.first().first();
493  nearInfo = boundaryTree.findNearest
494  (
495  sample,
496  patchBb.magSqr()
497  );
498 
499  if (!nearInfo.hit())
500  {
501  near.first().second().first() = Foam::sqr(GREAT);
502  near.first().second().second() = myRank;
503  near.second() = mySampleWorld;
504  }
505  else
506  {
507  const point& pt = nearInfo.point();
508 
509  near.first().second().first() = sample.distSqr(pt);
510  near.first().second().second() = myRank;
511  near.second() = mySampleWorld;
512  }
513  }
514  }
515  break;
516  }
517 
518  case NEARESTFACE:
519  {
520  if (samplePatch.size() && samplePatch != "none")
521  {
523  << "No need to supply a patch name when in "
524  << sampleModeNames_[mode] << " mode." << exit(FatalError);
525  }
526 
527  //- Note: face-diagonal decomposition
528  const meshSearchMeshObject& meshSearchEngine =
530 
531  forAll(samples, sampleI)
532  {
533  const point& sample = samples[sampleI];
534  nearInfoWorld& near = nearest[sampleI];
535 
536  label facei = meshSearchEngine.findNearestFace(sample);
537 
538  if (facei == -1)
539  {
540  near.first().second().first() = Foam::sqr(GREAT);
541  near.first().second().second() = myRank;
542  near.second() = mySampleWorld;
543  }
544  else
545  {
546  const point& fc = mesh.faceCentres()[facei];
547 
548  near.first().first() = pointIndexHit(true, fc, facei);
549  near.first().second().first() = sample.distSqr(fc);
550  near.first().second().second() = myRank;
551  near.second() = mySampleWorld;
552  }
553  }
554  break;
555  }
556 
557  case NEARESTPATCHFACEAMI:
558  {
559  // nothing to do here
560  return;
561  }
562 
563  default:
564  {
566  << "problem." << abort(FatalError);
567  }
568  }
569 }
570 
571 
573 (
574  const sampleMode mode,
575  const label myWorld,
576  const pointField& samples,
577  const labelList& wantedWorlds,
578  const labelList& origProcs,
579 
580  labelList& sampleProcs,
581  labelList& sampleIndices,
582  pointField& sampleLocations
583 ) const
584 {
585  DebugInFunction << nl;
586 
587  // Find the processor/cell containing the samples. Does not account
588  // for samples being found in two processors.
589 
590  const label myComm = getCommunicator(); // Get or create
591  const label myRank = Pstream::myProcNo(myComm);
592  const label nProcs = Pstream::nProcs(myComm);
593 
594  const label oldWarnComm(Pstream::warnComm);
595  Pstream::warnComm = myComm;
596 
597  wordList samplePatches(nProcs);
598  {
599  samplePatches[myRank] = samplePatch_;
600  Pstream::allGatherList(samplePatches, Pstream::msgType(), myComm);
601  }
602  wordList sampleRegions(nProcs);
603  {
604  sampleRegions[myRank] = sampleRegion_;
605  Pstream::allGatherList(sampleRegions, Pstream::msgType(), myComm);
606  }
607 
608 
609  // Find all the info for nearest
610  List<nearInfoWorld> nearest(samples.size());
611  forAll(nearest, samplei)
612  {
613  nearest[samplei].first() = nearInfo
614  (
615  pointIndexHit(),
616  Tuple2<scalar, label>(Foam::sqr(GREAT), -1)
617  );
618  nearest[samplei].second() = wantedWorlds[samplei];
619  }
620 
621 
622  // Extract samples to search for locally
623  {
624  DynamicList<label> localMap(samples.size());
625  forAll(wantedWorlds, samplei)
626  {
627  if (wantedWorlds[samplei] == myWorld)
628  {
629  localMap.append(samplei);
630  }
631  }
632 
633  if (localMap.size())
634  {
635  pointField localSamples(samples, localMap);
636  labelList localOrigProcs(origProcs, localMap);
637 
638  //Assume single patch to sample for now
639  const word localOrigPatch(samplePatches[localOrigProcs[0]]);
640  const word localOrigRegion(sampleRegions[localOrigProcs[0]]);
641  List<nearInfoWorld> localNearest(localSamples.size());
642 
643  if (debug)
644  {
645  Pout<< "Searching locally for " << localSamples.size()
646  << " samples on region:" << localOrigRegion
647  << " on patch:" << localOrigPatch << endl;
648  }
649  findLocalSamples
650  (
651  mode,
652  myWorld,
653  localOrigRegion,
654  localOrigPatch,
655  localSamples,
656  localNearest
657  );
658  UIndirectList<nearInfoWorld>(nearest, localMap) = localNearest;
659  }
660  }
661 
662 
663  // Find nearest - globally consistent
665  (
666  nearest,
669  myComm
670  );
671 
672  //if (debug)
673  //{
674  // Pout<< "** After combining:" << endl;
675  // forAll(nearest, samplei)
676  // {
677  // Pout<< " sample:" << samples[samplei]
678  // << " origating from proc:" << origProcs[samplei]
679  // << " to be found on world:" << wantedWorlds[samplei] << nl
680  // << " found on world:" << nearest[samplei].second() << nl
681  // << " found on proc:"
682  // << nearest[samplei].first().second().second() << nl
683  // << " found on patchfacei:"
684  // << nearest[samplei].first().first().index() << nl
685  // << " found at location:"
686  // << nearest[samplei].first().first().point() << nl;
687  // }
688  // Pout<< endl;
689  //}
690 
691  // Convert back into proc+local index
692  sampleProcs.setSize(samples.size());
693  sampleIndices.setSize(samples.size());
694  sampleLocations.setSize(samples.size());
695 
696  forAll(nearest, sampleI)
697  {
698  const nearInfo& ni = nearest[sampleI].first();
699 
700  if (!ni.first().hit())
701  {
702  sampleProcs[sampleI] = -1;
703  sampleIndices[sampleI] = -1;
704  sampleLocations[sampleI] = vector::max;
705  }
706  else
707  {
708  sampleProcs[sampleI] = ni.second().second();
709  sampleIndices[sampleI] = ni.first().index();
710  sampleLocations[sampleI] = ni.first().point();
711  }
712  }
713 
714  Pstream::warnComm = oldWarnComm;
715 }
716 
719 {
721 
722  static bool hasWarned = false;
723  if (mapPtr_)
724  {
726  << "Mapping already calculated" << exit(FatalError);
727  }
728 
729  DebugInFunction << nl;
730 
731 
732  // Early construction of tetBasePtIs since does parallel
733  // communication which might conflict with inter-world communicator
734  // (for multi-world operation)
735  (void)patch_.boundaryMesh().mesh().tetBasePtIs();
736 
737  const label myComm = getCommunicator(); // Get or create
738 
740  //if (sampleDatabase() && !sameWorld())
741  //{
742  // const word syncName("syncObjects");
743  // const polyMesh& mesh = patch_.boundaryMesh().mesh();
744  // const Time& runTime = mesh.time();
745  // const functionObjectList& fos = runTime.functionObjects();
746  //
747  // forAll(fos, i)
748  // {
749  // Pout<< "** FO:" << fos[i].name() << " tpye:" << fos[i].type()
750  // << endl;
751  // }
752  //
753  // if (fos.findObjectID(syncName) == -1)
754  // {
755  // //FatalErrorInFunction
756  // WarningInFunction
757  // << "Did not detect functionObject " << syncName
758  // << ". This is used to synchronise data inbetween worlds"
759  // << endl;
760  // }
761  //}
762 
763 
764  // Get points on face (since cannot use face-centres - might be off
765  // face-diagonal decomposed tets.
766  tmp<pointField> patchPoints(facePoints(patch_));
767 
768  // Get offsetted points
769  const pointField offsettedPoints(samplePoints(patchPoints()));
770 
771  // Do a sanity check - am I sampling my own patch?
772  // This only makes sense for a non-zero offset.
773  bool sampleMyself =
774  (
775  mode_ == NEARESTPATCHFACE
776  && sampleWorld() == UPstream::myWorld()
777  && sampleRegion() == patch_.boundaryMesh().mesh().name()
778  && samplePatch() == patch_.name()
779  );
780 
781  if (sampleMyself)
782  {
783  // Check offset
784  vectorField d(offsettedPoints-patchPoints());
785  bool coincident = (gAverage(mag(d)) <= ROOTVSMALL);
786 
787  if (sampleMyself && coincident)
788  {
790  << "Invalid offset " << d << endl
791  << "Offset is the vector added to the patch face centres to"
792  << " find the patch face supplying the data." << endl
793  << "Setting it to " << d
794  << " on the same patch, on the same region"
795  << " will find the faces themselves which does not make sense"
796  << " for anything but testing." << endl
797  << "patch:" << patch_.name() << endl
798  << "sampleRegion:" << sampleRegion() << endl
799  << "mode:" << sampleModeNames_[mode_] << endl
800  << "samplePatch:" << samplePatch() << endl
801  << "offsetMode:" << offsetModeNames_[offsetMode_] << endl;
802  }
803  }
804 
805  // Get local world and world-to-sample in index form
806  const label myWorld = Pstream::myWorldID();
807  const label mySampleWorld = Pstream::allWorlds().find(sampleWorld_);
808 
809 
810  // Get global list of all samples and the processor and face they come from.
811  pointField samples; // coordinates
812  labelList patchFaceWorlds; // world to sample
813  labelList patchFaceProcs; // originating processor
814  labelList patchFaces; // originating face on processor patch
815  pointField patchFc; // originating face centre
816  collectSamples
817  (
818  mySampleWorld, // world I want to sample
819  patchPoints,
820 
821  samples,
822  patchFaceWorlds,
823  patchFaceProcs,
824  patchFaces,
825  patchFc
826  );
827 
828  //if (debug)
829  //{
830  // forAll(samples, samplei)
831  // {
832  // Pout<< " sample:" << samples[samplei]
833  // << " origating from proc:" << patchFaceProcs[samplei]
834  // << " face:" << patchFaces[samplei]
835  // << " to be found on world:" << patchFaceWorlds[samplei] << nl;
836  // }
837  //}
838 
839  // Find processor and cell/face samples are in and actual location.
840  labelList sampleProcs;
841  labelList sampleIndices;
842  pointField sampleLocations;
843  findSamples
844  (
845  mode_,
846  myWorld, // my world (in index form)
847  samples,
848  patchFaceWorlds,
849  patchFaceProcs,
850 
851  sampleProcs,
852  sampleIndices,
853  sampleLocations
854  );
855 
856  // Check for samples that were not found. This will only happen for
857  // NEARESTCELL since finds cell containing a location
858  if (mode_ == NEARESTCELL)
859  {
860  label nNotFound = 0;
861  forAll(sampleProcs, sampleI)
862  {
863  if (sampleProcs[sampleI] == -1)
864  {
865  nNotFound++;
866  }
867  }
868  reduce(nNotFound, sumOp<label>(), Pstream::msgType(), myComm);
869 
870  if (nNotFound > 0)
871  {
872  if (!hasWarned)
873  {
875  << "Did not find " << nNotFound
876  << " out of " << sampleProcs.size() << " total samples."
877  << " Sampling these on owner cell centre instead." << endl
878  << "On patch " << patch_.name()
879  << " on region " << sampleRegion()
880  << " in mode " << sampleModeNames_[mode_] << endl
881  << "with offset mode " << offsetModeNames_[offsetMode_]
882  << ". Suppressing further warnings from " << type() << endl;
883 
884  hasWarned = true;
885  }
886 
887  // Collect the samples that cannot be found
888  DynamicList<label> subMap;
889  forAll(sampleProcs, sampleI)
890  {
891  if (sampleProcs[sampleI] == -1)
892  {
893  subMap.append(sampleI);
894  }
895  }
896 
897  // And re-search for pure nearest (should not fail)
898  labelList subSampleProcs;
899  labelList subSampleIndices;
900  pointField subSampleLocations;
901  findSamples
902  (
903  NEARESTONLYCELL,
904  myWorld, // my world (in index form)
905 
906  pointField(samples, subMap),
907  UIndirectList<label>(patchFaceWorlds, subMap)(),
908  UIndirectList<label>(patchFaceProcs, subMap)(),
909 
910  subSampleProcs,
911  subSampleIndices,
912  subSampleLocations
913  );
914 
915  // Insert
916  labelUIndList(sampleProcs, subMap) = subSampleProcs;
917  labelUIndList(sampleIndices, subMap) = subSampleIndices;
918  UIndirectList<point>(sampleLocations, subMap) = subSampleLocations;
919  }
920  }
921 
922  // Now we have all the data we need:
923  // - where sample originates from (so destination when mapping):
924  // patchFaces, patchFaceProcs.
925  // - cell/face sample is in (so source when mapping)
926  // sampleIndices, sampleProcs.
927 
928  if (Pstream::master(myComm))
929  {
930  forAll(samples, i)
931  {
932  if (sampleProcs[i] == -1)
933  {
934  FatalErrorInFunction << "did not find sample "
935  << patchFc[i] << " on patch " << patch_.name()
936  << " on region "
937  << patch_.boundaryMesh().mesh().name()
938  << " on processor " << patchFaceProcs[i]
939  << exit(FatalError);
940  }
941  }
942  }
943 
944 
945  if (debug && Pstream::master(myComm))
946  {
947  //forAll(samples, i)
948  //{
949  // Pout<< i << " need data in region "
950  // << patch_.boundaryMesh().mesh().name()
951  // << " for proc:" << patchFaceProcs[i]
952  // << " face:" << patchFaces[i]
953  // << " at:" << patchFc[i] << endl
954  // << "Found data in region " << sampleRegion()
955  // << " at proc:" << sampleProcs[i]
956  // << " face:" << sampleIndices[i]
957  // << " at:" << sampleLocations[i]
958  // << nl << endl;
959  //}
960 
961  OFstream str
962  (
963  patch_.boundaryMesh().mesh().time().path()
964  / patch_.name()
965  + "_mapped.obj"
966  );
967  Pout<< "Dumping mapping as lines from patch faceCentres to"
968  << " sampled cell/faceCentres/points to file " << str.name()
969  << endl;
970 
971  label vertI = 0;
972 
973  forAll(patchFc, i)
974  {
975  meshTools::writeOBJ(str, patchFc[i]);
976  vertI++;
977  meshTools::writeOBJ(str, sampleLocations[i]);
978  vertI++;
979  str << "l " << vertI-1 << ' ' << vertI << nl;
980  }
981  }
982 
983  // Determine schedule.
984  mapPtr_.reset(new mapDistribute(sampleProcs, patchFaceProcs, myComm));
985 
986  // Rework the schedule from indices into samples to cell data to send,
987  // face data to receive.
988 
989  labelListList& subMap = mapPtr_().subMap();
990  labelListList& constructMap = mapPtr_().constructMap();
991 
992  forAll(subMap, proci)
993  {
994  subMap[proci] = labelUIndList(sampleIndices, subMap[proci]);
995  constructMap[proci] = labelUIndList(patchFaces, constructMap[proci]);
996 
997  if (debug)
998  {
999  Pout<< "To proc:" << proci << " sending values of cells/faces:"
1000  << subMap[proci] << endl;
1001  Pout<< "From proc:" << proci
1002  << " receiving values of patch faces:"
1003  << constructMap[proci] << endl;
1004  }
1005  }
1006 
1007 
1008  // Redo constructSize
1009  mapPtr_().constructSize() = patch_.size();
1010 
1011  if (debug)
1012  {
1013  // Check that all elements get a value.
1014  bitSet used(patch_.size());
1015  forAll(constructMap, proci)
1016  {
1017  const labelList& map = constructMap[proci];
1018 
1019  forAll(map, i)
1020  {
1021  label facei = map[i];
1022 
1023  if (used.test(facei))
1024  {
1026  << "On patch " << patch_.name()
1027  << " patchface " << facei
1028  << " is assigned to more than once."
1029  << abort(FatalError);
1030  }
1031  else
1032  {
1033  used.set(facei);
1034  }
1035  }
1036  }
1037  forAll(used, facei)
1038  {
1039  if (!used.test(facei))
1040  {
1042  << "On patch " << patch_.name()
1043  << " patchface " << facei
1044  << " is never assigned to."
1045  << abort(FatalError);
1046  }
1047  }
1048  }
1049 
1050  updateMeshTime().setUpToDate();
1051  if (sameWorld())
1052  {
1053  updateSampleMeshTime().setUpToDate();
1054  }
1055 }
1056 
1059 const
1060 {
1061  const word surfType(surfDict_.getOrDefault<word>("type", "none"));
1062 
1063  if (!surfPtr_ && surfType != "none")
1064  {
1065  word surfName(surfDict_.getOrDefault("name", patch_.name()));
1066 
1067  const polyMesh& mesh = patch_.boundaryMesh().mesh();
1068 
1069  surfPtr_ =
1071  (
1072  surfType,
1073  IOobject
1074  (
1075  surfName,
1076  mesh.time().constant(),
1077  "triSurface",
1078  mesh,
1081  ),
1082  surfDict_
1083  );
1084  }
1085 
1086  return surfPtr_;
1087 }
1088 
1090 void Foam::mappedPatchBase::calcAMI() const
1091 {
1092  if (AMIPtr_->upToDate())
1093  {
1095  << "AMI already up-to-date"
1096  << endl;
1097 
1098  return;
1099  }
1100 
1101  DebugInFunction << nl;
1102 
1103  const label myComm = getCommunicator(); // Get or create
1104 
1105  // Pre-calculate surface (if any)
1106  const auto& surf = surfPtr();
1107 
1108  const label oldWorldComm(Pstream::worldComm);
1109  const label oldWarnComm(Pstream::warnComm);
1110 
1111  // Check if running locally
1112  if (sampleWorld_.empty() || sameWorld())
1113  {
1114  const polyPatch& nbr = samplePolyPatch();
1115 
1116  // Transform neighbour patch to local system
1117  const pointField nbrPoints(samplePoints(nbr.localPoints()));
1118 
1119  const primitivePatch nbrPatch0
1120  (
1122  (
1123  nbr.localFaces(),
1124  nbr.size()
1125  ),
1126  nbrPoints
1127  );
1128 
1129 
1130  if (debug)
1131  {
1132  OFstream os(patch_.name() + "_neighbourPatch-org.obj");
1133  meshTools::writeOBJ(os, samplePolyPatch().localFaces(), nbrPoints);
1134 
1135  OFstream osN(patch_.name() + "_neighbourPatch-trans.obj");
1136  meshTools::writeOBJ(osN, nbrPatch0, nbrPoints);
1137 
1138  OFstream osO(patch_.name() + "_ownerPatch.obj");
1139  meshTools::writeOBJ(osO, patch_.localFaces(), patch_.localPoints());
1140  }
1141 
1142  // Construct/apply AMI interpolation to determine addressing and
1143  // weights.
1144 
1145  // Change to use inter-world communicator
1146  Pstream::worldComm = myComm;
1148 
1149  AMIPtr_->calculate(patch_, nbrPatch0, surf);
1150 
1151  Pstream::warnComm = oldWarnComm;
1152  Pstream::worldComm = oldWorldComm;
1153  }
1154  else
1155  {
1156  faceList dummyFaces;
1157  pointField dummyPoints;
1158  const primitivePatch dummyPatch
1159  (
1160  SubList<face>(dummyFaces),
1161  dummyPoints
1162  );
1163 
1164  // Change to use inter-world communicator
1165  Pstream::worldComm = myComm;
1167 
1168  if (masterWorld())
1169  {
1170  // Construct/apply AMI interpolation to determine addressing
1171  // and weights. Have patch_ for src faces, 0 faces for the
1172  // target side
1173  AMIPtr_->calculate(patch_, dummyPatch, surf);
1174  }
1175  else
1176  {
1177  // Construct/apply AMI interpolation to determine addressing
1178  // and weights. Have 0 faces for src side, patch_ for the tgt
1179  // side
1180  AMIPtr_->calculate(dummyPatch, patch_, surf);
1181  }
1182  // Now the AMI addressing/weights will be from src side (on masterWorld
1183  // processors) to tgt side (on other processors)
1184 
1185  Pstream::warnComm = oldWarnComm;
1186  Pstream::worldComm = oldWorldComm;
1187  }
1188 }
1189 
1190 
1192 (
1193  const objectRegistry& obr,
1194  const wordList& names,
1195  const label index
1196 )
1197 {
1198  const objectRegistry& sub = obr.subRegistry(names[index], true);
1199  if (index == names.size()-1)
1200  {
1201  return sub;
1202  }
1203  else
1204  {
1205  return subRegistry(sub, names, index+1);
1206  }
1207 }
1208 
1209 
1210 // * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * * * * * //
1213 :
1214  patch_(pp),
1215  sampleWorld_(),
1216  sampleRegion_(patch_.boundaryMesh().mesh().name()),
1217  mode_(NEARESTPATCHFACE),
1218  samplePatch_(),
1219  coupleGroup_(),
1220  sampleDatabasePtr_(),
1221  offsetMode_(UNIFORM),
1222  offset_(Zero),
1223  offsets_(pp.size(), offset_),
1224  distance_(0),
1225  communicator_(-1), // Demand-driven (cached value)
1226  sameRegion_(true),
1227  mapPtr_(nullptr),
1228  AMIReverse_(false),
1229  AMIPtr_(new faceAreaWeightAMI(true, AMIReverse_)),
1230  surfPtr_(nullptr),
1231  surfDict_(fileName("surface"))
1232 {
1233  // NOTE: same region, no sample-world. Thus no world-world communication
1234 }
1235 
1236 
1238 (
1239  const polyPatch& pp,
1240  const word& sampleRegion,
1241  const sampleMode mode,
1242  const word& samplePatch,
1243  const vectorField& offsets
1244 )
1245 :
1247  (
1248  pp,
1249  sampleRegion,
1250  mode,
1251  samplePatch,
1252  scalar(0)
1253  )
1254 {
1256 }
1257 
1258 
1260 (
1261  const polyPatch& pp,
1262  const word& sampleRegion,
1263  const sampleMode mode,
1264  const word& samplePatch,
1265  const vector& uniformOffset
1266 )
1267 :
1269  (
1270  pp,
1271  sampleRegion,
1272  mode,
1273  samplePatch,
1274  scalar(0)
1275  )
1276 {
1277  mappedPatchBase::setOffset(uniformOffset);
1278 }
1279 
1280 
1282 (
1283  const polyPatch& pp,
1284  const word& sampleRegion,
1285  const sampleMode mode,
1286  const word& samplePatch,
1287  const scalar normalDistance
1288 )
1289 :
1290  patch_(pp),
1291  sampleWorld_(),
1292  sampleRegion_(sampleRegion),
1293  mode_(mode),
1294  samplePatch_(samplePatch),
1295  coupleGroup_(),
1296  sampleDatabasePtr_(),
1297  offsetMode_(NORMAL),
1298  offset_(Zero),
1299  offsets_(0),
1300  distance_(normalDistance),
1301  communicator_(-1), // Demand-driven (cached value)
1302  sameRegion_
1303  (
1304  sampleWorld_.empty()
1305  && sampleRegion_ == patch_.boundaryMesh().mesh().name()
1306  ),
1307  mapPtr_(nullptr),
1308  AMIReverse_(false),
1309  AMIPtr_(new faceAreaWeightAMI(true, AMIReverse_)),
1310  surfPtr_(nullptr),
1311  surfDict_(fileName("surface"))
1312 {
1314 }
1315 
1316 
1318 (
1319  const polyPatch& pp,
1320  const dictionary& dict
1321 )
1322 :
1323  patch_(pp),
1324  sampleWorld_(dict.getOrDefault<word>("sampleWorld", word::null)),
1325  sampleRegion_(dict.getOrDefault<word>("sampleRegion", word::null)),
1326  mode_(sampleModeNames_.get("sampleMode", dict)),
1327  samplePatch_(dict.getOrDefault<word>("samplePatch", word::null)),
1328  coupleGroup_(dict),
1329  sampleDatabasePtr_(readDatabase(dict)),
1330  offsetMode_(UNIFORM),
1331  offset_(Zero),
1332  offsets_(0),
1333  distance_(0),
1334  communicator_(-1), // Demand-driven (cached value)
1335  sameRegion_
1336  (
1337  sampleWorld_.empty()
1338  && sampleRegion_ == patch_.boundaryMesh().mesh().name()
1339  ),
1340  mapPtr_(nullptr),
1341  AMIReverse_
1342  (
1343  dict.getOrDefault
1344  (
1345  "reverseTarget", // AMIInterpolation uses this keyword
1346  dict.getOrDefault
1347  (
1348  "flipNormals",
1349  false
1350  )
1351  )
1352  ),
1353  AMIPtr_
1354  (
1356  (
1357  dict.getOrDefault("AMIMethod", faceAreaWeightAMI::typeName),
1358  dict,
1359  AMIReverse_
1360  )
1361  ),
1362  surfPtr_(nullptr),
1363  surfDict_(dict.subOrEmptyDict("surface"))
1364 {
1366 
1367  if (!coupleGroup_.valid())
1368  {
1369  if (sampleWorld_.empty() && sampleRegion_.empty())
1370  {
1371  // If no coupleGroup and no sampleRegion assume local region
1373  sameRegion_ = true;
1374  }
1375  }
1376 
1377  if (offsetModeNames_.readIfPresent("offsetMode", dict, offsetMode_))
1378  {
1379  switch (offsetMode_)
1380  {
1381  case UNIFORM:
1382  {
1383  dict.readEntry("offset", offset_);
1384  }
1385  break;
1386 
1387  case NONUNIFORM:
1388  {
1389  offsets_ = pointField("offsets", dict, patch_.size());
1390  }
1391  break;
1392 
1393  case NORMAL:
1394  {
1395  dict.readEntry("distance", distance_);
1396  }
1397  break;
1398  }
1399  }
1400  else if (dict.readIfPresent("offset", offset_))
1401  {
1402  offsetMode_ = UNIFORM;
1403  }
1404  else if (dict.found("offsets"))
1405  {
1407  offsets_ = pointField("offsets", dict, patch_.size());
1408  }
1409  else if (mode_ != NEARESTPATCHFACE && mode_ != NEARESTPATCHFACEAMI)
1410  {
1412  << "Please supply the offsetMode as one of "
1413  << offsetModeNames_
1414  << exit(FatalIOError);
1415  }
1416 }
1417 
1418 
1420 (
1421  const polyPatch& pp,
1422  const sampleMode mode,
1423  const dictionary& dict
1424 )
1425 :
1426  patch_(pp),
1427  sampleWorld_(dict.getOrDefault<word>("sampleWorld", word::null)),
1428  sampleRegion_(dict.getOrDefault<word>("sampleRegion", word::null)),
1429  mode_(mode),
1430  samplePatch_(dict.getOrDefault<word>("samplePatch", word::null)),
1431  coupleGroup_(dict), //dict.getOrDefault<word>("coupleGroup", word::null)),
1432  sampleDatabasePtr_(readDatabase(dict)),
1433  offsetMode_(UNIFORM),
1434  offset_(Zero),
1435  offsets_(0),
1436  distance_(0),
1437  communicator_(-1), // Demand-driven (cached value)
1438  sameRegion_
1439  (
1440  sampleWorld_.empty()
1441  && sampleRegion_ == patch_.boundaryMesh().mesh().name()
1442  ),
1443  mapPtr_(nullptr),
1444  AMIReverse_
1445  (
1446  dict.getOrDefault
1447  (
1448  "reverseTarget", // AMIInterpolation uses this keyword
1449  dict.getOrDefault
1450  (
1451  "flipNormals",
1452  false
1453  )
1454  )
1455  ),
1456  AMIPtr_
1457  (
1459  (
1460  dict.getOrDefault("AMIMethod", faceAreaWeightAMI::typeName),
1461  dict,
1462  AMIReverse_
1463  )
1464  ),
1465  surfPtr_(nullptr),
1466  surfDict_(dict.subOrEmptyDict("surface"))
1467 {
1469 
1471  {
1473  << "Construct from sampleMode and dictionary only applicable for "
1474  << " collocated patches in modes "
1477  << exit(FatalIOError);
1478  }
1479 
1480  if (!coupleGroup_.valid())
1481  {
1482  if (sampleWorld_.empty() && sampleRegion_.empty())
1483  {
1484  // If no coupleGroup and no sampleRegion assume local region
1486  sameRegion_ = true;
1487  }
1488  }
1489 }
1490 
1491 
1493 (
1494  const polyPatch& pp,
1495  const mappedPatchBase& mpb
1496 )
1497 :
1498  patch_(pp),
1499  sampleWorld_(mpb.sampleWorld_),
1500  sampleRegion_(mpb.sampleRegion_),
1501  mode_(mpb.mode_),
1502  samplePatch_(mpb.samplePatch_),
1503  coupleGroup_(mpb.coupleGroup_),
1504  sampleDatabasePtr_
1505  (
1506  mpb.sampleDatabasePtr_
1507  ? new fileName(mpb.sampleDatabasePtr_())
1508  : nullptr
1509  ),
1510  offsetMode_(mpb.offsetMode_),
1511  offset_(mpb.offset_),
1512  offsets_(mpb.offsets_),
1513  distance_(mpb.distance_),
1514  communicator_(mpb.communicator_),
1515  sameRegion_(mpb.sameRegion_),
1516  mapPtr_(nullptr),
1517  AMIReverse_(mpb.AMIReverse_),
1518  AMIPtr_(mpb.AMIPtr_->clone()),
1519  surfPtr_(nullptr),
1520  surfDict_(mpb.surfDict_)
1521 {}
1522 
1523 
1525 (
1526  const polyPatch& pp,
1527  const mappedPatchBase& mpb,
1528  const labelUList& mapAddressing
1529 )
1530 :
1531  patch_(pp),
1532  sampleWorld_(mpb.sampleWorld_),
1533  sampleRegion_(mpb.sampleRegion_),
1534  mode_(mpb.mode_),
1535  samplePatch_(mpb.samplePatch_),
1536  coupleGroup_(mpb.coupleGroup_),
1537  sampleDatabasePtr_
1538  (
1539  mpb.sampleDatabasePtr_
1540  ? new fileName(mpb.sampleDatabasePtr_())
1541  : nullptr
1542  ),
1543  offsetMode_(mpb.offsetMode_),
1544  offset_(mpb.offset_),
1545  offsets_
1546  (
1547  offsetMode_ == NONUNIFORM
1548  ? vectorField(mpb.offsets_, mapAddressing)
1549  : vectorField()
1550  ),
1551  distance_(mpb.distance_),
1552  communicator_(mpb.communicator_),
1553  sameRegion_(mpb.sameRegion_),
1554  mapPtr_(nullptr),
1555  AMIReverse_(mpb.AMIReverse_),
1556  AMIPtr_(mpb.AMIPtr_->clone()),
1557  surfPtr_(nullptr),
1558  surfDict_(mpb.surfDict_)
1559 {}
1560 
1561 
1562 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
1565 {
1566  clearOut();
1567 }
1568 
1571 {
1572  mapPtr_.reset(nullptr);
1573  surfPtr_.reset(nullptr);
1574  AMIPtr_->upToDate() = false;
1575 }
1576 
1577 
1578 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
1580 void Foam::mappedPatchBase::setOffset(const scalar normalDist)
1581 {
1582  clearOut();
1583  offsetMode_ = offsetMode::NORMAL;
1584  offset_ = Zero;
1585  offsets_.clear();
1586  distance_ = normalDist;
1587 }
1588 
1590 void Foam::mappedPatchBase::setOffset(const vector& uniformOffset)
1591 {
1592  clearOut();
1593  offsetMode_ = offsetMode::UNIFORM;
1594  offset_ = uniformOffset;
1595  offsets_.clear();
1596  distance_ = Zero;
1597 }
1598 
1600 void Foam::mappedPatchBase::setOffset(const vectorField& offsets)
1601 {
1602  clearOut();
1603  offsetMode_ = offsetMode::NONUNIFORM;
1604  offset_ = Zero;
1605  offsets_ = offsets;
1606  distance_ = Zero;
1607 }
1608 
1609 
1611 (
1612  const word& sampleRegion
1613 ) const
1614 {
1615  const polyMesh& thisMesh = patch_.boundaryMesh().mesh();
1616  return
1617  (
1618  sampleRegion.empty() || sampleRegion == thisMesh.name()
1619  ? thisMesh
1620  : thisMesh.time().lookupObject<polyMesh>(sampleRegion)
1621  );
1622 }
1623 
1624 
1626 (
1627  const word& sampleRegion,
1628  const word& samplePatch
1629 ) const
1630 {
1631  const polyMesh& nbrMesh = lookupMesh(sampleRegion);
1632 
1633  const label patchi = nbrMesh.boundaryMesh().findPatchID(samplePatch);
1634 
1635  if (patchi == -1)
1636  {
1638  << "Cannot find patch " << samplePatch
1639  << " in region " << sampleRegion_ << endl
1640  << exit(FatalError);
1641  }
1642  return nbrMesh.boundaryMesh()[patchi];
1643 }
1644 
1647 {
1648  if (UPstream::myWorld() != sampleWorld_)
1649  {
1651  << "sampleWorld : " << sampleWorld_
1652  << " is not the current world : " << UPstream::myWorld()
1653  << exit(FatalError);
1654  }
1655  return lookupMesh(sampleRegion());
1656 }
1657 
1660 {
1661  const polyMesh& nbrMesh = sampleMesh();
1662 
1663  const label patchi = nbrMesh.boundaryMesh().findPatchID(samplePatch());
1664 
1665  if (patchi == -1)
1666  {
1668  << "Cannot find patch " << samplePatch()
1669  << " in region " << sampleRegion_ << endl
1670  << "Valid patches are " << nbrMesh.boundaryMesh().names()
1671  << exit(FatalError);
1672  }
1673 
1674  return nbrMesh.boundaryMesh()[patchi];
1675 }
1676 
1677 
1679 (
1680  const pointField& fc
1681 ) const
1682 {
1683  auto tfld = tmp<pointField>::New(fc);
1684  auto& fld = tfld.ref();
1685 
1686  switch (offsetMode_)
1687  {
1688  case UNIFORM:
1689  {
1690  fld += offset_;
1691  break;
1692  }
1693  case NONUNIFORM:
1694  {
1695  fld += offsets_;
1696  break;
1697  }
1698  case NORMAL:
1699  {
1700  // Get outwards pointing normal
1701  vectorField n(patch_.faceAreas());
1702  n /= mag(n);
1703 
1704  fld += distance_*n;
1705  break;
1706  }
1707  }
1708 
1709  return tfld;
1710 }
1711 
1714 {
1715  return samplePoints(facePoints(patch_));
1716 }
1717 
1718 
1720 (
1721  const polyMesh& mesh,
1722  const label facei,
1723  const polyMesh::cellDecomposition decompMode
1724 )
1725 {
1726  const point& fc = mesh.faceCentres()[facei];
1727 
1728  switch (decompMode)
1729  {
1730  case polyMesh::FACE_PLANES:
1732  {
1733  // For both decompositions the face centre is guaranteed to be
1734  // on the face
1735  return pointIndexHit(true, fc, facei);
1736  }
1737  break;
1738 
1740  case polyMesh::CELL_TETS:
1741  {
1742  // Find the intersection of a ray from face centre to cell centre
1743  // Find intersection of (face-centre-decomposition) centre to
1744  // cell-centre with face-diagonal-decomposition triangles.
1745 
1746  const pointField& p = mesh.points();
1747  const face& f = mesh.faces()[facei];
1748 
1749  if (f.size() <= 3)
1750  {
1751  // Return centre of triangle.
1752  return pointIndexHit(true, fc, 0);
1753  }
1754 
1755  const label celli = mesh.faceOwner()[facei];
1756  const point& cc = mesh.cellCentres()[celli];
1757  vector d = fc-cc;
1758 
1759  const label fp0 = mesh.tetBasePtIs()[facei];
1760  const point& basePoint = p[f[fp0]];
1761 
1762  label fp = f.fcIndex(fp0);
1763  for (label i = 2; i < f.size(); i++)
1764  {
1765  const point& thisPoint = p[f[fp]];
1766  label nextFp = f.fcIndex(fp);
1767  const point& nextPoint = p[f[nextFp]];
1768 
1769  const triPointRef tri(basePoint, thisPoint, nextPoint);
1770  pointHit hitInfo = tri.intersection
1771  (
1772  cc,
1773  d,
1775  );
1776 
1777  if (hitInfo.hit() && hitInfo.distance() > 0)
1778  {
1779  return pointIndexHit(true, hitInfo.point(), i-2);
1780  }
1781 
1782  fp = nextFp;
1783  }
1784 
1785  // Fall-back
1786  return pointIndexHit(false, fc, -1);
1787  }
1788  break;
1789 
1790  default:
1791  {
1793  << "problem" << abort(FatalError);
1794  return pointIndexHit();
1795  }
1796  }
1797 }
1798 
1799 
1801 (
1802  const objectRegistry& obr,
1803  const fileName& path
1804 )
1805 {
1806  // Lookup (and create if non-existing) a registry using
1807  // '/' separated path. Like 'mkdir -p'
1808 
1809  fileName cleanedPath(path);
1810  cleanedPath.clean(); // Remove unneeded ".."
1811  const wordList names(cleanedPath.components());
1812 
1813  if (names.empty())
1814  {
1815  return obr;
1816  }
1817  else
1818  {
1819  return subRegistry(obr, names, 0);
1820  }
1821 }
1822 
1823 
1825 (
1826  const fileName& root,
1827  const label proci
1828 )
1829 {
1830  const word processorName("processor"+Foam::name(proci));
1831  return root/"send"/processorName;
1832 }
1833 
1835 Foam::fileName Foam::mappedPatchBase::sendPath(const label proci) const
1836 {
1837  return sendPath(sampleDatabasePath(), proci);
1838 }
1839 
1840 
1842 (
1843  const fileName& root,
1844  const label proci
1845 )
1846 {
1847  const word processorName("processor"+Foam::name(proci));
1848  return root/"receive"/processorName;
1849 }
1850 
1852 Foam::fileName Foam::mappedPatchBase::receivePath(const label proci) const
1853 {
1854  return receivePath(sampleDatabasePath(), proci);
1855 }
1856 
1857 
1859 (
1860  const objectRegistry& obr,
1861  dictionary& dict
1862 )
1863 {
1864  forAllIters(obr, iter)
1865  {
1866  regIOobject* objPtr = iter.val();
1867  const regIOobject& obj = *objPtr;
1868 
1869  if (isA<objectRegistry>(obj))
1870  {
1871  dictionary& d = dict.subDictOrAdd(obj.name());
1872  writeDict(dynamic_cast<const objectRegistry&>(obj), d);
1873  }
1874  else if
1875  (
1876  writeIOField<scalar>(obj, dict)
1877  || writeIOField<vector>(obj, dict)
1878  || writeIOField<sphericalTensor>(obj, dict)
1879  || writeIOField<symmTensor>(obj, dict)
1880  || writeIOField<tensor>(obj, dict)
1881  )
1882  {
1883  // IOField specialisation
1884  }
1885  else
1886  {
1887  // Not tested. No way of retrieving data anyway ...
1888  OTstream os;
1889  obj.writeData(os);
1890 
1891  primitiveEntry* pePtr = new primitiveEntry(obj.name(), os.tokens());
1892  dict.add(pePtr);
1893  }
1894  }
1895 }
1896 
1899 {
1900  // Reverse of writeDict - reads fields from dictionary into objectRegistry
1901  for (const entry& e : d)
1902  {
1903  if (e.isDict())
1904  {
1905  // Add sub registry
1906  objectRegistry& sub = const_cast<objectRegistry&>
1907  (
1908  obr.subRegistry(e.keyword(), true)
1909  );
1910 
1911  readDict(e.dict(), sub);
1912  }
1913  else
1914  {
1915  ITstream& is = e.stream();
1916  token tok(is);
1917 
1918  if
1919  (
1920  constructIOField<scalar>(e.keyword(), tok, is, obr)
1921  || constructIOField<vector>(e.keyword(), tok, is, obr)
1922  || constructIOField<sphericalTensor>(e.keyword(), tok, is, obr)
1923  || constructIOField<symmTensor>(e.keyword(), tok, is, obr)
1924  || constructIOField<tensor>(e.keyword(), tok, is, obr)
1925  )
1926  {
1927  // Done storing field
1928  }
1929  else
1930  {
1931  FatalErrorInFunction << "Unsupported type " << e.keyword()
1932  << exit(FatalError);
1933  }
1934  }
1935  }
1936 }
1937 
1940 {
1941  os.writeEntry("sampleMode", sampleModeNames_[mode_]);
1942  os.writeEntryIfDifferent<word>("sampleWorld", word::null, sampleWorld_);
1943  os.writeEntryIfDifferent<word>("sampleRegion", word::null, sampleRegion_);
1944  os.writeEntryIfDifferent<word>("samplePatch", word::null, samplePatch_);
1945 
1946  if (sampleDatabasePtr_)
1947  {
1948  os.writeEntry("sampleDatabase", Switch::name(true));
1949  // Write database path if differing
1951  (
1952  "sampleDatabasePath",
1954  sampleDatabasePtr_()
1955  );
1956  }
1957  coupleGroup_.write(os);
1958 
1959  if
1960  (
1961  offsetMode_ == UNIFORM
1962  && offset_ == vector::zero
1963  && (mode_ == NEARESTPATCHFACE || mode_ == NEARESTPATCHFACEAMI)
1964  )
1965  {
1966  // Collocated mode. No need to write offset data
1967  }
1968  else
1969  {
1970  os.writeEntry("offsetMode", offsetModeNames_[offsetMode_]);
1971 
1972  switch (offsetMode_)
1973  {
1974  case UNIFORM:
1975  {
1976  os.writeEntry("offset", offset_);
1977  break;
1978  }
1979  case NONUNIFORM:
1980  {
1981  offsets_.writeEntry("offsets", os);
1982  break;
1983  }
1984  case NORMAL:
1985  {
1986  os.writeEntry("distance", distance_);
1987  break;
1988  }
1989  }
1990  }
1991 
1992  if (mode_ == NEARESTPATCHFACEAMI)
1993  {
1994  if (AMIPtr_)
1995  {
1996  // Use AMI to write itself. Problem: outputs:
1997  // - restartUncoveredSourceFace
1998  // - reverseTarget (instead of flipNormals)
1999  AMIPtr_->write(os);
2000  }
2001  if (!surfDict_.empty())
2002  {
2003  surfDict_.writeEntry(surfDict_.dictName(), os);
2004  }
2005  }
2006 }
2007 
2008 
2009 // ************************************************************************* //
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.
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:584
static const word & myWorld()
My world.
Definition: UPstream.H:739
virtual const fileName & name() const
Get the name of the output serial stream. (eg, the name of the Fstream file name) ...
Definition: OSstream.H:128
use face normal + distance
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:51
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:235
dictionary dict
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:118
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:71
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:120
tmp< pointField > samplePoints() const
Get the sample points.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:578
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:55
const polyMesh & sampleMesh() const
Get the region mesh.
bool valid() const noexcept
Is a valid patchGroup (non-empty) name.
label find(const T &val, label pos=0) const
Find index of the first occurrence of the value.
Definition: UList.C:204
const labelIOList & tetBasePtIs() const
Return the tetBasePtIs.
Definition: polyMesh.C:893
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:150
bool test(const Key &key) const
Same as found() - return true if key exists in the set.
Definition: HashSet.H:213
static const meshSearchMeshObject & New(const polyMesh &mesh, Args &&... args)
Get existing or create a new MeshObject.
Definition: MeshObject.C:41
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)
const polyBoundaryMesh & boundaryMesh() const
Return boundaryMesh reference.
Definition: polyPatch.C:309
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:49
bool empty() const noexcept
True if the UList is empty (ie, size() is zero)
Definition: UListI.H:420
A token holds an item read from Istream.
Definition: token.H:64
static const fileName null
An empty fileName.
Definition: fileName.H:110
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:798
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:487
PointIndexHit< point > pointIndexHit
A PointIndexHit with a 3D point.
Definition: pointIndexHit.H:58
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:312
static int & msgType() noexcept
Message tag of standard messages.
Definition: UPstream.H:806
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:637
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)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:688
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:361
static label worldComm
Default world communicator (all processors). May differ from globalComm if local worlds are in use...
Definition: UPstream.H:361
const objectRegistry & subRegistry(const word &name, const bool forceCreate=false, const bool recursive=false) const
Lookup and return a const sub-objectRegistry.
cellDecomposition
Enumeration defining the decomposition of the cell for.
Definition: polyMesh.H:102
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.
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:493
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1066
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:413
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:53
Face area weighted Arbitrary Mesh Interface (AMI) method.
A keyword and a list of tokens comprise a primitiveEntry. A primitiveEntry can be read...
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 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:752
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) is 1 for serial run.
Definition: UPstream.H:656
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:289
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 INVALID.
Definition: exprTraits.C:52
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)
Definition: labelList.C:31
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:1104
#define forAllIters(container, iter)
Iterate across all elements in the container object.
Definition: stdFoam.H:328
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:203
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
void append(const T &val)
Copy append an element to the end of this list.
Definition: DynamicList.H:558
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1091
const vectorField & cellCentres() const
static label warnComm
Debugging: warn for use of any communicator differing from warnComm.
Definition: UPstream.H:366
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:327
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:55
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:89
wordList components(const char delim='/') const
Return path components as wordList.
Definition: fileName.C:530
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 void broadcasts(const label comm, Type &arg1, Args &&... args)
Broadcast multiple items to all processes in communicator.
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:927
sampleMode mode() const noexcept
What to sample.
UIndirectList< label > labelUIndList
UIndirectList of labels.
Definition: IndirectList.H:65
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...
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:607
Type gAverage(const FieldField< Field, Type > &f)
A List with indirect addressing. Like IndirectList but does not store addressing. ...
Definition: faMatrix.H:57
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)
Am I the master rank.
Definition: UPstream.H:672
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:441
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.
bool readIfPresent(const word &key, const dictionary &dict, EnumType &val) const
Find an entry if present, and assign to T val.
Definition: EnumI.H:125
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:69
label n
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:73
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:715
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:166
mode_t mode(const fileName &name, const bool followLink=true)
Return the file mode, normally following symbolic links.
Definition: POSIX.C:726
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:48
const vectorField & offsets() const noexcept
Offset vectors (from patch faces to destination mesh objects)
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:731
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:157