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-2015 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  (
151  "meshSearch::findNearestCellWalk(const point&, const label)"
152  ) << "illegal seedCell:" << seedCellI << exit(FatalError);
153  }
154 
155  // Walk in direction of face that decreases distance
156 
157  label curCellI = seedCellI;
158  scalar distanceSqr = magSqr(mesh_.cellCentres()[curCellI] - location);
159 
160  bool closer;
161 
162  do
163  {
164  // Try neighbours of curCellI
165  closer = findNearer
166  (
167  location,
168  mesh_.cellCentres(),
169  mesh_.cellCells()[curCellI],
170  curCellI,
171  distanceSqr
172  );
173  } while (closer);
174 
175  return curCellI;
176 }
177 
178 
179 // tree based searching
180 Foam::label Foam::meshSearch::findNearestFaceTree(const point& location) const
181 {
182  // Search nearest cell centre.
183  const indexedOctree<treeDataCell>& tree = cellTree();
184 
185  // Search with decent span
186  pointIndexHit info = tree.findNearest
187  (
188  location,
189  magSqr(tree.bb().max()-tree.bb().min())
190  );
191 
192  if (!info.hit())
193  {
194  // Search with desparate span
195  info = tree.findNearest(location, Foam::sqr(GREAT));
196  }
197 
198 
199  // Now check any of the faces of the nearest cell
200  const vectorField& centres = mesh_.faceCentres();
201  const cell& ownFaces = mesh_.cells()[info.index()];
202 
203  label nearestFaceI = ownFaces[0];
204  scalar minProximity = magSqr(centres[nearestFaceI] - location);
205 
206  findNearer
207  (
208  location,
209  centres,
210  ownFaces,
211  nearestFaceI,
212  minProximity
213  );
214 
215  return nearestFaceI;
216 }
217 
218 
219 // linear searching
220 Foam::label Foam::meshSearch::findNearestFaceLinear(const point& location) const
221 {
222  const vectorField& centres = mesh_.faceCentres();
223 
224  label nearestFaceI = 0;
225  scalar minProximity = magSqr(centres[nearestFaceI] - location);
226 
227  findNearer
228  (
229  location,
230  centres,
231  nearestFaceI,
232  minProximity
233  );
234 
235  return nearestFaceI;
236 }
237 
238 
239 // walking from seed
240 Foam::label Foam::meshSearch::findNearestFaceWalk
241 (
242  const point& location,
243  const label seedFaceI
244 ) const
245 {
246  if (seedFaceI < 0)
247  {
249  (
250  "meshSearch::findNearestFaceWalk(const point&, const label)"
251  ) << "illegal seedFace:" << seedFaceI << exit(FatalError);
252  }
253 
254  const vectorField& centres = mesh_.faceCentres();
255 
256 
257  // Walk in direction of face that decreases distance
258 
259  label curFaceI = seedFaceI;
260  scalar distanceSqr = magSqr(centres[curFaceI] - location);
261 
262  while (true)
263  {
264  label betterFaceI = curFaceI;
265 
266  findNearer
267  (
268  location,
269  centres,
270  mesh_.cells()[mesh_.faceOwner()[curFaceI]],
271  betterFaceI,
272  distanceSqr
273  );
274 
275  if (mesh_.isInternalFace(curFaceI))
276  {
277  findNearer
278  (
279  location,
280  centres,
281  mesh_.cells()[mesh_.faceNeighbour()[curFaceI]],
282  betterFaceI,
283  distanceSqr
284  );
285  }
286 
287  if (betterFaceI == curFaceI)
288  {
289  break;
290  }
291 
292  curFaceI = betterFaceI;
293  }
294 
295  return curFaceI;
296 }
297 
298 
299 Foam::label Foam::meshSearch::findCellLinear(const point& location) const
300 {
301  bool cellFound = false;
302  label n = 0;
303 
304  label cellI = -1;
305 
306  while ((!cellFound) && (n < mesh_.nCells()))
307  {
308  if (mesh_.pointInCell(location, n, cellDecompMode_))
309  {
310  cellFound = true;
311  cellI = n;
312  }
313  else
314  {
315  n++;
316  }
317  }
318  if (cellFound)
319  {
320  return cellI;
321  }
322  else
323  {
324  return -1;
325  }
326 }
327 
328 
329 // walking from seed
330 Foam::label Foam::meshSearch::findCellWalk
331 (
332  const point& location,
333  const label seedCellI
334 ) const
335 {
336  if (seedCellI < 0)
337  {
339  (
340  "meshSearch::findCellWalk(const point&, const label)"
341  ) << "illegal seedCell:" << seedCellI << exit(FatalError);
342  }
343 
344  if (mesh_.pointInCell(location, seedCellI, cellDecompMode_))
345  {
346  return seedCellI;
347  }
348 
349  // Walk in direction of face that decreases distance
350  label curCellI = seedCellI;
351  scalar nearestDistSqr = magSqr(mesh_.cellCentres()[curCellI] - location);
352 
353  while(true)
354  {
355  // Try neighbours of curCellI
356 
357  const cell& cFaces = mesh_.cells()[curCellI];
358 
359  label nearestCellI = -1;
360 
361  forAll(cFaces, i)
362  {
363  label faceI = cFaces[i];
364 
365  if (mesh_.isInternalFace(faceI))
366  {
367  label cellI = mesh_.faceOwner()[faceI];
368  if (cellI == curCellI)
369  {
370  cellI = mesh_.faceNeighbour()[faceI];
371  }
372 
373  // Check if this is the correct cell
374  if (mesh_.pointInCell(location, cellI, cellDecompMode_))
375  {
376  return cellI;
377  }
378 
379  // Also calculate the nearest cell
380  scalar distSqr = magSqr(mesh_.cellCentres()[cellI] - location);
381 
382  if (distSqr < nearestDistSqr)
383  {
384  nearestDistSqr = distSqr;
385  nearestCellI = cellI;
386  }
387  }
388  }
389 
390  if (nearestCellI == -1)
391  {
392  return -1;
393  }
394 
395  // Continue with the nearest cell
396  curCellI = nearestCellI;
397  }
398 
399  return -1;
400 }
401 
402 
403 Foam::label Foam::meshSearch::findNearestBoundaryFaceWalk
404 (
405  const point& location,
406  const label seedFaceI
407 ) const
408 {
409  if (seedFaceI < 0)
410  {
412  (
413  "meshSearch::findNearestBoundaryFaceWalk"
414  "(const point&, const label)"
415  ) << "illegal seedFace:" << seedFaceI << exit(FatalError);
416  }
417 
418  // Start off from seedFaceI
419 
420  label curFaceI = seedFaceI;
421 
422  const face& f = mesh_.faces()[curFaceI];
423 
424  scalar minDist = f.nearestPoint
425  (
426  location,
427  mesh_.points()
428  ).distance();
429 
430  bool closer;
431 
432  do
433  {
434  closer = false;
435 
436  // Search through all neighbouring boundary faces by going
437  // across edges
438 
439  label lastFaceI = curFaceI;
440 
441  const labelList& myEdges = mesh_.faceEdges()[curFaceI];
442 
443  forAll(myEdges, myEdgeI)
444  {
445  const labelList& neighbours = mesh_.edgeFaces()[myEdges[myEdgeI]];
446 
447  // Check any face which uses edge, is boundary face and
448  // is not curFaceI itself.
449 
450  forAll(neighbours, nI)
451  {
452  label faceI = neighbours[nI];
453 
454  if
455  (
456  (faceI >= mesh_.nInternalFaces())
457  && (faceI != lastFaceI)
458  )
459  {
460  const face& f = mesh_.faces()[faceI];
461 
462  pointHit curHit = f.nearestPoint
463  (
464  location,
465  mesh_.points()
466  );
467 
468  // If the face is closer, reset current face and distance
469  if (curHit.distance() < minDist)
470  {
471  minDist = curHit.distance();
472  curFaceI = faceI;
473  closer = true; // a closer neighbour has been found
474  }
475  }
476  }
477  }
478  } while (closer);
479 
480  return curFaceI;
481 }
482 
483 
484 Foam::vector Foam::meshSearch::offset
485 (
486  const point& bPoint,
487  const label bFaceI,
488  const vector& dir
489 ) const
490 {
491  // Get the neighbouring cell
492  label ownerCellI = mesh_.faceOwner()[bFaceI];
493 
494  const point& c = mesh_.cellCentres()[ownerCellI];
495 
496  // Typical dimension: distance from point on face to cell centre
497  scalar typDim = mag(c - bPoint);
498 
499  return tol_*typDim*dir;
500 }
501 
502 
503 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
504 
505 Foam::meshSearch::meshSearch
506 (
507  const polyMesh& mesh,
508  const polyMesh::cellDecomposition cellDecompMode
509 )
510 :
511  mesh_(mesh),
512  cellDecompMode_(cellDecompMode)
513 {
514  if (cellDecompMode_ == polyMesh::FACE_DIAG_TRIS)
515  {
516  // Force construction of face diagonals
517  (void)mesh.tetBasePtIs();
518  }
519 }
520 
521 
522 // Construct with a custom bounding box
523 Foam::meshSearch::meshSearch
524 (
525  const polyMesh& mesh,
526  const treeBoundBox& bb,
527  const polyMesh::cellDecomposition cellDecompMode
528 )
529 :
530  mesh_(mesh),
531  cellDecompMode_(cellDecompMode)
532 {
533  overallBbPtr_.reset(new treeBoundBox(bb));
534 
535  if (cellDecompMode_ == polyMesh::FACE_DIAG_TRIS)
536  {
537  // Force construction of face diagonals
538  (void)mesh.tetBasePtIs();
539  }
540 }
541 
542 
543 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
544 
546 {
547  clearOut();
548 }
549 
550 
551 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
552 
554  const
555 {
556  if (!boundaryTreePtr_.valid())
557  {
558  //
559  // Construct tree
560  //
561 
562  if (!overallBbPtr_.valid())
563  {
564  Random rndGen(261782);
565  overallBbPtr_.reset
566  (
567  new treeBoundBox(mesh_.points())
568  );
569 
570  treeBoundBox& overallBb = overallBbPtr_();
571  // Extend slightly and make 3D
572  overallBb = overallBb.extend(rndGen, 1e-4);
573  overallBb.min() -= point(ROOTVSMALL, ROOTVSMALL, ROOTVSMALL);
574  overallBb.max() += point(ROOTVSMALL, ROOTVSMALL, ROOTVSMALL);
575  }
576 
577  // all boundary faces (not just walls)
578  labelList bndFaces(mesh_.nFaces()-mesh_.nInternalFaces());
579  forAll(bndFaces, i)
580  {
581  bndFaces[i] = mesh_.nInternalFaces() + i;
582  }
583 
584  boundaryTreePtr_.reset
585  (
587  (
588  treeDataFace // all information needed to search faces
589  (
590  false, // do not cache bb
591  mesh_,
592  bndFaces // boundary faces only
593  ),
594  overallBbPtr_(), // overall search domain
595  8, // maxLevel
596  10, // leafsize
597  3.0 // duplicity
598  )
599  );
600  }
601 
602  return boundaryTreePtr_();
603 }
604 
605 
607 const
608 {
609  if (!cellTreePtr_.valid())
610  {
611  //
612  // Construct tree
613  //
614 
615  if (!overallBbPtr_.valid())
616  {
617  Random rndGen(261782);
618  overallBbPtr_.reset
619  (
620  new treeBoundBox(mesh_.points())
621  );
622 
623  treeBoundBox& overallBb = overallBbPtr_();
624  // Extend slightly and make 3D
625  overallBb = overallBb.extend(rndGen, 1e-4);
626  overallBb.min() -= point(ROOTVSMALL, ROOTVSMALL, ROOTVSMALL);
627  overallBb.max() += point(ROOTVSMALL, ROOTVSMALL, ROOTVSMALL);
628  }
629 
630  cellTreePtr_.reset
631  (
633  (
635  (
636  false, // not cache bb
637  mesh_,
638  cellDecompMode_ // cell decomposition mode for inside tests
639  ),
640  overallBbPtr_(),
641  8, // maxLevel
642  10, // leafsize
643  6.0 // duplicity
644  )
645  );
646  }
647 
648  return cellTreePtr_();
649 }
650 
651 
656 //bool Foam::meshSearch::pointInCell(const point& p, label cellI) const
657 //{
658 // if (faceDecomp_)
659 // {
660 // const point& ctr = mesh_.cellCentres()[cellI];
661 //
662 // vector dir(p - ctr);
663 // scalar magDir = mag(dir);
664 //
665 // // Check if any faces are hit by ray from cell centre to p.
666 // // If none -> p is in cell.
667 // const labelList& cFaces = mesh_.cells()[cellI];
668 //
669 // // Make sure half_ray does not pick up any faces on the wrong
670 // // side of the ray.
671 // scalar oldTol = intersection::setPlanarTol(0.0);
672 //
673 // forAll(cFaces, i)
674 // {
675 // label faceI = cFaces[i];
676 //
677 // pointHit inter = mesh_.faces()[faceI].ray
678 // (
679 // ctr,
680 // dir,
681 // mesh_.points(),
682 // intersection::HALF_RAY,
683 // intersection::VECTOR
684 // );
685 //
686 // if (inter.hit())
687 // {
688 // scalar dist = inter.distance();
689 //
690 // if (dist < magDir)
691 // {
692 // // Valid hit. Hit face so point is not in cell.
693 // intersection::setPlanarTol(oldTol);
694 //
695 // return false;
696 // }
697 // }
698 // }
699 //
700 // intersection::setPlanarTol(oldTol);
701 //
702 // // No face inbetween point and cell centre so point is inside.
703 // return true;
704 // }
705 // else
706 // {
707 // const labelList& f = mesh_.cells()[cellI];
708 // const labelList& owner = mesh_.faceOwner();
709 // const vectorField& cf = mesh_.faceCentres();
710 // const vectorField& Sf = mesh_.faceAreas();
711 //
712 // forAll(f, facei)
713 // {
714 // label nFace = f[facei];
715 // vector proj = p - cf[nFace];
716 // vector normal = Sf[nFace];
717 // if (owner[nFace] == cellI)
718 // {
719 // if ((normal & proj) > 0)
720 // {
721 // return false;
722 // }
723 // }
724 // else
725 // {
726 // if ((normal & proj) < 0)
727 // {
728 // return false;
729 // }
730 // }
731 // }
732 //
733 // return true;
734 // }
735 //}
736 
737 
739 (
740  const point& location,
741  const label seedCellI,
742  const bool useTreeSearch
743 ) const
744 {
745  if (seedCellI == -1)
746  {
747  if (useTreeSearch)
748  {
749  return findNearestCellTree(location);
750  }
751  else
752  {
753  return findNearestCellLinear(location);
754  }
755  }
756  else
757  {
758  return findNearestCellWalk(location, seedCellI);
759  }
760 }
761 
762 
764 (
765  const point& location,
766  const label seedFaceI,
767  const bool useTreeSearch
768 ) const
769 {
770  if (seedFaceI == -1)
771  {
772  if (useTreeSearch)
773  {
774  return findNearestFaceTree(location);
775  }
776  else
777  {
778  return findNearestFaceLinear(location);
779  }
780  }
781  else
782  {
783  return findNearestFaceWalk(location, seedFaceI);
784  }
785 }
786 
787 
789 (
790  const point& location,
791  const label seedCellI,
792  const bool useTreeSearch
793 ) const
794 {
795  // Find the nearest cell centre to this location
796  if (seedCellI == -1)
797  {
798  if (useTreeSearch)
799  {
800  return cellTree().findInside(location);
801  }
802  else
803  {
804  return findCellLinear(location);
805  }
806  }
807  else
808  {
809  return findCellWalk(location, seedCellI);
810  }
811 }
812 
813 
815 (
816  const point& location,
817  const label seedFaceI,
818  const bool useTreeSearch
819 ) const
820 {
821  if (seedFaceI == -1)
822  {
823  if (useTreeSearch)
824  {
825  const indexedOctree<treeDataFace>& tree = boundaryTree();
826 
827  pointIndexHit info = boundaryTree().findNearest
828  (
829  location,
830  magSqr(tree.bb().max()-tree.bb().min())
831  );
832 
833  if (!info.hit())
834  {
835  info = boundaryTree().findNearest
836  (
837  location,
838  Foam::sqr(GREAT)
839  );
840  }
841 
842  return tree.shapes().faceLabels()[info.index()];
843  }
844  else
845  {
846  scalar minDist = GREAT;
847 
848  label minFaceI = -1;
849 
850  for
851  (
852  label faceI = mesh_.nInternalFaces();
853  faceI < mesh_.nFaces();
854  faceI++
855  )
856  {
857  const face& f = mesh_.faces()[faceI];
858 
859  pointHit curHit =
860  f.nearestPoint
861  (
862  location,
863  mesh_.points()
864  );
865 
866  if (curHit.distance() < minDist)
867  {
868  minDist = curHit.distance();
869  minFaceI = faceI;
870  }
871  }
872  return minFaceI;
873  }
874  }
875  else
876  {
877  return findNearestBoundaryFaceWalk(location, seedFaceI);
878  }
879 }
880 
881 
883 (
884  const point& pStart,
885  const point& pEnd
886 ) const
887 {
888  pointIndexHit curHit = boundaryTree().findLine(pStart, pEnd);
889 
890  if (curHit.hit())
891  {
892  // Change index into octreeData into face label
893  curHit.setIndex(boundaryTree().shapes().faceLabels()[curHit.index()]);
894  }
895  return curHit;
896 }
897 
898 
900 (
901  const point& pStart,
902  const point& pEnd
903 ) const
904 {
906 
907  vector edgeVec = pEnd - pStart;
908  edgeVec /= mag(edgeVec);
909 
910  point pt = pStart;
911 
912  pointIndexHit bHit;
913  do
914  {
915  bHit = intersection(pt, pEnd);
916 
917  if (bHit.hit())
918  {
919  hits.append(bHit);
920 
921  const vector& area = mesh_.faceAreas()[bHit.index()];
922 
923  scalar typDim = Foam::sqrt(mag(area));
924 
925  if ((mag(bHit.hitPoint() - pEnd)/typDim) < SMALL)
926  {
927  break;
928  }
929 
930  // Restart from hitPoint shifted a little bit in the direction
931  // of the destination
932 
933  pt =
934  bHit.hitPoint()
935  + offset(bHit.hitPoint(), bHit.index(), edgeVec);
936  }
937 
938  } while (bHit.hit());
939 
940 
941  hits.shrink();
942 
943  return hits;
944 }
945 
946 
948 {
949  return (boundaryTree().getVolumeType(p) == volumeType::INSIDE);
950 }
951 
952 
953 // Delete all storage
955 {
956  boundaryTreePtr_.clear();
957  cellTreePtr_.clear();
958  overallBbPtr_.clear();
959 }
960 
961 
963 {
964  clearOut();
965 }
966 
967 
968 // ************************************************************************* //
dimensionedScalar sqrt(const dimensionedScalar &ds)
Simple random number generator.
Definition: Random.H:49
cachedRandom rndGen(label(0),-1)
const point & min() const
Minimum describing the bounding box.
Definition: boundBoxI.H:54
const indexedOctree< treeDataFace > & boundaryTree() const
Get (demand driven) reference to octree holding all.
Definition: meshSearch.C:553
This class describes the interaction of (usually) a face and a point. It carries the info of a succes...
Definition: PointIndexHit.H:53
vector point
Point is a vector.
Definition: point.H:41
List< pointIndexHit > intersections(const point &pStart, const point &pEnd) const
Find all intersections of boundary within segment pStart .. pEnd.
Definition: meshSearch.C:900
dimensioned< scalar > mag(const dimensioned< Type > &)
const labelList & tetBasePtIs() const
Return the tetBasePtIs.
Definition: polyMesh.C:864
dimensioned< scalar > magSqr(const dimensioned< Type > &)
const Type & shapes() const
Reference to shape.
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:739
const treeBoundBox & bb() const
Top bounding box.
bool isInside(const point &) const
Determine inside/outside status.
Definition: meshSearch.C:947
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
Standard boundBox + extra functionality for use in octree.
Definition: treeBoundBox.H:75
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Definition: DynamicListI.H:310
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
PointIndexHit< point > pointIndexHit
Definition: pointIndexHit.H:42
bool hit() const
Is there a hit.
label findNearestFace(const point &location, const label seedFaceI=-1, const bool useTreeSearch=true) const
Definition: meshSearch.C:764
Encapsulation of data needed to search for faces.
Definition: treeDataFace.H:58
void correct()
Correct for mesh geom/topo changes.
Definition: meshSearch.C:962
label index() const
Return index.
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
Namespace for OpenFOAM.
PointHit< point > pointHit
Definition: pointHit.H:41
const Point & hitPoint() const
Return hit point.
scalar distance(const vector &p1, const vector &p2)
Definition: curveTools.C:12
label n
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
const double e
Elementary charge.
Definition: doubleFloat.H:78
pointIndexHit intersection(const point &pStart, const point &pEnd) const
Find first intersection of boundary in segment [pStart, pEnd].
Definition: meshSearch.C:883
label findCell(const point &location, const label seedCellI=-1, const bool useTreeSearch=true) const
Find cell containing location.
Definition: meshSearch.C:789
treeBoundBox extend(Random &, const scalar s) const
Return slightly wider bounding box.
volScalarField & p
Definition: createFields.H:51
const labelList & faceLabels() const
Definition: treeDataFace.H:179
Field< vector > vectorField
Specialisation of Field<T> for vector.
#define forAll(list, i)
Definition: UList.H:421
DynamicList< T, SizeInc, SizeMult, SizeDiv > & shrink()
Shrink the allocated space to the number of elements used.
Definition: DynamicListI.H:258
Template functions to aid in the implementation of demand driven data.
const point & max() const
Maximum describing the bounding box.
Definition: boundBoxI.H:60
void setIndex(const label index)
#define FatalErrorIn(functionName)
Report an error message using Foam::FatalError.
Definition: error.H:314
scalar distance() const
Return distance to hit.
Definition: PointHit.H:139
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
const indexedOctree< treeDataCell > & cellTree() const
Get (demand driven) reference to octree holding all cells.
Definition: meshSearch.C:606
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects...
Definition: DynamicList.H:56
error FatalError
pointHit nearestPoint(const point &p, const pointField &) const
Return nearest point to face.
List< label > labelList
A List of labels.
Definition: labelList.H:56
label findNearestBoundaryFace(const point &location, const label seedFaceI=-1, const bool useTreeSearch=true) const
Find nearest boundary face.
Definition: meshSearch.C:815
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
void clearOut()
Delete all storage.
Definition: meshSearch.C:954
Foam::intersection.
Definition: intersection.H:49
static scalar tol_
Tolerance on linear dimensions.
Definition: meshSearch.H:163
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Encapsulation of data needed to search in/for cells. Used to find the cell containing a point (e...
Definition: treeDataCell.H:54
cellDecomposition
Enumeration defining the decomposition of the cell for.
Definition: polyMesh.H:98
~meshSearch()
Destructor.
Definition: meshSearch.C:545
defineTypeNameAndDebug(combustionModel, 0)