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