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