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-2020 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"
33#include "treeDataCell.H"
34#include "treeDataFace.H"
35
36// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37
38namespace Foam
39{
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().faceLabels()[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 =
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 = magSqr(points[pointi] - sample);
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 forAll(indices, i)
165 {
166 label pointi = indices[i];
167
168 scalar distSqr = magSqr(points[pointi] - sample);
169
170 if (distSqr < nearestDistSqr)
171 {
172 nearestDistSqr = distSqr;
173 nearestI = pointi;
174 nearer = true;
175 }
176 }
177
178 return nearer;
179}
180
181
182// tree based searching
183Foam::label Foam::meshSearch::findNearestCellTree(const point& location) const
184{
185 const indexedOctree<treeDataCell>& tree = cellTree();
186
187 pointIndexHit info = tree.findNearest
188 (
189 location,
190 magSqr(tree.bb().max()-tree.bb().min())
191 );
192
193 if (!info.hit())
194 {
195 info = tree.findNearest(location, Foam::sqr(GREAT));
196 }
197 return info.index();
198}
199
200
201Foam::label Foam::meshSearch::findNearestCellLinear(const point& location) const
202{
203 const vectorField& centres = mesh_.cellCentres();
204
205 label nearestIndex = 0;
206 scalar minProximity = magSqr(centres[nearestIndex] - location);
207
208 findNearer
209 (
210 location,
211 centres,
212 nearestIndex,
213 minProximity
214 );
215
216 return nearestIndex;
217}
218
219
220Foam::label Foam::meshSearch::findNearestCellWalk
221(
222 const point& location,
223 const label seedCelli
224) const
225{
226 if (seedCelli < 0)
227 {
229 << "illegal seedCell:" << seedCelli << exit(FatalError);
230 }
231
232 // Walk in direction of face that decreases distance
233
234 label curCelli = seedCelli;
235 scalar distanceSqr = magSqr(mesh_.cellCentres()[curCelli] - location);
236
237 bool closer;
238
239 do
240 {
241 // Try neighbours of curCelli
242 closer = findNearer
243 (
244 location,
245 mesh_.cellCentres(),
246 mesh_.cellCells()[curCelli],
247 curCelli,
248 distanceSqr
249 );
250 } while (closer);
251
252 return curCelli;
253}
254
255
256Foam::label Foam::meshSearch::findNearestFaceTree(const point& location) const
257{
258 // Search nearest cell centre.
259 const indexedOctree<treeDataCell>& tree = cellTree();
260
261 // Search with decent span
262 pointIndexHit info = tree.findNearest
263 (
264 location,
265 magSqr(tree.bb().max()-tree.bb().min())
266 );
267
268 if (!info.hit())
269 {
270 // Search with desperate span
271 info = tree.findNearest(location, Foam::sqr(GREAT));
272 }
273
274
275 // Now check any of the faces of the nearest cell
276 const vectorField& centres = mesh_.faceCentres();
277 const cell& ownFaces = mesh_.cells()[info.index()];
278
279 label nearestFacei = ownFaces[0];
280 scalar minProximity = magSqr(centres[nearestFacei] - location);
281
282 findNearer
283 (
284 location,
285 centres,
286 ownFaces,
287 nearestFacei,
288 minProximity
289 );
290
291 return nearestFacei;
292}
293
294
295Foam::label Foam::meshSearch::findNearestFaceLinear(const point& location) const
296{
297 const vectorField& centres = mesh_.faceCentres();
298
299 label nearestFacei = 0;
300 scalar minProximity = magSqr(centres[nearestFacei] - location);
301
302 findNearer
303 (
304 location,
305 centres,
306 nearestFacei,
307 minProximity
308 );
309
310 return nearestFacei;
311}
312
313
314Foam::label Foam::meshSearch::findNearestFaceWalk
315(
316 const point& location,
317 const label seedFacei
318) const
319{
320 if (seedFacei < 0)
321 {
323 << "illegal seedFace:" << seedFacei << exit(FatalError);
324 }
325
326 const vectorField& centres = mesh_.faceCentres();
327
328
329 // Walk in direction of face that decreases distance
330
331 label curFacei = seedFacei;
332 scalar distanceSqr = magSqr(centres[curFacei] - location);
333
334 while (true)
335 {
336 label betterFacei = curFacei;
337
338 findNearer
339 (
340 location,
341 centres,
342 mesh_.cells()[mesh_.faceOwner()[curFacei]],
343 betterFacei,
344 distanceSqr
345 );
346
347 if (mesh_.isInternalFace(curFacei))
348 {
349 findNearer
350 (
351 location,
352 centres,
353 mesh_.cells()[mesh_.faceNeighbour()[curFacei]],
354 betterFacei,
355 distanceSqr
356 );
357 }
358
359 if (betterFacei == curFacei)
360 {
361 break;
362 }
363
364 curFacei = betterFacei;
365 }
366
367 return curFacei;
368}
369
370
371Foam::label Foam::meshSearch::findCellLinear(const point& location) const
372{
373 bool cellFound = false;
374 label n = 0;
375
376 label celli = -1;
377
378 while ((!cellFound) && (n < mesh_.nCells()))
379 {
380 if (mesh_.pointInCell(location, n, cellDecompMode_))
381 {
382 cellFound = true;
383 celli = n;
384 }
385 else
386 {
387 n++;
388 }
389 }
390 if (cellFound)
391 {
392 return celli;
393 }
394
395 return -1;
396}
397
398
399Foam::label Foam::meshSearch::findCellWalk
400(
401 const point& location,
402 const label seedCelli
403) const
404{
405 if (seedCelli < 0)
406 {
408 << "illegal seedCell:" << seedCelli << exit(FatalError);
409 }
410
411 if (mesh_.pointInCell(location, seedCelli, cellDecompMode_))
412 {
413 return seedCelli;
414 }
415
416 // Walk in direction of face that decreases distance
417 label curCelli = seedCelli;
418 scalar nearestDistSqr = magSqr(mesh_.cellCentres()[curCelli] - location);
419
420 while(true)
421 {
422 // Try neighbours of curCelli
423
424 const cell& cFaces = mesh_.cells()[curCelli];
425
426 label nearestCelli = -1;
427
428 forAll(cFaces, i)
429 {
430 label facei = cFaces[i];
431
432 if (mesh_.isInternalFace(facei))
433 {
434 label celli = mesh_.faceOwner()[facei];
435 if (celli == curCelli)
436 {
437 celli = mesh_.faceNeighbour()[facei];
438 }
439
440 // Check if this is the correct cell
441 if (mesh_.pointInCell(location, celli, cellDecompMode_))
442 {
443 return celli;
444 }
445
446 // Also calculate the nearest cell
447 scalar distSqr = magSqr(mesh_.cellCentres()[celli] - location);
448
449 if (distSqr < nearestDistSqr)
450 {
451 nearestDistSqr = distSqr;
452 nearestCelli = celli;
453 }
454 }
455 }
456
457 if (nearestCelli == -1)
458 {
459 return -1;
460 }
461
462 // Continue with the nearest cell
463 curCelli = nearestCelli;
464 }
465
466 return -1;
467}
468
469
470Foam::label Foam::meshSearch::findNearestBoundaryFaceWalk
471(
472 const point& location,
473 const label seedFacei
474) const
475{
476 if (seedFacei < 0)
477 {
479 << "illegal seedFace:" << seedFacei << exit(FatalError);
480 }
481
482 // Start off from seedFacei
483
484 label curFacei = seedFacei;
485
486 const face& f = mesh_.faces()[curFacei];
487
488 scalar minDist = f.nearestPoint
489 (
490 location,
491 mesh_.points()
492 ).distance();
493
494 bool closer;
495
496 do
497 {
498 closer = false;
499
500 // Search through all neighbouring boundary faces by going
501 // across edges
502
503 label lastFacei = curFacei;
504
505 const labelList& myEdges = mesh_.faceEdges()[curFacei];
506
507 forAll(myEdges, myEdgeI)
508 {
509 const labelList& neighbours = mesh_.edgeFaces()[myEdges[myEdgeI]];
510
511 // Check any face which uses edge, is boundary face and
512 // is not curFacei itself.
513
514 forAll(neighbours, nI)
515 {
516 label facei = neighbours[nI];
517
518 if
519 (
520 (facei >= mesh_.nInternalFaces())
521 && (facei != lastFacei)
522 )
523 {
524 const face& f = mesh_.faces()[facei];
525
526 pointHit curHit = f.nearestPoint
527 (
528 location,
529 mesh_.points()
530 );
531
532 // If the face is closer, reset current face and distance
533 if (curHit.distance() < minDist)
534 {
535 minDist = curHit.distance();
536 curFacei = facei;
537 closer = true; // a closer neighbour has been found
538 }
539 }
540 }
541 }
542 } while (closer);
543
544 return curFacei;
545}
546
547
548// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
549
551(
552 const polyMesh& mesh,
553 const polyMesh::cellDecomposition cellDecompMode
554)
555:
556 mesh_(mesh),
557 cellDecompMode_(cellDecompMode)
558{
559 if
560 (
561 cellDecompMode_ == polyMesh::FACE_DIAG_TRIS
562 || cellDecompMode_ == polyMesh::CELL_TETS
563 )
564 {
565 // Force construction of face diagonals
566 (void)mesh.tetBasePtIs();
567 }
568}
569
570
572(
573 const polyMesh& mesh,
574 const treeBoundBox& bb,
575 const polyMesh::cellDecomposition cellDecompMode
576)
577:
578 mesh_(mesh),
579 cellDecompMode_(cellDecompMode)
580{
581 overallBbPtr_.reset(new treeBoundBox(bb));
582
583 if
584 (
585 cellDecompMode_ == polyMesh::FACE_DIAG_TRIS
586 || cellDecompMode_ == polyMesh::CELL_TETS
587 )
588 {
589 // Force construction of face diagonals
590 (void)mesh.tetBasePtIs();
591 }
592}
593
594
595// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
596
598{
599 clearOut();
600}
601
602
603// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
604
605const Foam::treeBoundBox& Foam::meshSearch::dataBoundBox() const
606{
607 if (!overallBbPtr_)
608 {
609 Random rndGen(261782);
610 overallBbPtr_.reset
611 (
612 new treeBoundBox(mesh_.points())
613 );
614
615 treeBoundBox& overallBb = overallBbPtr_();
616
617 // Extend slightly and make 3D
618 overallBb = overallBb.extend(rndGen, 1e-4);
619 overallBb.min() -= point::uniform(ROOTVSMALL);
620 overallBb.max() += point::uniform(ROOTVSMALL);
621 }
622
623 return *overallBbPtr_;
624}
625
626
629{
630 if (!boundaryTreePtr_)
631 {
632 // All boundary faces (not just walls)
633 labelList bndFaces
634 (
635 identity(mesh_.nBoundaryFaces(), mesh_.nInternalFaces())
636 );
637
638 boundaryTreePtr_.reset
639 (
641 (
642 treeDataFace // all information needed to search faces
643 (
644 false, // do not cache bb
645 mesh_,
646 bndFaces // boundary faces only
647 ),
648 dataBoundBox(), // overall search domain
649 8, // maxLevel
650 10, // leafsize
651 3.0 // duplicity
652 )
653 );
654 }
655
656 return *boundaryTreePtr_;
657}
658
659
660
663{
664 if (!nonCoupledBoundaryTreePtr_)
665 {
666 // All non-coupled boundary faces (not just walls)
667 const polyBoundaryMesh& patches = mesh_.boundaryMesh();
668
669 labelList bndFaces(mesh_.nBoundaryFaces());
670
671 label bndi = 0;
672 for (const polyPatch& pp : patches)
673 {
674 if (!pp.coupled())
675 {
676 forAll(pp, i)
677 {
678 bndFaces[bndi++] = pp.start()+i;
679 }
680 }
681 }
682 bndFaces.setSize(bndi);
683
684 nonCoupledBoundaryTreePtr_.reset
685 (
687 (
688 treeDataFace // all information needed to search faces
689 (
690 false, // do not cache bb
691 mesh_,
692 bndFaces // boundary faces only
693 ),
694 dataBoundBox(), // overall search domain
695 8, // maxLevel
696 10, // leafsize
697 3.0 // duplicity
698 )
699 );
700 }
701
702 return *nonCoupledBoundaryTreePtr_;
703}
704
705
708{
709 if (!cellTreePtr_)
710 {
711 cellTreePtr_.reset
712 (
714 (
716 (
717 false, // not cache bb
718 mesh_,
719 cellDecompMode_ // cell decomposition mode for inside tests
720 ),
721 dataBoundBox(), // overall search domain
722 8, // maxLevel
723 10, // leafsize
724 6.0 // duplicity
725 )
726 );
727 }
728
729 return *cellTreePtr_;
730}
731
732
734(
735 const point& location,
736 const label seedCelli,
737 const bool useTreeSearch
738) const
739{
740 if (seedCelli == -1)
741 {
742 if (useTreeSearch)
743 {
744 return findNearestCellTree(location);
745 }
746 else
747 {
748 return findNearestCellLinear(location);
749 }
750 }
751
752 return findNearestCellWalk(location, seedCelli);
753}
754
755
757(
758 const point& location,
759 const label seedFacei,
760 const bool useTreeSearch
761) const
762{
763 if (seedFacei == -1)
764 {
765 if (useTreeSearch)
766 {
767 return findNearestFaceTree(location);
768 }
769 else
770 {
771 return findNearestFaceLinear(location);
772 }
773 }
774
775 return findNearestFaceWalk(location, seedFacei);
776}
777
778
780(
781 const point& location,
782 const label seedCelli,
783 const bool useTreeSearch
784) const
785{
786 // Find the nearest cell centre to this location
787 if (seedCelli == -1)
788 {
789 if (useTreeSearch)
790 {
791 return cellTree().findInside(location);
792 }
793 else
794 {
795 return findCellLinear(location);
796 }
797 }
798
799 return findCellWalk(location, seedCelli);
800}
801
802
804(
805 const point& location,
806 const label seedFacei,
807 const bool useTreeSearch
808) const
809{
810 if (seedFacei == -1)
811 {
812 if (useTreeSearch)
813 {
814 const indexedOctree<treeDataFace>& tree = boundaryTree();
815
816 pointIndexHit info = boundaryTree().findNearest
817 (
818 location,
819 magSqr(tree.bb().max()-tree.bb().min())
820 );
821
822 if (!info.hit())
823 {
824 info = boundaryTree().findNearest
825 (
826 location,
827 Foam::sqr(GREAT)
828 );
829 }
830
831 return tree.shapes().faceLabels()[info.index()];
832 }
833 else
834 {
835 scalar minDist = GREAT;
836
837 label minFacei = -1;
838
839 for
840 (
841 label facei = mesh_.nInternalFaces();
842 facei < mesh_.nFaces();
843 facei++
844 )
845 {
846 const face& f = mesh_.faces()[facei];
847
848 pointHit curHit =
849 f.nearestPoint
850 (
851 location,
852 mesh_.points()
853 );
854
855 if (curHit.distance() < minDist)
856 {
857 minDist = curHit.distance();
858 minFacei = facei;
859 }
860 }
861 return minFacei;
862 }
863 }
864
865 return findNearestBoundaryFaceWalk(location, seedFacei);
866}
867
868
870(
871 const point& pStart,
872 const point& pEnd
873) const
874{
875 pointIndexHit curHit = boundaryTree().findLine(pStart, pEnd);
876
877 if (curHit.hit())
878 {
879 // Change index into octreeData into face label
880 curHit.setIndex(boundaryTree().shapes().faceLabels()[curHit.index()]);
881 }
882 return curHit;
883}
884
885
887(
888 const point& pStart,
889 const point& pEnd
890) const
891{
893 findUniqueIntersectOp iop(boundaryTree(), hits);
894
895 while (true)
896 {
897 // Get the next hit, or quit
898 pointIndexHit curHit = boundaryTree().findLine(pStart, pEnd, iop);
899 if (!curHit.hit()) break;
900
901 // Change index into octreeData into face label
902 curHit.setIndex(boundaryTree().shapes().faceLabels()[curHit.index()]);
903
904 hits.append(curHit);
905 }
906
907 hits.shrink();
908
909 return hits;
910}
911
912
914{
915 return (boundaryTree().getVolumeType(p) == volumeType::INSIDE);
916}
917
918
920{
921 boundaryTreePtr_.clear();
922 cellTreePtr_.clear();
923 overallBbPtr_.clear();
924}
925
926
928{
929 clearOut();
930}
931
932
933// ************************************************************************* //
label n
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:72
DynamicList< T, SizeMin > & shrink()
Shrink the allocated space to the number of elements used.
Definition: DynamicListI.H:434
void append(const T &val)
Copy append an element to the end of this list.
Definition: DynamicListI.H:503
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:77
void setSize(const label n)
Alias for resize()
Definition: List.H:218
Describes the interaction of a face and a point. It carries the info of a successful hit and (if succ...
Definition: PointHit.H:54
scalar distance() const noexcept
Return distance to hit.
Definition: PointHit.H:139
This class describes the interaction of (usually) a face and a point. It carries the info of a succes...
Definition: PointIndexHit.H:66
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:60
const point & min() const
Minimum describing the bounding box.
Definition: boundBoxI.H:91
const point & max() const
Maximum describing the bounding box.
Definition: boundBoxI.H:97
Particle-size distribution model wherein random samples are drawn from the doubly-truncated uniform p...
Definition: uniform.H:164
friend Ostream & operator(Ostream &, const faMatrix< Type > &)
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
const List< pointIndexHit > & hits_
Definition: meshSearch.C:54
findUniqueIntersectOp(const indexedOctree< treeDataFace > &tree, const List< pointIndexHit > &hits)
Construct from components.
Definition: meshSearch.C:60
const indexedOctree< treeDataFace > & tree_
Definition: meshSearch.C:52
Non-pointer based hierarchical recursive searching.
Definition: indexedOctree.H:74
const treeBoundBox & bb() const
Top bounding box.
const Type & shapes() const
Reference to shape.
Foam::intersection.
Definition: intersection.H:53
Various (local, not parallel) searches on polyMesh; uses (demand driven) octree to search.
Definition: meshSearch.H:61
label findNearestFace(const point &location, const label seedFacei=-1, const bool useTreeSearch=true) const
Definition: meshSearch.C:757
label findNearestBoundaryFace(const point &location, const label seedFacei=-1, const bool useTreeSearch=true) const
Find nearest boundary face.
Definition: meshSearch.C:804
const indexedOctree< treeDataFace > & boundaryTree() const
Demand-driven reference to octree holding all boundary faces.
Definition: meshSearch.C:628
const polyMesh & mesh() const
Definition: meshSearch.H:195
void correct()
Correct for mesh geom/topo changes.
Definition: meshSearch.C:927
~meshSearch()
Destructor.
Definition: meshSearch.C:597
bool isInside(const point &) const
Determine inside/outside status.
Definition: meshSearch.C:913
const indexedOctree< treeDataFace > & nonCoupledBoundaryTree() const
Definition: meshSearch.C:662
List< pointIndexHit > intersections(const point &pStart, const point &pEnd) const
Find all intersections of boundary within segment pStart .. pEnd.
Definition: meshSearch.C:887
label findCell(const point &location, const label seedCelli=-1, const bool useTreeSearch=true) const
Find cell containing location.
Definition: meshSearch.C:780
const indexedOctree< treeDataCell > & cellTree() const
Demand-driven reference to octree holding all cells.
Definition: meshSearch.C:707
static scalar tol_
Tolerance on linear dimensions.
Definition: meshSearch.H:163
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:734
void clearOut()
Delete all storage.
Definition: meshSearch.C:919
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
cellDecomposition
Enumeration defining the decomposition of the cell for.
Definition: polyMesh.H:101
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1108
const labelIOList & tetBasePtIs() const
Return the tetBasePtIs.
Definition: polyMesh.C:906
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:75
Cell-face mesh analysis engine.
Definition: primitiveMesh.H:79
const vectorField & faceAreas() const
Standard boundBox with extra functionality for use in octree.
Definition: treeBoundBox.H:89
treeBoundBox extend(Random &rndGen, const scalar s) const
Return slightly wider bounding box.
Encapsulation of data needed to search in/for cells. Used to find the cell containing a point (e....
Definition: treeDataCell.H:57
findIntersectOp(const indexedOctree< treeDataFace > &tree)
Definition: treeDataFace.C:141
Encapsulation of data needed to search for faces.
Definition: treeDataFace.H:60
@ INSIDE
A location inside the volume.
Definition: volumeType.H:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
volScalarField & p
const polyBoundaryMesh & patches
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
const pointField & points
Namespace for OpenFOAM.
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i)
Definition: labelList.C:38
List< label > labelList
A List of labels.
Definition: List.H:66
dimensionedSymmTensor sqr(const dimensionedVector &dv)
PointIndexHit< point > pointIndexHit
A PointIndexHit for 3D points.
Definition: pointIndexHit.H:46
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
vector point
Point is a vector.
Definition: point.H:43
PointHit< point > pointHit
A PointIndexHit for 3D points.
Definition: pointHit.H:44
Field< vector > vectorField
Specialisation of Field<T> for vector.
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition: HashSet.H:85
error FatalError
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
constexpr auto end(C &c) -> decltype(c.end())
Return iterator to the end of the container c.
Definition: stdFoam.H:158
labelList f(nPoints)
volScalarField & e
Definition: createFields.H:11
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
Random rndGen
Definition: createFields.H:23