Go to the documentation of this file.
34 #ifndef helperFunctionsGeometryQueries_H
35 #define helperFunctionsGeometryQueries_H
61 template<
class ListType>
62 bool isnan(
const ListType&);
65 template<
class ListType>
66 bool isinf(
const ListType&);
69 template<
class Face1,
class Face2>
78 template<
class Face1,
class Face2>
109 const point& edgePoint0,
110 const point& edgePoint1,
117 const point& edgePoint0,
118 const point& edgePoint1,
125 const point& edgePoint0,
126 const point& edgePoint1,
133 const point& edgePoint0,
134 const point& edgePoint1,
137 point& nearestOnEdge,
154 const tetrahedron<point, point>& tet
160 const triangle<point, point>& tria,
161 const point& lineStart,
162 const point& lineEnd,
206 const triangle<point, point>& tri,
234 FixedList<point, 2>& overlappingPart,
235 const scalar distTol = -1.0,
236 const scalar cosTol =
Foam::cos(5.0*(M_PI/180.0))
242 const triangle<point, point>& tri0,
243 const triangle<point, point>& tri1,
245 const scalar distTol = -1.0,
246 const scalar cosTol =
Foam::cos(5.0*(M_PI/180.0))
252 const triangle<point, point> &tri0,
253 const triangle<point, point> &tri1,
254 const scalar distTol = -1.0
259 const triangle<point, point>& tri0,
260 const triangle<point, point>& tri1,
262 const scalar distTol = -1.0
272 const scalar distTol = SMALL
281 const scalar distTol = SMALL
vectorField pointField
pointField is a vectorField.
label numberOfFaceGroups(const labelHashSet &containedElements, const point ¢re, const scalar range, const triSurf &surface)
find number of face groups within a given range
bool doFaceAndTriangleIntersect(const triSurf &surface, const label triI, const face &f, const pointField &facePoints)
check if the surface triangle and the face intersect
bool doTrianglesIntersect(const triangle< point, point > &tri0, const triangle< point, point > &tri1, const scalar distTol=-1.0)
check the existence of intersection between the two triangles
bool isVertexVisible(const point &p, const plane &pl)
check if the vertex is on the positive side of the face plane
bool vertexOnLine(const point &p, const edge &e, const pointField &ep)
check if the point p belongs to the edge e
bool planeIntersectsEdge(const point &start, const point &end, const plane &pl, point &intersection)
check if the edge intersects the plane
point nearestPointOnTheEdge(const point &edgePoint0, const point &edgePoint1, const point &p)
find the vertex on the line of the edge nearest to the point p
bool isnan(const ListType &)
check if a list has nan entries
scalar distanceOfPointFromTheEdge(const point &edgePoint0, const point &edgePoint1, const point &p)
find and return the distance between the edge and the point p
A dynamic list is a 1-D vector of objects of type T which resizes itself as necessary to accept the n...
point nearestPointOnTheTriangle(const triangle< point, point > &tri, const point &)
find the nearest point on the triangle to the given point
dimensionedScalar pMin("pMin", dimPressure, fluid)
bool lineFaceIntersection(const point &, const point &, const face &, const pointField &fp, point &intersection)
check if the line and the face intersect
bool nearestEdgePointToTheLine(const point &edgePoint0, const point &edgePoint1, const point &lp0, const point &lp1, point &nearestOnEdge, point &nearestOnLine)
find the nearest points on the edge and the line
faceList mergePatchFaces(const List< DynList< label > > &pfcs, const pointField &polyPoints)
bool isinf(const ListType &)
check if a list has inf entries
bool pointInsideFace(const point &p, const face &f, const vector &n, const pointField &fp, const scalar distTol=SMALL)
check if the point is inside or outside the face
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
point nearestPointOnTheEdgeExact(const point &edgePoint0, const point &edgePoint1, const point &p)
find the vertex on the edge nearest to the point p
bool isSharedEdgeConvex(const pointField &points, const Face1 &f1, const Face2 &f2)
check if the faces share a convex edge
bool isFaceConvexAndOk(const face &f, const pointField &fp, DynList< bool > &OkPoints)
check if the face is convex. Concave points are flagged false
bool triLineIntersection(const triangle< point, point > &tria, const point &lineStart, const point &lineEnd, point &intersection)
check if a line intersects the triangle, and return the intersection
const double e
Elementary charge.
Vector< scalar > vector
A scalar version of the templated Vector.
bool findMinimizerPoint(const DynList< point > &origins, const DynList< vector > &normals, point &pMin)
bool pointInTetrahedron(const point &p, const tetrahedron< point, point > &tet)
check if a vertex lies inside the tetrahedron
scalar angleBetweenFaces(const pointField &points, const Face1 &f1, const Face2 &f2)
angle between the two faces in radians
label numberOfEdgeGroups(const labelHashSet &containedEdges, const point ¢re, const scalar range, const triSurf &surface)
find the number of edge groups within the given range
bool boundBoxLineIntersection(const point &, const point &, const boundBox &)
check if the line intersects the bounding box
bool doTrianglesOverlap(const triangle< point, point > &tri0, const triangle< point, point > &tri1, DynList< point > &overlappingPolygon, const scalar distTol=-1.0, const scalar cosTol=Foam::cos(5.0 *(M_PI/180.0)))
check the existence of overlap between the two triangles
A class for triangulated surface used in the meshing process. It is derived from points and facets wi...
bool vertexInPlane(const point &p, const plane &pl)
check if the point p belongs to the plane
vector point
Point is a vector.
HashSet< label, Hash< label > > labelHashSet
A HashSet with label keys.
bool doEdgesOverlap(const point &e0p0, const point &e0p1, const point &e1p0, const point &e1p1, FixedList< point, 2 > &overlappingPart, const scalar distTol=-1.0, const scalar cosTol=Foam::cos(5.0 *(M_PI/180.0)))
check the existence of overlap between the two edges
dimensionedScalar cos(const dimensionedScalar &ds)