39 scalar treeDataFace::tolSqr =
sqr(1
e-6);
49template<
class ElementIds>
53 const ElementIds& elemIds
65 return treeBoundBox(mesh.points(), mesh.faces()[facei]);
74template<
class ElementIds>
77 const primitiveMesh&
mesh,
78 const ElementIds& elemIds
83 for (
const label facei : elemIds)
85 bb.add(
mesh.points(),
mesh.faces()[facei]);
151Foam::treeDataFace::usesFace(
const label facei)
const
153 return (!useSubset_ || isTreeFace_.test(facei));
157inline Foam::treeBoundBox
158Foam::treeDataFace::getBounds(
const label index)
const
160 const label facei = objectIndex(index);
161 return treeBoundBox(mesh_.points(), mesh_.faces()[facei]);
165void Foam::treeDataFace::update()
235 isTreeFace_(mesh_.nFaces(), faceLabels_),
252 isTreeFace_(mesh_.nFaces(), faceLabels_),
269 for (
const label index : indices)
271 const label facei = faceLabels_[index];
273 bb.
add(mesh_.points(), mesh_.faces()[facei]);
290 return mesh_.faceCentres();
296 const indexedOctree<treeDataFace>& oc,
314 const label index = info.
index();
319 <<
"Could not find " <<
sample <<
" in octree."
325 const label facei = objectIndex(index);
329 Pout<<
"getSampleType : sample:" << sample
330 <<
" nearest face:" << facei;
338 const face&
f = mesh_.faces()[facei];
339 const vector&
area = mesh_.faceAreas()[facei];
340 const point& fc = mesh_.faceCentres()[facei];
343 const point& curPt = curHit.point();
355 Pout<<
" -> face hit:" << curPt
356 <<
" comparing to face normal " <<
area <<
endl;
363 Pout<<
" -> face miss:" << curPt;
371 const scalar typDimSqr =
mag(area) + VSMALL;
373 for (
const label fp :
f)
375 const scalar relDistSqr = (
magSqr(
points[fp] - curPt)/typDimSqr);
377 if (relDistSqr < tolSqr)
386 for (
const label ptFacei : mesh_.pointFaces()[fp])
388 if (usesFace(ptFacei))
390 pointNormal +=
normalised(mesh_.faceAreas()[ptFacei]);
397 <<
" point normal:" << pointNormal
398 <<
" distance:" << relDistSqr <<
endl;
408 const scalar relDistSqr = (
magSqr(fc - curPt)/typDimSqr);
410 if (relDistSqr < tolSqr)
417 Pout<<
" -> centre hit:" << fc
418 <<
" distance:" << relDistSqr <<
endl;
429 for (
const label edgei : mesh_.faceEdges()[facei])
432 mesh_.edges()[edgei].line(
points).nearestDist(sample);
434 const scalar relDistSqr = edgeHit.point().distSqr(curPt)/typDimSqr;
436 if (relDistSqr < tolSqr)
445 for (
const label eFacei : mesh_.edgeFaces()[edgei])
447 if (usesFace(eFacei))
449 edgeNormal +=
normalised(mesh_.faceAreas()[eFacei]);
455 Pout<<
" -> real edge hit point:" << edgeHit.point()
456 <<
" comparing to edge normal:" << edgeNormal
479 const scalar relDistSqr = edgeHit.
point().
distSqr(curPt)/typDimSqr;
481 if (relDistSqr < tolSqr)
495 Pout<<
" -> internal edge hit point:" << edgeHit.point()
496 <<
" comparing to edge normal "
497 << 0.5*(nLeft + nRight)
505 0.5*(nLeft + nRight),
513 Pout<<
"Did not find sample " << sample
514 <<
" anywhere related to nearest face " << facei <<
endl
519 Pout<<
" vertex:" <<
f[fp] <<
" coord:" <<
points[
f[fp]]
545 : !searchBox.
overlaps(getBounds(index))
554 const label facei = objectIndex(index);
556 const face&
f = mesh_.faces()[facei];
573 const point& fc = mesh_.faceCentres()[facei];
597 const scalar radiusSqr
604 ? !bbs_[index].overlaps(centre, radiusSqr)
605 : !getBounds(index).overlaps(centre, radiusSqr)
611 const label facei = objectIndex(index);
617 if (
sqr(nearHit.distance()) < radiusSqr)
651 scalar& nearestDistSqr,
656 for (
const label index : indices)
658 const label facei = objectIndex(index);
665 if (distSqr < nearestDistSqr)
667 nearestDistSqr = distSqr;
669 nearestPoint = nearHit.
point();
675void Foam::treeDataFace::findNearestOp::operator()
680 scalar& nearestDistSqr,
685 tree_.shapes().findNearest
696void Foam::treeDataFace::findNearestOp::operator()
711bool Foam::treeDataFace::findIntersectOp::operator()
716 point& intersectionPoint
735 const vector dir(end - start);
746 if (inter.hit() && inter.distance() <= 1)
750 intersectionPoint = inter.point();
labelList faceLabels(nFaceLabels)
Minimal example by using system/controlDict.functions:
scalar distance() const noexcept
Return distance to hit.
bool hit() const noexcept
Is there a hit.
const point_type & point() const noexcept
Return the point, no checks.
label index() const noexcept
Return the hit index.
iterator begin() noexcept
Return an iterator to begin traversing the UList.
label rcIndex(const label i) const noexcept
The reverse circular index. The previous index in the list which returns to the last at the beginning...
label fcIndex(const label i) const noexcept
The forward circular index. The next index in the list which returns to the first at the end of the l...
scalar distSqr(const Vector< Cmpt > &v2) const
The L2-norm distance squared from another vector. The magSqr() of the difference.
void add(const boundBox &bb)
Extend to include the second box.
bool containsAny(const UList< point > &points) const
Contains any of the points? (inside or on edge).
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface....
A face is a list of labels corresponding to mesh vertices.
Non-pointer based hierarchical recursive searching.
static volumeType getSide(const vector &outsideNormal, const vector &vec)
Helper function to return the side. Returns outside if.
A range or interval of labels defined by a start and a size.
PointHit< Point > nearestDist(const Point &p) const
Return nearest distance to line from a given point.
virtual const faceList & faces() const
Return raw faces.
A patch is a list of labels that address the faces in the global face list.
Cell-face mesh analysis engine.
const vectorField & faceCentres() const
virtual const faceList & faces() const =0
Return faces.
virtual const pointField & points() const =0
Return mesh points.
A class for managing temporary objects.
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Standard boundBox with extra functionality for use in octree.
bool overlaps(const boundBox &bb) const
Overlaps with other bounding box, sphere etc?
bool intersects(const point &overallStart, const vector &overallVec, const point &start, const point &end, point &pt, direction &ptBits) const
Intersects segment; set point to intersection position and face,.
direction posBits(const point &pt) const
Position of point relative to bounding box.
findIntersectOp(const indexedOctree< treeDataFace > &tree)
findNearestOp(const indexedOctree< treeDataFace > &tree)
Encapsulation of data for searching on faces.
tmp< pointField > centres() const
Representative point cloud for contained shapes. One point per shape, corresponding to the face centr...
const labelList & faceLabels() const noexcept
The subset of face ids to use.
void findNearest(const labelUList &indices, const point &sample, scalar &nearestDistSqr, label &nearestIndex, point &nearestPoint) const
Calculates nearest (to sample) point in shape.
static treeBoundBox bounds(const primitiveMesh &mesh, const labelRange &range)
Return bounding box of specified range of faces.
const point & centre(const label index) const
Representative point (face centre) at shape index.
label objectIndex(const label index) const
Map from shape index to original (non-subset) face label.
treeDataFace(const bool cacheBb, const primitiveMesh &mesh)
Construct from mesh, using all faces.
volumeType getVolumeType(const indexedOctree< treeDataFace > &, const point &) const
Get type (inside,outside,mixed,unknown) of point w.r.t. surface.
static treeBoundBoxList boxes(const primitiveMesh &mesh)
Calculate and return bounding boxes for all mesh faces.
const primitiveMesh & mesh() const noexcept
Reference to the supporting mesh.
bool overlaps(const label index, const treeBoundBox &searchBox) const
Does (bb of) shape at index overlap searchBox.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
An enumeration wrapper for classification of a location as being inside/outside of a volume.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Namespace for handling debugging switches.
const wordList area
Standard area field types (scalar, vector, tensor, etc).
PointHit< point > pointHit
A PointHit with a 3D point.
List< label > labelList
A List of labels.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
quaternion normalised(const quaternion &q)
Return the normalised (unit) quaternion of the given quaternion.
static treeBoundBox boundsImpl(const primitiveMesh &mesh, const ElementIds &elemIds)
List< treeBoundBox > treeBoundBoxList
A List of treeBoundBox.
static treeBoundBoxList boxesImpl(const primitiveMesh &mesh, const ElementIds &elemIds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
line< point, const point & > linePointRef
A line using referred points.
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i), works like std::iota() but returning a...
errorManip< error > abort(error &err)
vector point
Point is a vector.
static constexpr const zero Zero
Global zero (0).
triangle< point, const point & > triPointRef
A triangle using referred points.
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
vectorField pointField
pointField is a vectorField.
PointIndexHit< point > pointIndexHit
A PointIndexHit with a 3D point.
UList< label > labelUList
A UList of labels.
bool ln(const fileName &src, const fileName &dst)
Create a softlink. dst should not exist. Returns true if successful.
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
Tree tree(triangles.begin(), triangles.end())
#define forAll(list, i)
Loop across all elements in list.