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 | Copyright 2015 OpenCFD Ltd.
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 
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 
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
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
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
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
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
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
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 
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
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 
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 
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 
497 (
498  const polyMesh& mesh,
499  const polyMesh::cellDecomposition cellDecompMode
500 )
501 :
502  mesh_(mesh),
503  cellDecompMode_(cellDecompMode)
504 {
505  if
506  (
507  cellDecompMode_ == polyMesh::FACE_DIAG_TRIS
508  || cellDecompMode_ == polyMesh::CELL_TETS
509  )
510  {
511  // Force construction of face diagonals
512  (void)mesh.tetBasePtIs();
513  }
514 }
515 
516 
517 // Construct with a custom bounding box
519 (
520  const polyMesh& mesh,
521  const treeBoundBox& bb,
522  const polyMesh::cellDecomposition cellDecompMode
523 )
524 :
525  mesh_(mesh),
526  cellDecompMode_(cellDecompMode)
527 {
528  overallBbPtr_.reset(new treeBoundBox(bb));
529 
530  if
531  (
532  cellDecompMode_ == polyMesh::FACE_DIAG_TRIS
533  || cellDecompMode_ == polyMesh::CELL_TETS
534  )
535  {
536  // Force construction of face diagonals
537  (void)mesh.tetBasePtIs();
538  }
539 }
540 
541 
542 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
543 
545 {
546  clearOut();
547 }
548 
549 
550 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
551 
553  const
554 {
555  if (!boundaryTreePtr_.valid())
556  {
557  //
558  // Construct tree
559  //
560 
561  if (!overallBbPtr_.valid())
562  {
563  Random rndGen(261782);
564  overallBbPtr_.reset
565  (
566  new treeBoundBox(mesh_.points())
567  );
568 
569  treeBoundBox& overallBb = overallBbPtr_();
570  // Extend slightly and make 3D
571  overallBb = overallBb.extend(rndGen, 1e-4);
572  overallBb.min() -= point(ROOTVSMALL, ROOTVSMALL, ROOTVSMALL);
573  overallBb.max() += point(ROOTVSMALL, ROOTVSMALL, ROOTVSMALL);
574  }
575 
576  // all boundary faces (not just walls)
577  labelList bndFaces(mesh_.nFaces()-mesh_.nInternalFaces());
578  forAll(bndFaces, i)
579  {
580  bndFaces[i] = mesh_.nInternalFaces() + i;
581  }
582 
583  boundaryTreePtr_.reset
584  (
586  (
587  treeDataFace // all information needed to search faces
588  (
589  false, // do not cache bb
590  mesh_,
591  bndFaces // boundary faces only
592  ),
593  overallBbPtr_(), // overall search domain
594  8, // maxLevel
595  10, // leafsize
596  3.0 // duplicity
597  )
598  );
599  }
600 
601  return boundaryTreePtr_();
602 }
603 
604 
606 const
607 {
608  if (!cellTreePtr_.valid())
609  {
610  //
611  // Construct tree
612  //
613 
614  if (!overallBbPtr_.valid())
615  {
616  Random rndGen(261782);
617  overallBbPtr_.reset
618  (
619  new treeBoundBox(mesh_.points())
620  );
621 
622  treeBoundBox& overallBb = overallBbPtr_();
623  // Extend slightly and make 3D
624  overallBb = overallBb.extend(rndGen, 1e-4);
625  overallBb.min() -= point(ROOTVSMALL, ROOTVSMALL, ROOTVSMALL);
626  overallBb.max() += point(ROOTVSMALL, ROOTVSMALL, ROOTVSMALL);
627  }
628 
629  cellTreePtr_.reset
630  (
632  (
634  (
635  false, // not cache bb
636  mesh_,
637  cellDecompMode_ // cell decomposition mode for inside tests
638  ),
639  overallBbPtr_(),
640  8, // maxLevel
641  10, // leafsize
642  6.0 // duplicity
643  )
644  );
645  }
646 
647  return cellTreePtr_();
648 }
649 
650 
655 //bool Foam::meshSearch::pointInCell(const point& p, label cellI) const
656 //{
657 // if (faceDecomp_)
658 // {
659 // const point& ctr = mesh_.cellCentres()[cellI];
660 //
661 // vector dir(p - ctr);
662 // scalar magDir = mag(dir);
663 //
664 // // Check if any faces are hit by ray from cell centre to p.
665 // // If none -> p is in cell.
666 // const labelList& cFaces = mesh_.cells()[cellI];
667 //
668 // // Make sure half_ray does not pick up any faces on the wrong
669 // // side of the ray.
670 // scalar oldTol = intersection::setPlanarTol(0.0);
671 //
672 // forAll(cFaces, i)
673 // {
674 // label faceI = cFaces[i];
675 //
676 // pointHit inter = mesh_.faces()[faceI].ray
677 // (
678 // ctr,
679 // dir,
680 // mesh_.points(),
681 // intersection::HALF_RAY,
682 // intersection::VECTOR
683 // );
684 //
685 // if (inter.hit())
686 // {
687 // scalar dist = inter.distance();
688 //
689 // if (dist < magDir)
690 // {
691 // // Valid hit. Hit face so point is not in cell.
692 // intersection::setPlanarTol(oldTol);
693 //
694 // return false;
695 // }
696 // }
697 // }
698 //
699 // intersection::setPlanarTol(oldTol);
700 //
701 // // No face inbetween point and cell centre so point is inside.
702 // return true;
703 // }
704 // else
705 // {
706 // const labelList& f = mesh_.cells()[cellI];
707 // const labelList& owner = mesh_.faceOwner();
708 // const vectorField& cf = mesh_.faceCentres();
709 // const vectorField& Sf = mesh_.faceAreas();
710 //
711 // forAll(f, facei)
712 // {
713 // label nFace = f[facei];
714 // vector proj = p - cf[nFace];
715 // vector normal = Sf[nFace];
716 // if (owner[nFace] == cellI)
717 // {
718 // if ((normal & proj) > 0)
719 // {
720 // return false;
721 // }
722 // }
723 // else
724 // {
725 // if ((normal & proj) < 0)
726 // {
727 // return false;
728 // }
729 // }
730 // }
731 //
732 // return true;
733 // }
734 //}
735 
736 
738 (
739  const point& location,
740  const label seedCellI,
741  const bool useTreeSearch
742 ) const
743 {
744  if (seedCellI == -1)
745  {
746  if (useTreeSearch)
747  {
748  return findNearestCellTree(location);
749  }
750  else
751  {
752  return findNearestCellLinear(location);
753  }
754  }
755  else
756  {
757  return findNearestCellWalk(location, seedCellI);
758  }
759 }
760 
761 
763 (
764  const point& location,
765  const label seedFaceI,
766  const bool useTreeSearch
767 ) const
768 {
769  if (seedFaceI == -1)
770  {
771  if (useTreeSearch)
772  {
773  return findNearestFaceTree(location);
774  }
775  else
776  {
777  return findNearestFaceLinear(location);
778  }
779  }
780  else
781  {
782  return findNearestFaceWalk(location, seedFaceI);
783  }
784 }
785 
786 
788 (
789  const point& location,
790  const label seedCellI,
791  const bool useTreeSearch
792 ) const
793 {
794  // Find the nearest cell centre to this location
795  if (seedCellI == -1)
796  {
797  if (useTreeSearch)
798  {
799  return cellTree().findInside(location);
800  }
801  else
802  {
803  return findCellLinear(location);
804  }
805  }
806  else
807  {
808  return findCellWalk(location, seedCellI);
809  }
810 }
811 
812 
814 (
815  const point& location,
816  const label seedFaceI,
817  const bool useTreeSearch
818 ) const
819 {
820  if (seedFaceI == -1)
821  {
822  if (useTreeSearch)
823  {
824  const indexedOctree<treeDataFace>& tree = boundaryTree();
825 
826  pointIndexHit info = boundaryTree().findNearest
827  (
828  location,
829  magSqr(tree.bb().max()-tree.bb().min())
830  );
831 
832  if (!info.hit())
833  {
834  info = boundaryTree().findNearest
835  (
836  location,
837  Foam::sqr(GREAT)
838  );
839  }
840 
841  return tree.shapes().faceLabels()[info.index()];
842  }
843  else
844  {
845  scalar minDist = GREAT;
846 
847  label minFaceI = -1;
848 
849  for
850  (
851  label faceI = mesh_.nInternalFaces();
852  faceI < mesh_.nFaces();
853  faceI++
854  )
855  {
856  const face& f = mesh_.faces()[faceI];
857 
858  pointHit curHit =
859  f.nearestPoint
860  (
861  location,
862  mesh_.points()
863  );
864 
865  if (curHit.distance() < minDist)
866  {
867  minDist = curHit.distance();
868  minFaceI = faceI;
869  }
870  }
871  return minFaceI;
872  }
873  }
874  else
875  {
876  return findNearestBoundaryFaceWalk(location, seedFaceI);
877  }
878 }
879 
880 
882 (
883  const point& pStart,
884  const point& pEnd
885 ) const
886 {
887  pointIndexHit curHit = boundaryTree().findLine(pStart, pEnd);
888 
889  if (curHit.hit())
890  {
891  // Change index into octreeData into face label
892  curHit.setIndex(boundaryTree().shapes().faceLabels()[curHit.index()]);
893  }
894  return curHit;
895 }
896 
897 
899 (
900  const point& pStart,
901  const point& pEnd
902 ) const
903 {
905 
906  vector edgeVec = pEnd - pStart;
907  edgeVec /= mag(edgeVec);
908 
909  point pt = pStart;
910 
911  pointIndexHit bHit;
912  do
913  {
914  bHit = intersection(pt, pEnd);
915 
916  if (bHit.hit())
917  {
918  hits.append(bHit);
919 
920  const vector& area = mesh_.faceAreas()[bHit.index()];
921 
922  scalar typDim = Foam::sqrt(mag(area));
923 
924  if ((mag(bHit.hitPoint() - pEnd)/typDim) < SMALL)
925  {
926  break;
927  }
928 
929  // Restart from hitPoint shifted a little bit in the direction
930  // of the destination
931 
932  pt =
933  bHit.hitPoint()
934  + offset(bHit.hitPoint(), bHit.index(), edgeVec);
935  }
936 
937  } while (bHit.hit());
938 
939 
940  hits.shrink();
941 
942  return hits;
943 }
944 
945 
947 {
948  return (boundaryTree().getVolumeType(p) == volumeType::INSIDE);
949 }
950 
951 
952 // Delete all storage
954 {
955  boundaryTreePtr_.clear();
956  cellTreePtr_.clear();
957  overallBbPtr_.clear();
958 }
959 
960 
962 {
963  clearOut();
964 }
965 
966 
967 // ************************************************************************* //
Foam::Random
Simple random number generator.
Definition: Random.H:49
Foam::volumeType::INSIDE
@ INSIDE
Definition: volumeType.H:63
Foam::PointIndexHit::index
label index() const
Return index.
Definition: PointIndexHit.H:111
Foam::polyMesh::FACE_DIAG_TRIS
@ FACE_DIAG_TRIS
Definition: polyMesh.H:105
Foam::polyMesh::cellDecomposition
cellDecomposition
Enumeration defining the decomposition of the cell for.
Definition: polyMesh.H:98
p
p
Definition: pEqn.H:62
Foam::PointIndexHit::setIndex
void setIndex(const label index)
Definition: PointIndexHit.H:168
Foam::boundBox::max
const point & max() const
Maximum describing the bounding box.
Definition: boundBoxI.H:60
Foam::meshSearch::findNearestFaceLinear
label findNearestFaceLinear(const point &) const
Definition: meshSearch.C:218
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::PointHit
This class describes the interaction of a face and a point. It carries the info of a successful hit a...
Definition: PointHit.H:51
Foam::DynamicList
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:56
Foam::polyMesh::tetBasePtIs
const labelList & tetBasePtIs() const
Return the tetBasePtIs.
Definition: polyMesh.C:818
demandDrivenData.H
Template functions to aid in the implementation of demand driven data.
Foam::treeBoundBox
Standard boundBox + extra functionality for use in octree.
Definition: treeBoundBox.H:75
indexedOctree.H
Foam::PointIndexHit::hit
bool hit() const
Is there a hit.
Definition: PointIndexHit.H:105
Foam::indexedOctree::shapes
const Type & shapes() const
Reference to shape.
Definition: indexedOctree.H:449
Foam::meshSearch::findNearestBoundaryFace
label findNearestBoundaryFace(const point &location, const label seedFaceI=-1, const bool useTreeSearch=true) const
Find nearest boundary face.
Definition: meshSearch.C:814
Foam::meshSearch::findNearer
static bool findNearer(const point &sample, const pointField &points, label &nearestI, scalar &nearestDistSqr)
Updates nearestI, nearestDistSqr from any closer ones.
Definition: meshSearch.C:47
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
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:899
Foam::meshSearch::findCellWalk
label findCellWalk(const point &, const label) const
Walk from seed. Does not 'go around' boundary, just returns.
Definition: meshSearch.C:327
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
treeDataFace.H
Foam::meshSearch::cellTree
const indexedOctree< treeDataCell > & cellTree() const
Get (demand driven) reference to octree holding all cells.
Definition: meshSearch.C:605
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::PointIndexHit
This class describes the interaction of (usually) a face and a point. It carries the info of a succes...
Definition: PointIndexHit.H:53
Foam::meshSearch::clearOut
void clearOut()
Delete all storage.
Definition: meshSearch.C:953
Foam::meshSearch::correct
void correct()
Correct for mesh geom/topo changes.
Definition: meshSearch.C:961
Foam::PointHit::distance
scalar distance() const
Return distance to hit.
Definition: PointHit.H:139
Foam::meshSearch::findNearestFace
label findNearestFace(const point &location, const label seedFaceI=-1, const bool useTreeSearch=true) const
Definition: meshSearch.C:763
Foam::treeDataCell
Encapsulation of data needed to search in/for cells. Used to find the cell containing a point (e....
Definition: treeDataCell.H:54
Foam::label
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
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::meshSearch::findCellLinear
label findCellLinear(const point &) const
Cell containing location. Linear search.
Definition: meshSearch.C:295
Foam::meshSearch::findNearestCellWalk
label findNearestCellWalk(const point &, const label) const
Walk from seed. Does not 'go around' boundary, just returns.
Definition: meshSearch.C:142
Foam::treeBoundBox::extend
treeBoundBox extend(Random &, const scalar s) const
Return slightly wider bounding box.
Definition: treeBoundBoxI.H:317
Foam::PointIndexHit::hitPoint
const Point & hitPoint() const
Return hit point.
Definition: PointIndexHit.H:117
Foam::indexedOctree
Non-pointer based hierarchical recursive searching.
Definition: treeDataTriSurface.H:47
Foam::DynamicList::shrink
DynamicList< T, SizeInc, SizeMult, SizeDiv > & shrink()
Shrink the allocated space to the number of elements used.
Definition: DynamicListI.H:258
Foam::boundBox::min
const point & min() const
Minimum describing the bounding box.
Definition: boundBoxI.H:54
Foam::meshSearch::findNearestBoundaryFaceWalk
label findNearestBoundaryFaceWalk(const point &location, const label seedFaceI) const
Walk from seed to find nearest boundary face. Gets stuck in.
Definition: meshSearch.C:398
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:738
Foam::meshSearch::isInside
bool isInside(const point &) const
Determine inside/outside status.
Definition: meshSearch.C:946
Foam::FatalError
error FatalError
Foam::meshSearch::offset
vector offset(const point &bPoint, const label bFaceI, const vector &dir) const
Calculate offset vector in direction dir with as length a.
Definition: meshSearch.C:476
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
treeDataCell.H
Foam::e
const double e
Elementary charge.
Definition: doubleFloat.H:94
Foam::meshSearch::findNearestFaceTree
label findNearestFaceTree(const point &) const
Definition: meshSearch.C:178
Foam::DynamicList::append
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Foam::indexedOctree::bb
const treeBoundBox & bb() const
Top bounding box.
Definition: indexedOctree.H:468
Foam::meshSearch::findNearestCellLinear
label findNearestCellLinear(const point &) const
Nearest cell centre going through all cells.
Definition: meshSearch.C:121
Foam::polyMesh::CELL_TETS
@ CELL_TETS
Definition: polyMesh.H:107
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:49
meshSearch.H
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: HashTable.H:59
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:142
Foam::meshSearch::findCell
label findCell(const point &location, const label seedCellI=-1, const bool useTreeSearch=true) const
Find cell containing location.
Definition: meshSearch.C:788
Foam::treeDataFace
Encapsulation of data needed to search for faces.
Definition: treeDataFace.H:58
Foam::meshSearch::~meshSearch
~meshSearch()
Destructor.
Definition: meshSearch.C:544
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::intersection
Foam::intersection.
Definition: intersection.H:49
Foam::constant::universal::c
const dimensionedScalar c
Speed of light in a vacuum.
Foam::indexedOctree::findNearest
void findNearest(const label nodeI, const linePointRef &ln, treeBoundBox &tightest, label &nearestShapeI, point &linePoint, point &nearestPoint, const FindNearestOp &fnOp) const
Find nearest point to line.
Definition: indexedOctree.C:590
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
DynamicList.H
Foam::point
vector point
Point is a vector.
Definition: point.H:41
Foam::meshSearch::meshSearch
meshSearch(const meshSearch &)
Disallow default bitwise copy construct.
Foam::meshSearch::findNearestFaceWalk
label findNearestFaceWalk(const point &, const label) const
Definition: meshSearch.C:239
Foam::meshSearch::intersection
pointIndexHit intersection(const point &pStart, const point &pEnd) const
Find first intersection of boundary in segment [pStart, pEnd].
Definition: meshSearch.C:882
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::cell
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:56
Foam::meshSearch::boundaryTree
const indexedOctree< treeDataFace > & boundaryTree() const
Get (demand driven) reference to octree holding all.
Definition: meshSearch.C:552
Foam::magSqr
dimensioned< scalar > magSqr(const dimensioned< Type > &)
rndGen
cachedRandom rndGen(label(0), -1)
Foam::meshSearch::findNearestCellTree
label findNearestCellTree(const point &) const
Nearest cell centre using octree.
Definition: meshSearch.C:102
Foam::meshSearch::tol_
static scalar tol_
Tolerance on linear dimensions.
Definition: meshSearch.H:163