1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright held by original author
7 -------------------------------------------------------------------------------
9 This file is part of OpenFOAM.
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 the
13 Free Software Foundation; either version 2 of the License, or (at your
14 option) any later version.
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
21 You should have received a copy of the GNU General Public License
22 along with OpenFOAM; if not, write to the Free Software Foundation,
23 Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 \*---------------------------------------------------------------------------*/
27 #include "meshSearch.H"
29 #include "indexedOctree.H"
30 #include "DynamicList.H"
31 #include "demandDrivenData.H"
32 #include "treeDataCell.H"
33 #include "treeDataFace.H"
34 #include "treeDataPoint.H"
36 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38 defineTypeNameAndDebug(Foam::meshSearch, 0);
40 Foam::scalar Foam::meshSearch::tol_ = 1E-3;
43 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
45 bool Foam::meshSearch::findNearer
48 const pointField& points,
50 scalar& nearestDistSqr
55 forAll(points, pointI)
57 scalar distSqr = magSqr(points[pointI] - sample);
59 if (distSqr < nearestDistSqr)
61 nearestDistSqr = distSqr;
71 bool Foam::meshSearch::findNearer
74 const pointField& points,
75 const labelList& indices,
77 scalar& nearestDistSqr
84 label pointI = indices[i];
86 scalar distSqr = magSqr(points[pointI] - sample);
88 if (distSqr < nearestDistSqr)
90 nearestDistSqr = distSqr;
100 // tree based searching
101 Foam::label Foam::meshSearch::findNearestCellTree(const point& location) const
103 const indexedOctree<treeDataPoint>& tree = cellCentreTree();
105 scalar span = tree.bb().mag();
107 pointIndexHit info = tree.findNearest(location, Foam::sqr(span));
111 info = tree.findNearest(location, Foam::sqr(GREAT));
118 Foam::label Foam::meshSearch::findNearestCellLinear(const point& location) const
120 const vectorField& centres = mesh_.cellCentres();
122 label nearestIndex = 0;
123 scalar minProximity = magSqr(centres[nearestIndex] - location);
138 Foam::label Foam::meshSearch::findNearestCellWalk
140 const point& location,
141 const label seedCellI
148 "meshSearch::findNearestCellWalk(const point&, const label)"
149 ) << "illegal seedCell:" << seedCellI << exit(FatalError);
152 // Walk in direction of face that decreases distance
154 label curCellI = seedCellI;
155 scalar distanceSqr = magSqr(mesh_.cellCentres()[curCellI] - location);
161 // Try neighbours of curCellI
166 mesh_.cellCells()[curCellI],
176 // tree based searching
177 Foam::label Foam::meshSearch::findNearestFaceTree(const point& location) const
179 // Search nearest cell centre.
180 const indexedOctree<treeDataPoint>& tree = cellCentreTree();
182 scalar span = tree.bb().mag();
184 // Search with decent span
185 pointIndexHit info = tree.findNearest(location, Foam::sqr(span));
189 // Search with desparate span
190 info = tree.findNearest(location, Foam::sqr(GREAT));
194 // Now check any of the faces of the nearest cell
195 const vectorField& centres = mesh_.faceCentres();
196 const cell& ownFaces = mesh_.cells()[info.index()];
198 label nearestFaceI = ownFaces[0];
199 scalar minProximity = magSqr(centres[nearestFaceI] - location);
215 Foam::label Foam::meshSearch::findNearestFaceLinear(const point& location) const
217 const vectorField& centres = mesh_.faceCentres();
219 label nearestFaceI = 0;
220 scalar minProximity = magSqr(centres[nearestFaceI] - location);
235 Foam::label Foam::meshSearch::findNearestFaceWalk
237 const point& location,
238 const label seedFaceI
245 "meshSearch::findNearestFaceWalk(const point&, const label)"
246 ) << "illegal seedFace:" << seedFaceI << exit(FatalError);
249 const vectorField& centres = mesh_.faceCentres();
252 // Walk in direction of face that decreases distance
254 label curFaceI = seedFaceI;
255 scalar distanceSqr = magSqr(centres[curFaceI] - location);
259 label betterFaceI = curFaceI;
265 mesh_.cells()[mesh_.faceOwner()[curFaceI]],
270 if (mesh_.isInternalFace(curFaceI))
276 mesh_.cells()[mesh_.faceNeighbour()[curFaceI]],
282 if (betterFaceI == curFaceI)
287 curFaceI = betterFaceI;
294 Foam::label Foam::meshSearch::findCellLinear(const point& location) const
296 bool cellFound = false;
301 while ((!cellFound) && (n < mesh_.nCells()))
303 if (pointInCell(location, n))
324 Foam::label Foam::meshSearch::findNearestBoundaryFaceWalk
326 const point& location,
327 const label seedFaceI
334 "meshSearch::findNearestBoundaryFaceWalk"
335 "(const point&, const label)"
336 ) << "illegal seedFace:" << seedFaceI << exit(FatalError);
339 // Start off from seedFaceI
341 label curFaceI = seedFaceI;
343 const face& f = mesh_.faces()[curFaceI];
345 scalar minDist = f.nearestPoint
357 // Search through all neighbouring boundary faces by going
360 label lastFaceI = curFaceI;
362 const labelList& myEdges = mesh_.faceEdges()[curFaceI];
364 forAll (myEdges, myEdgeI)
366 const labelList& neighbours = mesh_.edgeFaces()[myEdges[myEdgeI]];
368 // Check any face which uses edge, is boundary face and
369 // is not curFaceI itself.
371 forAll(neighbours, nI)
373 label faceI = neighbours[nI];
377 (faceI >= mesh_.nInternalFaces())
378 && (faceI != lastFaceI)
381 const face& f = mesh_.faces()[faceI];
383 pointHit curHit = f.nearestPoint
389 // If the face is closer, reset current face and distance
390 if (curHit.distance() < minDist)
392 minDist = curHit.distance();
394 closer = true; // a closer neighbour has been found
405 Foam::vector Foam::meshSearch::offset
412 // Get the neighbouring cell
413 label ownerCellI = mesh_.faceOwner()[bFaceI];
415 const point& c = mesh_.cellCentres()[ownerCellI];
417 // Typical dimension: distance from point on face to cell centre
418 scalar typDim = mag(c - bPoint);
420 return tol_*typDim*dir;
424 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
426 // Construct from components
427 Foam::meshSearch::meshSearch(const polyMesh& mesh, const bool faceDecomp)
430 faceDecomp_(faceDecomp),
431 cloud_(mesh_, IDLList<passiveParticle>()),
432 boundaryTreePtr_(NULL),
434 cellCentreTreePtr_(NULL)
438 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
440 Foam::meshSearch::~meshSearch()
446 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
448 const Foam::indexedOctree<Foam::treeDataFace>& Foam::meshSearch::boundaryTree()
451 if (!boundaryTreePtr_)
457 // all boundary faces (not just walls)
458 labelList bndFaces(mesh_.nFaces()-mesh_.nInternalFaces());
461 bndFaces[i] = mesh_.nInternalFaces() + i;
464 treeBoundBox overallBb(mesh_.points());
465 Random rndGen(123456);
466 overallBb = overallBb.extend(rndGen, 1E-4);
467 overallBb.min() -= point(ROOTVSMALL, ROOTVSMALL, ROOTVSMALL);
468 overallBb.max() += point(ROOTVSMALL, ROOTVSMALL, ROOTVSMALL);
470 boundaryTreePtr_ = new indexedOctree<treeDataFace>
472 treeDataFace // all information needed to search faces
474 false, // do not cache bb
476 bndFaces // boundary faces only
478 overallBb, // overall search domain
485 return *boundaryTreePtr_;
489 const Foam::indexedOctree<Foam::treeDataCell>& Foam::meshSearch::cellTree()
498 treeBoundBox overallBb(mesh_.points());
499 Random rndGen(123456);
500 overallBb = overallBb.extend(rndGen, 1E-4);
501 overallBb.min() -= point(ROOTVSMALL, ROOTVSMALL, ROOTVSMALL);
502 overallBb.max() += point(ROOTVSMALL, ROOTVSMALL, ROOTVSMALL);
504 cellTreePtr_ = new indexedOctree<treeDataCell>
508 false, // not cache bb
511 overallBb, // overall search domain
518 return *cellTreePtr_;
523 const Foam::indexedOctree<Foam::treeDataPoint>&
524 Foam::meshSearch::cellCentreTree() const
526 if (!cellCentreTreePtr_)
532 treeBoundBox overallBb(mesh_.cellCentres());
533 Random rndGen(123456);
534 overallBb = overallBb.extend(rndGen, 1E-4);
535 overallBb.min() -= point(ROOTVSMALL, ROOTVSMALL, ROOTVSMALL);
536 overallBb.max() += point(ROOTVSMALL, ROOTVSMALL, ROOTVSMALL);
538 cellCentreTreePtr_ = new indexedOctree<treeDataPoint>
540 treeDataPoint(mesh_.cellCentres()),
541 overallBb, // overall search domain
543 10.0, // maximum ratio of cubes v.s. cells
544 100.0 // max. duplicity; n/a since no bounding boxes.
548 return *cellCentreTreePtr_;
552 // Is the point in the cell
553 // Works by checking if there is a face inbetween the point and the cell
555 // Check for internal uses proper face decomposition or just average normal.
556 bool Foam::meshSearch::pointInCell(const point& p, label cellI) const
560 const point& ctr = mesh_.cellCentres()[cellI];
563 scalar magDir = mag(dir);
565 // Check if any faces are hit by ray from cell centre to p.
566 // If none -> p is in cell.
567 const labelList& cFaces = mesh_.cells()[cellI];
569 // Make sure half_ray does not pick up any faces on the wrong
571 scalar oldTol = intersection::setPlanarTol(0.0);
575 label faceI = cFaces[i];
577 const face& f = mesh_.faces()[faceI];
581 pointHit inter = f.ray
586 intersection::HALF_RAY,
592 scalar dist = inter.distance();
596 // Valid hit. Hit face so point is not in cell.
597 intersection::setPlanarTol(oldTol);
605 intersection::setPlanarTol(oldTol);
607 // No face inbetween point and cell centre so point is inside.
612 const labelList& f = mesh_.cells()[cellI];
613 const labelList& owner = mesh_.faceOwner();
614 const vectorField& cf = mesh_.faceCentres();
615 const vectorField& Sf = mesh_.faceAreas();
619 label nFace = f[facei];
620 vector proj = p - cf[nFace];
621 vector normal = Sf[nFace];
622 if (owner[nFace] == cellI)
624 if ((normal & proj) > 0)
631 if ((normal & proj) < 0)
643 Foam::label Foam::meshSearch::findNearestCell
645 const point& location,
646 const label seedCellI,
647 const bool useTreeSearch
654 return findNearestCellTree(location);
658 return findNearestCellLinear(location);
663 return findNearestCellWalk(location, seedCellI);
668 Foam::label Foam::meshSearch::findNearestFace
670 const point& location,
671 const label seedFaceI,
672 const bool useTreeSearch
679 return findNearestFaceTree(location);
683 return findNearestFaceLinear(location);
688 return findNearestFaceWalk(location, seedFaceI);
693 Foam::label Foam::meshSearch::findCell
695 const point& location,
696 const label seedCellI,
697 const bool useTreeSearch
700 // Find the nearest cell centre to this location
701 label nearCellI = findNearestCell(location, seedCellI, useTreeSearch);
705 Pout<< "findCell : nearCellI:" << nearCellI
706 << " ctr:" << mesh_.cellCentres()[nearCellI]
710 if (pointInCell(location, nearCellI))
718 // Start searching from cell centre of nearCell
719 point curPoint = mesh_.cellCentres()[nearCellI];
721 vector edgeVec = location - curPoint;
722 edgeVec /= mag(edgeVec);
726 // Walk neighbours (using tracking) to get closer
727 passiveParticle tracker(cloud_, curPoint, nearCellI);
731 Pout<< "findCell : tracked from curPoint:" << curPoint
732 << " nearCellI:" << nearCellI;
735 tracker.track(location);
739 Pout<< " to " << tracker.position()
740 << " need:" << location
741 << " onB:" << tracker.onBoundary()
742 << " cell:" << tracker.cell()
743 << " face:" << tracker.face() << endl;
746 if (!tracker.onBoundary())
748 // stopped not on boundary -> reached location
749 return tracker.cell();
752 // stopped on boundary face. Compare positions
753 scalar typDim = sqrt(mag(mesh_.faceAreas()[tracker.face()]));
755 if ((mag(tracker.position() - location)/typDim) < SMALL)
757 return tracker.cell();
760 // tracking stopped at boundary. Find next boundary
761 // intersection from current point (shifted out a little bit)
762 // in the direction of the destination
766 + offset(tracker.position(), tracker.face(), edgeVec);
770 Pout<< "Searching for next boundary from curPoint:"
772 << " to location:" << location << endl;
774 pointIndexHit curHit = intersection(curPoint, location);
777 Pout<< "Returned from line search with ";
788 nearCellI = mesh_.faceOwner()[curHit.index()];
791 + offset(curHit.hitPoint(), curHit.index(), edgeVec);
798 return findCellLinear(location);
805 Foam::label Foam::meshSearch::findNearestBoundaryFace
807 const point& location,
808 const label seedFaceI,
809 const bool useTreeSearch
816 const indexedOctree<treeDataFace>& tree = boundaryTree();
818 scalar span = tree.bb().mag();
820 pointIndexHit info = boundaryTree().findNearest
828 info = boundaryTree().findNearest
835 return tree.shapes().faceLabels()[info.index()];
839 scalar minDist = GREAT;
845 label faceI = mesh_.nInternalFaces();
846 faceI < mesh_.nFaces();
850 const face& f = mesh_.faces()[faceI];
859 if (curHit.distance() < minDist)
861 minDist = curHit.distance();
870 return findNearestBoundaryFaceWalk(location, seedFaceI);
875 Foam::pointIndexHit Foam::meshSearch::intersection
881 pointIndexHit curHit = boundaryTree().findLine(pStart, pEnd);
885 // Change index into octreeData into face label
886 curHit.setIndex(boundaryTree().shapes().faceLabels()[curHit.index()]);
892 Foam::List<Foam::pointIndexHit> Foam::meshSearch::intersections
898 DynamicList<pointIndexHit> hits;
900 vector edgeVec = pEnd - pStart;
901 edgeVec /= mag(edgeVec);
908 bHit = intersection(pt, pEnd);
914 const vector& area = mesh_.faceAreas()[bHit.index()];
916 scalar typDim = Foam::sqrt(mag(area));
918 if ((mag(bHit.hitPoint() - pEnd)/typDim) < SMALL)
923 // Restart from hitPoint shifted a little bit in the direction
924 // of the destination
928 + offset(bHit.hitPoint(), bHit.index(), edgeVec);
940 bool Foam::meshSearch::isInside(const point& p) const
943 boundaryTree().getVolumeType(p)
944 == indexedOctree<treeDataFace>::INSIDE;
948 // Delete all storage
949 void Foam::meshSearch::clearOut()
951 deleteDemandDrivenData(boundaryTreePtr_);
952 deleteDemandDrivenData(cellTreePtr_);
953 deleteDemandDrivenData(cellCentreTreePtr_);
957 void Foam::meshSearch::correct()
963 // ************************************************************************* //