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 | Copyright (C) 2011-2016 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "mappedPatchBase.H"
28 #include "ListListOps.H"
29 #include "meshSearchMeshObject.H"
30 #include "meshTools.H"
31 #include "OFstream.H"
32 #include "Random.H"
33 #include "treeDataFace.H"
34 #include "treeDataPoint.H"
35 #include "indexedOctree.H"
36 #include "polyMesh.H"
37 #include "polyPatch.H"
38 #include "Time.H"
39 #include "mapDistribute.H"
40 #include "SubField.H"
41 #include "triPointRef.H"
42 #include "syncTools.H"
43 #include "treeDataCell.H"
44 #include "DynamicField.H"
45 
46 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
47 
48 namespace Foam
49 {
50  defineTypeNameAndDebug(mappedPatchBase, 0);
51 
52  template<>
53  const char* Foam::NamedEnum
54  <
56  6
57  >::names[] =
58  {
59  "nearestCell",
60  "nearestPatchFace",
61  "nearestPatchFaceAMI",
62  "nearestPatchPoint",
63  "nearestFace",
64  "nearestOnlyCell"
65  };
66 
67  template<>
68  const char* Foam::NamedEnum
69  <
71  3
72  >::names[] =
73  {
74  "uniform",
75  "nonuniform",
76  "normal"
77  };
78 }
79 
80 
83 
86 
87 
88 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
89 
91 (
92  const polyPatch& pp
93 ) const
94 {
95  const polyMesh& mesh = pp.boundaryMesh().mesh();
96 
97  // Force construction of min-tet decomp
98  (void)mesh.tetBasePtIs();
99 
100  // Initialise to face-centre
101  tmp<pointField> tfacePoints(new pointField(patch_.size()));
102  pointField& facePoints = tfacePoints.ref();
103 
104  forAll(pp, facei)
105  {
106  facePoints[facei] = facePoint
107  (
108  mesh,
109  pp.start()+facei,
111  ).rawPoint();
112  }
113 
114  return tfacePoints;
115 }
116 
117 
119 (
120  const pointField& facePoints,
121  pointField& samples,
122  labelList& patchFaceProcs,
123  labelList& patchFaces,
124  pointField& patchFc
125 ) const
126 {
127  // Collect all sample points and the faces they come from.
128  {
129  List<pointField> globalFc(Pstream::nProcs());
130  globalFc[Pstream::myProcNo()] = facePoints;
131  Pstream::gatherList(globalFc);
132  Pstream::scatterList(globalFc);
133  // Rework into straight list
134  patchFc = ListListOps::combine<pointField>
135  (
136  globalFc,
138  );
139  }
140 
141  {
142  List<pointField> globalSamples(Pstream::nProcs());
143  globalSamples[Pstream::myProcNo()] = samplePoints(facePoints);
144  Pstream::gatherList(globalSamples);
145  Pstream::scatterList(globalSamples);
146  // Rework into straight list
147  samples = ListListOps::combine<pointField>
148  (
149  globalSamples,
151  );
152  }
153 
154  {
155  labelListList globalFaces(Pstream::nProcs());
156  globalFaces[Pstream::myProcNo()] = identity(patch_.size());
157  // Distribute to all processors
158  Pstream::gatherList(globalFaces);
159  Pstream::scatterList(globalFaces);
160 
161  patchFaces = ListListOps::combine<labelList>
162  (
163  globalFaces,
165  );
166  }
167 
168  {
169  labelList nPerProc(Pstream::nProcs());
170  nPerProc[Pstream::myProcNo()] = patch_.size();
171  Pstream::gatherList(nPerProc);
172  Pstream::scatterList(nPerProc);
173 
174  patchFaceProcs.setSize(patchFaces.size());
175 
176  label sampleI = 0;
177  forAll(nPerProc, proci)
178  {
179  for (label i = 0; i < nPerProc[proci]; i++)
180  {
181  patchFaceProcs[sampleI++] = proci;
182  }
183  }
184  }
185 }
186 
187 
188 // Find the processor/cell containing the samples. Does not account
189 // for samples being found in two processors.
191 (
192  const sampleMode mode,
193  const pointField& samples,
194  labelList& sampleProcs,
195  labelList& sampleIndices,
196  pointField& sampleLocations
197 ) const
198 {
199  // Lookup the correct region
200  const polyMesh& mesh = sampleMesh();
201 
202  // All the info for nearest. Construct to miss
203  List<nearInfo> nearest(samples.size());
204 
205  switch (mode)
206  {
207  case NEARESTCELL:
208  {
209  if (samplePatch_.size() && samplePatch_ != "none")
210  {
212  << "No need to supply a patch name when in "
213  << sampleModeNames_[mode] << " mode." << exit(FatalError);
214  }
215 
216  //- Note: face-diagonal decomposition
217  const indexedOctree<Foam::treeDataCell>& tree = mesh.cellTree();
218 
219  forAll(samples, sampleI)
220  {
221  const point& sample = samples[sampleI];
222 
223  label celli = tree.findInside(sample);
224 
225  if (celli == -1)
226  {
227  nearest[sampleI].second().first() = Foam::sqr(GREAT);
228  nearest[sampleI].second().second() = Pstream::myProcNo();
229  }
230  else
231  {
232  const point& cc = mesh.cellCentres()[celli];
233 
234  nearest[sampleI].first() = pointIndexHit
235  (
236  true,
237  cc,
238  celli
239  );
240  nearest[sampleI].second().first() = magSqr(cc-sample);
241  nearest[sampleI].second().second() = Pstream::myProcNo();
242  }
243  }
244  break;
245  }
246 
247  case NEARESTONLYCELL:
248  {
249  if (samplePatch_.size() && samplePatch_ != "none")
250  {
252  << "No need to supply a patch name when in "
253  << sampleModeNames_[mode] << " mode." << exit(FatalError);
254  }
255 
256  //- Note: face-diagonal decomposition
257  const indexedOctree<Foam::treeDataCell>& tree = mesh.cellTree();
258 
259  forAll(samples, sampleI)
260  {
261  const point& sample = samples[sampleI];
262 
263  nearest[sampleI].first() = tree.findNearest(sample, sqr(GREAT));
264  nearest[sampleI].second().first() = magSqr
265  (
266  nearest[sampleI].first().hitPoint()
267  -sample
268  );
269  nearest[sampleI].second().second() = Pstream::myProcNo();
270  }
271  break;
272  }
273 
274  case NEARESTPATCHFACE:
275  {
276  Random rndGen(123456);
277 
278  const polyPatch& pp = samplePolyPatch();
279 
280  if (pp.empty())
281  {
282  forAll(samples, sampleI)
283  {
284  nearest[sampleI].second().first() = Foam::sqr(GREAT);
285  nearest[sampleI].second().second() = Pstream::myProcNo();
286  }
287  }
288  else
289  {
290  // patch faces
291  const labelList patchFaces(identity(pp.size()) + pp.start());
292 
293  treeBoundBox patchBb
294  (
295  treeBoundBox(pp.points(), pp.meshPoints()).extend
296  (
297  rndGen,
298  1e-4
299  )
300  );
301  patchBb.min() -= point(ROOTVSMALL, ROOTVSMALL, ROOTVSMALL);
302  patchBb.max() += point(ROOTVSMALL, ROOTVSMALL, ROOTVSMALL);
303 
304  indexedOctree<treeDataFace> boundaryTree
305  (
306  treeDataFace // all information needed to search faces
307  (
308  false, // do not cache bb
309  mesh,
310  patchFaces // boundary faces only
311  ),
312  patchBb, // overall search domain
313  8, // maxLevel
314  10, // leafsize
315  3.0 // duplicity
316  );
317 
318  forAll(samples, sampleI)
319  {
320  const point& sample = samples[sampleI];
321 
322  pointIndexHit& nearInfo = nearest[sampleI].first();
323  nearInfo = boundaryTree.findNearest
324  (
325  sample,
326  magSqr(patchBb.span())
327  );
328 
329  if (!nearInfo.hit())
330  {
331  nearest[sampleI].second().first() = Foam::sqr(GREAT);
332  nearest[sampleI].second().second() =
334  }
335  else
336  {
337  point fc(pp[nearInfo.index()].centre(pp.points()));
338  nearInfo.setPoint(fc);
339  nearest[sampleI].second().first() = magSqr(fc-sample);
340  nearest[sampleI].second().second() =
342  }
343  }
344  }
345  break;
346  }
347 
348  case NEARESTPATCHPOINT:
349  {
350  Random rndGen(123456);
351 
352  const polyPatch& pp = samplePolyPatch();
353 
354  if (pp.empty())
355  {
356  forAll(samples, sampleI)
357  {
358  nearest[sampleI].second().first() = Foam::sqr(GREAT);
359  nearest[sampleI].second().second() = Pstream::myProcNo();
360  }
361  }
362  else
363  {
364  // patch (local) points
365  treeBoundBox patchBb
366  (
367  treeBoundBox(pp.points(), pp.meshPoints()).extend
368  (
369  rndGen,
370  1e-4
371  )
372  );
373  patchBb.min() -= point(ROOTVSMALL, ROOTVSMALL, ROOTVSMALL);
374  patchBb.max() += point(ROOTVSMALL, ROOTVSMALL, ROOTVSMALL);
375 
376  indexedOctree<treeDataPoint> boundaryTree
377  (
378  treeDataPoint // all information needed to search faces
379  (
380  mesh.points(),
381  pp.meshPoints() // selection of points to search on
382  ),
383  patchBb, // overall search domain
384  8, // maxLevel
385  10, // leafsize
386  3.0 // duplicity
387  );
388 
389  forAll(samples, sampleI)
390  {
391  const point& sample = samples[sampleI];
392 
393  pointIndexHit& nearInfo = nearest[sampleI].first();
394  nearInfo = boundaryTree.findNearest
395  (
396  sample,
397  magSqr(patchBb.span())
398  );
399 
400  if (!nearInfo.hit())
401  {
402  nearest[sampleI].second().first() = Foam::sqr(GREAT);
403  nearest[sampleI].second().second() =
405  }
406  else
407  {
408  const point& pt = nearInfo.hitPoint();
409 
410  nearest[sampleI].second().first() = magSqr(pt-sample);
411  nearest[sampleI].second().second() =
413  }
414  }
415  }
416  break;
417  }
418 
419  case NEARESTFACE:
420  {
421  if (samplePatch().size() && samplePatch() != "none")
422  {
424  << "No need to supply a patch name when in "
425  << sampleModeNames_[mode] << " mode." << exit(FatalError);
426  }
427 
428  //- Note: face-diagonal decomposition
429  const meshSearchMeshObject& meshSearchEngine =
431 
432  forAll(samples, sampleI)
433  {
434  const point& sample = samples[sampleI];
435 
436  label facei = meshSearchEngine.findNearestFace(sample);
437 
438  if (facei == -1)
439  {
440  nearest[sampleI].second().first() = Foam::sqr(GREAT);
441  nearest[sampleI].second().second() = Pstream::myProcNo();
442  }
443  else
444  {
445  const point& fc = mesh.faceCentres()[facei];
446 
447  nearest[sampleI].first() = pointIndexHit
448  (
449  true,
450  fc,
451  facei
452  );
453  nearest[sampleI].second().first() = magSqr(fc-sample);
454  nearest[sampleI].second().second() = Pstream::myProcNo();
455  }
456  }
457  break;
458  }
459 
460  case NEARESTPATCHFACEAMI:
461  {
462  // nothing to do here
463  return;
464  }
465 
466  default:
467  {
469  << "problem." << abort(FatalError);
470  }
471  }
472 
473 
474  // Find nearest. Combine on master.
477 
478 
479  if (debug)
480  {
482  << "mesh " << sampleRegion() << " : " << endl;
483 
484  forAll(nearest, sampleI)
485  {
486  label proci = nearest[sampleI].second().second();
487  label localI = nearest[sampleI].first().index();
488 
489  Info<< " " << sampleI << " coord:"<< samples[sampleI]
490  << " found on processor:" << proci
491  << " in local cell/face/point:" << localI
492  << " with location:" << nearest[sampleI].first().rawPoint()
493  << endl;
494  }
495  }
496 
497  // Convert back into proc+local index
498  sampleProcs.setSize(samples.size());
499  sampleIndices.setSize(samples.size());
500  sampleLocations.setSize(samples.size());
501 
502  forAll(nearest, sampleI)
503  {
504  if (!nearest[sampleI].first().hit())
505  {
506  sampleProcs[sampleI] = -1;
507  sampleIndices[sampleI] = -1;
508  sampleLocations[sampleI] = vector::max;
509  }
510  else
511  {
512  sampleProcs[sampleI] = nearest[sampleI].second().second();
513  sampleIndices[sampleI] = nearest[sampleI].first().index();
514  sampleLocations[sampleI] = nearest[sampleI].first().hitPoint();
515  }
516  }
517 }
518 
519 
521 {
522  static bool hasWarned = false;
523  if (mapPtr_.valid())
524  {
526  << "Mapping already calculated" << exit(FatalError);
527  }
528 
529  // Get points on face (since cannot use face-centres - might be off
530  // face-diagonal decomposed tets.
531  tmp<pointField> patchPoints(facePoints(patch_));
532 
533  // Get offsetted points
534  const pointField offsettedPoints(samplePoints(patchPoints()));
535 
536  // Do a sanity check - am I sampling my own patch?
537  // This only makes sense for a non-zero offset.
538  bool sampleMyself =
539  (
540  mode_ == NEARESTPATCHFACE
541  && sampleRegion() == patch_.boundaryMesh().mesh().name()
542  && samplePatch() == patch_.name()
543  );
544 
545  // Check offset
546  vectorField d(offsettedPoints-patchPoints());
547  bool coincident = (gAverage(mag(d)) <= ROOTVSMALL);
548 
549  if (sampleMyself && coincident)
550  {
552  << "Invalid offset " << d << endl
553  << "Offset is the vector added to the patch face centres to"
554  << " find the patch face supplying the data." << endl
555  << "Setting it to " << d
556  << " on the same patch, on the same region"
557  << " will find the faces themselves which does not make sense"
558  << " for anything but testing." << endl
559  << "patch_:" << patch_.name() << endl
560  << "sampleRegion_:" << sampleRegion() << endl
561  << "mode_:" << sampleModeNames_[mode_] << endl
562  << "samplePatch_:" << samplePatch() << endl
563  << "offsetMode_:" << offsetModeNames_[offsetMode_] << endl;
564  }
565 
566  // Get global list of all samples and the processor and face they come from.
568  labelList patchFaceProcs;
569  labelList patchFaces;
570  pointField patchFc;
571  collectSamples
572  (
573  patchPoints,
574  samples,
575  patchFaceProcs,
576  patchFaces,
577  patchFc
578  );
579 
580  // Find processor and cell/face samples are in and actual location.
581  labelList sampleProcs;
582  labelList sampleIndices;
583  pointField sampleLocations;
584  findSamples(mode_, samples, sampleProcs, sampleIndices, sampleLocations);
585 
586  // Check for samples that were not found. This will only happen for
587  // NEARESTCELL since finds cell containing a location
588  if (mode_ == NEARESTCELL)
589  {
590  label nNotFound = 0;
591  forAll(sampleProcs, sampleI)
592  {
593  if (sampleProcs[sampleI] == -1)
594  {
595  nNotFound++;
596  }
597  }
598  reduce(nNotFound, sumOp<label>());
599 
600  if (nNotFound > 0)
601  {
602  if (!hasWarned)
603  {
605  << "Did not find " << nNotFound
606  << " out of " << sampleProcs.size() << " total samples."
607  << " Sampling these on owner cell centre instead." << endl
608  << "On patch " << patch_.name()
609  << " on region " << sampleRegion()
610  << " in mode " << sampleModeNames_[mode_] << endl
611  << "with offset mode " << offsetModeNames_[offsetMode_]
612  << ". Suppressing further warnings from " << type() << endl;
613 
614  hasWarned = true;
615  }
616 
617  // Collect the samples that cannot be found
618  DynamicList<label> subMap;
619  DynamicField<point> subSamples;
620 
621  forAll(sampleProcs, sampleI)
622  {
623  if (sampleProcs[sampleI] == -1)
624  {
625  subMap.append(sampleI);
626  subSamples.append(samples[sampleI]);
627  }
628  }
629 
630  // And re-search for pure nearest (should not fail)
631  labelList subSampleProcs;
632  labelList subSampleIndices;
633  pointField subSampleLocations;
634  findSamples
635  (
636  NEARESTONLYCELL,
637  subSamples,
638  subSampleProcs,
639  subSampleIndices,
640  subSampleLocations
641  );
642 
643  // Insert
644  UIndirectList<label>(sampleProcs, subMap) = subSampleProcs;
645  UIndirectList<label>(sampleIndices, subMap) = subSampleIndices;
646  UIndirectList<point>(sampleLocations, subMap) = subSampleLocations;
647  }
648  }
649 
650  // Now we have all the data we need:
651  // - where sample originates from (so destination when mapping):
652  // patchFaces, patchFaceProcs.
653  // - cell/face sample is in (so source when mapping)
654  // sampleIndices, sampleProcs.
655 
656  //forAll(samples, i)
657  //{
658  // Info<< i << " need data in region "
659  // << patch_.boundaryMesh().mesh().name()
660  // << " for proc:" << patchFaceProcs[i]
661  // << " face:" << patchFaces[i]
662  // << " at:" << patchFc[i] << endl
663  // << "Found data in region " << sampleRegion()
664  // << " at proc:" << sampleProcs[i]
665  // << " face:" << sampleIndices[i]
666  // << " at:" << sampleLocations[i]
667  // << nl << endl;
668  //}
669 
670 
671 
672  if (debug && Pstream::master())
673  {
674  OFstream str
675  (
676  patch_.boundaryMesh().mesh().time().path()
677  / patch_.name()
678  + "_mapped.obj"
679  );
680  Pout<< "Dumping mapping as lines from patch faceCentres to"
681  << " sampled cell/faceCentres/points to file " << str.name()
682  << endl;
683 
684  label vertI = 0;
685 
686  forAll(patchFc, i)
687  {
688  meshTools::writeOBJ(str, patchFc[i]);
689  vertI++;
690  meshTools::writeOBJ(str, sampleLocations[i]);
691  vertI++;
692  str << "l " << vertI-1 << ' ' << vertI << nl;
693  }
694  }
695 
696  // Determine schedule.
697  mapPtr_.reset(new mapDistribute(sampleProcs, patchFaceProcs));
698 
699  // Rework the schedule from indices into samples to cell data to send,
700  // face data to receive.
701 
702  labelListList& subMap = mapPtr_().subMap();
703  labelListList& constructMap = mapPtr_().constructMap();
704 
705  forAll(subMap, proci)
706  {
707  subMap[proci] = UIndirectList<label>
708  (
709  sampleIndices,
710  subMap[proci]
711  );
712  constructMap[proci] = UIndirectList<label>
713  (
714  patchFaces,
715  constructMap[proci]
716  );
717 
718  //if (debug)
719  //{
720  // Pout<< "To proc:" << proci << " sending values of cells/faces:"
721  // << subMap[proci] << endl;
722  // Pout<< "From proc:" << proci
723  // << " receiving values of patch faces:"
724  // << constructMap[proci] << endl;
725  //}
726  }
727 
728  // Redo constructSize
729  mapPtr_().constructSize() = patch_.size();
730 
731  if (debug)
732  {
733  // Check that all elements get a value.
734  PackedBoolList used(patch_.size());
735  forAll(constructMap, proci)
736  {
737  const labelList& map = constructMap[proci];
738 
739  forAll(map, i)
740  {
741  label facei = map[i];
742 
743  if (used[facei] == 0)
744  {
745  used[facei] = 1;
746  }
747  else
748  {
750  << "On patch " << patch_.name()
751  << " patchface " << facei
752  << " is assigned to more than once."
753  << abort(FatalError);
754  }
755  }
756  }
757  forAll(used, facei)
758  {
759  if (used[facei] == 0)
760  {
762  << "On patch " << patch_.name()
763  << " patchface " << facei
764  << " is never assigned to."
765  << abort(FatalError);
766  }
767  }
768  }
769 }
770 
771 
773 const
774 {
775  const word surfType(surfDict_.lookupOrDefault<word>("type", "none"));
776 
777  if (!surfPtr_.valid() && surfType != "none")
778  {
779  word surfName(surfDict_.lookupOrDefault("name", patch_.name()));
780 
781  const polyMesh& mesh = patch_.boundaryMesh().mesh();
782 
783  surfPtr_ =
785  (
786  surfType,
787  IOobject
788  (
789  surfName,
790  mesh.time().constant(),
791  "triSurface",
792  mesh,
795  ),
796  surfDict_
797  );
798  }
799 
800  return surfPtr_;
801 }
802 
803 
805 {
806  if (AMIPtr_.valid())
807  {
809  << "AMI already calculated" << exit(FatalError);
810  }
811 
812  AMIPtr_.clear();
813 
814  if (debug)
815  {
816  const polyPatch& nbr = samplePolyPatch();
817 
818  pointField nbrPoints(nbr.localPoints());
819 
820  OFstream os(patch_.name() + "_neighbourPatch-org.obj");
821  meshTools::writeOBJ(os, samplePolyPatch().localFaces(), nbrPoints);
822 
823  // transform neighbour patch to local system
824  primitivePatch nbrPatch0
825  (
827  (
828  nbr.localFaces(),
829  nbr.size()
830  ),
831  nbrPoints
832  );
833 
834  OFstream osN(patch_.name() + "_neighbourPatch-trans.obj");
835  meshTools::writeOBJ(osN, nbrPatch0, nbrPoints);
836 
837  OFstream osO(patch_.name() + "_ownerPatch.obj");
838  meshTools::writeOBJ(osO, patch_.localFaces(), patch_.localPoints());
839  }
840 
841  // Construct/apply AMI interpolation to determine addressing and weights
842  AMIPtr_.reset
843  (
845  (
846  patch_,
847  samplePolyPatch(), // nbrPatch0,
848  surfPtr(),
850  true,
852  -1,
853  AMIReverse_
854  )
855  );
856 }
857 
858 
859 // Hack to read old (List-based) format. See Field.C. The difference
860 // is only that in case of falling back to old format it expects a non-uniform
861 // list instead of a single vector.
863 (
864  const word& keyword,
865  const dictionary& dict,
866  const label size
867 )
868 {
869  tmp<pointField> tfld(new pointField());
870  pointField& fld = tfld.ref();
871 
872  if (size)
873  {
874  ITstream& is = dict.lookup(keyword);
875 
876  // Read first token
877  token firstToken(is);
878 
879  if (firstToken.isWord())
880  {
881  if (firstToken.wordToken() == "uniform")
882  {
883  fld.setSize(size);
884  fld = pTraits<vector>(is);
885  }
886  else if (firstToken.wordToken() == "nonuniform")
887  {
888  is >> static_cast<List<vector>&>(fld);
889  if (fld.size() != size)
890  {
892  (
893  dict
894  ) << "size " << fld.size()
895  << " is not equal to the given value of " << size
896  << exit(FatalIOError);
897  }
898  }
899  else
900  {
902  (
903  dict
904  ) << "expected keyword 'uniform' or 'nonuniform', found "
905  << firstToken.wordToken()
906  << exit(FatalIOError);
907  }
908  }
909  else
910  {
911  if (is.version() == 2.0)
912  {
914  (
915  dict
916  ) << "expected keyword 'uniform' or 'nonuniform', "
917  "assuming List format for backwards compatibility."
918  "Foam version 2.0." << endl;
919 
920  is.putBack(firstToken);
921  is >> static_cast<List<vector>&>(fld);
922  }
923  }
924  }
925  return tfld;
926 }
927 
928 
929 // * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * * * * * //
930 
932 (
933  const polyPatch& pp
934 )
935 :
936  patch_(pp),
937  sampleRegion_(patch_.boundaryMesh().mesh().name()),
938  mode_(NEARESTPATCHFACE),
939  samplePatch_(""),
940  coupleGroup_(),
941  offsetMode_(UNIFORM),
942  offset_(Zero),
943  offsets_(pp.size(), offset_),
944  distance_(0),
945  sameRegion_(sampleRegion_ == patch_.boundaryMesh().mesh().name()),
946  mapPtr_(nullptr),
947  AMIPtr_(nullptr),
948  AMIReverse_(false),
949  surfPtr_(nullptr),
950  surfDict_(fileName("surface"))
951 {}
952 
953 
955 (
956  const polyPatch& pp,
957  const word& sampleRegion,
958  const sampleMode mode,
959  const word& samplePatch,
960  const vectorField& offsets
961 )
962 :
963  patch_(pp),
964  sampleRegion_(sampleRegion),
965  mode_(mode),
966  samplePatch_(samplePatch),
967  coupleGroup_(),
968  offsetMode_(NONUNIFORM),
969  offset_(Zero),
970  offsets_(offsets),
971  distance_(0),
972  sameRegion_(sampleRegion_ == patch_.boundaryMesh().mesh().name()),
973  mapPtr_(nullptr),
974  AMIPtr_(nullptr),
975  AMIReverse_(false),
976  surfPtr_(nullptr),
977  surfDict_(fileName("surface"))
978 {}
979 
980 
982 (
983  const polyPatch& pp,
984  const word& sampleRegion,
985  const sampleMode mode,
986  const word& samplePatch,
987  const vector& offset
988 )
989 :
990  patch_(pp),
991  sampleRegion_(sampleRegion),
992  mode_(mode),
993  samplePatch_(samplePatch),
994  coupleGroup_(),
995  offsetMode_(UNIFORM),
996  offset_(offset),
997  offsets_(0),
998  distance_(0),
999  sameRegion_(sampleRegion_ == patch_.boundaryMesh().mesh().name()),
1000  mapPtr_(nullptr),
1001  AMIPtr_(nullptr),
1002  AMIReverse_(false),
1003  surfPtr_(nullptr),
1004  surfDict_(fileName("surface"))
1005 {}
1006 
1007 
1010  const polyPatch& pp,
1011  const word& sampleRegion,
1012  const sampleMode mode,
1013  const word& samplePatch,
1014  const scalar distance
1015 )
1016 :
1017  patch_(pp),
1018  sampleRegion_(sampleRegion),
1019  mode_(mode),
1020  samplePatch_(samplePatch),
1021  coupleGroup_(),
1022  offsetMode_(NORMAL),
1023  offset_(Zero),
1024  offsets_(0),
1025  distance_(distance),
1026  sameRegion_(sampleRegion_ == patch_.boundaryMesh().mesh().name()),
1027  mapPtr_(nullptr),
1028  AMIPtr_(nullptr),
1029  AMIReverse_(false),
1030  surfPtr_(nullptr),
1031  surfDict_(fileName("surface"))
1032 {}
1033 
1034 
1037  const polyPatch& pp,
1038  const dictionary& dict
1039 )
1040 :
1041  patch_(pp),
1042  sampleRegion_(dict.lookupOrDefault<word>("sampleRegion", "")),
1043  mode_(sampleModeNames_.read(dict.lookup("sampleMode"))),
1044  samplePatch_(dict.lookupOrDefault<word>("samplePatch", "")),
1045  coupleGroup_(dict),
1046  offsetMode_(UNIFORM),
1047  offset_(Zero),
1048  offsets_(0),
1049  distance_(0.0),
1050  sameRegion_(sampleRegion_ == patch_.boundaryMesh().mesh().name()),
1051  mapPtr_(nullptr),
1052  AMIPtr_(nullptr),
1053  AMIReverse_(dict.lookupOrDefault<bool>("flipNormals", false)),
1054  surfPtr_(nullptr),
1055  surfDict_(dict.subOrEmptyDict("surface"))
1056 {
1057  if (!coupleGroup_.valid())
1058  {
1059  if (sampleRegion_.empty())
1060  {
1061  // If no coupleGroup and no sampleRegion assume local region
1062  sampleRegion_ = patch_.boundaryMesh().mesh().name();
1063  sameRegion_ = true;
1064  }
1065  }
1066 
1067  if (dict.found("offsetMode"))
1068  {
1069  offsetMode_ = offsetModeNames_.read(dict.lookup("offsetMode"));
1070 
1071  switch(offsetMode_)
1072  {
1073  case UNIFORM:
1074  {
1075  offset_ = point(dict.lookup("offset"));
1076  }
1077  break;
1078 
1079  case NONUNIFORM:
1080  {
1081  //offsets_ = pointField(dict.lookup("offsets"));
1082  offsets_ = readListOrField("offsets", dict, patch_.size());
1083  }
1084  break;
1085 
1086  case NORMAL:
1087  {
1088  distance_ = readScalar(dict.lookup("distance"));
1089  }
1090  break;
1091  }
1092  }
1093  else if (dict.found("offset"))
1094  {
1095  offsetMode_ = UNIFORM;
1096  offset_ = point(dict.lookup("offset"));
1097  }
1098  else if (dict.found("offsets"))
1099  {
1100  offsetMode_ = NONUNIFORM;
1101  //offsets_ = pointField(dict.lookup("offsets"));
1102  offsets_ = readListOrField("offsets", dict, patch_.size());
1103  }
1104  else if (mode_ != NEARESTPATCHFACE && mode_ != NEARESTPATCHFACEAMI)
1105  {
1107  (
1108  dict
1109  ) << "Please supply the offsetMode as one of "
1111  << exit(FatalIOError);
1112  }
1113 }
1114 
1115 
1118  const polyPatch& pp,
1119  const sampleMode mode,
1120  const dictionary& dict
1121 )
1122 :
1123  patch_(pp),
1124  sampleRegion_(dict.lookupOrDefault<word>("sampleRegion", "")),
1125  mode_(mode),
1126  samplePatch_(dict.lookupOrDefault<word>("samplePatch", "")),
1127  coupleGroup_(dict), //dict.lookupOrDefault<word>("coupleGroup", "")),
1128  offsetMode_(UNIFORM),
1129  offset_(Zero),
1130  offsets_(0),
1131  distance_(0.0),
1132  sameRegion_(sampleRegion_ == patch_.boundaryMesh().mesh().name()),
1133  mapPtr_(nullptr),
1134  AMIPtr_(nullptr),
1135  AMIReverse_(dict.lookupOrDefault<bool>("flipNormals", false)),
1136  surfPtr_(nullptr),
1137  surfDict_(dict.subOrEmptyDict("surface"))
1138 {
1139  if (mode != NEARESTPATCHFACE && mode != NEARESTPATCHFACEAMI)
1140  {
1142  (
1143  dict
1144  ) << "Construct from sampleMode and dictionary only applicable for "
1145  << " collocated patches in modes "
1146  << sampleModeNames_[NEARESTPATCHFACE] << ','
1147  << sampleModeNames_[NEARESTPATCHFACEAMI]
1148  << exit(FatalIOError);
1149  }
1150 
1151 
1152  if (!coupleGroup_.valid())
1153  {
1154  if (sampleRegion_.empty())
1155  {
1156  // If no coupleGroup and no sampleRegion assume local region
1157  sampleRegion_ = patch_.boundaryMesh().mesh().name();
1158  sameRegion_ = true;
1159  }
1160  }
1161 }
1162 
1163 
1166  const polyPatch& pp,
1167  const mappedPatchBase& mpb
1168 )
1169 :
1170  patch_(pp),
1171  sampleRegion_(mpb.sampleRegion_),
1172  mode_(mpb.mode_),
1173  samplePatch_(mpb.samplePatch_),
1174  coupleGroup_(mpb.coupleGroup_),
1175  offsetMode_(mpb.offsetMode_),
1176  offset_(mpb.offset_),
1177  offsets_(mpb.offsets_),
1178  distance_(mpb.distance_),
1179  sameRegion_(mpb.sameRegion_),
1180  mapPtr_(nullptr),
1181  AMIPtr_(nullptr),
1182  AMIReverse_(mpb.AMIReverse_),
1183  surfPtr_(nullptr),
1184  surfDict_(mpb.surfDict_)
1185 {}
1186 
1187 
1190  const polyPatch& pp,
1191  const mappedPatchBase& mpb,
1192  const labelUList& mapAddressing
1193 )
1194 :
1195  patch_(pp),
1196  sampleRegion_(mpb.sampleRegion_),
1197  mode_(mpb.mode_),
1198  samplePatch_(mpb.samplePatch_),
1199  coupleGroup_(mpb.coupleGroup_),
1200  offsetMode_(mpb.offsetMode_),
1201  offset_(mpb.offset_),
1202  offsets_
1203  (
1204  offsetMode_ == NONUNIFORM
1205  ? vectorField(mpb.offsets_, mapAddressing)
1206  : vectorField(0)
1207  ),
1208  distance_(mpb.distance_),
1209  sameRegion_(mpb.sameRegion_),
1210  mapPtr_(nullptr),
1211  AMIPtr_(nullptr),
1212  AMIReverse_(mpb.AMIReverse_),
1213  surfPtr_(nullptr),
1214  surfDict_(mpb.surfDict_)
1215 {}
1216 
1217 
1218 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
1219 
1221 {
1222  clearOut();
1223 }
1224 
1225 
1227 {
1228  mapPtr_.clear();
1229  AMIPtr_.clear();
1230  surfPtr_.clear();
1231 }
1232 
1233 
1234 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
1235 
1237 {
1238  return patch_.boundaryMesh().mesh().time().lookupObject<polyMesh>
1239  (
1240  sampleRegion()
1241  );
1242 }
1243 
1244 
1246 {
1247  const polyMesh& nbrMesh = sampleMesh();
1248 
1249  const label patchi = nbrMesh.boundaryMesh().findPatchID(samplePatch());
1250 
1251  if (patchi == -1)
1252  {
1254  << "Cannot find patch " << samplePatch()
1255  << " in region " << sampleRegion_ << endl
1256  << "Valid patches are " << nbrMesh.boundaryMesh().names()
1257  << exit(FatalError);
1258  }
1259 
1260  return nbrMesh.boundaryMesh()[patchi];
1261 }
1262 
1263 
1266  const pointField& fc
1267 ) const
1268 {
1269  tmp<pointField> tfld(new pointField(fc));
1270  pointField& fld = tfld.ref();
1271 
1272  switch (offsetMode_)
1273  {
1274  case UNIFORM:
1275  {
1276  fld += offset_;
1277  break;
1278  }
1279  case NONUNIFORM:
1280  {
1281  fld += offsets_;
1282  break;
1283  }
1284  case NORMAL:
1285  {
1286  // Get outwards pointing normal
1287  vectorField n(patch_.faceAreas());
1288  n /= mag(n);
1289 
1290  fld += distance_*n;
1291  break;
1292  }
1293  }
1294 
1295  return tfld;
1296 }
1297 
1298 
1300 {
1301  return samplePoints(facePoints(patch_));
1302 }
1303 
1304 
1307  const polyMesh& mesh,
1308  const label facei,
1309  const polyMesh::cellDecomposition decompMode
1310 )
1311 {
1312  const point& fc = mesh.faceCentres()[facei];
1313 
1314  switch (decompMode)
1315  {
1316  case polyMesh::FACE_PLANES:
1318  {
1319  // For both decompositions the face centre is guaranteed to be
1320  // on the face
1321  return pointIndexHit(true, fc, facei);
1322  }
1323  break;
1324 
1326  case polyMesh::CELL_TETS:
1327  {
1328  // Find the intersection of a ray from face centre to cell centre
1329  // Find intersection of (face-centre-decomposition) centre to
1330  // cell-centre with face-diagonal-decomposition triangles.
1331 
1332  const pointField& p = mesh.points();
1333  const face& f = mesh.faces()[facei];
1334 
1335  if (f.size() <= 3)
1336  {
1337  // Return centre of triangle.
1338  return pointIndexHit(true, fc, 0);
1339  }
1340 
1341  label celli = mesh.faceOwner()[facei];
1342  const point& cc = mesh.cellCentres()[celli];
1343  vector d = fc-cc;
1344 
1345  const label fp0 = mesh.tetBasePtIs()[facei];
1346  const point& basePoint = p[f[fp0]];
1347 
1348  label fp = f.fcIndex(fp0);
1349  for (label i = 2; i < f.size(); i++)
1350  {
1351  const point& thisPoint = p[f[fp]];
1352  label nextFp = f.fcIndex(fp);
1353  const point& nextPoint = p[f[nextFp]];
1354 
1355  const triPointRef tri(basePoint, thisPoint, nextPoint);
1356  pointHit hitInfo = tri.intersection
1357  (
1358  cc,
1359  d,
1361  );
1362 
1363  if (hitInfo.hit() && hitInfo.distance() > 0)
1364  {
1365  return pointIndexHit(true, hitInfo.hitPoint(), i-2);
1366  }
1367 
1368  fp = nextFp;
1369  }
1370 
1371  // Fall-back
1372  return pointIndexHit(false, fc, -1);
1373  }
1374  break;
1375 
1376  default:
1377  {
1379  << "problem" << abort(FatalError);
1380  return pointIndexHit();
1381  }
1382  }
1383 }
1384 
1385 
1387 {
1388  os.writeKeyword("sampleMode") << sampleModeNames_[mode_]
1389  << token::END_STATEMENT << nl;
1390  if (!sampleRegion_.empty())
1391  {
1392  os.writeKeyword("sampleRegion") << sampleRegion_
1393  << token::END_STATEMENT << nl;
1394  }
1395  if (!samplePatch_.empty())
1396  {
1397  os.writeKeyword("samplePatch") << samplePatch_
1398  << token::END_STATEMENT << nl;
1399  }
1400  coupleGroup_.write(os);
1401 
1402  if
1403  (
1404  offsetMode_ == UNIFORM
1405  && offset_ == vector::zero
1406  && (mode_ == NEARESTPATCHFACE || mode_ == NEARESTPATCHFACEAMI)
1407  )
1408  {
1409  // Collocated mode. No need to write offset data
1410  }
1411  else
1412  {
1413  os.writeKeyword("offsetMode") << offsetModeNames_[offsetMode_]
1414  << token::END_STATEMENT << nl;
1415 
1416  switch (offsetMode_)
1417  {
1418  case UNIFORM:
1419  {
1420  os.writeKeyword("offset") << offset_ << token::END_STATEMENT
1421  << nl;
1422  break;
1423  }
1424  case NONUNIFORM:
1425  {
1426  offsets_.writeEntry("offsets", os);
1427  break;
1428  }
1429  case NORMAL:
1430  {
1431  os.writeKeyword("distance") << distance_ << token::END_STATEMENT
1432  << nl;
1433  break;
1434  }
1435  }
1436 
1437  if (mode_ == NEARESTPATCHFACEAMI)
1438  {
1439  if (AMIReverse_)
1440  {
1441  os.writeKeyword("flipNormals") << AMIReverse_
1442  << token::END_STATEMENT << nl;
1443  }
1444 
1445  if (!surfDict_.empty())
1446  {
1447  os.writeKeyword(surfDict_.dictName());
1448  os << surfDict_;
1449  }
1450  }
1451  }
1452 }
1453 
1454 
1455 // ************************************************************************* //
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:424
word samplePatch_
Patch (if in sampleMode NEARESTPATCH*)
virtual const fileName & name() const
Return the name of the stream.
Definition: OSstream.H:88
static wordList words()
The set of names as a list of words.
Definition: NamedEnum.C:105
static void scatterList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Scatter data. Reverse of gatherList.
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:431
static void listCombineScatter(const List< commsStruct > &comms, List< T > &Value, const int tag, const label comm)
Scatter data. Reverse of combineGather.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
bool isWord() const
Definition: tokenI.H:213
A triangle primitive used to calculate face normals and swept volumes.
Definition: triangle.H:58
vector offset_
Offset vector (uniform)
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
A class for handling file names.
Definition: fileName.H:69
Holds (reference to) pointField. Encapsulation of data needed for octree searches. Used for searching for nearest point. No bounding boxes around points. Only overlaps and calcNearest are implemented, rest makes little sense.
Definition: treeDataPoint.H:59
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
word sampleRegion_
Region to sample.
const double e
Elementary charge.
Definition: doubleFloat.H:78
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
MeshObject wrapper around meshSearch(mesh).
tmp< pointField > samplePoints() const
Get the sample points.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
A 2-tuple for storing two objects of different types.
Definition: HashTable.H:66
const polyMesh & sampleMesh() const
Get the region mesh.
const labelIOList & tetBasePtIs() const
Return the tetBasePtIs.
Definition: polyMesh.C:841
Encapsulation of data needed to search for faces.
Definition: treeDataFace.H:58
bool sameRegion_
Same region.
const Type1 & first() const
Return first.
Definition: Tuple2.H:94
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const polyBoundaryMesh & boundaryMesh() const
Return boundaryMesh reference.
Definition: polyPatch.C:278
virtual void write(Ostream &) const
Write as a dictionary.
const word & wordToken() const
Definition: tokenI.H:218
Output to file stream.
Definition: OFstream.H:82
const coupleGroupIdentifier coupleGroup_
PatchGroup (if in sampleMode NEARESTPATCH*)
label findNearestFace(const point &location, const label seedFacei=-1, const bool useTreeSearch=true) const
Definition: meshSearch.C:669
virtual ~mappedPatchBase()
Destructor.
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:163
A token holds items read from Istream.
Definition: token.H:69
label fcIndex(const label i) const
Return the forward circular index, i.e. the next index.
Definition: UListI.H:58
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:418
void putBack(const token &)
Put back token.
Definition: Istream.C:30
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
Traits class for primitives.
Definition: pTraits.H:50
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:412
label findPatchID(const word &patchName) const
Find patch index given a name.
labelList identity(const label len)
Create identity map (map[i] == i) of given length.
Definition: ListOps.C:104
This class describes the interaction of (usually) a face and a point. It carries the info of a succes...
Definition: PointIndexHit.H:53
PointIndexHit< point > pointIndexHit
Definition: pointIndexHit.H:42
void writeOBJ(Ostream &os, const point &pt)
Write obj representation of point.
Definition: meshTools.C:203
versionNumber version() const
Return the stream version.
Definition: IOstream.H:399
T & first()
Return the first element of the list.
Definition: UListI.H:114
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
Initialise the NamedEnum HashTable from the static list of names.
Definition: NamedEnum.H:52
static void listCombineGather(const List< commsStruct > &comms, List< T > &Value, const CombineOp &cop, const int tag, const label comm)
cellDecomposition
Enumeration defining the decomposition of the cell for.
Definition: polyMesh.H:98
vectorField offsets_
Offset vector (nonuniform)
void setPoint(const Point &p)
Macros for easy insertion into run-time selection tables.
const labelList & meshPoints() const
Return labelList of mesh points in patch. They are constructed.
scalar distance_
Offset distance (normal)
void collectSamples(const pointField &facePoints, pointField &, labelList &patchFaceProcs, labelList &patchFaces, pointField &patchFc) const
Collect single list of samples and originating processor+face.
static const meshSearchMeshObject & New(const polyMesh &mesh)
Definition: MeshObject.C:44
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1011
const List< Face > & localFaces() const
Return patch faces addressing into local point list.
scalarField samples(nIntervals, 0)
A list of faces which address into the list of points.
const Point & hitPoint() const
Return hit point.
A List obtained as a section of another List.
Definition: SubList.H:53
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
dynamicFvMesh & mesh
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< ' ';}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< ' ';}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< ' ';}gmvFile<< nl;forAll(lagrangianScalarNames, i){ const word &name=lagrangianScalarNames[i];IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
A class for handling words, derived from string.
Definition: word.H:59
label size() const
Return the number of elements in the list.
virtual const fileName & name() const
Return the name of the stream.
Definition: IOstream.H:297
const Point & hitPoint() const
Return hit point.
Definition: PointHit.H:126
wordList names() const
Return a list of patch names.
mode_t mode(const fileName &, const bool followLink=true)
Return the file mode.
Definition: POSIX.C:467
const word & constant() const
Return constant name.
Definition: TimePaths.H:124
cachedRandom rndGen(label(0), -1)
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1049
static tmp< pointField > readListOrField(const word &keyword, const dictionary &dict, const label size)
Helper to read field or non-uniform list from dictionary.
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Definition: DynamicListI.H:292
bool hit() const
Is there a hit.
const Field< PointType > & points() const
Return reference to global points.
void calcMapping() const
Calculate mapping.
const polyMesh & mesh() const
Return the mesh reference.
Determines a mapping between patch face centres and mesh cell or face centres and processors they&#39;re ...
bool hit() const
Is there a hit.
Definition: PointHit.H:120
Dynamically sized Field.
Definition: DynamicField.H:49
static const zero Zero
Definition: zero.H:91
offsetMode offsetMode_
How to obtain samples.
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if succesful.
Definition: doubleScalar.H:63
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1036
const vectorField & cellCentres() const
errorManip< error > abort(error &err)
Definition: errorManip.H:131
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: HashTable.H:61
Simple random number generator.
Definition: Random.H:49
tmp< pointField > facePoints(const polyPatch &) const
Get the points from face-centre-decomposition face centres.
const autoPtr< Foam::searchableSurface > & surfPtr() const
Return a pointer to the AMI projection surface.
dimensioned< scalar > magSqr(const dimensioned< Type > &)
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
const Field< PointType > & localPoints() const
Return pointField of points in patch.
static const char nl
Definition: Ostream.H:262
label findInside(const point &) const
Find shape containing point. Only implemented for certain.
label facePoint(const int facei, const block &block, const label i, const label j)
const Time & time() const
Return time.
defineTypeNameAndDebug(combustionModel, 0)
static pointIndexHit facePoint(const polyMesh &, const label facei, const polyMesh::cellDecomposition)
Get a point on the face given a face decomposition method:
Ostream & writeKeyword(const keyType &)
Write the keyword followed by an appropriate indentation.
Definition: Ostream.C:54
labelList f(nPoints)
offsetMode
How to project face centres.
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
static const NamedEnum< sampleMode, 6 > sampleModeNames_
fileName::Type type(const fileName &, const bool followLink=true)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:485
void findSamples(const sampleMode mode, const pointField &, labelList &sampleProcs, labelList &sampleIndices, pointField &sampleLocations) const
Find cells/faces containing samples.
dictionary surfDict_
Dictionary storing projection surface description.
const vectorField & faceCentres() const
void setSize(const label)
Reset size of List.
Definition: List.C:281
const indexedOctree< treeDataCell > & cellTree() const
Return the cell search tree.
Definition: polyMesh.C:875
A bit-packed bool list.
static label nProcs(const label communicator=0)
Number of processes in parallel run.
Definition: UPstream.H:400
label patchi
Class containing processor-to-processor mapping information.
void calcAMI() const
Calculate AMI interpolator.
vector point
Point is a vector.
Definition: point.H:41
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
#define WarningInFunction
Report a warning using Foam::Warning.
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:331
Type gAverage(const FieldField< Field, Type > &f)
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
A List with indirect addressing.
Definition: fvMatrix.H:106
sampleMode
Mesh items to sample.
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:300
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 bool AMIReverse_
Flag to indicate that slave patch should be reversed for AMI.
Standard boundBox + extra functionality for use in octree.
Definition: treeBoundBox.H:87
messageStream Info
dimensioned< scalar > mag(const dimensioned< Type > &)
label n
Field< vector > vectorField
Specialisation of Field<T> for vector.
#define IOWarningInFunction(ios)
Report an IO warning using Foam::Warning.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
scalar distance() const
Return distance to hit.
Definition: PointHit.H:139
const point & min() const
Minimum describing the bounding box.
Definition: boundBoxI.H:54
static const NamedEnum< offsetMode, 3 > offsetModeNames_
volScalarField & p
A class for managing temporary objects.
Definition: PtrList.H:53
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:427
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
static void gatherList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Gather data but keep individual values separate.
const polyPatch & samplePolyPatch() const
Get the patch on the region.
label index() const
Return index.
mappedPatchBase(const polyPatch &)
Construct from patch.
DynamicField< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Input token stream.
Definition: ITstream.H:49
dictionary subOrEmptyDict(const word &, const bool mustRead=false) const
Find and return a sub-dictionary as a copy, or.
Definition: dictionary.C:727
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:576
const sampleMode mode_
What to sample.
IOerror FatalIOError
#define InfoInFunction
Report an information message using Foam::Info.