helperFunctionsGeometryQueries.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | cfMesh: A library for mesh generation
4  \\ / O peration |
5  \\ / A nd | Author: Franjo Juretic (franjo.juretic@c-fields.com)
6  \\/ M anipulation | Copyright (C) Creative Fields, Ltd.
7 -------------------------------------------------------------------------------
8 License
9  This file is part of cfMesh.
10 
11  cfMesh 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 3 of the License, or (at your
14  option) any later version.
15 
16  cfMesh 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 cfMesh. If not, see <http://www.gnu.org/licenses/>.
23 
24 InNamespace
25  Foam::help
26 
27 Description
28  Geometry queries useful for mesh generation
29 
30 SourceFiles
31 
32 \*---------------------------------------------------------------------------*/
33 
34 #ifndef helperFunctionsGeometryQueries_H
35 #define helperFunctionsGeometryQueries_H
36 
37 #include "DynList.H"
38 #include "plane.H"
39 #include "face.H"
40 
41 #include "triSurf.H"
42 #include "triangle.H"
43 #include "tetrahedron.H"
44 #include "boolList.H"
45 #include "HashSet.H"
46 
47 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
48 
49 class boundBox;
50 
51 namespace Foam
52 {
53 
54 /*---------------------------------------------------------------------------*\
55  Namespace help functions Declaration
56 \*---------------------------------------------------------------------------*/
57 
58 namespace help
59 {
60  //- check if a list has nan entries
61  template<class ListType>
62  bool isnan(const ListType&);
63 
64  //- check if a list has inf entries
65  template<class ListType>
66  bool isinf(const ListType&);
67 
68  //- check if the faces share a convex edge
69  template<class Face1, class Face2>
70  inline bool isSharedEdgeConvex
71  (
72  const pointField& points,
73  const Face1& f1,
74  const Face2& f2
75  );
76 
77  //- angle between the two faces in radians
78  template<class Face1, class Face2>
79  inline scalar angleBetweenFaces
80  (
81  const pointField& points,
82  const Face1& f1,
83  const Face2& f2
84  );
85 
86  // merge faces which are in the same patch
87  // This is need to merge the triagles which are generated for cells
88  // cut by more than on boundary regions
90  (
91  const List< DynList<label> >& pfcs,
92  const pointField& polyPoints
93  );
94 
95  //- check if the point p belongs to the edge e
96  inline bool vertexOnLine
97  (
98  const point& p,
99  const edge& e,
100  const pointField& ep
101  );
102 
103  //- check if the point p belongs to the plane
104  inline bool vertexInPlane(const point& p, const plane& pl);
105 
106  //- find the vertex on the line of the edge nearest to the point p
108  (
109  const point& edgePoint0,
110  const point& edgePoint1,
111  const point& p
112  );
113 
114  //- find the vertex on the edge nearest to the point p
116  (
117  const point& edgePoint0,
118  const point& edgePoint1,
119  const point& p
120  );
121 
122  //- find and return the distance between the edge and the point p
123  inline scalar distanceOfPointFromTheEdge
124  (
125  const point& edgePoint0,
126  const point& edgePoint1,
127  const point& p
128  );
129 
130  //- find the nearest points on the edge and the line
131  inline bool nearestEdgePointToTheLine
132  (
133  const point& edgePoint0,
134  const point& edgePoint1,
135  const point& lp0,
136  const point& lp1,
137  point& nearestOnEdge,
138  point& nearestOnLine
139  );
140 
141  //- check if the edge intersects the plane
142  inline bool planeIntersectsEdge
143  (
144  const point& start,
145  const point& end,
146  const plane& pl,
147  point& intersection
148  );
149 
150  //- check if a vertex lies inside the tetrahedron
151  inline bool pointInTetrahedron
152  (
153  const point& p,
154  const tetrahedron<point, point>& tet
155  );
156 
157  //- check if a line intersects the triangle, and return the intersection
158  inline bool triLineIntersection
159  (
160  const triangle<point, point>& tria,
161  const point& lineStart,
162  const point& lineEnd,
163  point& intersection
164  );
165 
166  //- check if a line intersects the triangle and return the intersection
167  inline bool triLineIntersection
168  (
169  const triSurf&,
170  const label,
171  const point&,
172  const point&,
173  point&
174  );
175 
176  //- check if the line intersects the bounding box
177  inline bool boundBoxLineIntersection
178  (
179  const point&,
180  const point&,
181  const boundBox&
182  );
183 
184  //- check if the line and the face intersect
185  inline bool lineFaceIntersection
186  (
187  const point&,
188  const point&,
189  const face&,
190  const pointField& fp,
191  point& intersection
192  );
193 
194  //- check if the surface triangle and the face intersect
195  inline bool doFaceAndTriangleIntersect
196  (
197  const triSurf& surface,
198  const label triI,
199  const face& f,
200  const pointField& facePoints
201  );
202 
203  //- find the nearest point on the triangle to the given point
205  (
206  const triangle<point, point>& tri,
207  const point&
208  );
209 
210  //- find the nearest vertex on the surface triangle to the given point
212  (
213  const label,
214  const triSurf&,
215  const point&
216  );
217 
218  //- find the minimiser point from a point and the given planes
219  //- returns true if the minimizer exists
220  inline bool findMinimizerPoint
221  (
222  const DynList<point>& origins,
223  const DynList<vector>& normals,
224  point& pMin
225  );
226 
227  //- check the existence of overlap between the two edges
228  inline bool doEdgesOverlap
229  (
230  const point& e0p0,
231  const point& e0p1,
232  const point& e1p0,
233  const point& e1p1,
234  FixedList<point, 2>& overlappingPart,
235  const scalar distTol = -1.0,
236  const scalar cosTol = Foam::cos(5.0*(M_PI/180.0)) // cosine tolerance
237  );
238 
239  //- check the existence of overlap between the two triangles
240  inline bool doTrianglesOverlap
241  (
242  const triangle<point, point>& tri0,
243  const triangle<point, point>& tri1,
244  DynList<point>& overlappingPolygon,
245  const scalar distTol = -1.0,
246  const scalar cosTol = Foam::cos(5.0*(M_PI/180.0)) // cosine tolerance
247  );
248 
249  //- check the existence of intersection between the two triangles
250  inline bool doTrianglesIntersect
251  (
252  const triangle<point, point> &tri0,
253  const triangle<point, point> &tri1,
254  const scalar distTol = -1.0
255  );
256 
257  inline bool doTrianglesIntersect
258  (
259  const triangle<point, point>& tri0,
260  const triangle<point, point>& tri1,
261  DynList<point>& intersectionPoints,
262  const scalar distTol = -1.0
263  );
264 
265  //- check if the point is inside or outside the face
266  inline bool pointInsideFace
267  (
268  const point& p,
269  const face& f,
270  const vector& n,
271  const pointField& fp,
272  const scalar distTol = SMALL
273  );
274 
275  //- check if the point is inside or outside the face
276  inline bool pointInsideFace
277  (
278  const point& p,
279  const face& f,
280  const pointField& fp,
281  const scalar distTol = SMALL
282  );
283 
284  //- check if the face is convex. Concave points are flagged false
285  inline bool isFaceConvexAndOk
286  (
287  const face& f,
288  const pointField& fp,
289  DynList<bool>& OkPoints
290  );
291 
292  //- check if the vertex is on the positive side of the face plane
293  inline bool isVertexVisible(const point& p, const plane& pl);
294 
295  //- find number of face groups within a given range
297  (
298  const labelHashSet& containedElements,
299  const point& centre,
300  const scalar range,
301  const triSurf& surface
302  );
303 
304  //- find the number of edge groups within the given range
306  (
307  const labelHashSet& containedEdges,
308  const point& centre,
309  const scalar range,
310  const triSurf& surface
311  );
312 
313 } // End namespace help
314 
315 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
316 
317 } // End namespace Foam
318 
319 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
320 
322 
323 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
324 
325 #endif
326 
327 // ************************************************************************* //
Foam::pointField
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
Foam::help::numberOfFaceGroups
label numberOfFaceGroups(const labelHashSet &containedElements, const point &centre, const scalar range, const triSurf &surface)
find number of face groups within a given range
Definition: helperFunctionsGeometryQueriesI.H:1670
Foam::help::doFaceAndTriangleIntersect
bool doFaceAndTriangleIntersect(const triSurf &surface, const label triI, const face &f, const pointField &facePoints)
check if the surface triangle and the face intersect
Definition: helperFunctionsGeometryQueriesI.H:805
triSurf.H
boolList.H
Foam::help::doTrianglesIntersect
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
Definition: helperFunctionsGeometryQueriesI.H:1149
p
p
Definition: pEqn.H:62
Foam::help::isVertexVisible
bool isVertexVisible(const point &p, const plane &pl)
check if the vertex is on the positive side of the face plane
Foam::help::vertexOnLine
bool vertexOnLine(const point &p, const edge &e, const pointField &ep)
check if the point p belongs to the edge e
Definition: helperFunctionsGeometryQueriesI.H:308
Foam::help::planeIntersectsEdge
bool planeIntersectsEdge(const point &start, const point &end, const plane &pl, point &intersection)
check if the edge intersects the plane
Definition: helperFunctionsGeometryQueriesI.H:338
Foam::help::nearestPointOnTheEdge
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
Definition: helperFunctionsGeometryQueriesI.H:1616
Foam::help::isnan
bool isnan(const ListType &)
check if a list has nan entries
Definition: helperFunctionsGeometryQueriesI.H:53
triangle.H
Foam::help::distanceOfPointFromTheEdge
scalar distanceOfPointFromTheEdge(const point &edgePoint0, const point &edgePoint1, const point &p)
find and return the distance between the edge and the point p
Definition: helperFunctionsGeometryQueriesI.H:1660
DynList
A dynamic list is a 1-D vector of objects of type T which resizes itself as necessary to accept the n...
Foam::help::nearestPointOnTheTriangle
point nearestPointOnTheTriangle(const triangle< point, point > &tri, const point &)
find the nearest point on the triangle to the given point
Definition: helperFunctionsGeometryQueriesI.H:475
pMin
dimensionedScalar pMin("pMin", dimPressure, fluid)
face.H
Foam::help::lineFaceIntersection
bool lineFaceIntersection(const point &, const point &, const face &, const pointField &fp, point &intersection)
check if the line and the face intersect
Definition: helperFunctionsGeometryQueriesI.H:778
Foam::help::nearestEdgePointToTheLine
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
Definition: helperFunctionsGeometryQueriesI.H:417
Foam::help::mergePatchFaces
faceList mergePatchFaces(const List< DynList< label > > &pfcs, const pointField &polyPoints)
Definition: helperFunctionsGeometryQueriesI.H:232
Foam::help::isinf
bool isinf(const ListType &)
check if a list has inf entries
Definition: helperFunctionsGeometryQueriesI.H:63
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::help::pointInsideFace
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
Definition: helperFunctionsGeometryQueriesI.H:1491
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
plane.H
f1
scalar f1
Definition: createFields.H:28
Foam::help::nearestPointOnTheEdgeExact
point nearestPointOnTheEdgeExact(const point &edgePoint0, const point &edgePoint1, const point &p)
find the vertex on the edge nearest to the point p
Definition: helperFunctionsGeometryQueriesI.H:1633
HashSet.H
Foam::help::isSharedEdgeConvex
bool isSharedEdgeConvex(const pointField &points, const Face1 &f1, const Face2 &f2)
check if the faces share a convex edge
Definition: helperFunctionsGeometryQueriesI.H:74
Foam::help::isFaceConvexAndOk
bool isFaceConvexAndOk(const face &f, const pointField &fp, DynList< bool > &OkPoints)
check if the face is convex. Concave points are flagged false
Definition: helperFunctionsGeometryQueriesI.H:1553
Foam::help::triLineIntersection
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
Definition: helperFunctionsGeometryQueriesI.H:653
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::e
const double e
Elementary charge.
Definition: doubleFloat.H:94
Foam::vector
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
Foam::help::findMinimizerPoint
bool findMinimizerPoint(const DynList< point > &origins, const DynList< vector > &normals, point &pMin)
Definition: helperFunctionsGeometryQueriesI.H:614
helperFunctionsGeometryQueriesI.H
Foam::help::pointInTetrahedron
bool pointInTetrahedron(const point &p, const tetrahedron< point, point > &tet)
check if a vertex lies inside the tetrahedron
Definition: helperFunctionsGeometryQueriesI.H:365
Foam::help::angleBetweenFaces
scalar angleBetweenFaces(const pointField &points, const Face1 &f1, const Face2 &f2)
angle between the two faces in radians
Definition: helperFunctionsGeometryQueriesI.H:132
Foam::help::numberOfEdgeGroups
label numberOfEdgeGroups(const labelHashSet &containedEdges, const point &centre, const scalar range, const triSurf &surface)
find the number of edge groups within the given range
Definition: helperFunctionsGeometryQueriesI.H:1758
Foam::help::boundBoxLineIntersection
bool boundBoxLineIntersection(const point &, const point &, const boundBox &)
check if the line intersects the bounding box
Definition: helperFunctionsGeometryQueriesI.H:723
Foam::help::doTrianglesOverlap
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
Definition: helperFunctionsGeometryQueriesI.H:957
range
scalar range
Definition: LISASMDCalcMethod1.H:12
f
labelList f(nPoints)
Foam::faceList
List< face > faceList
Definition: faceListFwd.H:43
triSurf
A class for triangulated surface used in the meshing process. It is derived from points and facets wi...
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::help::vertexInPlane
bool vertexInPlane(const point &p, const plane &pl)
check if the point p belongs to the plane
Definition: helperFunctionsGeometryQueriesI.H:322
List
Definition: Test.C:19
Foam::point
vector point
Point is a vector.
Definition: point.H:41
Foam::labelHashSet
HashSet< label, Hash< label > > labelHashSet
A HashSet with label keys.
Definition: HashSet.H:210
Foam::help::doEdgesOverlap
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
Definition: helperFunctionsGeometryQueriesI.H:872
DynList.H
tetrahedron.H
Foam::cos
dimensionedScalar cos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:256