collection_searchableSurface.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2025 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 
27 #include "SortableList.H"
28 #include "Time.H"
29 #include "ListOps.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36  namespace searchableSurfaces
37  {
39 
41  (
43  collection,
45  );
46 
48  (
50  collection,
51  dictionary,
52  searchableSurfaceCollection,
53  "searchableSurfaceCollection"
54  );
55  }
56 }
57 
58 
59 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
60 
61 void Foam::searchableSurfaces::collection::findNearest
62 (
63  const pointField& samples,
64  scalarField& minDistSqr,
65  List<pointIndexHit>& nearestInfo,
66  labelList& nearestSurf
67 ) const
68 {
69  // Initialise
70  nearestInfo.setSize(samples.size());
71  nearestInfo = pointIndexHit();
72  nearestSurf.setSize(samples.size());
73  nearestSurf = -1;
74 
75  List<pointIndexHit> hitInfo(samples.size());
76 
77  const scalarField localMinDistSqr(samples.size(), great);
78 
79  forAll(subGeom_, surfI)
80  {
81  subGeom_[surfI].findNearest
82  (
83  cmptDivide // Transform then divide
84  (
85  transform_[surfI].localPosition(samples),
86  scale_[surfI]
87  ),
88  localMinDistSqr,
89  hitInfo
90  );
91 
92  forAll(hitInfo, pointi)
93  {
94  if (hitInfo[pointi].hit())
95  {
96  // Rework back into global coordinate sys. Multiply then
97  // transform
98  point globalPt = transform_[surfI].globalPosition
99  (
101  (
102  hitInfo[pointi].rawPoint(),
103  scale_[surfI]
104  )
105  );
106 
107  scalar distSqr = magSqr(globalPt - samples[pointi]);
108 
109  if (distSqr < minDistSqr[pointi])
110  {
111  minDistSqr[pointi] = distSqr;
112  nearestInfo[pointi].setPoint(globalPt);
113  nearestInfo[pointi].setHit();
114  nearestInfo[pointi].setIndex
115  (
116  hitInfo[pointi].index()
117  + indexOffset_[surfI]
118  );
119  nearestSurf[pointi] = surfI;
120  }
121  }
122  }
123  }
124 }
125 
126 
127 // Sort hits into per-surface bins. Misses are rejected. Maintains map back
128 // to position
129 void Foam::searchableSurfaces::collection::sortHits
130 (
131  const List<pointIndexHit>& info,
132  List<List<pointIndexHit>>& surfInfo,
133  labelListList& infoMap
134 ) const
135 {
136  // Count hits per surface.
137  labelList nHits(subGeom_.size(), 0);
138 
139  forAll(info, pointi)
140  {
141  if (info[pointi].hit())
142  {
143  label index = info[pointi].index();
144  label surfI = findLower(indexOffset_, index+1);
145  nHits[surfI]++;
146  }
147  }
148 
149  // Per surface the hit
150  surfInfo.setSize(subGeom_.size());
151  // Per surface the original position
152  infoMap.setSize(subGeom_.size());
153 
154  forAll(surfInfo, surfI)
155  {
156  surfInfo[surfI].setSize(nHits[surfI]);
157  infoMap[surfI].setSize(nHits[surfI]);
158  }
159  nHits = 0;
160 
161  forAll(info, pointi)
162  {
163  if (info[pointi].hit())
164  {
165  label index = info[pointi].index();
166  label surfI = findLower(indexOffset_, index+1);
167 
168  // Store for correct surface and adapt indices back to local
169  // ones
170  label localI = nHits[surfI]++;
171  surfInfo[surfI][localI] = pointIndexHit
172  (
173  info[pointi].hit(),
174  info[pointi].rawPoint(),
175  index-indexOffset_[surfI]
176  );
177  infoMap[surfI][localI] = pointi;
178  }
179  }
180 }
181 
182 
183 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
184 
186 (
187  const IOobject& io,
188  const dictionary& dict
189 )
190 :
191  searchableSurface(io),
192  instance_(dict.size()),
193  scale_(dict.size()),
194  transform_(dict.size()),
195  subGeom_(dict.size()),
196  mergeSubRegions_(dict.lookup("mergeSubRegions")),
197  indexOffset_(dict.size()+1)
198 {
199  Info<< "SearchableCollection : " << name() << endl;
200 
201  label surfI = 0;
202  label startIndex = 0;
204  {
205  if (dict.isDict(iter().keyword()))
206  {
207  instance_[surfI] = iter().keyword();
208 
209  const dictionary& subDict = dict.subDict(instance_[surfI]);
210 
211  scale_[surfI] = subDict.lookup<vector>("scale", dimless);
212  transform_.set
213  (
214  surfI,
216  (
217  io.db(),
218  subDict.subDict("transform")
219  )
220  );
221 
222  const word subGeomName(subDict.lookup("surface"));
223  // Pout<< "Trying to find " << subGeomName << endl;
224 
225  const searchableSurface& s =
226  io.db().lookupObject<searchableSurface>(subGeomName);
227 
228  // I don't know yet how to handle the globalSize combined with
229  // regionOffset. Would cause non-consecutive indices locally
230  // if all indices offset by globalSize() of the local region...
231  if (s.size() != s.globalSize())
232  {
234  << "Cannot use a distributed surface in a collection."
235  << exit(FatalError);
236  }
237 
238  subGeom_.set(surfI, &const_cast<searchableSurface&>(s));
239 
240  indexOffset_[surfI] = startIndex;
241  startIndex += subGeom_[surfI].size();
242 
243  Info<< " instance : " << instance_[surfI] << endl;
244  Info<< " surface : " << s.name() << endl;
245  Info<< " scale : " << scale_[surfI] << endl;
246  Info<< " coordsys : " << transform_[surfI] << endl;
247 
248  surfI++;
249  }
250  }
251  indexOffset_[surfI] = startIndex;
252 
253  instance_.setSize(surfI);
254  scale_.setSize(surfI);
255  transform_.setSize(surfI);
256  subGeom_.setSize(surfI);
257  indexOffset_.setSize(surfI+1);
258 
259  // Bounds is the overall bounds
261 
262  forAll(subGeom_, surfI)
263  {
264  const boundBox& surfBb = subGeom_[surfI].bounds();
265 
266  // Transform back to global coordinate sys.
267  const point surfBbMin = transform_[surfI].globalPosition
268  (
270  (
271  surfBb.min(),
272  scale_[surfI]
273  )
274  );
275  const point surfBbMax = transform_[surfI].globalPosition
276  (
278  (
279  surfBb.max(),
280  scale_[surfI]
281  )
282  );
283 
284  bounds().min() = min(bounds().min(), surfBbMin);
285  bounds().max() = max(bounds().max(), surfBbMax);
286  }
287 }
288 
289 
290 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
291 
293 {}
294 
295 
296 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
297 
299 {
300  if (regions_.size() == 0)
301  {
302  regionOffset_.setSize(subGeom_.size());
303 
304  DynamicList<word> allRegions;
305  forAll(subGeom_, surfI)
306  {
307  regionOffset_[surfI] = allRegions.size();
308 
309  if (mergeSubRegions_)
310  {
311  // Single name regardless how many regions subsurface has
312  allRegions.append(instance_[surfI] + "_" + Foam::name(surfI));
313  }
314  else
315  {
316  const wordList& subRegions = subGeom_[surfI].regions();
317 
318  forAll(subRegions, i)
319  {
320  allRegions.append(instance_[surfI] + "_" + subRegions[i]);
321  }
322  }
323  }
324  regions_.transfer(allRegions.shrink());
325  }
326  return regions_;
327 }
328 
329 
331 {
332  return indexOffset_.last();
333 }
334 
335 
338 {
339  tmp<pointField> tCtrs = tmp<pointField>(new pointField(size()));
340  pointField& ctrs = tCtrs.ref();
341 
342  // Append individual coordinates
343  label coordI = 0;
344 
345  forAll(subGeom_, surfI)
346  {
347  const pointField subCoords(subGeom_[surfI].coordinates());
348 
349  forAll(subCoords, i)
350  {
351  ctrs[coordI++] = transform_[surfI].globalPosition
352  (
354  (
355  subCoords[i],
356  scale_[surfI]
357  )
358  );
359  }
360  }
361 
362  return tCtrs;
363 }
364 
365 
367 (
368  pointField& centres,
369  scalarField& radiusSqr
370 ) const
371 {
372  centres.setSize(size());
373  radiusSqr.setSize(centres.size());
374 
375  // Append individual coordinates
376  label coordI = 0;
377 
378  forAll(subGeom_, surfI)
379  {
380  scalar maxScale = cmptMax(scale_[surfI]);
381 
382  pointField subCentres;
383  scalarField subRadiusSqr;
384  subGeom_[surfI].boundingSpheres(subCentres, subRadiusSqr);
385 
386  forAll(subCentres, i)
387  {
388  centres[coordI] = transform_[surfI].globalPosition
389  (
391  (
392  subCentres[i],
393  scale_[surfI]
394  )
395  );
396  radiusSqr[coordI] = maxScale*subRadiusSqr[i];
397  coordI++;
398  }
399  }
400 }
401 
402 
405 {
406  // Get overall size
407  label nPoints = 0;
408 
409  forAll(subGeom_, surfI)
410  {
411  nPoints += subGeom_[surfI].points()().size();
412  }
413 
415  pointField& pts = tPts.ref();
416 
417  // Append individual coordinates
418  nPoints = 0;
419 
420  forAll(subGeom_, surfI)
421  {
422  const pointField subCoords(subGeom_[surfI].points());
423 
424  forAll(subCoords, i)
425  {
426  pts[nPoints++] = transform_[surfI].globalPosition
427  (
429  (
430  subCoords[i],
431  scale_[surfI]
432  )
433  );
434  }
435  }
436 
437  return tPts;
438 }
439 
440 
441 void Foam::searchableSurfaces::collection::findNearest
442 (
443  const pointField& samples,
444  const scalarField& nearestDistSqr,
445  List<pointIndexHit>& nearestInfo
446 ) const
447 {
448  // How to scale distance?
449  scalarField minDistSqr(nearestDistSqr);
450 
451  labelList nearestSurf;
452  findNearest
453  (
454  samples,
455  minDistSqr,
456  nearestInfo,
457  nearestSurf
458  );
459 }
460 
461 
463 (
464  const pointField& start,
465  const pointField& end,
466  List<pointIndexHit>& info
467 ) const
468 {
469  info.setSize(start.size());
470  info = pointIndexHit();
471 
472  // Current nearest (to start) intersection
473  pointField nearest(end);
474 
475  List<pointIndexHit> hitInfo(start.size());
476 
477  forAll(subGeom_, surfI)
478  {
479  // Starting point
481  (
482  transform_[surfI].localPosition
483  (
484  start
485  ),
486  scale_[surfI]
487  );
488 
489  // Current best end point
491  (
492  transform_[surfI].localPosition
493  (
494  nearest
495  ),
496  scale_[surfI]
497  );
498 
499  subGeom_[surfI].findLine(e0, e1, hitInfo);
500 
501  forAll(hitInfo, pointi)
502  {
503  if (hitInfo[pointi].hit())
504  {
505  // Transform back to global coordinate sys.
506  nearest[pointi] = transform_[surfI].globalPosition
507  (
509  (
510  hitInfo[pointi].rawPoint(),
511  scale_[surfI]
512  )
513  );
514  info[pointi] = hitInfo[pointi];
515  info[pointi].rawPoint() = nearest[pointi];
516  info[pointi].setIndex
517  (
518  hitInfo[pointi].index()
519  + indexOffset_[surfI]
520  );
521  }
522  }
523  }
524 
525 
526  // Debug check
527  if (false)
528  {
529  forAll(info, pointi)
530  {
531  if (info[pointi].hit())
532  {
533  vector n(end[pointi] - start[pointi]);
534  scalar magN = mag(n);
535 
536  if (magN > small)
537  {
538  n /= mag(n);
539 
540  scalar s = ((info[pointi].rawPoint()-start[pointi])&n);
541 
542  if (s < 0 || s > 1)
543  {
545  << "point:" << info[pointi]
546  << " s:" << s
547  << " outside vector "
548  << " start:" << start[pointi]
549  << " end:" << end[pointi]
550  << abort(FatalError);
551  }
552  }
553  }
554  }
555  }
556 }
557 
558 
560 (
561  const pointField& start,
562  const pointField& end,
563  List<pointIndexHit>& info
564 ) const
565 {
566  // To be done ...
567  findLine(start, end, info);
568 }
569 
570 
572 (
573  const pointField& start,
574  const pointField& end,
576 ) const
577 {
578  // To be done. Assume for now only one intersection.
579  List<pointIndexHit> nearestInfo;
580  findLine(start, end, nearestInfo);
581 
582  info.setSize(start.size());
583  forAll(info, pointi)
584  {
585  if (nearestInfo[pointi].hit())
586  {
587  info[pointi].setSize(1);
588  info[pointi][0] = nearestInfo[pointi];
589  }
590  else
591  {
592  info[pointi].clear();
593  }
594  }
595 }
596 
597 
599 (
600  const List<pointIndexHit>& info,
601  labelList& region
602 ) const
603 {
604  if (subGeom_.size() == 0)
605  {}
606  else if (subGeom_.size() == 1)
607  {
608  if (mergeSubRegions_)
609  {
610  region.setSize(info.size());
611  region = regionOffset_[0];
612  }
613  else
614  {
615  subGeom_[0].getRegion(info, region);
616  }
617  }
618  else
619  {
620  // Multiple surfaces. Sort by surface.
621 
622  // Per surface the hit
623  List<List<pointIndexHit>> surfInfo;
624  // Per surface the original position
625  List<List<label>> infoMap;
626  sortHits(info, surfInfo, infoMap);
627 
628  region.setSize(info.size());
629  region = -1;
630 
631  // Do region tests
632 
633  if (mergeSubRegions_)
634  {
635  // Actually no need for surfInfo. Just take region for surface.
636  forAll(infoMap, surfI)
637  {
638  const labelList& map = infoMap[surfI];
639  forAll(map, i)
640  {
641  region[map[i]] = regionOffset_[surfI];
642  }
643  }
644  }
645  else
646  {
647  forAll(infoMap, surfI)
648  {
649  labelList surfRegion;
650  subGeom_[surfI].getRegion(surfInfo[surfI], surfRegion);
651 
652  const labelList& map = infoMap[surfI];
653  forAll(map, i)
654  {
655  region[map[i]] = regionOffset_[surfI] + surfRegion[i];
656  }
657  }
658  }
659  }
660 }
661 
662 
664 (
665  const List<pointIndexHit>& info,
666  vectorField& normal
667 ) const
668 {
669  if (subGeom_.size() == 0)
670  {}
671  else if (subGeom_.size() == 1)
672  {
673  subGeom_[0].getNormal(info, normal);
674  }
675  else
676  {
677  // Multiple surfaces. Sort by surface.
678 
679  // Per surface the hit
680  List<List<pointIndexHit>> surfInfo;
681  // Per surface the original position
682  List<List<label>> infoMap;
683  sortHits(info, surfInfo, infoMap);
684 
685  normal.setSize(info.size());
686 
687  // Do region tests
688  forAll(surfInfo, surfI)
689  {
690  vectorField surfNormal;
691  subGeom_[surfI].getNormal(surfInfo[surfI], surfNormal);
692 
693  // Transform back to global coordinate sys.
694  surfNormal = transform_[surfI].globalVector(surfNormal);
695 
696  const labelList& map = infoMap[surfI];
697  forAll(map, i)
698  {
699  normal[map[i]] = surfNormal[i];
700  }
701  }
702  }
703 }
704 
705 
707 (
708  const pointField& points,
709  List<volumeType>& volType
710 ) const
711 {
713  << "Volume type not supported for collection."
714  << exit(FatalError);
715 }
716 
717 
719 (
720  const List<treeBoundBox>& bbs,
721  const bool keepNonLocal,
723  autoPtr<distributionMap>& pointMap
724 )
725 {
726  forAll(subGeom_, surfI)
727  {
728  // Note:Transform the bounding boxes? Something like
729  // pointField bbPoints =
730  // cmptDivide
731  // (
732  // transform_[surfI].localPosition
733  // (
734  // bbs[i].points()
735  // ),
736  // scale_[surfI]
737  // );
738  // treeBoundBox newBb(bbPoints);
739 
740  // Note: what to do with faceMap, pointMap from multiple surfaces?
741  subGeom_[surfI].distribute
742  (
743  bbs,
744  keepNonLocal,
745  faceMap,
746  pointMap
747  );
748  }
749 }
750 
751 
753 {
754  forAll(subGeom_, surfI)
755  {
756  subGeom_[surfI].setField
757  (
758  static_cast<const labelList&>
759  (
761  (
762  values,
763  subGeom_[surfI].size(),
764  indexOffset_[surfI]
765  )
766  )
767  );
768  }
769 }
770 
771 
773 (
774  const List<pointIndexHit>& info,
775  labelList& values
776 ) const
777 {
778  if (subGeom_.size() == 0)
779  {}
780  else if (subGeom_.size() == 1)
781  {
782  subGeom_[0].getField(info, values);
783  }
784  else
785  {
786  // Multiple surfaces. Sort by surface.
787 
788  // Per surface the hit
789  List<List<pointIndexHit>> surfInfo;
790  // Per surface the original position
791  List<List<label>> infoMap;
792  sortHits(info, surfInfo, infoMap);
793 
794  // Do surface tests
795  forAll(surfInfo, surfI)
796  {
797  labelList surfValues;
798  subGeom_[surfI].getField(surfInfo[surfI], surfValues);
799 
800  if (surfValues.size())
801  {
802  // Size values only when we have a surface that supports it.
803  values.setSize(info.size());
804 
805  const labelList& map = infoMap[surfI];
806  forAll(map, i)
807  {
808  values[map[i]] = surfValues[i];
809  }
810  }
811  }
812  }
813 }
814 
815 
816 // ************************************************************************* //
Various functions to operate on Lists.
label n
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:476
Macros for easy insertion into run-time selection tables.
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Definition: DynamicListI.H:296
DynamicList< T, SizeInc, SizeMult, SizeDiv > & shrink()
Shrink the allocated space to the number of elements used.
Definition: DynamicListI.H:252
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
const objectRegistry & db() const
Return the local objectRegistry.
Definition: IOobject.C:309
const word & name() const
Return name.
Definition: IOobject.H:307
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
void clear()
Clear the list, i.e. set size to zero.
Definition: ListI.H:125
void setSize(const label)
Reset size of List.
Definition: List.C:281
A List obtained as a section of another List.
Definition: SubList.H:56
static const Form max
Definition: VectorSpace.H:120
static const Form min
Definition: VectorSpace.H:121
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
A bounding box defined in terms of the points at its extremities.
Definition: boundBox.H:59
const point & min() const
Minimum point defining the bounding box.
Definition: boundBoxI.H:60
const point & max() const
Maximum point defining the bounding box.
Definition: boundBoxI.H:66
static autoPtr< coordinateSystem > New(const objectRegistry &obr, const dictionary &dict)
Select constructed from dictionary and objectRegistry.
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:740
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:849
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type and name.
Base class of (analytical or triangulated) surface. Encapsulates all the search routines....
const boundBox & bounds() const
Return const reference to boundBox.
Makes a collection of surface geometries by copying from an existing defined surface geometry....
virtual void getVolumeType(const pointField &, List< volumeType > &) const
Determine type (inside/outside/mixed) for point. unknown if.
virtual void findLine(const pointField &start, const pointField &end, List< pointIndexHit > &) const
Find first intersection on segment from start to end.
virtual void findLineAny(const pointField &start, const pointField &end, List< pointIndexHit > &) const
Return any intersection on segment from start to end.
virtual label size() const
Range of local indices that can be returned.
virtual void boundingSpheres(pointField &centres, scalarField &radiusSqr) const
Get bounding spheres (centre and radius squared), one per element.
virtual void setField(const labelList &values)
WIP. Store element-wise field.
virtual void getField(const List< pointIndexHit > &, labelList &) const
WIP. From a set of hits (points and.
virtual void getNormal(const List< pointIndexHit > &, vectorField &normal) const
From a set of points and indices get the normal.
collection(const IOobject &io, const dictionary &dict)
Construct from dictionary (used by searchableSurface)
virtual void getRegion(const List< pointIndexHit > &, labelList &region) const
From a set of points and indices get the region.
virtual const wordList & regions() const
Names of regions.
virtual tmp< pointField > coordinates() const
Get representative set of element coordinates.
virtual void findLineAll(const pointField &start, const pointField &end, List< List< pointIndexHit >> &) const
Get all intersections in order from start to end.
virtual void distribute(const List< treeBoundBox > &, const bool keepNonLocal, autoPtr< distributionMap > &faceMap, autoPtr< distributionMap > &pointMap)
Set bounds of surface. Bounds currently set as list of.
virtual tmp< pointField > points() const
Get the points that define the surface.
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:197
A class for handling words, derived from string.
Definition: word.H:62
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
volScalarField scalarField(fieldObject, mesh)
const pointField & points
label nPoints
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(lagrangian::Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(lagrangian::Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(lagrangian::Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.name(), lagrangian::cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
addToRunTimeSelectionTable(searchableSurface, box, dictionary)
addBackwardCompatibleToRunTimeSelectionTable(searchableSurface, box, dictionary, searchableBox, "searchableBox")
barycentric coordinates(const polyMesh &mesh, const point &position, const label celli, const label facei, const label faceTrii, const scalar stepFraction)
Return the coordinates given the position and tet topology.
Definition: tracking.C:1259
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
void cmptMultiply(LagrangianPatchField< Type > &f, const LagrangianPatchField< Type > &f1, const LagrangianPatchField< Type > &f2)
List< label > labelList
A List of labels.
Definition: labelList.H:56
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
PointIndexHit< point > pointIndexHit
Definition: pointIndexHit.H:42
void cmptDivide(LagrangianPatchField< Type > &f, const LagrangianPatchField< Type > &f1, const LagrangianPatchField< Type > &f2)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:258
label findLower(const ListType &, typename ListType::const_reference, const label stary, const BinaryOp &bop)
Find last element < given value in sorted list and return index,.
errorManip< error > abort(error &err)
Definition: errorManip.H:131
const dimensionSet dimless
messageStream Info
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
void mag(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:57
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
void cmptMax(FieldField< Field, typename FieldField< Field, Type >::cmptType > &cf, const FieldField< Field, Type > &f)
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
error FatalError
void magSqr(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
dictionary dict
scalarField samples(nIntervals, 0)