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 | Copyright (C) 2011-2016 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "meshSearch.H"
27 #include "polyMesh.H"
28 #include "indexedOctree.H"
29 #include "DynamicList.H"
30 #include "demandDrivenData.H"
31 #include "treeDataCell.H"
32 #include "treeDataFace.H"
33 
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38 defineTypeNameAndDebug(meshSearch, 0);
39 
40 scalar meshSearch::tol_ = 1e-3;
41 }
42 
43 
44 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
45 
46 bool Foam::meshSearch::findNearer
47 (
48  const point& sample,
49  const pointField& points,
50  label& nearestI,
51  scalar& nearestDistSqr
52 )
53 {
54  bool nearer = false;
55 
56  forAll(points, pointi)
57  {
58  scalar distSqr = magSqr(points[pointi] - sample);
59 
60  if (distSqr < nearestDistSqr)
61  {
62  nearestDistSqr = distSqr;
63  nearestI = pointi;
64  nearer = true;
65  }
66  }
67 
68  return nearer;
69 }
70 
71 
72 bool Foam::meshSearch::findNearer
73 (
74  const point& sample,
75  const pointField& points,
76  const labelList& indices,
77  label& nearestI,
78  scalar& nearestDistSqr
79 )
80 {
81  bool nearer = false;
82 
83  forAll(indices, i)
84  {
85  label pointi = indices[i];
86 
87  scalar distSqr = magSqr(points[pointi] - sample);
88 
89  if (distSqr < nearestDistSqr)
90  {
91  nearestDistSqr = distSqr;
92  nearestI = pointi;
93  nearer = true;
94  }
95  }
96 
97  return nearer;
98 }
99 
100 
101 // tree based searching
102 Foam::label Foam::meshSearch::findNearestCellTree(const point& location) const
103 {
104  const indexedOctree<treeDataCell>& tree = cellTree();
105 
106  pointIndexHit info = tree.findNearest
107  (
108  location,
109  magSqr(tree.bb().max()-tree.bb().min())
110  );
111 
112  if (!info.hit())
113  {
114  info = tree.findNearest(location, Foam::sqr(GREAT));
115  }
116  return info.index();
117 }
118 
119 
120 // linear searching
121 Foam::label Foam::meshSearch::findNearestCellLinear(const point& location) const
122 {
123  const vectorField& centres = mesh_.cellCentres();
124 
125  label nearestIndex = 0;
126  scalar minProximity = magSqr(centres[nearestIndex] - location);
127 
128  findNearer
129  (
130  location,
131  centres,
132  nearestIndex,
133  minProximity
134  );
135 
136  return nearestIndex;
137 }
138 
139 
140 // walking from seed
141 Foam::label Foam::meshSearch::findNearestCellWalk
142 (
143  const point& location,
144  const label seedCelli
145 ) const
146 {
147  if (seedCelli < 0)
148  {
150  << "illegal seedCell:" << seedCelli << exit(FatalError);
151  }
152 
153  // Walk in direction of face that decreases distance
154 
155  label curCelli = seedCelli;
156  scalar distanceSqr = magSqr(mesh_.cellCentres()[curCelli] - location);
157 
158  bool closer;
159 
160  do
161  {
162  // Try neighbours of curCelli
163  closer = findNearer
164  (
165  location,
166  mesh_.cellCentres(),
167  mesh_.cellCells()[curCelli],
168  curCelli,
169  distanceSqr
170  );
171  } while (closer);
172 
173  return curCelli;
174 }
175 
176 
177 // tree based searching
178 Foam::label Foam::meshSearch::findNearestFaceTree(const point& location) const
179 {
180  // Search nearest cell centre.
181  const indexedOctree<treeDataCell>& tree = cellTree();
182 
183  // Search with decent span
184  pointIndexHit info = tree.findNearest
185  (
186  location,
187  magSqr(tree.bb().max()-tree.bb().min())
188  );
189 
190  if (!info.hit())
191  {
192  // Search with desparate span
193  info = tree.findNearest(location, Foam::sqr(GREAT));
194  }
195 
196 
197  // Now check any of the faces of the nearest cell
198  const vectorField& centres = mesh_.faceCentres();
199  const cell& ownFaces = mesh_.cells()[info.index()];
200 
201  label nearestFacei = ownFaces[0];
202  scalar minProximity = magSqr(centres[nearestFacei] - location);
203 
204  findNearer
205  (
206  location,
207  centres,
208  ownFaces,
209  nearestFacei,
210  minProximity
211  );
212 
213  return nearestFacei;
214 }
215 
216 
217 // linear searching
218 Foam::label Foam::meshSearch::findNearestFaceLinear(const point& location) const
219 {
220  const vectorField& centres = mesh_.faceCentres();
221 
222  label nearestFacei = 0;
223  scalar minProximity = magSqr(centres[nearestFacei] - location);
224 
225  findNearer
226  (
227  location,
228  centres,
229  nearestFacei,
230  minProximity
231  );
232 
233  return nearestFacei;
234 }
235 
236 
237 // walking from seed
238 Foam::label Foam::meshSearch::findNearestFaceWalk
239 (
240  const point& location,
241  const label seedFacei
242 ) const
243 {
244  if (seedFacei < 0)
245  {
247  << "illegal seedFace:" << seedFacei << exit(FatalError);
248  }
249 
250  const vectorField& centres = mesh_.faceCentres();
251 
252 
253  // Walk in direction of face that decreases distance
254 
255  label curFacei = seedFacei;
256  scalar distanceSqr = magSqr(centres[curFacei] - location);
257 
258  while (true)
259  {
260  label betterFacei = curFacei;
261 
262  findNearer
263  (
264  location,
265  centres,
266  mesh_.cells()[mesh_.faceOwner()[curFacei]],
267  betterFacei,
268  distanceSqr
269  );
270 
271  if (mesh_.isInternalFace(curFacei))
272  {
273  findNearer
274  (
275  location,
276  centres,
277  mesh_.cells()[mesh_.faceNeighbour()[curFacei]],
278  betterFacei,
279  distanceSqr
280  );
281  }
282 
283  if (betterFacei == curFacei)
284  {
285  break;
286  }
287 
288  curFacei = betterFacei;
289  }
290 
291  return curFacei;
292 }
293 
294 
295 Foam::label Foam::meshSearch::findCellLinear(const point& location) const
296 {
297  bool cellFound = false;
298  label n = 0;
299 
300  label celli = -1;
301 
302  while ((!cellFound) && (n < mesh_.nCells()))
303  {
304  if (mesh_.pointInCell(location, n, cellDecompMode_))
305  {
306  cellFound = true;
307  celli = n;
308  }
309  else
310  {
311  n++;
312  }
313  }
314  if (cellFound)
315  {
316  return celli;
317  }
318  else
319  {
320  return -1;
321  }
322 }
323 
324 
325 // walking from seed
326 Foam::label Foam::meshSearch::findCellWalk
327 (
328  const point& location,
329  const label seedCelli
330 ) const
331 {
332  if (seedCelli < 0)
333  {
335  << "illegal seedCell:" << seedCelli << exit(FatalError);
336  }
337 
338  if (mesh_.pointInCell(location, seedCelli, cellDecompMode_))
339  {
340  return seedCelli;
341  }
342 
343  // Walk in direction of face that decreases distance
344  label curCelli = seedCelli;
345  scalar nearestDistSqr = magSqr(mesh_.cellCentres()[curCelli] - location);
346 
347  while(true)
348  {
349  // Try neighbours of curCelli
350 
351  const cell& cFaces = mesh_.cells()[curCelli];
352 
353  label nearestCelli = -1;
354 
355  forAll(cFaces, i)
356  {
357  label facei = cFaces[i];
358 
359  if (mesh_.isInternalFace(facei))
360  {
361  label celli = mesh_.faceOwner()[facei];
362  if (celli == curCelli)
363  {
364  celli = mesh_.faceNeighbour()[facei];
365  }
366 
367  // Check if this is the correct cell
368  if (mesh_.pointInCell(location, celli, cellDecompMode_))
369  {
370  return celli;
371  }
372 
373  // Also calculate the nearest cell
374  scalar distSqr = magSqr(mesh_.cellCentres()[celli] - location);
375 
376  if (distSqr < nearestDistSqr)
377  {
378  nearestDistSqr = distSqr;
379  nearestCelli = celli;
380  }
381  }
382  }
383 
384  if (nearestCelli == -1)
385  {
386  return -1;
387  }
388 
389  // Continue with the nearest cell
390  curCelli = nearestCelli;
391  }
392 
393  return -1;
394 }
395 
396 
397 Foam::label Foam::meshSearch::findNearestBoundaryFaceWalk
398 (
399  const point& location,
400  const label seedFacei
401 ) const
402 {
403  if (seedFacei < 0)
404  {
406  << "illegal seedFace:" << seedFacei << exit(FatalError);
407  }
408 
409  // Start off from seedFacei
410 
411  label curFacei = seedFacei;
412 
413  const face& f = mesh_.faces()[curFacei];
414 
415  scalar minDist = f.nearestPoint
416  (
417  location,
418  mesh_.points()
419  ).distance();
420 
421  bool closer;
422 
423  do
424  {
425  closer = false;
426 
427  // Search through all neighbouring boundary faces by going
428  // across edges
429 
430  label lastFacei = curFacei;
431 
432  const labelList& myEdges = mesh_.faceEdges()[curFacei];
433 
434  forAll(myEdges, myEdgeI)
435  {
436  const labelList& neighbours = mesh_.edgeFaces()[myEdges[myEdgeI]];
437 
438  // Check any face which uses edge, is boundary face and
439  // is not curFacei itself.
440 
441  forAll(neighbours, nI)
442  {
443  label facei = neighbours[nI];
444 
445  if
446  (
447  (facei >= mesh_.nInternalFaces())
448  && (facei != lastFacei)
449  )
450  {
451  const face& f = mesh_.faces()[facei];
452 
453  pointHit curHit = f.nearestPoint
454  (
455  location,
456  mesh_.points()
457  );
458 
459  // If the face is closer, reset current face and distance
460  if (curHit.distance() < minDist)
461  {
462  minDist = curHit.distance();
463  curFacei = facei;
464  closer = true; // a closer neighbour has been found
465  }
466  }
467  }
468  }
469  } while (closer);
470 
471  return curFacei;
472 }
473 
474 
475 Foam::vector Foam::meshSearch::offset
476 (
477  const point& bPoint,
478  const label bFacei,
479  const vector& dir
480 ) const
481 {
482  // Get the neighbouring cell
483  label ownerCelli = mesh_.faceOwner()[bFacei];
484 
485  const point& c = mesh_.cellCentres()[ownerCelli];
486 
487  // Typical dimension: distance from point on face to cell centre
488  scalar typDim = mag(c - bPoint);
489 
490  return tol_*typDim*dir;
491 }
492 
493 
494 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
495 
496 Foam::meshSearch::meshSearch
497 (
498  const polyMesh& mesh,
499  const polyMesh::cellDecomposition cellDecompMode
500 )
501 :
502  mesh_(mesh),
503  cellDecompMode_(cellDecompMode)
504 {
505  if (cellDecompMode_ == polyMesh::FACE_DIAG_TRIS)
506  {
507  // Force construction of face diagonals
508  (void)mesh.tetBasePtIs();
509  }
510 }
511 
512 
513 // Construct with a custom bounding box
514 Foam::meshSearch::meshSearch
515 (
516  const polyMesh& mesh,
517  const treeBoundBox& bb,
518  const polyMesh::cellDecomposition cellDecompMode
519 )
520 :
521  mesh_(mesh),
522  cellDecompMode_(cellDecompMode)
523 {
524  overallBbPtr_.reset(new treeBoundBox(bb));
525 
526  if (cellDecompMode_ == polyMesh::FACE_DIAG_TRIS)
527  {
528  // Force construction of face diagonals
529  (void)mesh.tetBasePtIs();
530  }
531 }
532 
533 
534 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
535 
537 {
538  clearOut();
539 }
540 
541 
542 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
543 
545  const
546 {
547  if (!boundaryTreePtr_.valid())
548  {
549  //
550  // Construct tree
551  //
552 
553  if (!overallBbPtr_.valid())
554  {
555  Random rndGen(261782);
556  overallBbPtr_.reset
557  (
558  new treeBoundBox(mesh_.points())
559  );
560 
561  treeBoundBox& overallBb = overallBbPtr_();
562  // Extend slightly and make 3D
563  overallBb = overallBb.extend(rndGen, 1e-4);
564  overallBb.min() -= point(ROOTVSMALL, ROOTVSMALL, ROOTVSMALL);
565  overallBb.max() += point(ROOTVSMALL, ROOTVSMALL, ROOTVSMALL);
566  }
567 
568  // all boundary faces (not just walls)
569  labelList bndFaces(mesh_.nFaces()-mesh_.nInternalFaces());
570  forAll(bndFaces, i)
571  {
572  bndFaces[i] = mesh_.nInternalFaces() + i;
573  }
574 
575  boundaryTreePtr_.reset
576  (
578  (
579  treeDataFace // all information needed to search faces
580  (
581  false, // do not cache bb
582  mesh_,
583  bndFaces // boundary faces only
584  ),
585  overallBbPtr_(), // overall search domain
586  8, // maxLevel
587  10, // leafsize
588  3.0 // duplicity
589  )
590  );
591  }
592 
593  return boundaryTreePtr_();
594 }
595 
596 
598 const
599 {
600  if (!cellTreePtr_.valid())
601  {
602  //
603  // Construct tree
604  //
605 
606  if (!overallBbPtr_.valid())
607  {
608  Random rndGen(261782);
609  overallBbPtr_.reset
610  (
611  new treeBoundBox(mesh_.points())
612  );
613 
614  treeBoundBox& overallBb = overallBbPtr_();
615  // Extend slightly and make 3D
616  overallBb = overallBb.extend(rndGen, 1e-4);
617  overallBb.min() -= point(ROOTVSMALL, ROOTVSMALL, ROOTVSMALL);
618  overallBb.max() += point(ROOTVSMALL, ROOTVSMALL, ROOTVSMALL);
619  }
620 
621  cellTreePtr_.reset
622  (
624  (
626  (
627  false, // not cache bb
628  mesh_,
629  cellDecompMode_ // cell decomposition mode for inside tests
630  ),
631  overallBbPtr_(),
632  8, // maxLevel
633  10, // leafsize
634  6.0 // duplicity
635  )
636  );
637  }
638 
639  return cellTreePtr_();
640 }
641 
642 
647 //bool Foam::meshSearch::pointInCell(const point& p, label celli) const
648 //{
649 // if (faceDecomp_)
650 // {
651 // const point& ctr = mesh_.cellCentres()[celli];
652 //
653 // vector dir(p - ctr);
654 // scalar magDir = mag(dir);
655 //
656 // // Check if any faces are hit by ray from cell centre to p.
657 // // If none -> p is in cell.
658 // const labelList& cFaces = mesh_.cells()[celli];
659 //
660 // // Make sure half_ray does not pick up any faces on the wrong
661 // // side of the ray.
662 // scalar oldTol = intersection::setPlanarTol(0.0);
663 //
664 // forAll(cFaces, i)
665 // {
666 // label facei = cFaces[i];
667 //
668 // pointHit inter = mesh_.faces()[facei].ray
669 // (
670 // ctr,
671 // dir,
672 // mesh_.points(),
673 // intersection::HALF_RAY,
674 // intersection::VECTOR
675 // );
676 //
677 // if (inter.hit())
678 // {
679 // scalar dist = inter.distance();
680 //
681 // if (dist < magDir)
682 // {
683 // // Valid hit. Hit face so point is not in cell.
684 // intersection::setPlanarTol(oldTol);
685 //
686 // return false;
687 // }
688 // }
689 // }
690 //
691 // intersection::setPlanarTol(oldTol);
692 //
693 // // No face inbetween point and cell centre so point is inside.
694 // return true;
695 // }
696 // else
697 // {
698 // const labelList& f = mesh_.cells()[celli];
699 // const labelList& owner = mesh_.faceOwner();
700 // const vectorField& cf = mesh_.faceCentres();
701 // const vectorField& Sf = mesh_.faceAreas();
702 //
703 // forAll(f, facei)
704 // {
705 // label nFace = f[facei];
706 // vector proj = p - cf[nFace];
707 // vector normal = Sf[nFace];
708 // if (owner[nFace] == celli)
709 // {
710 // if ((normal & proj) > 0)
711 // {
712 // return false;
713 // }
714 // }
715 // else
716 // {
717 // if ((normal & proj) < 0)
718 // {
719 // return false;
720 // }
721 // }
722 // }
723 //
724 // return true;
725 // }
726 //}
727 
728 
730 (
731  const point& location,
732  const label seedCelli,
733  const bool useTreeSearch
734 ) const
735 {
736  if (seedCelli == -1)
737  {
738  if (useTreeSearch)
739  {
740  return findNearestCellTree(location);
741  }
742  else
743  {
744  return findNearestCellLinear(location);
745  }
746  }
747  else
748  {
749  return findNearestCellWalk(location, seedCelli);
750  }
751 }
752 
753 
755 (
756  const point& location,
757  const label seedFacei,
758  const bool useTreeSearch
759 ) const
760 {
761  if (seedFacei == -1)
762  {
763  if (useTreeSearch)
764  {
765  return findNearestFaceTree(location);
766  }
767  else
768  {
769  return findNearestFaceLinear(location);
770  }
771  }
772  else
773  {
774  return findNearestFaceWalk(location, seedFacei);
775  }
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  else
799  {
800  return findCellWalk(location, seedCelli);
801  }
802 }
803 
804 
806 (
807  const point& location,
808  const label seedFacei,
809  const bool useTreeSearch
810 ) const
811 {
812  if (seedFacei == -1)
813  {
814  if (useTreeSearch)
815  {
816  const indexedOctree<treeDataFace>& tree = boundaryTree();
817 
818  pointIndexHit info = boundaryTree().findNearest
819  (
820  location,
821  magSqr(tree.bb().max()-tree.bb().min())
822  );
823 
824  if (!info.hit())
825  {
826  info = boundaryTree().findNearest
827  (
828  location,
829  Foam::sqr(GREAT)
830  );
831  }
832 
833  return tree.shapes().faceLabels()[info.index()];
834  }
835  else
836  {
837  scalar minDist = GREAT;
838 
839  label minFacei = -1;
840 
841  for
842  (
843  label facei = mesh_.nInternalFaces();
844  facei < mesh_.nFaces();
845  facei++
846  )
847  {
848  const face& f = mesh_.faces()[facei];
849 
850  pointHit curHit =
851  f.nearestPoint
852  (
853  location,
854  mesh_.points()
855  );
856 
857  if (curHit.distance() < minDist)
858  {
859  minDist = curHit.distance();
860  minFacei = facei;
861  }
862  }
863  return minFacei;
864  }
865  }
866  else
867  {
868  return findNearestBoundaryFaceWalk(location, seedFacei);
869  }
870 }
871 
872 
874 (
875  const point& pStart,
876  const point& pEnd
877 ) const
878 {
879  pointIndexHit curHit = boundaryTree().findLine(pStart, pEnd);
880 
881  if (curHit.hit())
882  {
883  // Change index into octreeData into face label
884  curHit.setIndex(boundaryTree().shapes().faceLabels()[curHit.index()]);
885  }
886  return curHit;
887 }
888 
889 
891 (
892  const point& pStart,
893  const point& pEnd
894 ) const
895 {
897 
898  vector edgeVec = pEnd - pStart;
899  edgeVec /= mag(edgeVec);
900 
901  point pt = pStart;
902 
903  pointIndexHit bHit;
904  do
905  {
906  bHit = intersection(pt, pEnd);
907 
908  if (bHit.hit())
909  {
910  hits.append(bHit);
911 
912  const vector& area = mesh_.faceAreas()[bHit.index()];
913 
914  scalar typDim = Foam::sqrt(mag(area));
915 
916  if ((mag(bHit.hitPoint() - pEnd)/typDim) < SMALL)
917  {
918  break;
919  }
920 
921  // Restart from hitPoint shifted a little bit in the direction
922  // of the destination
923 
924  pt =
925  bHit.hitPoint()
926  + offset(bHit.hitPoint(), bHit.index(), edgeVec);
927  }
928 
929  } while (bHit.hit());
930 
931 
932  hits.shrink();
933 
934  return hits;
935 }
936 
937 
939 {
940  return (boundaryTree().getVolumeType(p) == volumeType::INSIDE);
941 }
942 
943 
944 // Delete all storage
946 {
947  boundaryTreePtr_.clear();
948  cellTreePtr_.clear();
949  overallBbPtr_.clear();
950 }
951 
952 
954 {
955  clearOut();
956 }
957 
958 
959 // ************************************************************************* //
cachedRandom rndGen(label(0),-1)
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
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
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
List< pointIndexHit > intersections(const point &pStart, const point &pEnd) const
Find all intersections of boundary within segment pStart .. pEnd.
Definition: meshSearch.C:891
const double e
Elementary charge.
Definition: doubleFloat.H:78
const point & min() const
Minimum describing the bounding box.
Definition: boundBoxI.H:54
error FatalError
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:730
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
Encapsulation of data needed to search for faces.
Definition: treeDataFace.H:58
bool hit() const
Is there a hit.
static scalar tol_
Tolerance on linear dimensions.
Definition: meshSearch.H:163
dimensionedSymmTensor sqr(const dimensionedVector &dv)
void setIndex(const label index)
label findNearestFace(const point &location, const label seedFacei=-1, const bool useTreeSearch=true) const
Definition: meshSearch.C:755
treeBoundBox extend(Random &, const scalar s) const
Return slightly wider bounding box.
dimensionedScalar sqrt(const dimensionedScalar &ds)
const point & max() const
Maximum describing the bounding box.
Definition: boundBoxI.H:60
scalar distance(const vector &p1, const vector &p2)
Definition: curveTools.C:12
This class describes the interaction of (usually) a face and a point. It carries the info of a succes...
Definition: PointIndexHit.H:53
PointIndexHit< point > pointIndexHit
Definition: pointIndexHit.H:42
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
pointIndexHit intersection(const point &pStart, const point &pEnd) const
Find first intersection of boundary in segment [pStart, pEnd].
Definition: meshSearch.C:874
label findCell(const point &location, const label seedCelli=-1, const bool useTreeSearch=true) const
Find cell containing location.
Definition: meshSearch.C:780
cellDecomposition
Enumeration defining the decomposition of the cell for.
Definition: polyMesh.H:98
bool isInside(const point &) const
Determine inside/outside status.
Definition: meshSearch.C:938
Foam::intersection.
Definition: intersection.H:49
scalar distance() const
Return distance to hit.
Definition: PointHit.H:139
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
void correct()
Correct for mesh geom/topo changes.
Definition: meshSearch.C:953
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects...
Definition: DynamicList.H:56
const indexedOctree< treeDataFace > & boundaryTree() const
Get (demand driven) reference to octree holding all.
Definition: meshSearch.C:544
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Definition: DynamicListI.H:292
List< label > labelList
A List of labels.
Definition: labelList.H:56
Simple random number generator.
Definition: Random.H:49
Templated 3D Vector derived from VectorSpace adding construction from 3 components, element access using x(), y() and z() member functions and the inner-product (dot-product) and cross product operators.
Definition: Vector.H:57
dimensioned< scalar > magSqr(const dimensioned< Type > &)
DynamicList< T, SizeInc, SizeMult, SizeDiv > & shrink()
Shrink the allocated space to the number of elements used.
Definition: DynamicListI.H:240
const treeBoundBox & bb() const
Top bounding box.
defineTypeNameAndDebug(combustionModel, 0)
const labelList & tetBasePtIs() const
Return the tetBasePtIs.
Definition: polyMesh.C:818
pointHit nearestPoint(const point &p, const pointField &) const
Return nearest point to face.
Encapsulation of data needed to search in/for cells. Used to find the cell containing a point (e...
Definition: treeDataCell.H:54
void clearOut()
Delete all storage.
Definition: meshSearch.C:945
Template functions to aid in the implementation of demand driven data.
vector point
Point is a vector.
Definition: point.H:41
~meshSearch()
Destructor.
Definition: meshSearch.C:536
const Type & shapes() const
Reference to shape.
label index() const
Return index.
const Point & hitPoint() const
Return hit point.
Standard boundBox + extra functionality for use in octree.
Definition: treeBoundBox.H:87
dimensioned< scalar > mag(const dimensioned< Type > &)
label n
Field< vector > vectorField
Specialisation of Field<T> for vector.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
PointHit< point > pointHit
Definition: pointHit.H:41
label findNearestBoundaryFace(const point &location, const label seedFacei=-1, const bool useTreeSearch=true) const
Find nearest boundary face.
Definition: meshSearch.C:806
volScalarField & p
const indexedOctree< treeDataCell > & cellTree() const
Get (demand driven) reference to octree holding all cells.
Definition: meshSearch.C:597
const labelList & faceLabels() const
Definition: treeDataFace.H:179
Namespace for OpenFOAM.