triSurfaceTools.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2020 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
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
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
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
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 Class
25  Foam::triSurfaceTools
26 
27 Description
28  A collection of tools for triSurface.
29 
30 SourceFiles
31  triSurfaceTools.C
32 
33 \*---------------------------------------------------------------------------*/
34 
35 #ifndef triSurfaceTools_H
36 #define triSurfaceTools_H
37 
38 #include "boolList.H"
39 #include "pointField.H"
40 #include "DynamicList.H"
41 #include "HashSet.H"
42 #include "FixedList.H"
43 #include "vector2D.H"
44 #include "triPointRef.H"
45 #include "surfaceLocation.H"
46 
47 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
48 
49 namespace Foam
50 {
51 
52 // Forward declaration of classes
53 class triSurface;
54 class edge;
55 class labelledTri;
56 class polyBoundaryMesh;
57 class plane;
58 class boundBox;
59 
60 /*---------------------------------------------------------------------------*\
61  Class triSurfaceTools Declaration
62 \*---------------------------------------------------------------------------*/
63 
64 class triSurfaceTools
65 {
66  // Private Member Functions
67 
68  // Refinement
69 
70  enum refineType
71  {
72  NONE,
73  RED,
74  GREEN
75  };
76  static void calcRefineStatus
77  (
78  const triSurface& surf,
79  const label facei,
80  List<refineType>& refine
81  );
82  static void greenRefine
83  (
84  const triSurface& surf,
85  const label facei,
86  const label edgeI,
87  const label newPointi,
88  DynamicList<labelledTri>& newFaces
89  );
90  static triSurface doRefine
91  (
92  const triSurface& surf,
93  const List<refineType>& refineStatus
94  );
95 
96 
97  // Coarsening
98 
99  static scalar faceCosAngle
100  (
101  const point& pStart,
102  const point& pEnd,
103  const point& pLeft,
104  const point& pRight
105  );
106 
107  static void protectNeighbours
108  (
109  const triSurface& surf,
110  const label vertI,
111  labelList& faceStatus
112  );
113 
114  //- Faces to collapse because of edge collapse
115  static labelHashSet getCollapsedFaces
116  (
117  const triSurface& surf,
118  label edgeI
119  );
120 
121  // Return value of faceUsed for faces using vertI (local numbering).
122  // Used internally.
123  static label vertexUsesFace
124  (
125  const triSurface& surf,
126  const labelHashSet& faceUsed,
127  const label vertI
128  );
129 
130  // Get new connections between faces (because of edge collapse)
131  // in form of tables:
132  // - given edge get other edge
133  // - given edge get other face
134  // A face using point v1 on edge will get connected to a face using
135  // point v2 if they share a common vertex
136  // (but not a common edge since then the triangles collapse to
137  // nothing)
138  static void getMergedEdges
139  (
140  const triSurface& surf,
141  const label edgeI,
142  const labelHashSet& collapsedFaces,
143  HashTable<label, label, Hash<label>>& edgeToEdge,
144  HashTable<label, label, Hash<label>>& edgeToFace
145  );
146 
147  //- Calculates (cos of) angle across edgeI of facei,
148  // taking into account updated addressing (resulting from edge
149  // collapse)
150  static scalar edgeCosAngle
151  (
152  const triSurface& surf,
153  const label v1,
154  const point& pt,
155  const labelHashSet& collapsedFaces,
156  const HashTable<label, label, Hash<label>>& edgeToEdge,
157  const HashTable<label, label, Hash<label>>& edgeToFace,
158  const label facei,
159  const label edgeI
160  );
161 
162  //- Calculate minimum (cos of) edge angle using addressing from
163  // collapsing
164  // edge to v1 at pt. Returns 1 if v1 is on edge without neighbours
165  // (and hence no edge angle can be defined)
166  static scalar collapseMinCosAngle
167  (
168  const triSurface& surf,
169  const label v1,
170  const point& pt,
171  const labelHashSet& collapsedFaces,
172  const HashTable<label, label, Hash<label>>& edgeToEdge,
173  const HashTable<label, label, Hash<label>>& edgeToFace
174  );
175 
176  //- Like collapseMinCosAngle but return true for value < minCos
177  bool collapseCreatesFold
178  (
179  const triSurface& surf,
180  const label v1,
181  const point& pt,
182  const labelHashSet& collapsedFaces,
183  const HashTable<label, label, Hash<label>>& edgeToEdge,
184  const HashTable<label, label, Hash<label>>& edgeToFace,
185  const scalar minCos
186  );
187 
190  // static bool collapseCreatesDuplicates
191  //(
192  // const triSurface& surf,
193  // const label edgeI,
194  // const HashTable<bool, label, Hash<label>>& collapsedFaces
195  //);
196 
197  // Tracking
198 
199  //- Finds the triangle edge/point cut by the plane between
200  // a point inside/on edge of a triangle and a point outside.
201  // Returns
202  // - location on edge/point and hit()
203  // - or miss() if no intersection found
204  static surfaceLocation cutEdge
205  (
206  const triSurface& s,
207  const label triI,
208  const label excludeEdgeI,
209  const label excludePointi,
210  const point& triPoint,
211  const plane& cutPlane,
212  const point& toPoint
213  );
214 
215  //- Checks if current is on the same triangle as the endpoint
216  // and shifts it there. If so updates current and sets a hit.
217  static void snapToEnd
218  (
219  const triSurface& s,
220  const surfaceLocation& endInfo,
221  surfaceLocation& current
222  );
223 
224  //- Visits faces eFaces around start. Does not visit triangle
225  // start.triangle() nor edge excludeEdgeI.
226  // Returns edge, triangle (if more than one choice) which gets
227  // us nearer endpoint.
228  // Returns
229  // - hit() if triangle contains endpoint
230  // - triangle()=-1 if no triangle found
231  // - nearest triangle/edge otherwise
232  static surfaceLocation visitFaces
233  (
234  const triSurface& s,
235  const labelList& eFaces,
236  const surfaceLocation& start,
237  const label excludeEdgeI,
238  const label excludePointi,
239  const surfaceLocation& end,
240  const plane& cutPlane
241  );
242 
243 
244 public:
245 
246  // OBJ writing
247 
248  //- Write pointField to OBJ format file
249  static void writeOBJ
250  (
251  const fileName& fName,
252  const pointField& pts
253  );
254 
255  //- Write vertex subset to OBJ format file
256  static void writeOBJ
257  (
258  const triSurface& surf,
259  const fileName& fName,
260  const boolList& markedVerts
261  );
262 
263 
264  // Additional addressing
265 
266  //- Get all triangles using edge endpoint
267  static void getVertexTriangles
268  (
269  const triSurface& surf,
270  const label edgeI,
271  labelList& edgeTris
272  );
273 
274  //- Get all vertices (local numbering) connected to vertices of edge
276  (
277  const triSurface& surf,
278  const edge& e
279  );
280 
282  // static void orderVertices
283  //(
284  // const labelledTri& f,
285  // const label v1,
286  // const label v2,
287  // label& vA,
288  // label& vB
289  //);
290 
291  //- Get face connected to edge not facei
292  static label otherFace
293  (
294  const triSurface& surf,
295  const label facei,
296  const label edgeI
297  );
298 
299  //- Get the two edges on facei counterclockwise after edgeI
300  static void otherEdges
301  (
302  const triSurface& surf,
303  const label facei,
304  const label edgeI,
305  label& e1,
306  label& e2
307  );
308 
309  //- Get the two vertices (local numbering) on facei counterclockwise
310  // vertI
311  static void otherVertices
312  (
313  const triSurface& surf,
314  const label facei,
315  const label vertI,
316  label& vert1I,
317  label& vert2I
318  );
319 
320  //- Get edge opposite vertex (local numbering)
321  static label oppositeEdge
322  (
323  const triSurface& surf,
324  const label facei,
325  const label vertI
326  );
327 
328  //- Get vertex (local numbering) opposite edge
329  static label oppositeVertex
330  (
331  const triSurface& surf,
332  const label facei,
333  const label edgeI
334  );
335 
336  //- Returns edge label connecting v1, v2 (local numbering)
337  static label getEdge
338  (
339  const triSurface& surf,
340  const label vert1I,
341  const label vert2I
342  );
343 
344  //- Return index of triangle (or -1) using all three edges
345  static label getTriangle
346  (
347  const triSurface& surf,
348  const label e0I,
349  const label e1I,
350  const label e2I
351  );
352 
353  // Coarsening
354 
355  //- Create new triSurface by collapsing edges to edge mids.
357  (
358  const triSurface& surf,
359  const labelList& collapsableEdges
360  );
361 
362 
363  //- Face collapse status.
364  // anyEdge: any edge can be collapsed
365  // noEdge: no edge can be collapsed
366  // collapsed: already collapsed
367  // >0: edge label that can be collapsed
368  static const label ANYEDGE;
369  static const label NOEDGE;
370  static const label COLLAPSED;
371 
372  //- Create new triSurface by collapsing edges to specified
373  // positions. faceStatus allows
374  // explicit control over which faces need to be protected (see above).
375  // faceStatus gets updated to protect collapsing already collapsed
376  // faces.
378  (
379  const triSurface& surf,
380  const labelList& collapsableEdges,
381  const pointField& edgeMids,
382  labelList& faceStatus
383  );
384 
385 
386  // Refinement
387 
388  //- Refine edges by splitting to opposite vertex
389  static triSurface greenRefine
390  (
391  const triSurface& surf,
392  const labelList& refineEdges
393  );
394 
395  //- Refine face by splitting all edges. Neighbouring face is
396  // greenRefine'd.
398  (
399  const triSurface& surf,
400  const labelList& refineFaces
401  );
402 
403 
404  // Geometric
405 
406  //- Returns element in edgeIndices with minimum length
407  static label minEdge
408  (
409  const triSurface& surf,
410  const labelList& edgeIndices
411  );
412 
413  //- Returns element in edgeIndices with minimum length
414  static label maxEdge
415  (
416  const triSurface& surf,
417  const labelList& edgeIndices
418  );
419 
420  //- Merge points within distance
421  static triSurface mergePoints
422  (
423  const triSurface& surf,
424  const scalar mergeTol
425  );
426 
427  //- Triangle (unit) normal. If nearest point to triangle on edge use
428  // edge normal (calculated on the fly); if on vertex use vertex normal.
429  // Uses planarTol.
430  static vector surfaceNormal
431  (
432  const triSurface& surf,
433  const label nearestFacei,
434  const point& nearestPt
435  );
436 
437  //- On which side of surface
438  enum sideType
439  {
440  UNKNOWN, // cannot be determined (e.g. non-manifold)
441  INSIDE, // inside
442  OUTSIDE // outside
443  };
444 
445  //- If nearest point is on edgeI, determine on which side of surface
446  // sample is.
447  static sideType edgeSide
448  (
449  const triSurface& surf,
450  const point& sample,
451  const point& nearestPoint,
452  const label edgeI
453  );
454 
455  //- Given nearest point (to sample) on surface determines which side
456  // sample is. Uses either face normal, edge normal or point normal
457  // (non-trivial). Uses triangle::classify.
458  static sideType surfaceSide
459  (
460  const triSurface& surf,
461  const point& sample,
462  const label nearestFacei
463  );
464 
465  // Triangulation of faces
466 
467  //- Simple triangulation of (selected patches of) boundaryMesh. Needs
468  // polyMesh (or polyBoundaryMesh) since only at this level are the
469  // triangles on neighbouring patches connected.
470  static triSurface triangulate
471  (
472  const polyBoundaryMesh& mBesh,
474  const bool verbose = false
475  );
476 
477 
478  static triSurface triangulate
479  (
480  const polyBoundaryMesh& bMesh,
482  const boundBox& bBox,
483  const bool verbose = false
484  );
485 
486 
487  //- Face-centre triangulation of (selected patches of) boundaryMesh.
488  // Needs
489  // polyMesh (or polyBoundaryMesh) since only at this level are the
490  // triangles on neighbouring patches connected.
492  (
493  const polyBoundaryMesh& mBesh,
495  const bool verbose = false
496  );
497 
498 
499  // Triangulation and interpolation
500 
501  //- Calculate linear interpolation weights for point (guaranteed to be
502  // inside triangle)
503  static void calcInterpolationWeights
504  (
505  const triPointRef&,
506  const point&,
507  FixedList<scalar, 3>& weights
508  );
509 
510  // Calculate weighting factors from samplePts to triangle it is in.
511  // Uses linear search to find triangle.
512  // Vertices are:
513  // (a b c) : vertices of the triangle abc the point is in
514  // or if the point is outside all triangles:
515  // (a b -1) : the edge ab the point is nearest to.
516  // (a -1 -1) : the vertex a the point is nearest to
517  static void calcInterpolationWeights
518  (
519  const triSurface& s,
520  const pointField& samplePts,
521  List<FixedList<label, 3>>& verts,
522  List<FixedList<scalar, 3>>& weights
523  );
524 
525  //- Do unconstrained Delaunay of points. Returns triSurface with 3D
526  // points with z=0. All triangles in region 0.
527  static triSurface delaunay2D(const List<vector2D>&);
528 
529 
530  // Tracking
531 
532  //- Test point on plane of triangle to see if on edge or point or inside
534  (
535  const triSurface&,
536  const label triI,
537  const point& trianglePoint
538  );
539 
540  //- Track on surface to get closer to point.
541  // Possible situations:
542  // - 1. reached endpoint
543  // - 2. reached edge (normal situation)
544  // - 3. reached end of surface (edge on single face)
545  // Input:
546  // - starting position+triangle/edge/point (so has to be on surface!)
547  // - (optional) previous position+triangle/edge/point to prevent
548  // going back. Set index (of triangle/edge/point) to -1 if not
549  // used.
550  // - end position+triangle/edge/point (so has to be on surface!)
551  // - plane to follow. Has to go through end point!
552  // Returns:
553  // - true if end point reached (situation 1)
554  // - new position+triangle/edge/point
555  // Caller has to check for situation 3 by checking that triangle()
556  // is not set.
558  (
559  const triSurface&,
560  const surfaceLocation& start,
561  const surfaceLocation& end,
562  const plane& cutPlane
563  );
564 
565  //- Track from edge to edge across surface. Uses trackToEdge.
566  // Not really useful by itself, more example of how to use trackToEdge.
567  // endInfo should be location on surface.
568  // hitInfo should be initialised to starting location (on surface as
569  // well). Upon return is set to end location.
570  static void track
571  (
572  const triSurface&,
573  const surfaceLocation& endInfo,
574  const plane& cutPlane,
575  surfaceLocation& hitInfo
576  );
577 };
578 
579 
580 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
581 
582 } // End namespace Foam
583 
584 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
585 
586 #endif
587 
588 // ************************************************************************* //
static labelList getVertexVertices(const triSurface &surf, const edge &e)
Get all vertices (local numbering) connected to vertices of edge.
A triangle primitive used to calculate face areas and swept volumes.
Definition: triangle.H:57
A class for handling file names.
Definition: fileName.H:79
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
static void track(const triSurface &, const surfaceLocation &endInfo, const plane &cutPlane, surfaceLocation &hitInfo)
Track from edge to edge across surface. Uses trackToEdge.
A bounding box defined in terms of the points at its extremities.
Definition: boundBox.H:58
static triSurface redGreenRefine(const triSurface &surf, const labelList &refineFaces)
Refine face by splitting all edges. Neighbouring face is.
static label maxEdge(const triSurface &surf, const labelList &edgeIndices)
Returns element in edgeIndices with minimum length.
Geometric class that creates a 2D plane and can return the intersection point between a line and the ...
Definition: plane.H:60
static vector surfaceNormal(const triSurface &surf, const label nearestFacei, const point &nearestPt)
Triangle (unit) normal. If nearest point to triangle on edge use.
static const label NOEDGE
static triSurface delaunay2D(const List< vector2D > &)
Do unconstrained Delaunay of points. Returns triSurface with 3D.
static triSurface collapseEdges(const triSurface &surf, const labelList &collapsableEdges)
Create new triSurface by collapsing edges to edge mids.
static label otherFace(const triSurface &surf, const label facei, const label edgeI)
Get face connected to edge not facei.
triSurface triangulateFaceCentre(const polyBoundaryMesh &mBesh, const labelHashSet &includePatches, const bool verbose=false)
Face-centre triangulation of (selected patches of) boundaryMesh.
static label oppositeVertex(const triSurface &surf, const label facei, const label edgeI)
Get vertex (local numbering) opposite edge.
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
static triSurface mergePoints(const triSurface &surf, const scalar mergeTol)
Merge points within distance.
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects...
Definition: DynamicList.H:56
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:58
static sideType edgeSide(const triSurface &surf, const point &sample, const point &nearestPoint, const label edgeI)
If nearest point is on edgeI, determine on which side of surface.
static label minEdge(const triSurface &surf, const labelList &edgeIndices)
Returns element in edgeIndices with minimum length.
static void otherVertices(const triSurface &surf, const label facei, const label vertI, label &vert1I, label &vert2I)
Get the two vertices (local numbering) on facei counterclockwise.
PointFrompoint toPoint(const Foam::point &p)
static const label COLLAPSED
An STL-conforming hash table.
Definition: HashTable.H:61
static surfaceLocation trackToEdge(const triSurface &, const surfaceLocation &start, const surfaceLocation &end, const plane &cutPlane)
Track on surface to get closer to point.
Foam::polyBoundaryMesh.
A collection of tools for triSurface.
sideType
On which side of surface.
labelHashSet includePatches
static void calcInterpolationWeights(const triPointRef &, const point &, FixedList< scalar, 3 > &weights)
Calculate linear interpolation weights for point (guaranteed to be.
static label getTriangle(const triSurface &surf, const label e0I, const label e1I, const label e2I)
Return index of triangle (or -1) using all three edges.
static sideType surfaceSide(const triSurface &surf, const point &sample, const label nearestFacei)
Given nearest point (to sample) on surface determines which side.
label newPointi
Definition: readKivaGrid.H:501
static void otherEdges(const triSurface &surf, const label facei, const label edgeI, label &e1, label &e2)
Get the two edges on facei counterclockwise after edgeI.
const polyBoundaryMesh & bMesh
Hash function class for primitives. All non-primitives used to hash entries on hash tables likely nee...
Definition: Hash.H:52
const doubleScalar e
Elementary charge.
Definition: doubleScalar.H:105
static triSurface triangulate(const polyBoundaryMesh &mBesh, const labelHashSet &includePatches, const bool verbose=false)
Simple triangulation of (selected patches of) boundaryMesh. Needs.
Contains information about location on a triSurface.
static void writeOBJ(const fileName &fName, const pointField &pts)
Write pointField to OBJ format file.
Triangulated surface description with patch information.
Definition: triSurface.H:66
static const label ANYEDGE
Face collapse status.
static label getEdge(const triSurface &surf, const label vert1I, const label vert2I)
Returns edge label connecting v1, v2 (local numbering)
static void getVertexTriangles(const triSurface &surf, const label edgeI, labelList &edgeTris)
Get all triangles using edge endpoint.
static surfaceLocation classify(const triSurface &, const label triI, const point &trianglePoint)
Test point on plane of triangle to see if on edge or point or inside.
static label oppositeEdge(const triSurface &surf, const label facei, const label vertI)
Get edge opposite vertex (local numbering)
Namespace for OpenFOAM.