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-2021 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 "triPointRef.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  autoPtr<fileName> dbNamePtr_;
95 
96  if (dict.found("sampleDatabase"))
97  {
98  const bool useDb = dict.get<bool>("sampleDatabase");
99  if (useDb)
100  {
101  dbNamePtr_.set
102  (
103  new fileName
104  (
105  dict.lookupOrDefault<fileName>
106  (
107  "sampleDatabasePath",
108  fileName::null
109  )
110  )
111  );
112  }
113  }
114  else if (dict.found("sampleDatabasePath"))
115  {
116  dbNamePtr_.set(new fileName(dict.get<fileName>("sampleDatabasePath")));
117  }
118 
119  return dbNamePtr_;
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::worldComm;
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  ).rawPoint();
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(Pstream::warnComm);
197  Pstream::warnComm = myComm;
198 
199  if (debug & 2)
200  {
201  Perr<< "patch: " << patch_.name()
202  << "[rank=" << myRank << " procs=" << nProcs
203  << " comm=" << myComm << "] collect samples" << endl;
204  }
205 
206  // Collect all sample points and the faces they come from.
207  {
208  List<pointField> globalFc(nProcs);
209  globalFc[myRank] = facePoints;
210 
211  Pstream::gatherList(globalFc, Pstream::msgType(), myComm);
212  Pstream::scatterList(globalFc, Pstream::msgType(), myComm);
213 
214  // Rework into straight list
215  patchFc = ListListOps::combine<pointField>
216  (
217  globalFc,
219  );
220  }
221 
222  {
223  List<pointField> globalSamples(nProcs);
224  globalSamples[myRank] = samplePoints(facePoints);
225  Pstream::gatherList(globalSamples, Pstream::msgType(), myComm);
226  Pstream::scatterList(globalSamples, Pstream::msgType(), myComm);
227  // Rework into straight list
228  samples = ListListOps::combine<pointField>
229  (
230  globalSamples,
232  );
233  }
234 
235  {
236  labelListList globalFaces(nProcs);
237  globalFaces[myRank] = identity(patch_.size());
238  // Distribute to all processors
239  Pstream::gatherList(globalFaces, Pstream::msgType(), myComm);
240  Pstream::scatterList(globalFaces, Pstream::msgType(), myComm);
241 
242  patchFaces = ListListOps::combine<labelList>
243  (
244  globalFaces,
246  );
247  }
248 
249  {
250  labelList procToWorldIndex(nProcs);
251  procToWorldIndex[myRank] = mySampleWorld;
252  Pstream::gatherList(procToWorldIndex, Pstream::msgType(), myComm);
253  Pstream::scatterList(procToWorldIndex, Pstream::msgType(), myComm);
254 
255  labelList nPerProc(nProcs);
256  nPerProc[myRank] = patch_.size();
257  Pstream::gatherList(nPerProc, Pstream::msgType(), myComm);
258  Pstream::scatterList(nPerProc, Pstream::msgType(), myComm);
259 
260  patchFaceWorlds.setSize(patchFaces.size());
261  patchFaceProcs.setSize(patchFaces.size());
262 
263  label sampleI = 0;
264  forAll(nPerProc, proci)
265  {
266  for (label i = 0; i < nPerProc[proci]; i++)
267  {
268  patchFaceWorlds[sampleI] = procToWorldIndex[proci];
269  patchFaceProcs[sampleI] = proci;
270  sampleI++;
271  }
272  }
273  }
274 
275  Pstream::warnComm = oldWarnComm;
276 }
277 
278 
280 (
281  const sampleMode mode,
282 
283  const label mySampleWorld, // local world to sample == my own world
284  const word& sampleRegion, // local region to sample
285  const word& samplePatch, // local patch to sample
286 
287  const pointField& samples,
288  List<nearInfoWorld>& nearest
289 ) const
290 {
291  DebugInFunction << nl;
292 
293  const label myComm = getCommunicator(); // Get or create
294 
295  // Find the local cell containing the samples
296  const label myRank = Pstream::myProcNo(myComm);
297 
298  // Lookup the correct region
299  const polyMesh& mesh = lookupMesh(sampleRegion);
300 
301  // All the info for nearest. Construct to miss
302  nearest.setSize(samples.size());
303  nearInfoWorld miss;
304  {
305  miss.first().second() = Tuple2<scalar, label>(Foam::sqr(GREAT), -1);
306  miss.second() = -1; // set world to be ignored
307  }
308  nearest = miss;
309 
310  switch (mode)
311  {
312  case NEARESTCELL:
313  {
314  if (samplePatch.size() && samplePatch != "none")
315  {
317  << "No need to supply a patch name when in "
318  << sampleModeNames_[mode] << " mode." << exit(FatalError);
319  }
320 
321  //- Note: face-diagonal decomposition
323 
324  forAll(samples, sampleI)
325  {
326  const point& sample = samples[sampleI];
327  nearInfoWorld& near = nearest[sampleI];
328 
329  label celli = tree.findInside(sample);
330 
331  if (celli == -1)
332  {
333  near.first().second().first() = Foam::sqr(GREAT);
334  near.first().second().second() = myRank;
335  near.second() = mySampleWorld;
336  }
337  else
338  {
339  const point& cc = mesh.cellCentres()[celli];
340 
341  near.first().first() = pointIndexHit
342  (
343  true,
344  cc,
345  celli
346  );
347  near.first().second().first() = magSqr(cc-sample);
348  near.first().second().second() = myRank;
349  near.second() = mySampleWorld;
350  }
351  }
352  break;
353  }
354 
355  case NEARESTONLYCELL:
356  {
357  if (samplePatch.size() && samplePatch != "none")
358  {
360  << "No need to supply a patch name when in "
361  << sampleModeNames_[mode] << " mode." << exit(FatalError);
362  }
363 
364  //- Note: face-diagonal decomposition
366 
367  forAll(samples, sampleI)
368  {
369  const point& sample = samples[sampleI];
370  nearInfoWorld& near = nearest[sampleI];
371 
372  near.first().first() = tree.findNearest(sample, sqr(GREAT));
373  near.first().second().first() = magSqr
374  (
375  near.first().first().hitPoint()
376  -sample
377  );
378  near.first().second().second() = myRank;
379  near.second() = mySampleWorld;
380  }
381  break;
382  }
383 
384  case NEARESTPATCHFACE:
385  {
386  Random rndGen(123456);
387 
388  const polyPatch& pp = lookupPatch(sampleRegion, samplePatch);
389 
390  if (pp.empty())
391  {
392  forAll(samples, sampleI)
393  {
394  nearInfoWorld& near = nearest[sampleI];
395  near.first().second().first() = Foam::sqr(GREAT);
396  near.first().second().second() = myRank;
397  near.second() = mySampleWorld;
398  }
399  }
400  else
401  {
402  treeBoundBox patchBb
403  (
404  treeBoundBox(pp.points(), pp.meshPoints()).extend
405  (
406  rndGen,
407  1e-4
408  )
409  );
410  patchBb.min() -= point::uniform(ROOTVSMALL);
411  patchBb.max() += point::uniform(ROOTVSMALL);
412 
413  indexedOctree<treeDataFace> boundaryTree
414  (
415  treeDataFace // all information needed to search faces
416  (
417  false, // do not cache bb
418  mesh,
419  identity(pp.range()) // boundary faces only
420  ),
421  patchBb, // overall search domain
422  8, // maxLevel
423  10, // leafsize
424  3.0 // duplicity
425  );
426 
427  forAll(samples, sampleI)
428  {
429  const point& sample = samples[sampleI];
430 
431  nearInfoWorld& near = nearest[sampleI];
432  pointIndexHit& nearInfo = near.first().first();
433  nearInfo = boundaryTree.findNearest
434  (
435  sample,
436  magSqr(patchBb.span())
437  );
438 
439  if (!nearInfo.hit())
440  {
441  near.first().second().first() = Foam::sqr(GREAT);
442  near.first().second().second() = myRank;
443  near.second() = mySampleWorld;
444  }
445  else
446  {
447  point fc(pp[nearInfo.index()].centre(pp.points()));
448  nearInfo.setPoint(fc);
449  near.first().second().first() = magSqr(fc-sample);
450  near.first().second().second() = myRank;
451  near.second() = mySampleWorld;
452  }
453  }
454  }
455  break;
456  }
457 
458  case NEARESTPATCHPOINT:
459  {
460  Random rndGen(123456);
461 
462  const polyPatch& pp = lookupPatch(sampleRegion, samplePatch);
463 
464  if (pp.empty())
465  {
466  forAll(samples, sampleI)
467  {
468  nearInfoWorld& near = nearest[sampleI];
469  near.first().second().first() = Foam::sqr(GREAT);
470  near.first().second().second() = myRank;
471  near.second() = mySampleWorld;
472  }
473  }
474  else
475  {
476  // patch (local) points
477  treeBoundBox patchBb
478  (
479  treeBoundBox(pp.points(), pp.meshPoints()).extend
480  (
481  rndGen,
482  1e-4
483  )
484  );
485  patchBb.min() -= point::uniform(ROOTVSMALL);
486  patchBb.max() += point::uniform(ROOTVSMALL);
487 
488  indexedOctree<treeDataPoint> boundaryTree
489  (
490  treeDataPoint // all information needed to search faces
491  (
492  mesh.points(),
493  pp.meshPoints() // selection of points to search on
494  ),
495  patchBb, // overall search domain
496  8, // maxLevel
497  10, // leafsize
498  3.0 // duplicity
499  );
500 
501  forAll(samples, sampleI)
502  {
503  const point& sample = samples[sampleI];
504 
505  nearInfoWorld& near = nearest[sampleI];
506  pointIndexHit& nearInfo = near.first().first();
507  nearInfo = boundaryTree.findNearest
508  (
509  sample,
510  magSqr(patchBb.span())
511  );
512 
513  if (!nearInfo.hit())
514  {
515  near.first().second().first() = Foam::sqr(GREAT);
516  near.first().second().second() = myRank;
517  near.second() = mySampleWorld;
518  }
519  else
520  {
521  const point& pt = nearInfo.hitPoint();
522 
523  near.first().second().first() = magSqr(pt-sample);
524  near.first().second().second() = myRank;
525  near.second() = mySampleWorld;
526  }
527  }
528  }
529  break;
530  }
531 
532  case NEARESTFACE:
533  {
534  if (samplePatch.size() && samplePatch != "none")
535  {
537  << "No need to supply a patch name when in "
538  << sampleModeNames_[mode] << " mode." << exit(FatalError);
539  }
540 
541  //- Note: face-diagonal decomposition
542  const meshSearchMeshObject& meshSearchEngine =
544 
545  forAll(samples, sampleI)
546  {
547  const point& sample = samples[sampleI];
548  nearInfoWorld& near = nearest[sampleI];
549 
550  label facei = meshSearchEngine.findNearestFace(sample);
551 
552  if (facei == -1)
553  {
554  near.first().second().first() = Foam::sqr(GREAT);
555  near.first().second().second() = myRank;
556  near.second() = mySampleWorld;
557  }
558  else
559  {
560  const point& fc = mesh.faceCentres()[facei];
561 
562  near.first().first() = pointIndexHit(true, fc, facei);
563  near.first().second().first() = magSqr(fc-sample);
564  near.first().second().second() = myRank;
565  near.second() = mySampleWorld;
566  }
567  }
568  break;
569  }
570 
571  case NEARESTPATCHFACEAMI:
572  {
573  // nothing to do here
574  return;
575  }
576 
577  default:
578  {
580  << "problem." << abort(FatalError);
581  }
582  }
583 }
584 
585 
587 (
588  const sampleMode mode,
589  const label myWorld,
590  const pointField& samples,
591  const labelList& wantedWorlds,
592  const labelList& origProcs,
593 
594  labelList& sampleProcs,
595  labelList& sampleIndices,
596  pointField& sampleLocations
597 ) const
598 {
599  DebugInFunction << nl;
600 
601  // Find the processor/cell containing the samples. Does not account
602  // for samples being found in two processors.
603 
604  const label myComm = getCommunicator(); // Get or create
605  const label myRank = Pstream::myProcNo(myComm);
606  const label nProcs = Pstream::nProcs(myComm);
607 
608  const label oldWarnComm(Pstream::warnComm);
609  Pstream::warnComm = myComm;
610 
611  wordList samplePatches(nProcs);
612  {
613  samplePatches[myRank] = samplePatch_;
614  Pstream::gatherList(samplePatches, Pstream::msgType(), myComm);
615  Pstream::scatterList(samplePatches, Pstream::msgType(), myComm);
616  }
617  wordList sampleRegions(nProcs);
618  {
619  sampleRegions[myRank] = sampleRegion_;
620  Pstream::gatherList(sampleRegions, Pstream::msgType(), myComm);
621  Pstream::scatterList(sampleRegions, Pstream::msgType(), myComm);
622  }
623 
624 
625  // Find all the info for nearest
626  List<nearInfoWorld> nearest(samples.size());
627  forAll(nearest, samplei)
628  {
629  nearest[samplei].first() = nearInfo
630  (
631  pointIndexHit(),
632  Tuple2<scalar, label>(Foam::sqr(GREAT), -1)
633  );
634  nearest[samplei].second() = wantedWorlds[samplei];
635  }
636 
637 
638  // Extract samples to search for locally
639  {
640  DynamicList<label> localMap(samples.size());
641  forAll(wantedWorlds, samplei)
642  {
643  if (wantedWorlds[samplei] == myWorld)
644  {
645  localMap.append(samplei);
646  }
647  }
648 
649  if (localMap.size())
650  {
651  pointField localSamples(samples, localMap);
652  labelList localOrigProcs(origProcs, localMap);
653 
654  //Assume single patch to sample for now
655  const word localOrigPatch(samplePatches[localOrigProcs[0]]);
656  const word localOrigRegion(sampleRegions[localOrigProcs[0]]);
657  List<nearInfoWorld> localNearest(localSamples.size());
658 
659  if (debug)
660  {
661  Pout<< "Searching locally for " << localSamples.size()
662  << " samples on region:" << localOrigRegion
663  << " on patch:" << localOrigPatch << endl;
664  }
665  findLocalSamples
666  (
667  mode,
668  myWorld,
669  localOrigRegion,
670  localOrigPatch,
671  localSamples,
672  localNearest
673  );
674  UIndirectList<nearInfoWorld>(nearest, localMap) = localNearest;
675  }
676  }
677 
678 
679  // Find nearest. Combine on master.
681  (
682  nearest,
685  myComm
686  );
687  Pstream::listCombineScatter(nearest, Pstream::msgType(), myComm);
688 
689  //if (debug)
690  //{
691  // Pout<< "** After combining:" << endl;
692  // forAll(nearest, samplei)
693  // {
694  // Pout<< " sample:" << samples[samplei]
695  // << " origating from proc:" << origProcs[samplei]
696  // << " to be found on world:" << wantedWorlds[samplei] << nl
697  // << " found on world:" << nearest[samplei].second() << nl
698  // << " found on proc:"
699  // << nearest[samplei].first().second().second() << nl
700  // << " found on patchfacei:"
701  // << nearest[samplei].first().first().index() << nl
702  // << " found at location:"
703  // << nearest[samplei].first().first().rawPoint() << nl;
704  // }
705  // Pout<< endl;
706  //}
707 
708  // Convert back into proc+local index
709  sampleProcs.setSize(samples.size());
710  sampleIndices.setSize(samples.size());
711  sampleLocations.setSize(samples.size());
712 
713  forAll(nearest, sampleI)
714  {
715  const nearInfo& ni = nearest[sampleI].first();
716 
717  if (!ni.first().hit())
718  {
719  sampleProcs[sampleI] = -1;
720  sampleIndices[sampleI] = -1;
721  sampleLocations[sampleI] = vector::max;
722  }
723  else
724  {
725  sampleProcs[sampleI] = ni.second().second();
726  sampleIndices[sampleI] = ni.first().index();
727  sampleLocations[sampleI] = ni.first().hitPoint();
728  }
729  }
730 
731  Pstream::warnComm = oldWarnComm;
732 }
733 
734 
736 {
737  static bool hasWarned = false;
738  if (mapPtr_)
739  {
741  << "Mapping already calculated" << exit(FatalError);
742  }
743 
744  DebugInFunction << nl;
745 
746  const label myComm = getCommunicator(); // Get or create
747 
749  //if (sampleDatabase() && !sameWorld())
750  //{
751  // const word syncName("syncObjects");
752  // const polyMesh& mesh = patch_.boundaryMesh().mesh();
753  // const Time& runTime = mesh.time();
754  // const functionObjectList& fos = runTime.functionObjects();
755  //
756  // forAll(fos, i)
757  // {
758  // Pout<< "** FO:" << fos[i].name() << " tpye:" << fos[i].type()
759  // << endl;
760  // }
761  //
762  // if (fos.findObjectID(syncName) == -1)
763  // {
764  // //FatalErrorInFunction
765  // WarningInFunction
766  // << "Did not detect functionObject " << syncName
767  // << ". This is used to synchronise data inbetween worlds"
768  // << endl;
769  // }
770  //}
771 
772 
773  // Get points on face (since cannot use face-centres - might be off
774  // face-diagonal decomposed tets.
775  tmp<pointField> patchPoints(facePoints(patch_));
776 
777  // Get offsetted points
778  const pointField offsettedPoints(samplePoints(patchPoints()));
779 
780  // Do a sanity check - am I sampling my own patch?
781  // This only makes sense for a non-zero offset.
782  bool sampleMyself =
783  (
784  mode_ == NEARESTPATCHFACE
785  && sampleWorld() == UPstream::myWorld()
786  && sampleRegion() == patch_.boundaryMesh().mesh().name()
787  && samplePatch() == patch_.name()
788  );
789 
790  if (sampleMyself)
791  {
792  // Check offset
793  vectorField d(offsettedPoints-patchPoints());
794  bool coincident = (gAverage(mag(d)) <= ROOTVSMALL);
795 
796  if (sampleMyself && coincident)
797  {
799  << "Invalid offset " << d << endl
800  << "Offset is the vector added to the patch face centres to"
801  << " find the patch face supplying the data." << endl
802  << "Setting it to " << d
803  << " on the same patch, on the same region"
804  << " will find the faces themselves which does not make sense"
805  << " for anything but testing." << endl
806  << "patch:" << patch_.name() << endl
807  << "sampleRegion:" << sampleRegion() << endl
808  << "mode:" << sampleModeNames_[mode_] << endl
809  << "samplePatch:" << samplePatch() << endl
810  << "offsetMode:" << offsetModeNames_[offsetMode_] << endl;
811  }
812  }
813 
814  // Get local world and world-to-sample in index form
815  const label myWorld = Pstream::myWorldID();
816  const label mySampleWorld = Pstream::allWorlds().find(sampleWorld_);
817 
818 
819  // Get global list of all samples and the processor and face they come from.
820  pointField samples; // coordinates
821  labelList patchFaceWorlds; // world to sample
822  labelList patchFaceProcs; // originating processor
823  labelList patchFaces; // originating face on processor patch
824  pointField patchFc; // originating face centre
825  collectSamples
826  (
827  mySampleWorld, // world I want to sample
828  patchPoints,
829 
830  samples,
831  patchFaceWorlds,
832  patchFaceProcs,
833  patchFaces,
834  patchFc
835  );
836 
837  //if (debug)
838  //{
839  // forAll(samples, samplei)
840  // {
841  // Pout<< " sample:" << samples[samplei]
842  // << " origating from proc:" << patchFaceProcs[samplei]
843  // << " face:" << patchFaces[samplei]
844  // << " to be found on world:" << patchFaceWorlds[samplei] << nl;
845  // }
846  //}
847 
848  // Find processor and cell/face samples are in and actual location.
849  labelList sampleProcs;
850  labelList sampleIndices;
851  pointField sampleLocations;
852  findSamples
853  (
854  mode_,
855  myWorld, // my world (in index form)
856  samples,
857  patchFaceWorlds,
858  patchFaceProcs,
859 
860  sampleProcs,
861  sampleIndices,
862  sampleLocations
863  );
864 
865  // Check for samples that were not found. This will only happen for
866  // NEARESTCELL since finds cell containing a location
867  if (mode_ == NEARESTCELL)
868  {
869  label nNotFound = 0;
870  forAll(sampleProcs, sampleI)
871  {
872  if (sampleProcs[sampleI] == -1)
873  {
874  nNotFound++;
875  }
876  }
877  reduce(nNotFound, sumOp<label>(), Pstream::msgType(), myComm);
878 
879  if (nNotFound > 0)
880  {
881  if (!hasWarned)
882  {
884  << "Did not find " << nNotFound
885  << " out of " << sampleProcs.size() << " total samples."
886  << " Sampling these on owner cell centre instead." << endl
887  << "On patch " << patch_.name()
888  << " on region " << sampleRegion()
889  << " in mode " << sampleModeNames_[mode_] << endl
890  << "with offset mode " << offsetModeNames_[offsetMode_]
891  << ". Suppressing further warnings from " << type() << endl;
892 
893  hasWarned = true;
894  }
895 
896  // Collect the samples that cannot be found
897  DynamicList<label> subMap;
898  forAll(sampleProcs, sampleI)
899  {
900  if (sampleProcs[sampleI] == -1)
901  {
902  subMap.append(sampleI);
903  }
904  }
905 
906  // And re-search for pure nearest (should not fail)
907  labelList subSampleProcs;
908  labelList subSampleIndices;
909  pointField subSampleLocations;
910  findSamples
911  (
912  NEARESTONLYCELL,
913  myWorld, // my world (in index form)
914 
915  pointField(samples, subMap),
916  UIndirectList<label>(patchFaceWorlds, subMap)(),
917  UIndirectList<label>(patchFaceProcs, subMap)(),
918 
919  subSampleProcs,
920  subSampleIndices,
921  subSampleLocations
922  );
923 
924  // Insert
925  labelUIndList(sampleProcs, subMap) = subSampleProcs;
926  labelUIndList(sampleIndices, subMap) = subSampleIndices;
927  UIndirectList<point>(sampleLocations, subMap) = subSampleLocations;
928  }
929  }
930 
931  // Now we have all the data we need:
932  // - where sample originates from (so destination when mapping):
933  // patchFaces, patchFaceProcs.
934  // - cell/face sample is in (so source when mapping)
935  // sampleIndices, sampleProcs.
936 
937  if (Pstream::master(myComm))
938  {
939  forAll(samples, i)
940  {
941  if (sampleProcs[i] == -1)
942  {
943  FatalErrorInFunction << "did not find sample "
944  << patchFc[i] << " on patch " << patch_.name()
945  << " on region "
946  << patch_.boundaryMesh().mesh().name()
947  << " on processor " << patchFaceProcs[i]
948  << exit(FatalError);
949  }
950  }
951  }
952 
953 
954  if (debug && Pstream::master(myComm))
955  {
956  //forAll(samples, i)
957  //{
958  // Pout<< i << " need data in region "
959  // << patch_.boundaryMesh().mesh().name()
960  // << " for proc:" << patchFaceProcs[i]
961  // << " face:" << patchFaces[i]
962  // << " at:" << patchFc[i] << endl
963  // << "Found data in region " << sampleRegion()
964  // << " at proc:" << sampleProcs[i]
965  // << " face:" << sampleIndices[i]
966  // << " at:" << sampleLocations[i]
967  // << nl << endl;
968  //}
969 
970  OFstream str
971  (
972  patch_.boundaryMesh().mesh().time().path()
973  / patch_.name()
974  + "_mapped.obj"
975  );
976  Pout<< "Dumping mapping as lines from patch faceCentres to"
977  << " sampled cell/faceCentres/points to file " << str.name()
978  << endl;
979 
980  label vertI = 0;
981 
982  forAll(patchFc, i)
983  {
984  meshTools::writeOBJ(str, patchFc[i]);
985  vertI++;
986  meshTools::writeOBJ(str, sampleLocations[i]);
987  vertI++;
988  str << "l " << vertI-1 << ' ' << vertI << nl;
989  }
990  }
991 
992  // Determine schedule.
993  mapPtr_.reset(new mapDistribute(sampleProcs, patchFaceProcs, myComm));
994 
995  // Rework the schedule from indices into samples to cell data to send,
996  // face data to receive.
997 
998  labelListList& subMap = mapPtr_().subMap();
999  labelListList& constructMap = mapPtr_().constructMap();
1000 
1001  forAll(subMap, proci)
1002  {
1003  subMap[proci] = labelUIndList(sampleIndices, subMap[proci]);
1004  constructMap[proci] = labelUIndList(patchFaces, constructMap[proci]);
1005 
1006  if (debug)
1007  {
1008  Pout<< "To proc:" << proci << " sending values of cells/faces:"
1009  << subMap[proci] << endl;
1010  Pout<< "From proc:" << proci
1011  << " receiving values of patch faces:"
1012  << constructMap[proci] << endl;
1013  }
1014  }
1015 
1016 
1017  // Redo constructSize
1018  mapPtr_().constructSize() = patch_.size();
1019 
1020  if (debug)
1021  {
1022  // Check that all elements get a value.
1023  bitSet used(patch_.size());
1024  forAll(constructMap, proci)
1025  {
1026  const labelList& map = constructMap[proci];
1027 
1028  forAll(map, i)
1029  {
1030  label facei = map[i];
1031 
1032  if (used.test(facei))
1033  {
1035  << "On patch " << patch_.name()
1036  << " patchface " << facei
1037  << " is assigned to more than once."
1038  << abort(FatalError);
1039  }
1040  else
1041  {
1042  used.set(facei);
1043  }
1044  }
1045  }
1046  forAll(used, facei)
1047  {
1048  if (!used.test(facei))
1049  {
1051  << "On patch " << patch_.name()
1052  << " patchface " << facei
1053  << " is never assigned to."
1054  << abort(FatalError);
1055  }
1056  }
1057  }
1058 }
1059 
1060 
1062 const
1063 {
1064  const word surfType(surfDict_.getOrDefault<word>("type", "none"));
1065 
1066  if (!surfPtr_ && surfType != "none")
1067  {
1068  word surfName(surfDict_.getOrDefault("name", patch_.name()));
1069 
1070  const polyMesh& mesh = patch_.boundaryMesh().mesh();
1071 
1072  surfPtr_ =
1074  (
1075  surfType,
1076  IOobject
1077  (
1078  surfName,
1079  mesh.time().constant(),
1080  "triSurface",
1081  mesh,
1084  ),
1085  surfDict_
1086  );
1087  }
1088 
1089  return surfPtr_;
1090 }
1091 
1092 
1094 {
1095  if (AMIPtr_->upToDate())
1096  {
1098  << "AMI already up-to-date"
1099  << endl;
1100 
1101  return;
1102  }
1103 
1104  DebugInFunction << nl;
1105 
1106  const label myComm = getCommunicator(); // Get or create
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, surfPtr());
1150 
1151  Pstream::warnComm = oldWarnComm;
1152  Pstream::worldComm = oldWorldComm;
1153  }
1154  else
1155  {
1156  faceList dummyFaces;
1157  pointField dummyPoints;
1158  const primitivePatch dummyPatch
1159  (
1161  (
1162  dummyFaces
1163  ),
1164  dummyPoints
1165  );
1166 
1167  // Change to use inter-world communicator
1168  Pstream::worldComm = myComm;
1170 
1171  if (masterWorld())
1172  {
1173  // Construct/apply AMI interpolation to determine addressing
1174  // and weights. Have patch_ for src faces, 0 faces for the
1175  // target side
1176  AMIPtr_->calculate(patch_, dummyPatch, surfPtr());
1177  }
1178  else
1179  {
1180  // Construct/apply AMI interpolation to determine addressing
1181  // and weights. Have 0 faces for src side, patch_ for the tgt
1182  // side
1183  AMIPtr_->calculate(dummyPatch, patch_, surfPtr());
1184  }
1185  // Now the AMI addressing/weights will be from src side (on masterWorld
1186  // processors) to tgt side (on other processors)
1187 
1188  Pstream::warnComm = oldWarnComm;
1189  Pstream::worldComm = oldWorldComm;
1190  }
1191 }
1192 
1193 
1196  const objectRegistry& obr,
1197  const wordList& names,
1198  const label index
1199 )
1200 {
1201  const objectRegistry& sub = obr.subRegistry(names[index], true);
1202  if (index == names.size()-1)
1203  {
1204  return sub;
1205  }
1206  else
1207  {
1208  return subRegistry(sub, names, index+1);
1209  }
1210 }
1211 
1212 
1213 // * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * * * * * //
1214 
1216 :
1217  patch_(pp),
1218  sampleWorld_(),
1219  sampleRegion_(patch_.boundaryMesh().mesh().name()),
1220  mode_(NEARESTPATCHFACE),
1221  samplePatch_(),
1222  coupleGroup_(),
1223  sampleDatabasePtr_(),
1224  offsetMode_(UNIFORM),
1225  offset_(Zero),
1226  offsets_(pp.size(), offset_),
1227  distance_(0),
1228  communicator_(-1), // Demand-driven (cached value)
1229  sameRegion_(true),
1230  mapPtr_(nullptr),
1231  AMIReverse_(false),
1232  AMIPtr_(new faceAreaWeightAMI(true, AMIReverse_)),
1233  surfPtr_(nullptr),
1234  surfDict_(fileName("surface"))
1235 {
1236  // NOTE: same region, no sample-world. Thus no world-world communication
1237 }
1238 
1239 
1242  const polyPatch& pp,
1243  const word& sampleRegion,
1244  const sampleMode mode,
1245  const word& samplePatch,
1246  const vectorField& offsets
1247 )
1248 :
1250  (
1251  pp,
1252  sampleRegion,
1253  mode,
1254  samplePatch,
1255  scalar(0)
1256  )
1257 {
1258  mappedPatchBase::setOffset(offsets);
1259 }
1260 
1261 
1264  const polyPatch& pp,
1265  const word& sampleRegion,
1266  const sampleMode mode,
1267  const word& samplePatch,
1268  const vector& uniformOffset
1269 )
1270 :
1272  (
1273  pp,
1274  sampleRegion,
1275  mode,
1276  samplePatch,
1277  scalar(0)
1278  )
1279 {
1280  mappedPatchBase::setOffset(uniformOffset);
1281 }
1282 
1283 
1286  const polyPatch& pp,
1287  const word& sampleRegion,
1288  const sampleMode mode,
1289  const word& samplePatch,
1290  const scalar normalDistance
1291 )
1292 :
1293  patch_(pp),
1294  sampleWorld_(),
1295  sampleRegion_(sampleRegion),
1296  mode_(mode),
1297  samplePatch_(samplePatch),
1298  coupleGroup_(),
1299  sampleDatabasePtr_(),
1300  offsetMode_(NORMAL),
1301  offset_(Zero),
1302  offsets_(0),
1303  distance_(normalDistance),
1304  communicator_(-1), // Demand-driven (cached value)
1305  sameRegion_
1306  (
1307  sampleWorld_.empty()
1308  && sampleRegion_ == patch_.boundaryMesh().mesh().name()
1309  ),
1310  mapPtr_(nullptr),
1311  AMIReverse_(false),
1312  AMIPtr_(new faceAreaWeightAMI(true, AMIReverse_)),
1313  surfPtr_(nullptr),
1314  surfDict_(fileName("surface"))
1315 {
1316  addWorldConnection();
1317 }
1318 
1319 
1322  const polyPatch& pp,
1323  const dictionary& dict
1324 )
1325 :
1326  patch_(pp),
1327  sampleWorld_(dict.getOrDefault<word>("sampleWorld", word::null)),
1328  sampleRegion_(dict.getOrDefault<word>("sampleRegion", word::null)),
1329  mode_(sampleModeNames_.get("sampleMode", dict)),
1330  samplePatch_(dict.getOrDefault<word>("samplePatch", word::null)),
1331  coupleGroup_(dict),
1332  sampleDatabasePtr_(readDatabase(dict)),
1333  offsetMode_(UNIFORM),
1334  offset_(Zero),
1335  offsets_(0),
1336  distance_(0),
1337  communicator_(-1), // Demand-driven (cached value)
1338  sameRegion_
1339  (
1340  sampleWorld_.empty()
1341  && sampleRegion_ == patch_.boundaryMesh().mesh().name()
1342  ),
1343  mapPtr_(nullptr),
1344  AMIReverse_
1345  (
1347  (
1348  "reverseTarget", // AMIInterpolation uses this keyword
1350  (
1351  "flipNormals",
1352  false
1353  )
1354  )
1355  ),
1356  AMIPtr_
1357  (
1359  (
1360  dict.getOrDefault("AMIMethod", faceAreaWeightAMI::typeName),
1361  dict,
1362  AMIReverse_
1363  )
1364  ),
1365  surfPtr_(nullptr),
1366  surfDict_(dict.subOrEmptyDict("surface"))
1367 {
1368  addWorldConnection();
1369 
1370  if (!coupleGroup_.valid())
1371  {
1372  if (sampleWorld_.empty() && sampleRegion_.empty())
1373  {
1374  // If no coupleGroup and no sampleRegion assume local region
1375  sampleRegion_ = patch_.boundaryMesh().mesh().name();
1376  sameRegion_ = true;
1377  }
1378  }
1379 
1380  if (offsetModeNames_.readIfPresent("offsetMode", dict, offsetMode_))
1381  {
1382  switch (offsetMode_)
1383  {
1384  case UNIFORM:
1385  {
1386  dict.readEntry("offset", offset_);
1387  }
1388  break;
1389 
1390  case NONUNIFORM:
1391  {
1392  offsets_ = pointField("offsets", dict, patch_.size());
1393  }
1394  break;
1395 
1396  case NORMAL:
1397  {
1398  dict.readEntry("distance", distance_);
1399  }
1400  break;
1401  }
1402  }
1403  else if (dict.readIfPresent("offset", offset_))
1404  {
1405  offsetMode_ = UNIFORM;
1406  }
1407  else if (dict.found("offsets"))
1408  {
1409  offsetMode_ = NONUNIFORM;
1410  offsets_ = pointField("offsets", dict, patch_.size());
1411  }
1412  else if (mode_ != NEARESTPATCHFACE && mode_ != NEARESTPATCHFACEAMI)
1413  {
1415  << "Please supply the offsetMode as one of "
1416  << offsetModeNames_
1417  << exit(FatalIOError);
1418  }
1419 }
1420 
1421 
1424  const polyPatch& pp,
1425  const sampleMode mode,
1426  const dictionary& dict
1427 )
1428 :
1429  patch_(pp),
1430  sampleWorld_(dict.getOrDefault<word>("sampleWorld", word::null)),
1431  sampleRegion_(dict.getOrDefault<word>("sampleRegion", word::null)),
1432  mode_(mode),
1433  samplePatch_(dict.getOrDefault<word>("samplePatch", word::null)),
1434  coupleGroup_(dict), //dict.getOrDefault<word>("coupleGroup", word::null)),
1435  sampleDatabasePtr_(readDatabase(dict)),
1436  offsetMode_(UNIFORM),
1437  offset_(Zero),
1438  offsets_(0),
1439  distance_(0),
1440  communicator_(-1), // Demand-driven (cached value)
1441  sameRegion_
1442  (
1443  sampleWorld_.empty()
1444  && sampleRegion_ == patch_.boundaryMesh().mesh().name()
1445  ),
1446  mapPtr_(nullptr),
1447  AMIReverse_
1448  (
1450  (
1451  "reverseTarget", // AMIInterpolation uses this keyword
1453  (
1454  "flipNormals",
1455  false
1456  )
1457  )
1458  ),
1459  AMIPtr_
1460  (
1462  (
1463  dict.getOrDefault("AMIMethod", faceAreaWeightAMI::typeName),
1464  dict,
1465  AMIReverse_
1466  )
1467  ),
1468  surfPtr_(nullptr),
1469  surfDict_(dict.subOrEmptyDict("surface"))
1470 {
1471  addWorldConnection();
1472 
1473  if (mode != NEARESTPATCHFACE && mode != NEARESTPATCHFACEAMI)
1474  {
1476  << "Construct from sampleMode and dictionary only applicable for "
1477  << " collocated patches in modes "
1478  << sampleModeNames_[NEARESTPATCHFACE] << ','
1479  << sampleModeNames_[NEARESTPATCHFACEAMI]
1480  << exit(FatalIOError);
1481  }
1482 
1483  if (!coupleGroup_.valid())
1484  {
1485  if (sampleWorld_.empty() && sampleRegion_.empty())
1486  {
1487  // If no coupleGroup and no sampleRegion assume local region
1488  sampleRegion_ = patch_.boundaryMesh().mesh().name();
1489  sameRegion_ = true;
1490  }
1491  }
1492 }
1493 
1494 
1497  const polyPatch& pp,
1498  const mappedPatchBase& mpb
1499 )
1500 :
1501  patch_(pp),
1502  sampleWorld_(mpb.sampleWorld_),
1503  sampleRegion_(mpb.sampleRegion_),
1504  mode_(mpb.mode_),
1505  samplePatch_(mpb.samplePatch_),
1506  coupleGroup_(mpb.coupleGroup_),
1507  sampleDatabasePtr_
1508  (
1509  mpb.sampleDatabasePtr_
1510  ? new fileName(mpb.sampleDatabasePtr_())
1511  : nullptr
1512  ),
1513  offsetMode_(mpb.offsetMode_),
1514  offset_(mpb.offset_),
1515  offsets_(mpb.offsets_),
1516  distance_(mpb.distance_),
1517  communicator_(mpb.communicator_),
1518  sameRegion_(mpb.sameRegion_),
1519  mapPtr_(nullptr),
1520  AMIReverse_(mpb.AMIReverse_),
1521  AMIPtr_(mpb.AMIPtr_->clone()),
1522  surfPtr_(nullptr),
1523  surfDict_(mpb.surfDict_)
1524 {}
1525 
1526 
1529  const polyPatch& pp,
1530  const mappedPatchBase& mpb,
1531  const labelUList& mapAddressing
1532 )
1533 :
1534  patch_(pp),
1535  sampleWorld_(mpb.sampleWorld_),
1536  sampleRegion_(mpb.sampleRegion_),
1537  mode_(mpb.mode_),
1538  samplePatch_(mpb.samplePatch_),
1539  coupleGroup_(mpb.coupleGroup_),
1540  sampleDatabasePtr_
1541  (
1542  mpb.sampleDatabasePtr_
1543  ? new fileName(mpb.sampleDatabasePtr_())
1544  : nullptr
1545  ),
1546  offsetMode_(mpb.offsetMode_),
1547  offset_(mpb.offset_),
1548  offsets_
1549  (
1550  offsetMode_ == NONUNIFORM
1551  ? vectorField(mpb.offsets_, mapAddressing)
1552  : vectorField()
1553  ),
1554  distance_(mpb.distance_),
1555  communicator_(mpb.communicator_),
1556  sameRegion_(mpb.sameRegion_),
1557  mapPtr_(nullptr),
1558  AMIReverse_(mpb.AMIReverse_),
1559  AMIPtr_(mpb.AMIPtr_->clone()),
1560  surfPtr_(nullptr),
1561  surfDict_(mpb.surfDict_)
1562 {}
1563 
1564 
1565 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
1566 
1568 {
1569  clearOut();
1570 }
1571 
1572 
1574 {
1575  mapPtr_.reset(nullptr);
1576  surfPtr_.reset(nullptr);
1577  AMIPtr_->upToDate() = false;
1578 }
1579 
1580 
1581 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
1582 
1583 void Foam::mappedPatchBase::setOffset(const scalar normalDist)
1584 {
1585  clearOut();
1586  offsetMode_ = offsetMode::NORMAL;
1587  offset_ = Zero;
1588  offsets_.clear();
1589  distance_ = normalDist;
1590 }
1591 
1592 
1593 void Foam::mappedPatchBase::setOffset(const vector& uniformOffset)
1594 {
1595  clearOut();
1596  offsetMode_ = offsetMode::UNIFORM;
1597  offset_ = uniformOffset;
1598  offsets_.clear();
1599  distance_ = Zero;
1600 }
1601 
1602 
1604 {
1605  clearOut();
1606  offsetMode_ = offsetMode::NONUNIFORM;
1607  offset_ = Zero;
1608  offsets_ = offsets;
1609  distance_ = Zero;
1610 }
1611 
1612 
1615  const word& sampleRegion
1616 ) const
1617 {
1618  const polyMesh& thisMesh = patch_.boundaryMesh().mesh();
1619  return
1620  (
1621  sampleRegion.empty() || sampleRegion == thisMesh.name()
1622  ? thisMesh
1623  : thisMesh.time().lookupObject<polyMesh>(sampleRegion)
1624  );
1625 }
1626 
1627 
1630  const word& sampleRegion,
1631  const word& samplePatch
1632 ) const
1633 {
1634  const polyMesh& nbrMesh = lookupMesh(sampleRegion);
1635 
1636  const label patchi = nbrMesh.boundaryMesh().findPatchID(samplePatch);
1637 
1638  if (patchi == -1)
1639  {
1641  << "Cannot find patch " << samplePatch
1642  << " in region " << sampleRegion_ << endl
1643  << exit(FatalError);
1644  }
1645  return nbrMesh.boundaryMesh()[patchi];
1646 }
1647 
1648 
1650 {
1651  if (UPstream::myWorld() != sampleWorld_)
1652  {
1654  << "sampleWorld : " << sampleWorld_
1655  << " is not the current world : " << UPstream::myWorld()
1656  << exit(FatalError);
1657  }
1658  return lookupMesh(sampleRegion());
1659 }
1660 
1661 
1663 {
1664  const polyMesh& nbrMesh = sampleMesh();
1665 
1666  const label patchi = nbrMesh.boundaryMesh().findPatchID(samplePatch());
1667 
1668  if (patchi == -1)
1669  {
1671  << "Cannot find patch " << samplePatch()
1672  << " in region " << sampleRegion_ << endl
1673  << "Valid patches are " << nbrMesh.boundaryMesh().names()
1674  << exit(FatalError);
1675  }
1676 
1677  return nbrMesh.boundaryMesh()[patchi];
1678 }
1679 
1680 
1683  const pointField& fc
1684 ) const
1685 {
1686  auto tfld = tmp<pointField>::New(fc);
1687  auto& fld = tfld.ref();
1688 
1689  switch (offsetMode_)
1690  {
1691  case UNIFORM:
1692  {
1693  fld += offset_;
1694  break;
1695  }
1696  case NONUNIFORM:
1697  {
1698  fld += offsets_;
1699  break;
1700  }
1701  case NORMAL:
1702  {
1703  // Get outwards pointing normal
1704  vectorField n(patch_.faceAreas());
1705  n /= mag(n);
1706 
1707  fld += distance_*n;
1708  break;
1709  }
1710  }
1711 
1712  return tfld;
1713 }
1714 
1715 
1717 {
1718  return samplePoints(facePoints(patch_));
1719 }
1720 
1721 
1724  const polyMesh& mesh,
1725  const label facei,
1726  const polyMesh::cellDecomposition decompMode
1727 )
1728 {
1729  const point& fc = mesh.faceCentres()[facei];
1730 
1731  switch (decompMode)
1732  {
1733  case polyMesh::FACE_PLANES:
1735  {
1736  // For both decompositions the face centre is guaranteed to be
1737  // on the face
1738  return pointIndexHit(true, fc, facei);
1739  }
1740  break;
1741 
1743  case polyMesh::CELL_TETS:
1744  {
1745  // Find the intersection of a ray from face centre to cell centre
1746  // Find intersection of (face-centre-decomposition) centre to
1747  // cell-centre with face-diagonal-decomposition triangles.
1748 
1749  const pointField& p = mesh.points();
1750  const face& f = mesh.faces()[facei];
1751 
1752  if (f.size() <= 3)
1753  {
1754  // Return centre of triangle.
1755  return pointIndexHit(true, fc, 0);
1756  }
1757 
1758  const label celli = mesh.faceOwner()[facei];
1759  const point& cc = mesh.cellCentres()[celli];
1760  vector d = fc-cc;
1761 
1762  const label fp0 = mesh.tetBasePtIs()[facei];
1763  const point& basePoint = p[f[fp0]];
1764 
1765  label fp = f.fcIndex(fp0);
1766  for (label i = 2; i < f.size(); i++)
1767  {
1768  const point& thisPoint = p[f[fp]];
1769  label nextFp = f.fcIndex(fp);
1770  const point& nextPoint = p[f[nextFp]];
1771 
1772  const triPointRef tri(basePoint, thisPoint, nextPoint);
1773  pointHit hitInfo = tri.intersection
1774  (
1775  cc,
1776  d,
1778  );
1779 
1780  if (hitInfo.hit() && hitInfo.distance() > 0)
1781  {
1782  return pointIndexHit(true, hitInfo.hitPoint(), i-2);
1783  }
1784 
1785  fp = nextFp;
1786  }
1787 
1788  // Fall-back
1789  return pointIndexHit(false, fc, -1);
1790  }
1791  break;
1792 
1793  default:
1794  {
1796  << "problem" << abort(FatalError);
1797  return pointIndexHit();
1798  }
1799  }
1800 }
1801 
1802 
1805  const objectRegistry& obr,
1806  const fileName& path
1807 )
1808 {
1809  // Lookup (and create if non-existing) a registry using
1810  // '/' separated path. Like 'mkdir -p'
1811 
1812  fileName cleanedPath(path);
1813  cleanedPath.clean(); // Remove unneeded ".."
1814  const wordList names(cleanedPath.components());
1815 
1816  if (names.empty())
1817  {
1818  return obr;
1819  }
1820  else
1821  {
1822  return subRegistry(obr, names, 0);
1823  }
1824 }
1825 
1826 
1829  const fileName& root,
1830  const label proci
1831 )
1832 {
1833  const word processorName("processor"+Foam::name(proci));
1834  return root/"send"/processorName;
1835 }
1836 
1837 
1839 {
1840  return sendPath(sampleDatabasePath(), proci);
1841 }
1842 
1843 
1846  const fileName& root,
1847  const label proci
1848 )
1849 {
1850  const word processorName("processor"+Foam::name(proci));
1851  return root/"receive"/processorName;
1852 }
1853 
1854 
1856 {
1857  return receivePath(sampleDatabasePath(), proci);
1858 }
1859 
1860 
1863  const objectRegistry& obr,
1864  dictionary& dict
1865 )
1866 {
1867  forAllIters(obr, iter)
1868  {
1869  regIOobject* objPtr = iter.val();
1870  const regIOobject& obj = *objPtr;
1871 
1872  if (isA<objectRegistry>(obj))
1873  {
1874  dictionary& d = dict.subDictOrAdd(obj.name());
1875  writeDict(dynamic_cast<const objectRegistry&>(obj), d);
1876  }
1877  else if
1878  (
1879  writeIOField<scalar>(obj, dict)
1880  || writeIOField<vector>(obj, dict)
1881  || writeIOField<sphericalTensor>(obj, dict)
1882  || writeIOField<symmTensor>(obj, dict)
1883  || writeIOField<tensor>(obj, dict)
1884  )
1885  {
1886  // IOField specialisation
1887  }
1888  else
1889  {
1890  // Not tested. No way of retrieving data anyway ...
1891  OTstream os;
1892  obj.writeData(os);
1893 
1894  primitiveEntry* pePtr = new primitiveEntry(obj.name(), os.tokens());
1895  dict.add(pePtr);
1896  }
1897  }
1898 }
1899 
1900 
1902 {
1903  // Reverse of writeDict - reads fields from dictionary into objectRegistry
1904  for (const entry& e : d)
1905  {
1906  if (e.isDict())
1907  {
1908  // Add sub registry
1909  objectRegistry& sub = const_cast<objectRegistry&>
1910  (
1911  obr.subRegistry(e.keyword(), true)
1912  );
1913 
1914  readDict(e.dict(), sub);
1915  }
1916  else
1917  {
1918  ITstream& is = e.stream();
1919  token tok(is);
1920 
1921  if
1922  (
1923  constructIOField<scalar>(e.keyword(), tok, is, obr)
1924  || constructIOField<vector>(e.keyword(), tok, is, obr)
1925  || constructIOField<sphericalTensor>(e.keyword(), tok, is, obr)
1926  || constructIOField<symmTensor>(e.keyword(), tok, is, obr)
1927  || constructIOField<tensor>(e.keyword(), tok, is, obr)
1928  )
1929  {
1930  // Done storing field
1931  }
1932  else
1933  {
1934  FatalErrorInFunction << "Unsupported type " << e.keyword()
1935  << exit(FatalError);
1936  }
1937  }
1938  }
1939 }
1940 
1941 
1943 {
1944  os.writeEntry("sampleMode", sampleModeNames_[mode_]);
1945  os.writeEntryIfDifferent<word>("sampleWorld", word::null, sampleWorld_);
1946  os.writeEntryIfDifferent<word>("sampleRegion", word::null, sampleRegion_);
1947  os.writeEntryIfDifferent<word>("samplePatch", word::null, samplePatch_);
1948 
1949  if (sampleDatabasePtr_)
1950  {
1951  os.writeEntry("sampleDatabase", Switch::name(true));
1952  // Write database path if differing
1954  (
1955  "sampleDatabasePath",
1957  sampleDatabasePtr_()
1958  );
1959  }
1960  coupleGroup_.write(os);
1961 
1962  if
1963  (
1964  offsetMode_ == UNIFORM
1965  && offset_ == vector::zero
1966  && (mode_ == NEARESTPATCHFACE || mode_ == NEARESTPATCHFACEAMI)
1967  )
1968  {
1969  // Collocated mode. No need to write offset data
1970  }
1971  else
1972  {
1973  os.writeEntry("offsetMode", offsetModeNames_[offsetMode_]);
1974 
1975  switch (offsetMode_)
1976  {
1977  case UNIFORM:
1978  {
1979  os.writeEntry("offset", offset_);
1980  break;
1981  }
1982  case NONUNIFORM:
1983  {
1984  offsets_.writeEntry("offsets", os);
1985  break;
1986  }
1987  case NORMAL:
1988  {
1989  os.writeEntry("distance", distance_);
1990  break;
1991  }
1992  }
1993  }
1994 
1995  if (mode_ == NEARESTPATCHFACEAMI)
1996  {
1997  if (AMIPtr_)
1998  {
1999  // Use AMI to write itself. Problem: outputs:
2000  // - restartUncoveredSourceFace
2001  // - reverseTarget (instead of flipNormals)
2002  AMIPtr_->write(os);
2003  }
2004  if (!surfDict_.empty())
2005  {
2006  surfDict_.writeEntry(surfDict_.dictName(), os);
2007  }
2008  }
2009 }
2010 
2011 
2012 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::PrimitivePatch::points
const Field< point_type > & points() const noexcept
Return reference to global points.
Definition: PrimitivePatch.H:299
Foam::entry
A keyword and a list of tokens is an 'entry'.
Definition: entry.H:67
Foam::IOobject::NO_WRITE
Definition: IOobject.H:195
Foam::UPstream::warnComm
static label warnComm
Debugging: warn for use of any communicator differing from warnComm.
Definition: UPstream.H:296
Foam::pointField
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
runTime
engineTime & runTime
Definition: createEngineTime.H:13
Foam::Random
Random number generator.
Definition: Random.H:59
Foam::polyMesh::points
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1069
Foam::Ostream::writeEntryIfDifferent
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:248
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
meshTools.H
Foam::mappedPatchBase::clearOut
void clearOut()
Definition: mappedPatchBase.C:1573
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::Enum
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition: IOstreamOption.H:57
Foam::mappedPatchBase::communicator_
label communicator_
Communicator.
Definition: mappedPatchBase.H:259
Foam::polyMesh::FACE_DIAG_TRIS
Definition: polyMesh.H:107
Foam::mappedPatchBase::sendPath
static fileName sendPath(const fileName &root, const label proci)
Helper: return path to store data to be sent to processor i.
Definition: mappedPatchBase.C:1828
Foam::Time
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:73
Foam::mappedPatchBase::readDict
static void readDict(const dictionary &d, objectRegistry &obr)
(recursively) construct and register IOFields from dictionary
Definition: mappedPatchBase.C:1901
Foam::polyMesh::cellDecomposition
cellDecomposition
Enumeration defining the decomposition of the cell for.
Definition: polyMesh.H:100
Foam::accessOp
Definition: UList.H:668
Foam::mappedPatchBase::AMIReverse_
const bool AMIReverse_
Flag to indicate that slave patch should be reversed for AMI.
Definition: mappedPatchBase.H:278
Foam::primitiveEntry
A keyword and a list of tokens comprise a primitiveEntry. A primitiveEntry can be read,...
Definition: primitiveEntry.H:63
Foam::mappedPatchBase::subRegistry
static const objectRegistry & subRegistry(const objectRegistry &obr, const wordList &names, const label index)
Definition: mappedPatchBase.C:1195
Foam::fileName::components
wordList components(const char delim='/') const
Return path components as wordList.
Definition: fileName.C:461
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::VectorSpace< Vector< Cmpt >, Cmpt, 3 >::uniform
static Vector< Cmpt > uniform(const Cmpt &s)
Return a VectorSpace with all elements = s.
Definition: VectorSpaceI.H:164
Foam::fileName
A class for handling file names.
Definition: fileName.H:73
SubField.H
Foam::intersection::HALF_RAY
Definition: intersection.H:75
Foam::mappedPatchBase::setOffset
void setOffset(const scalar normalDist)
Change to normal offset with given distance.
Definition: mappedPatchBase.C:1583
Foam::bitSet
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition: bitSet.H:63
polyPatch.H
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::PointHit
Describes the interaction of a face and a point. It carries the info of a successful hit and (if succ...
Definition: PointHit.H:53
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
meshSearchMeshObject.H
Foam::DynamicList< label >
Foam::dictionary::found
bool found(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Search for an entry (const access) with the given keyword.
Definition: dictionaryI.H:87
Foam::gAverage
Type gAverage(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:604
Foam::meshTools::writeOBJ
void writeOBJ(Ostream &os, const point &pt)
Write obj representation of a point.
Definition: meshTools.C:203
ListListOps.H
Foam::Pstream::scatterList
static void scatterList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Scatter data. Reverse of gatherList.
Definition: gatherScatterList.C:215
Foam::regIOobject::writeData
virtual bool writeData(Ostream &) const =0
Pure virtual writeData function.
Foam::SubList
A List obtained as a section of another List.
Definition: SubList.H:54
Foam::mappedPatchBase::sampleWorld_
word sampleWorld_
World to sample.
Definition: mappedPatchBase.H:229
Foam::mappedPatchBase::samplePolyPatch
const polyPatch & samplePolyPatch() const
Get the patch on the region.
Definition: mappedPatchBase.C:1662
Foam::treeBoundBox
Standard boundBox with extra functionality for use in octree.
Definition: treeBoundBox.H:86
Foam::MeshObject< polyMesh, GeometricMeshObject, meshSearchMeshObject >::New
static const meshSearchMeshObject & New(const polyMesh &mesh, Args &&... args)
Get existing or create a new MeshObject.
Definition: MeshObject.C:48
indexedOctree.H
Foam::mappedPatchBase::mappedPatchBase
mappedPatchBase(const polyPatch &)
Construct from patch.
Definition: mappedPatchBase.C:1215
Foam::mappedPatchBase::lookupMesh
const polyMesh & lookupMesh(const word &region) const
Lookup mesh.
Definition: mappedPatchBase.C:1614
Foam::mappedPatchBase::offset_
vector offset_
Offset vector (uniform)
Definition: mappedPatchBase.H:250
Foam::mappedPatchBase::findLocalSamples
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.
Definition: mappedPatchBase.C:280
Foam::PointHit::distance
scalar distance() const noexcept
Return distance to hit.
Definition: PointHit.H:139
Foam::UPstream::master
static bool master(const label communicator=worldComm)
Am I the master process.
Definition: UPstream.H:457
Foam::mappedPatchBase::mode_
const sampleMode mode_
What to sample.
Definition: mappedPatchBase.H:235
Foam::HashSet::test
bool test(const Key &key) const noexcept
Same as found() - return true if key exists in the set.
Definition: HashSet.H:180
Foam::PointHit::hitPoint
const point_type & hitPoint() const
Return the hit point. Fatal if not hit.
Definition: PointHit.H:145
Foam::mappedPatchBase::offsetModeNames_
static const Enum< offsetMode > offsetModeNames_
Definition: mappedPatchBase.H:139
Foam::FatalIOError
IOerror FatalIOError
Foam::mappedPatchBase::writeDict
static void writeDict(const objectRegistry &obr, dictionary &dict)
Convert objectRegistry contents into dictionary.
Definition: mappedPatchBase.C:1862
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:444
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::token
A token holds an item read from Istream.
Definition: token.H:68
Foam::mappedPatchBase::sameWorld
bool sameWorld() const
Is sample world the local world?
Definition: mappedPatchBaseI.H:169
Foam::mappedPatchBase::findSamples
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.
Definition: mappedPatchBase.C:587
Foam::Pout
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
Foam::mappedPatchBase::samplePatch_
word samplePatch_
Patch (if in sampleMode NEARESTPATCH*)
Definition: mappedPatchBase.H:238
triPointRef.H
Foam::AMIInterpolation::New
static autoPtr< AMIInterpolation > New(const word &modelName, const dictionary &dict, const bool reverseTarget=false)
Selector for dictionary.
Definition: AMIInterpolationNew.C:33
Foam::mappedPatchBase
Determines a mapping between patch face centres and mesh cell or face centres and processors they're ...
Definition: mappedPatchBase.H:112
polyMesh.H
syncTools.H
Foam::dictionary::subDictOrAdd
dictionary & subDictOrAdd(const word &keyword, enum keyType::option matchOpt=keyType::REGEX)
Find and return a sub-dictionary for manipulation.
Definition: dictionary.C:500
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:77
Foam::mappedPatchBase::coupleGroup_
const coupleGroupIdentifier coupleGroup_
PatchGroup (if in sampleMode NEARESTPATCH*)
Definition: mappedPatchBase.H:241
Foam::polyBoundaryMesh::names
wordList names() const
Return a list of patch names.
Definition: polyBoundaryMesh.C:555
Foam::boundBox::min
const point & min() const
Minimum describing the bounding box.
Definition: boundBoxI.H:91
Foam::sumOp
Definition: ops.H:213
Foam::mappedPatchBase::~mappedPatchBase
virtual ~mappedPatchBase()
Destructor.
Definition: mappedPatchBase.C:1567
multiWorldConnectionsObject.H
Foam::treeDataPoint
Holds (reference to) pointField. Encapsulation of data needed for octree searches....
Definition: treeDataPoint.H:62
Foam::mappedPatchBase::sampleModeNames_
static const Enum< sampleMode > sampleModeNames_
Definition: mappedPatchBase.H:137
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::mappedPatchBase::surfPtr
const autoPtr< Foam::searchableSurface > & surfPtr() const
Return a pointer to the AMI projection surface.
Definition: mappedPatchBase.C:1061
OFstream.H
Foam::faceAreaWeightAMI
Face area weighted Arbitrary Mesh Interface (AMI) method.
Definition: faceAreaWeightAMI.H:53
Foam::vectorField
Field< vector > vectorField
Specialisation of Field<T> for vector.
Definition: primitiveFieldsFwd.H:54
Foam::meshSearchMeshObject
MeshObject wrapper around meshSearch(mesh).
Definition: meshSearchMeshObject.H:51
Foam::magSqr
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
Foam::UPstream::allWorlds
static const wordList & allWorlds() noexcept
All worlds.
Definition: UPstream.H:483
treeDataFace.H
Foam::facePoint
label facePoint(const int facei, const block &block, const label i, const label j)
Definition: blockMeshMergeTopological.C:212
Foam::objectRegistry
Registry of regIOobjects.
Definition: objectRegistry.H:60
Foam::mode
mode_t mode(const fileName &name, const bool followLink=true)
Return the file mode, normally following symbolic links.
Definition: MSwindows.C:564
Foam::mappedPatchBase::sampleRegion_
word sampleRegion_
Region to sample.
Definition: mappedPatchBase.H:232
Foam::multiWorldConnections::New
static const multiWorldConnections & New(const Time &runTime)
Access mesh object.
Definition: multiWorldConnectionsObject.C:191
Foam::mappedPatchBase::nearestWorldEqOp
Definition: mappedPatchBase.H:194
Foam::mappedPatchBase::samplePoints
tmp< pointField > samplePoints() const
Get the sample points.
Definition: mappedPatchBase.C:1716
Foam::triangle
A triangle primitive used to calculate face normals and swept volumes.
Definition: triangle.H:59
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::mappedPatchBase::sampleMesh
const polyMesh & sampleMesh() const
Get the region mesh.
Definition: mappedPatchBase.C:1649
Foam::meshSearch::findNearestFace
label findNearestFace(const point &location, const label seedFacei=-1, const bool useTreeSearch=true) const
Definition: meshSearch.C:757
Foam::reduce
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
Definition: PstreamReduceOps.H:51
Foam::PointIndexHit
This class describes the interaction of (usually) a face and a point. It carries the info of a succes...
Definition: PointIndexHit.H:52
Foam::autoPtr::set
void set(T *p) noexcept
Deprecated(2018-02) Identical to reset().
Definition: autoPtr.H:288
Foam::polyBoundaryMesh::mesh
const polyMesh & mesh() const noexcept
Return the mesh reference.
Definition: polyBoundaryMesh.H:152
Foam::UPstream::myWorldID
static label myWorldID()
My worldID.
Definition: UPstream.H:495
Foam::Field< vector >
treeDataPoint.H
Foam::polyPatch::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundaryMesh reference.
Definition: polyPatch.C:315
Foam::ITstream
An input stream of tokens.
Definition: ITstream.H:52
Foam::polyMesh::FACE_CENTRE_TRIS
Definition: polyMesh.H:104
DebugInFunction
#define DebugInFunction
Report an information message using Foam::Info.
Definition: messageStream.H:388
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:68
Foam::OSstream::name
virtual const fileName & name() const
Get the name of the stream.
Definition: OSstream.H:107
Foam::polyPatch::range
labelRange range() const
Return start/size range of this patch.
Definition: polyPatch.H:367
Foam::List::setSize
void setSize(const label n)
Alias for resize()
Definition: List.H:222
Foam::DynamicList::append
DynamicList< T, SizeMin > & append(const T &val)
Append an element to the end of this list.
Definition: DynamicListI.H:511
Foam::mapDistribute
Class containing processor-to-processor mapping information.
Definition: mapDistribute.H:163
Foam::multiWorldConnections::getCommByName
label getCommByName(const word &otherWorld) const
Get communicator for myWorld to other world connection by NAME.
Definition: multiWorldConnectionsObject.C:385
Foam::polyMesh::faceOwner
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1107
Foam::indexedOctree< Foam::treeDataCell >
Foam::dictionary::readEntry
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, bool mandatory=true) const
Definition: dictionaryTemplates.C:302
samples
scalarField samples(nIntervals, Zero)
Foam::IOstream::name
virtual const fileName & name() const
Return the name of the stream.
Definition: IOstream.C:40
Foam::Switch::name
static const char * name(const bool b) noexcept
A string representation of bool as "false" / "true".
Definition: Switch.C:145
Foam::objectRegistry::lookupObject
const Type & lookupObject(const word &name, const bool recursive=false) const
Definition: objectRegistryTemplates.C:434
Foam::polyBoundaryMesh::findPatchID
label findPatchID(const word &patchName, const bool allowNotFound=true) const
Find patch index given a name, return -1 if not found.
Definition: polyBoundaryMesh.C:765
fld
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))
Definition: gmvOutputLagrangian.H:23
Foam::indexedOctree::findInside
label findInside(const point &) const
Find shape containing point. Only implemented for certain.
Definition: indexedOctree.C:2629
forAllIters
#define forAllIters(container, iter)
Iterate across all elements in the container object.
Definition: stdFoam.H:223
Foam::UPstream::myWorld
static const word & myWorld()
My world.
Definition: UPstream.H:501
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::mappedPatchBase::offsets_
vectorField offsets_
Offset vector (nonuniform)
Definition: mappedPatchBase.H:253
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
os
OBJstream os(runTime.globalPath()/outputName)
Foam::mappedPatchBase::receivePath
static fileName receivePath(const fileName &root, const label proci)
Definition: mappedPatchBase.C:1845
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::mappedPatchBase::distance_
scalar distance_
Offset distance (normal)
Definition: mappedPatchBase.H:256
Foam::PrimitivePatch::localFaces
const List< face_type > & localFaces() const
Return patch faces addressing into local point list.
Definition: PrimitivePatch.C:317
Foam::mappedPatchBase::calcAMI
void calcAMI() const
Calculate AMI interpolator.
Definition: mappedPatchBase.C:1093
Foam::Ostream::write
virtual bool write(const token &tok)=0
Write token to stream or otherwise handle it.
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:144
treeDataCell.H
Foam::mappedPatchBase::collectSamples
void collectSamples(const label mySampleWorld, const pointField &facePoints, pointField &samples, labelList &patchFaceWorlds, labelList &patchFaceProcs, labelList &patchFaces, pointField &patchFc) const
Definition: mappedPatchBase.C:180
Foam::searchableSurface::New
static autoPtr< searchableSurface > New(const word &surfaceType, const IOobject &io, const dictionary &dict)
Return a reference to the selected searchableSurface.
Definition: searchableSurface.C:43
Foam::polyPatch::start
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:361
Foam::dictionary::subOrEmptyDict
dictionary subOrEmptyDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX, const bool mandatory=false) const
Definition: dictionary.C:540
Random.H
Foam::fileName::null
static const fileName null
An empty fileName.
Definition: fileName.H:101
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::PrimitivePatch::localPoints
const Field< point_type > & localPoints() const
Return pointField of points in patch.
Definition: PrimitivePatch.C:359
Foam::OFstream
Output to file stream, using an OSstream.
Definition: OFstream.H:53
Foam::Pstream::listCombineGather
static void listCombineGather(const List< commsStruct > &comms, List< T > &Value, const CombineOp &cop, const int tag, const label comm)
Definition: combineGatherScatter.C:290
Foam::Perr
prefixOSstream Perr
OSstream wrapped stderr (std::cerr) with parallel prefix.
Foam::mappedPatchBase::calcMapping
void calcMapping() const
Calculate mapping.
Definition: mappedPatchBase.C:735
Foam::mappedPatchBase::addWorldConnection
bool addWorldConnection()
Add a world-world connection.
Definition: mappedPatchBase.C:123
Foam::IOobject::name
const word & name() const noexcept
Return name.
Definition: IOobjectI.H:65
Time.H
Foam::autoPtr< Foam::fileName >
Foam::polyMesh::cellTree
const indexedOctree< treeDataCell > & cellTree() const
Return the cell search tree.
Definition: polyMesh.C:926
Foam::regIOobject
regIOobject is an abstract class derived from IOobject to handle automatic object registration with t...
Definition: regIOobject.H:73
Foam::primitiveMesh::cellCentres
const vectorField & cellCentres() const
Definition: primitiveMeshCellCentresAndVols.C:84
Foam::multiWorldConnections
Centralized handling of multi-world MPI connections.
Definition: multiWorldConnectionsObject.H:55
Foam::polyMesh::CELL_TETS
Definition: polyMesh.H:109
Foam::polyMesh::faces
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1094
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::Pstream::gatherList
static void gatherList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Gather data but keep individual values separate.
Definition: gatherScatterList.C:52
DynamicField.H
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::pointIndexHit
PointIndexHit< point > pointIndexHit
A PointIndexHit for 3D points.
Definition: pointIndexHit.H:46
Foam::UPstream::msgType
static int & msgType() noexcept
Message tag of standard messages.
Definition: UPstream.H:540
Foam::UPstream::myProcNo
static int myProcNo(const label communicator=worldComm)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:463
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Foam::OTstream
A simple output token stream that can be used to build token lists. Always UNCOMPRESSED.
Definition: OTstream.H:56
Foam::VectorSpace< Vector< scalar >, scalar, 3 >::max
static const Vector< scalar > max
Definition: VectorSpace.H:117
mapDistribute.H
f
labelList f(nPoints)
Foam::objectRegistry::subRegistry
const objectRegistry & subRegistry(const word &name, const bool forceCreate=false, const bool recursive=false) const
Lookup and return a const sub-objectRegistry.
Definition: objectRegistry.C:174
faceAreaWeightAMI.H
Foam::mappedPatchBase::sampleDatabasePtr_
const autoPtr< fileName > sampleDatabasePtr_
Empty or location of database.
Definition: mappedPatchBase.H:244
Foam::Vector< scalar >
Foam::UPstream::worldComm
static label worldComm
Default communicator (all processors)
Definition: UPstream.H:293
Foam::List< label >
Foam::HashSetOps::used
labelHashSet used(const bitSet &select)
Convert a bitset to a labelHashSet of the indices used.
Definition: HashOps.C:35
Foam::treeDataFace
Encapsulation of data needed to search for faces.
Definition: treeDataFace.H:59
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::primitiveMesh::faceCentres
const vectorField & faceCentres() const
Definition: primitiveMeshFaceCentresAndAreas.C:77
Foam::type
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:590
Foam::UList< label >
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::mappedPatchBase::write
virtual void write(Ostream &os) const
Write as a dictionary.
Definition: mappedPatchBase.C:1942
Foam::identity
labelList identity(const label len, label start=0)
Create identity map of the given length with (map[i] == i)
Definition: labelList.C:38
Foam::tmp::New
static tmp< T > New(Args &&... args)
Construct tmp of T with forwarding arguments.
path
fileName path(UMean.rootPath()/UMean.caseName()/"graphs"/UMean.instance())
Foam::word::null
static const word null
An empty word.
Definition: word.H:80
Foam::mappedPatchBase::offsetMode_
offsetMode offsetMode_
How to obtain samples.
Definition: mappedPatchBase.H:247
Foam::mappedPatchBase::readDatabase
static autoPtr< fileName > readDatabase(const dictionary &dict)
Read optional database name from dictionary.
Definition: mappedPatchBase.C:90
Foam::triangle::intersection
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:439
Foam::Ostream::writeEntry
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:236
Foam::Tuple2::second
const T2 & second() const noexcept
Return second.
Definition: Tuple2.H:130
Foam::mappedPatchBase::offsetMode
offsetMode
How to project face centres.
Definition: mappedPatchBase.H:130
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
Foam::fileName::clean
static bool clean(std::string &str)
Definition: fileName.C:199
rndGen
Random rndGen
Definition: createFields.H:23
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:280
Foam::VectorSpace< Vector< scalar >, scalar, 3 >::zero
static const Vector< scalar > zero
Definition: VectorSpace.H:115
Foam::UIndirectList
A List with indirect addressing.
Definition: faMatrix.H:60
Foam::boundaryMesh
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface....
Definition: boundaryMesh.H:62
Foam::mappedPatchBase::facePoints
tmp< pointField > facePoints(const polyPatch &) const
Definition: mappedPatchBase.C:152
Foam::polyMesh::FACE_PLANES
Definition: polyMesh.H:102
Foam::mappedPatchBase::getWorldCommunicator
label getWorldCommunicator() const
Get the communicator for the world-world connection.
Definition: mappedPatchBase.C:138
Foam::mappedPatchBase::surfDict_
dictionary surfDict_
Dictionary storing projection surface description.
Definition: mappedPatchBase.H:287
Foam::dictionary::add
entry * add(entry *entryPtr, bool mergeEntry=false)
Add a new entry.
Definition: dictionary.C:640
Foam::mappedPatchBase::AMIPtr_
autoPtr< AMIPatchToPatchInterpolation > AMIPtr_
Pointer to AMI interpolator.
Definition: mappedPatchBase.H:281
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:72
FatalIOErrorInFunction
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:473
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::mappedPatchBase::patch_
const polyPatch & patch_
Patch to sample.
Definition: mappedPatchBase.H:226
Foam::mappedPatchBase::sameRegion_
bool sameRegion_
Same region.
Definition: mappedPatchBase.H:262
Foam::Tuple2
A 2-tuple for storing two objects of dissimilar types. The container is similar in purpose to std::pa...
Definition: stringOps.H:60
Foam::Pstream::listCombineScatter
static void listCombineScatter(const List< commsStruct > &comms, List< T > &Value, const int tag, const label comm)
Scatter data. Reverse of combineGather.
Definition: combineGatherScatter.C:432
Foam::mappedPatchBase::lookupPatch
const polyPatch & lookupPatch(const word &sampleRegion, const word &samplePatch) const
Lookup patch.
Definition: mappedPatchBase.C:1629
Foam::PrimitivePatch::meshPoints
const labelList & meshPoints() const
Return labelList of mesh points in patch.
Definition: PrimitivePatch.C:330
Foam::TimePaths::constant
const word & constant() const
Return constant name.
Definition: TimePathsI.H:96
Foam::dictionary::getOrDefault
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:148
Foam::PtrListOps::names
List< word > names(const UPtrList< T > &list, const UnaryMatchPredicate &matcher)
OTstream.H
Foam::mappedPatchBase::facePoint
static pointIndexHit facePoint(const polyMesh &, const label facei, const polyMesh::cellDecomposition)
Get a point on the face given a face decomposition method:
Definition: mappedPatchBase.C:1723
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::Tuple2::first
const T1 & first() const noexcept
Return first.
Definition: Tuple2.H:118
Foam::PointHit::hit
bool hit() const noexcept
Is there a hit.
Definition: PointHit.H:121
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:328
Foam::mappedPatchBase::sampleMode
sampleMode
Mesh items to sample.
Definition: mappedPatchBase.H:119
Foam::objectRegistry::time
const Time & time() const noexcept
Return time registry.
Definition: objectRegistry.H:178
Foam::UPstream::nProcs
static label nProcs(const label communicator=worldComm)
Number of processes in parallel run, and 1 for serial run.
Definition: UPstream.H:445
Foam::dictionary::readIfPresent
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:405
Foam::HashSet::set
bool set(const Key &key)
Same as insert (no value to overwrite)
Definition: HashSet.H:197
Foam::PrimitivePatch
A list of faces which address into the list of points.
Definition: PrimitivePatch.H:79
sample
Minimal example by using system/controlDict.functions:
Foam::IOobject::MUST_READ
Definition: IOobject.H:185
Foam::labelUIndList
UIndirectList< label > labelUIndList
UIndirectList of labels.
Definition: UIndirectList.H:58
mappedPatchBase.H
Foam::polyMesh::tetBasePtIs
const labelIOList & tetBasePtIs() const
Return the tetBasePtIs.
Definition: polyMesh.C:892