Loading...
Searching...
No Matches
meshSearch.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2018 OpenFOAM Foundation
9 Copyright (C) 2015-2022 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "meshSearch.H"
30#include "polyMesh.H"
31#include "indexedOctree.H"
32#include "DynamicList.H"
34#include "treeDataFace.H"
35
36// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37
38namespace Foam
41
42 scalar meshSearch::tol_ = 1e-3;
43
44 // Intersection operation that checks previous successful hits so that they
45 // are not duplicated
47 :
49 {
50 public:
51
53
55
56 public:
57
58 //- Construct from components
60 (
62 const List<pointIndexHit>& hits
63 )
64 :
66 tree_(tree),
67 hits_(hits)
68 {}
69
70 //- Calculate intersection of triangle with ray. Sets result
71 // accordingly
72 bool operator()
73 (
74 const label index,
75 const point& start,
76 const point& end,
77 point& intersectionPoint
78 ) const
79 {
80 const primitiveMesh& mesh = tree_.shapes().mesh();
81
82 // Check whether this hit has already happened. If the new face
83 // index is the same as an existing hit then return no new hit. If
84 // the new face shares a point with an existing hit face and the
85 // line passes through both faces in the same direction, then this
86 // is also assumed to be a duplicate hit.
87 const label newFacei = tree_.shapes().objectIndex(index);
88 const face& newFace = mesh.faces()[newFacei];
89 const scalar newDot = mesh.faceAreas()[newFacei] & (end - start);
90 forAll(hits_, hiti)
91 {
92 const label oldFacei = hits_[hiti].index();
93 const face& oldFace = mesh.faces()[oldFacei];
94 const scalar oldDot =
95 mesh.faceAreas()[oldFacei] & (end - start);
96
97 if
98 (
99 hits_[hiti].index() == newFacei
100 || (
101 newDot*oldDot > 0
102 && (labelHashSet(newFace) & labelHashSet(oldFace)).size()
103 )
104 )
105 {
106 return false;
107 }
108 }
109
110 const bool hit =
111 treeDataFace::findIntersectOp::operator()
112 (
113 index,
114 start,
115 end,
116 intersectionPoint
117 );
118
119 return hit;
120 }
121 };
122}
123
124
125// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
126
127bool Foam::meshSearch::findNearer
128(
129 const point& sample,
130 const pointField& points,
131 label& nearestI,
132 scalar& nearestDistSqr
133)
134{
135 bool nearer = false;
136
137 forAll(points, pointi)
138 {
139 scalar distSqr = sample.distSqr(points[pointi]);
140
141 if (distSqr < nearestDistSqr)
142 {
143 nearestDistSqr = distSqr;
144 nearestI = pointi;
145 nearer = true;
146 }
147 }
148
149 return nearer;
150}
151
152
153bool Foam::meshSearch::findNearer
154(
155 const point& sample,
156 const pointField& points,
157 const labelList& indices,
158 label& nearestI,
159 scalar& nearestDistSqr
160)
161{
162 bool nearer = false;
163
164 for (const label pointi : indices)
165 {
166 scalar distSqr = sample.distSqr(points[pointi]);
167
168 if (distSqr < nearestDistSqr)
169 {
170 nearestDistSqr = distSqr;
171 nearestI = pointi;
172 nearer = true;
173 }
174 }
175
176 return nearer;
177}
178
179
180// tree based searching
181Foam::label Foam::meshSearch::findNearestCellTree(const point& location) const
182{
183 const indexedOctree<treeDataCell>& tree = cellTree();
184
185 pointIndexHit info = tree.findNearest(location, tree.bb().magSqr());
186
187 if (!info.hit())
188 {
189 info = tree.findNearest(location, Foam::sqr(GREAT));
190 }
191 return info.index();
192}
193
194
195Foam::label Foam::meshSearch::findNearestCellLinear(const point& location) const
196{
197 const vectorField& centres = mesh_.cellCentres();
198
199 label nearestIndex = 0;
200 scalar minProximity = location.distSqr(centres[nearestIndex]);
201
202 findNearer
203 (
204 location,
205 centres,
206 nearestIndex,
207 minProximity
208 );
209
210 return nearestIndex;
211}
212
213
214Foam::label Foam::meshSearch::findNearestCellWalk
215(
216 const point& location,
217 const label seedCelli
218) const
219{
220 if (seedCelli < 0)
221 {
223 << "illegal seedCell:" << seedCelli << exit(FatalError);
224 }
225
226 // Walk in direction of face that decreases distance
227
228 label curCelli = seedCelli;
229 scalar distanceSqr = location.distSqr(mesh_.cellCentres()[curCelli]);
230
231 bool closer;
232
233 do
234 {
235 // Try neighbours of curCelli
236 closer = findNearer
237 (
238 location,
239 mesh_.cellCentres(),
240 mesh_.cellCells()[curCelli],
241 curCelli,
242 distanceSqr
243 );
244 } while (closer);
245
246 return curCelli;
247}
248
249
250Foam::label Foam::meshSearch::findNearestFaceTree(const point& location) const
251{
252 // Search nearest cell centre.
253 const indexedOctree<treeDataCell>& tree = cellTree();
254
255 // Search with decent span
256 pointIndexHit info = tree.findNearest(location, tree.bb().magSqr());
257
258 if (!info.hit())
259 {
260 // Search with desperate span
261 info = tree.findNearest(location, Foam::sqr(GREAT));
262 }
263
264
265 // Now check any of the faces of the nearest cell
266 const vectorField& centres = mesh_.faceCentres();
267 const cell& ownFaces = mesh_.cells()[info.index()];
268
269 label nearestFacei = ownFaces[0];
270 scalar minProximity = location.distSqr(centres[nearestFacei]);
271
272 findNearer
273 (
274 location,
275 centres,
276 ownFaces,
277 nearestFacei,
278 minProximity
279 );
280
281 return nearestFacei;
282}
283
284
285Foam::label Foam::meshSearch::findNearestFaceLinear(const point& location) const
286{
287 const vectorField& centres = mesh_.faceCentres();
288
289 label nearestFacei = 0;
290 scalar minProximity = location.distSqr(centres[nearestFacei]);
291
292 findNearer
293 (
294 location,
295 centres,
296 nearestFacei,
297 minProximity
298 );
299
300 return nearestFacei;
301}
302
303
304Foam::label Foam::meshSearch::findNearestFaceWalk
305(
306 const point& location,
307 const label seedFacei
308) const
309{
310 if (seedFacei < 0)
311 {
313 << "illegal seedFace:" << seedFacei << exit(FatalError);
314 }
315
316 const vectorField& centres = mesh_.faceCentres();
317
318
319 // Walk in direction of face that decreases distance
320
321 label curFacei = seedFacei;
322 scalar distanceSqr = location.distSqr(centres[curFacei]);
323
324 while (true)
325 {
326 label betterFacei = curFacei;
327
328 findNearer
329 (
330 location,
331 centres,
332 mesh_.cells()[mesh_.faceOwner()[curFacei]],
333 betterFacei,
334 distanceSqr
335 );
336
337 if (mesh_.isInternalFace(curFacei))
338 {
339 findNearer
340 (
341 location,
342 centres,
343 mesh_.cells()[mesh_.faceNeighbour()[curFacei]],
344 betterFacei,
345 distanceSqr
346 );
347 }
348
349 if (betterFacei == curFacei)
350 {
351 break;
352 }
353
354 curFacei = betterFacei;
355 }
356
357 return curFacei;
358}
359
360
361Foam::label Foam::meshSearch::findCellLinear(const point& location) const
362{
363 bool cellFound = false;
364 label n = 0;
365
366 label celli = -1;
367
368 while ((!cellFound) && (n < mesh_.nCells()))
369 {
370 if (mesh_.pointInCell(location, n, cellDecompMode_))
371 {
372 cellFound = true;
373 celli = n;
374 }
375 else
376 {
377 n++;
378 }
379 }
380 if (cellFound)
381 {
382 return celli;
383 }
384
385 return -1;
386}
387
388
389Foam::label Foam::meshSearch::findCellWalk
390(
391 const point& location,
392 const label seedCelli
393) const
394{
395 if (seedCelli < 0)
396 {
398 << "illegal seedCell:" << seedCelli << exit(FatalError);
399 }
400
401 if (mesh_.pointInCell(location, seedCelli, cellDecompMode_))
402 {
403 return seedCelli;
404 }
405
406 // Walk in direction of face that decreases distance
407 label curCelli = seedCelli;
408 scalar nearestDistSqr = location.distSqr(mesh_.cellCentres()[curCelli]);
409
410 while(true)
411 {
412 // Try neighbours of curCelli
413
414 const cell& cFaces = mesh_.cells()[curCelli];
415
416 label nearestCelli = -1;
417
418 forAll(cFaces, i)
419 {
420 label facei = cFaces[i];
421
422 if (mesh_.isInternalFace(facei))
423 {
424 label celli = mesh_.faceOwner()[facei];
425 if (celli == curCelli)
426 {
427 celli = mesh_.faceNeighbour()[facei];
428 }
429
430 // Check if this is the correct cell
431 if (mesh_.pointInCell(location, celli, cellDecompMode_))
432 {
433 return celli;
434 }
435
436 // Also calculate the nearest cell
437 scalar distSqr = location.distSqr(mesh_.cellCentres()[celli]);
438
439 if (distSqr < nearestDistSqr)
440 {
441 nearestDistSqr = distSqr;
442 nearestCelli = celli;
443 }
444 }
445 }
446
447 if (nearestCelli == -1)
448 {
449 return -1;
450 }
451
452 // Continue with the nearest cell
453 curCelli = nearestCelli;
454 }
455
456 return -1;
457}
458
459
460Foam::label Foam::meshSearch::findNearestBoundaryFaceWalk
461(
462 const point& location,
463 const label seedFacei
464) const
465{
466 if (seedFacei < 0)
467 {
469 << "illegal seedFace:" << seedFacei << exit(FatalError);
470 }
471
472 // Start off from seedFacei
473
474 label curFacei = seedFacei;
475
476 const face& f = mesh_.faces()[curFacei];
477
478 scalar minDist = f.nearestPoint
479 (
480 location,
481 mesh_.points()
482 ).distance();
483
484 bool closer;
485
486 do
487 {
488 closer = false;
489
490 // Search through all neighbouring boundary faces by going
491 // across edges
492
493 label lastFacei = curFacei;
494
495 const labelList& myEdges = mesh_.faceEdges()[curFacei];
496
497 forAll(myEdges, myEdgeI)
498 {
499 const labelList& neighbours = mesh_.edgeFaces()[myEdges[myEdgeI]];
500
501 // Check any face which uses edge, is boundary face and
502 // is not curFacei itself.
503
504 forAll(neighbours, nI)
505 {
506 label facei = neighbours[nI];
507
508 if
509 (
510 (facei >= mesh_.nInternalFaces())
511 && (facei != lastFacei)
512 )
513 {
514 const face& f = mesh_.faces()[facei];
515
516 pointHit curHit = f.nearestPoint
517 (
518 location,
519 mesh_.points()
520 );
521
522 // If the face is closer, reset current face and distance
523 if (curHit.distance() < minDist)
524 {
525 minDist = curHit.distance();
526 curFacei = facei;
527 closer = true; // a closer neighbour has been found
528 }
529 }
530 }
531 }
532 } while (closer);
533
534 return curFacei;
535}
536
537
538// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
539
540Foam::meshSearch::meshSearch
541(
542 const polyMesh& mesh,
543 const polyMesh::cellDecomposition cellDecompMode
544)
545:
546 mesh_(mesh),
547 cellDecompMode_(cellDecompMode)
548{
549 if
550 (
551 cellDecompMode_ == polyMesh::FACE_DIAG_TRIS
552 || cellDecompMode_ == polyMesh::CELL_TETS
553 )
554 {
555 // Force construction of face diagonals
556 (void)mesh.tetBasePtIs();
557 }
558}
560
561Foam::meshSearch::meshSearch
562(
563 const polyMesh& mesh,
564 const treeBoundBox& bb,
565 const polyMesh::cellDecomposition cellDecompMode
566)
567:
568 mesh_(mesh),
569 cellDecompMode_(cellDecompMode)
570{
571 overallBbPtr_.reset(new treeBoundBox(bb));
572
573 if
574 (
575 cellDecompMode_ == polyMesh::FACE_DIAG_TRIS
576 || cellDecompMode_ == polyMesh::CELL_TETS
577 )
578 {
579 // Force construction of face diagonals
580 (void)mesh.tetBasePtIs();
581 }
582}
583
584
585// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
586
588{
589 clearOut();
590}
591
592
593// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
594
595const Foam::treeBoundBox& Foam::meshSearch::dataBoundBox() const
596{
597 if (!overallBbPtr_)
598 {
599 Random rndGen(261782);
600 overallBbPtr_.reset
601 (
602 new treeBoundBox(mesh_.points())
603 );
604
605 treeBoundBox& overallBb = overallBbPtr_();
606
607 // Extend slightly and make 3D
608 overallBb.inflate(rndGen, 1e-4, ROOTVSMALL);
609 }
610
611 return *overallBbPtr_;
612}
613
617{
618 if (!boundaryTreePtr_)
619 {
620 // All boundary faces (not just walls)
621 const labelRange bndFaces
622 (
623 mesh_.nInternalFaces(),
624 mesh_.nBoundaryFaces()
625 );
626
627 boundaryTreePtr_.reset
628 (
630 (
631 treeDataFace(mesh_, bndFaces), // Boundary faces
632
633 dataBoundBox(), // overall search domain
634 8, // maxLevel
635 10, // leafsize
636 3.0 // duplicity
637 )
638 );
639 }
640
641 return *boundaryTreePtr_;
642}
643
644
648{
649 if (!nonCoupledBoundaryTreePtr_)
650 {
651 // All non-coupled boundary faces (not just walls)
652 const polyBoundaryMesh& patches = mesh_.boundaryMesh();
653
654 labelList bndFaces(mesh_.nBoundaryFaces());
655
656 label bndi = 0;
657 for (const polyPatch& pp : patches)
658 {
659 if (!pp.coupled())
660 {
661 forAll(pp, i)
662 {
663 bndFaces[bndi++] = pp.start()+i;
664 }
665 }
666 }
667 bndFaces.resize(bndi);
668
669 nonCoupledBoundaryTreePtr_.reset
670 (
672 (
673 treeDataFace(mesh_, bndFaces), // Boundary faces
674
675 dataBoundBox(), // overall search domain
676 8, // maxLevel
677 10, // leafsize
678 3.0 // duplicity
679 )
680 );
681 }
682
683 return *nonCoupledBoundaryTreePtr_;
684}
685
689{
690 if (!cellTreePtr_)
691 {
692 cellTreePtr_.reset
693 (
695 (
697 (
698 false, // not cache bb
699 mesh_,
700 cellDecompMode_ // cell decomposition mode for inside tests
701 ),
702 dataBoundBox(), // overall search domain
703 8, // maxLevel
704 10, // leafsize
705 6.0 // duplicity
706 )
707 );
708 }
709
710 return *cellTreePtr_;
711}
713
715(
716 const point& location,
717 const label seedCelli,
718 const bool useTreeSearch
719) const
720{
721 if (seedCelli == -1)
722 {
723 if (useTreeSearch)
724 {
725 return findNearestCellTree(location);
726 }
727 else
728 {
729 return findNearestCellLinear(location);
730 }
731 }
732
733 return findNearestCellWalk(location, seedCelli);
734}
736
738(
739 const point& location,
740 const label seedFacei,
741 const bool useTreeSearch
742) const
743{
744 if (seedFacei == -1)
745 {
746 if (useTreeSearch)
747 {
748 return findNearestFaceTree(location);
749 }
750 else
751 {
752 return findNearestFaceLinear(location);
753 }
754 }
755
756 return findNearestFaceWalk(location, seedFacei);
757}
759
761(
762 const point& location,
763 const label seedCelli,
764 const bool useTreeSearch
765) const
766{
767 // Find the nearest cell centre to this location
768 if (seedCelli == -1)
769 {
770 if (useTreeSearch)
771 {
772 return cellTree().findInside(location);
773 }
774 else
775 {
776 return findCellLinear(location);
777 }
778 }
779
780 return findCellWalk(location, seedCelli);
781}
783
785(
786 const point& location,
787 const label seedFacei,
788 const bool useTreeSearch
789) const
790{
791 if (seedFacei == -1)
792 {
793 if (useTreeSearch)
794 {
796
797 pointIndexHit info = boundaryTree().findNearest
798 (
799 location,
800 tree.bb().magSqr()
801 );
802
803 if (!info.hit())
804 {
805 info = boundaryTree().findNearest
806 (
807 location,
808 Foam::sqr(GREAT)
809 );
810 }
811
812 return tree.shapes().objectIndex(info.index());
813 }
814 else
815 {
816 scalar minDist = GREAT;
817
818 label minFacei = -1;
819
820 for
821 (
822 label facei = mesh_.nInternalFaces();
823 facei < mesh_.nFaces();
824 facei++
825 )
826 {
827 const face& f = mesh_.faces()[facei];
828
829 pointHit curHit =
830 f.nearestPoint
831 (
832 location,
833 mesh_.points()
834 );
835
836 if (curHit.distance() < minDist)
837 {
838 minDist = curHit.distance();
839 minFacei = facei;
840 }
841 }
842 return minFacei;
843 }
844 }
845
846 return findNearestBoundaryFaceWalk(location, seedFacei);
847}
849
851(
852 const point& pStart,
853 const point& pEnd
854) const
855{
856 pointIndexHit curHit = boundaryTree().findLine(pStart, pEnd);
857
858 if (curHit.hit())
859 {
860 // Change index into octreeData into face label
861 curHit.setIndex(boundaryTree().shapes().objectIndex(curHit.index()));
862 }
863 return curHit;
864}
866
868(
869 const point& pStart,
870 const point& pEnd
871) const
872{
875
876 while (true)
877 {
878 // Get the next hit, or quit
879 pointIndexHit curHit = boundaryTree().findLine(pStart, pEnd, iop);
880 if (!curHit.hit()) break;
881
882 // Change index into octreeData into face label
883 curHit.setIndex(boundaryTree().shapes().objectIndex(curHit.index()));
884
885 hits.append(curHit);
886 }
887
888 hits.shrink();
889
890 return hits;
891}
893
894bool Foam::meshSearch::isInside(const point& p) const
895{
896 return (boundaryTree().getVolumeType(p) == volumeType::INSIDE);
897}
899
901{
902 boundaryTreePtr_.clear();
903 cellTreePtr_.clear();
904 overallBbPtr_.clear();
905}
907
909{
910 clearOut();
911}
912
913
914// ************************************************************************* //
label n
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition DynamicList.H:68
void append(const T &val)
Copy append an element to the end of this list.
DynamicList< T, SizeMin > & shrink()
Calls shrink_to_fit() and returns a reference to the DynamicList.
Minimal example by using system/controlDict.functions:
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
void resize(const label len)
Adjust allocated size of list.
Definition ListI.H:153
scalar distance() const noexcept
Return distance to hit.
Definition pointHit.H:169
void setIndex(const label index) noexcept
Set the index.
label index() const noexcept
Return the hit index.
bool hit() const noexcept
Is there a hit?
Random number generator.
Definition Random.H:56
void inflate(const scalar factor)
Expand box by factor*mag(span) in all dimensions.
Definition boundBoxI.H:381
A cell is defined as a list of faces with extra functionality.
Definition cell.H:56
A face is a list of labels corresponding to mesh vertices.
Definition face.H:71
const List< pointIndexHit > & hits_
Definition meshSearch.C:47
findUniqueIntersectOp(const indexedOctree< treeDataFace > &tree, const List< pointIndexHit > &hits)
Construct from components.
Definition meshSearch.C:55
const indexedOctree< treeDataFace > & tree_
Definition meshSearch.C:45
Non-pointer based hierarchical recursive searching.
A range or interval of labels defined by a start and a size.
Definition labelRange.H:66
Various (local, not parallel) searches on polyMesh; uses (demand driven) octree to search.
Definition meshSearch.H:57
label findNearestFace(const point &location, const label seedFacei=-1, const bool useTreeSearch=true) const
Definition meshSearch.C:736
label findNearestBoundaryFace(const point &location, const label seedFacei=-1, const bool useTreeSearch=true) const
Find nearest boundary face.
Definition meshSearch.C:783
pointIndexHit intersection(const point &pStart, const point &pEnd) const
Find first intersection of boundary in segment [pStart, pEnd].
Definition meshSearch.C:849
const indexedOctree< treeDataFace > & boundaryTree() const
Demand-driven reference to octree holding all boundary faces.
Definition meshSearch.C:614
const polyMesh & mesh() const
Definition meshSearch.H:236
void correct()
Correct for mesh geom/topo changes.
Definition meshSearch.C:906
~meshSearch()
Destructor.
Definition meshSearch.C:585
bool isInside(const point &) const
Determine inside/outside status.
Definition meshSearch.C:892
const indexedOctree< treeDataFace > & nonCoupledBoundaryTree() const
Demand-driven reference to octree holding all non-coupled boundary faces.
Definition meshSearch.C:645
List< pointIndexHit > intersections(const point &pStart, const point &pEnd) const
Find all intersections of boundary within segment pStart .. pEnd.
Definition meshSearch.C:866
label findCell(const point &location, const label seedCelli=-1, const bool useTreeSearch=true) const
Find cell containing location.
Definition meshSearch.C:759
const indexedOctree< treeDataCell > & cellTree() const
Demand-driven reference to octree holding all cells.
Definition meshSearch.C:686
static scalar tol_
Tolerance on linear dimensions.
Definition meshSearch.H:196
label findNearestCell(const point &location, const label seedCelli=-1, const bool useTreeSearch=true) const
Find nearest cell in terms of cell centre.
Definition meshSearch.C:713
void clearOut()
Delete all storage.
Definition meshSearch.C:898
A polyBoundaryMesh is a polyPatch list with registered IO, a reference to the associated polyMesh,...
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
cellDecomposition
Enumeration defining the decomposition of the cell for.
Definition polyMesh.H:105
A patch is a list of labels that address the faces in the global face list.
Definition polyPatch.H:73
Cell-face mesh analysis engine.
Standard boundBox with extra functionality for use in octree.
Encapsulation of data needed to search in/for cells. Used to find the cell containing a point (e....
findIntersectOp(const indexedOctree< treeDataFace > &tree)
Encapsulation of data for searching on faces.
@ INSIDE
A location inside the volume.
Definition volumeType.H:65
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
volScalarField & p
const polyBoundaryMesh & patches
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
const pointField & points
Namespace for OpenFOAM.
PointHit< point > pointHit
A PointHit with a 3D point.
Definition pointHit.H:51
List< label > labelList
A List of labels.
Definition List.H:62
dimensionedSymmTensor sqr(const dimensionedVector &dv)
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition HashSet.H:85
Field< vector > vectorField
Specialisation of Field<T> for vector.
vector point
Point is a vector.
Definition point.H:37
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
vectorField pointField
pointField is a vectorField.
PointIndexHit< point > pointIndexHit
A PointIndexHit with a 3D point.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
labelList f(nPoints)
Tree tree(triangles.begin(), triangles.end())
volScalarField & e
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
Random rndGen