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