33 bool Foam::triSurfaceSearch::checkUniqueHit
36 const DynamicList<pointIndexHit, 1, 1>& hits,
44 const labelledTri&
f =
surface()[currHit.index()];
46 f.nearestPointClassify
58 const label nearPointi =
f[nearLabel];
65 const label pointFacei = pointFaces[pI];
67 if (pointFacei != currHit.index())
73 if (hit.index() == pointFacei)
88 const label edgeI = fEdges[nearLabel];
94 const label edgeFacei = edgeFaces[fi];
96 if (edgeFacei != currHit.index())
102 if (hit.index() == edgeFacei)
105 const vector currHitNormal =
108 const vector existingHitNormal =
111 const label signCurrHit =
112 pos0(currHitNormal & lineVec);
114 const label signExistingHit =
115 pos0(existingHitNormal & lineVec);
117 if (signCurrHit == signExistingHit)
154 if (
dict.readIfPresent(
"tolerance", tolerance_) && tolerance_ > 0)
156 Info<<
" using intersection tolerance " << tolerance_ <<
endl;
160 if (
dict.readIfPresent(
"maxTreeDepth", maxTreeDepth_) && maxTreeDepth_ > 0)
162 Info<<
" using maximum tree depth " << maxTreeDepth_ <<
endl;
170 const scalar tolerance,
171 const label maxTreeDepth
175 tolerance_(tolerance),
176 maxTreeDepth_(maxTreeDepth),
200 if (treePtr_.empty())
205 if (surface().size())
213 <<
"Surface does not have compact point numbering."
214 <<
" Of " << surface().points().size()
216 <<
" are used. This might give problems in some routines."
259 if (!tree().bb().contains(sample))
261 inside[sampleI] =
false;
265 inside[sampleI] =
true;
269 inside[sampleI] =
false;
292 info[i] = octree.findNearest
311 const scalar nearestDistSqr = 0.25*
magSqr(span);
313 return tree().findNearest(pt, nearestDistSqr);
333 info[i] = octree.findLine
380 info.setSize(start.
size());
409 vector lineVec = end[i] - start[i];
410 lineVec /=
mag(lineVec) + vSmall;
412 if (checkUniqueHit(inter, hits, lineVec))
425 info[i].transfer(hits);
#define forAll(list, i)
Loop across all elements in list.
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
void clear()
Clear the addressed list, i.e. set the size to zero.
void size(const label)
Override size to be inconsistent with allocated storage.
void setSize(const label)
Reset size of List.
This class describes the interaction of (usually) a face and a point. It carries the info of a succes...
label index() const
Return index.
bool hit() const
Is there a hit.
const Map< label > & meshPointMap() const
Mesh point map. Given the global point index find its.
const labelListList & pointFaces() const
Return point-face addressing.
const labelListList & edgeFaces() const
Return edge-face addressing.
const labelListList & faceEdges() const
Return face-edge addressing.
const Field< PointType > & faceNormals() const
Return face normals for patch.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Non-pointer based hierarchical recursive searching.
static scalar & perturbTol()
Get the perturbation tolerance.
pointIndexHit findLineAny(const point &start, const point &end) const
Find any intersection of line between start and end.
Standard boundBox + extra functionality for use in octree.
treeBoundBox extend(const scalar s) const
Return asymmetrically extended bounding box, with guaranteed.
Encapsulation of data needed to search on PrimitivePatches.
pointIndexHit nearest(const point &, const vector &span) const
Calculate nearest point on surface for single searchPoint. Returns.
boolList calcInside(const pointField &searchPoints) const
Calculate for each searchPoint inside/outside status.
void findLineAll(const pointField &start, const pointField &end, List< List< pointIndexHit >> &info) const
Calculate all intersections from start to end.
const indexedOctree< treeDataTriSurface > & tree() const
Demand driven construction of the octree.
~triSurfaceSearch()
Destructor.
void findLine(const pointField &start, const pointField &end, List< pointIndexHit > &info) const
void findLineAny(const pointField &start, const pointField &end, List< pointIndexHit > &info) const
const triSurface & surface() const
Return reference to the surface.
triSurfaceSearch(const triSurface &)
Construct from surface. Holds reference to surface!
void clearOut()
Clear storage.
void findNearest(const pointField &samples, const scalarField &nearestDistSqr, List< pointIndexHit > &info) const
Triangulated surface description with patch information.
#define WarningInFunction
Report a warning using Foam::Warning.
treeDataPrimitivePatch< triSurface > treeDataTriSurface
List< label > labelList
A List of labels.
dimensionedScalar pos0(const dimensionedScalar &ds)
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
PointIndexHit< point > pointIndexHit
Ostream & endl(Ostream &os)
Add newline and flush stream.
Vector< scalar > vector
A scalar version of the templated Vector.
dimensioned< scalar > mag(const dimensioned< Type > &)
dimensioned< scalar > magSqr(const dimensioned< Type > &)
scalarField samples(nIntervals, 0)