triSurfaceTools.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2018 OpenFOAM Foundation
9  Copyright (C) 2019 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 Class
28  Foam::triSurfaceTools
29 
30 Description
31  A collection of tools for triSurface.
32 
33 Note
34  The curvature calculation is an implementation of the algorithm from:
35  \verbatim
36  "Estimating Curvatures and their Derivatives on Triangle Meshes"
37  by S. Rusinkiewicz
38  3DPVT'04 Proceedings of the 3D Data Processing,
39  Visualization, and Transmission, 2nd International Symposium
40  Pages 486-493
41  http://gfx.cs.princeton.edu/pubs/_2004_ECA/curvpaper.pdf
42  \endverbatim
43 
44 SourceFiles
45  triSurfaceTools.C
46  triSurfaceCloseness.C
47  triSurfaceCurvature.C
48 
49 \*---------------------------------------------------------------------------*/
50 
51 #ifndef triSurfaceTools_H
52 #define triSurfaceTools_H
53 
54 #include "boolList.H"
55 #include "pointField.H"
56 #include "vectorField.H"
57 #include "triadFieldFwd.H"
58 #include "DynamicList.H"
59 #include "HashSet.H"
60 #include "Map.H"
61 #include "FixedList.H"
62 #include "Pair.H"
63 #include "vector2D.H"
64 #include "triPointRef.H"
65 #include "surfaceLocation.H"
66 
67 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
68 
69 namespace Foam
70 {
71 
72 // Forward declaration of classes
73 class boundBox;
74 class edge;
75 class labelledTri;
76 class polyBoundaryMesh;
77 class plane;
78 class triSurface;
79 class face;
80 class Time;
81 template<class Face> class MeshedSurface;
82 
83 /*---------------------------------------------------------------------------*\
84  Class triSurfaceTools Declaration
85 \*---------------------------------------------------------------------------*/
86 
87 class triSurfaceTools
88 {
89  // Private Member Functions
90 
91  // Refinement
92 
93  enum refineType
94  {
95  NONE,
96  RED,
97  GREEN
98  };
99  static void calcRefineStatus
100  (
101  const triSurface& surf,
102  const label facei,
103  List<refineType>& refine
104  );
105  static void greenRefine
106  (
107  const triSurface& surf,
108  const label facei,
109  const label edgeI,
110  const label newPointi,
111  DynamicList<labelledTri>& newFaces
112  );
113  static triSurface doRefine
114  (
115  const triSurface& surf,
116  const List<refineType>& refineStatus
117  );
118 
119 
120  // Coarsening
121 
122  static scalar faceCosAngle
123  (
124  const point& pStart,
125  const point& pEnd,
126  const point& pLeft,
127  const point& pRight
128  );
129 
130  static void protectNeighbours
131  (
132  const triSurface& surf,
133  const label vertI,
134  labelList& faceStatus
135  );
136 
137  //- Faces to collapse because of edge collapse
138  static labelHashSet getCollapsedFaces
139  (
140  const triSurface& surf,
141  label edgeI
142  );
143 
144  // Return value of faceUsed for faces using vertI (local numbering).
145  // Used internally.
146  static label vertexUsesFace
147  (
148  const triSurface& surf,
149  const labelHashSet& faceUsed,
150  const label vertI
151  );
152 
153  // Get new connections between faces (because of edge collapse)
154  // in form of tables:
155  // - given edge get other edge
156  // - given edge get other face
157  // A face using point v1 on edge will get connected to a face using
158  // point v2 if they share a common vertex
159  // (but not a common edge since then the triangles collapse to
160  // nothing)
161  static void getMergedEdges
162  (
163  const triSurface& surf,
164  const label edgeI,
165  const labelHashSet& collapsedFaces,
166  Map<label>& edgeToEdge,
167  Map<label>& edgeToFace
168  );
169 
170  //- Calculates (cos of) angle across edgeI of facei,
171  // taking into account updated addressing (resulting from edge
172  // collapse)
173  static scalar edgeCosAngle
174  (
175  const triSurface& surf,
176  const label v1,
177  const point& pt,
178  const labelHashSet& collapsedFaces,
179  const Map<label>& edgeToEdge,
180  const Map<label>& edgeToFace,
181  const label facei,
182  const label edgeI
183  );
184 
185  //- Calculate minimum (cos of) edge angle using addressing from
186  // collapsing
187  // edge to v1 at pt. Returns 1 if v1 is on edge without neighbours
188  // (and hence no edge angle can be defined)
189  static scalar collapseMinCosAngle
190  (
191  const triSurface& surf,
192  const label v1,
193  const point& pt,
194  const labelHashSet& collapsedFaces,
195  const Map<label>& edgeToEdge,
196  const Map<label>& edgeToFace
197  );
198 
199  //- Like collapseMinCosAngle but return true for value < minCos
200  bool collapseCreatesFold
201  (
202  const triSurface& surf,
203  const label v1,
204  const point& pt,
205  const labelHashSet& collapsedFaces,
206  const Map<label>& edgeToEdge,
207  const Map<label>& edgeToFace,
208  const scalar minCos
209  );
210 
213  //static bool collapseCreatesDuplicates
214  //(
215  // const triSurface& surf,
216  // const label edgeI,
217  // const Map<bool>& collapsedFaces
218  //);
219 
220  // Tracking
221 
222  //- Finds the triangle edge/point cut by the plane between
223  // a point inside/on edge of a triangle and a point outside.
224  // Returns
225  // - location on edge/point and hit()
226  // - or miss() if no intersection found
227  static surfaceLocation cutEdge
228  (
229  const triSurface& s,
230  const label triI,
231  const label excludeEdgeI,
232  const label excludePointi,
233  const point& triPoint,
234  const plane& cutPlane,
235  const point& toPoint
236  );
237 
238  //- Checks if current is on the same triangle as the endpoint
239  // and shifts it there. If so updates current and sets a hit.
240  static void snapToEnd
241  (
242  const triSurface& s,
243  const surfaceLocation& endInfo,
244  surfaceLocation& current
245  );
246 
247  //- Visits faces eFaces around start. Does not visit triangle
248  // start.triangle() nor edge excludeEdgeI.
249  // Returns edge, triangle (if more than one choice) which gets
250  // us nearer endpoint.
251  // Returns
252  // - hit() if triangle contains endpoint
253  // - triangle()=-1 if no triangle found
254  // - nearest triangle/edge otherwise
255  static surfaceLocation visitFaces
256  (
257  const triSurface& s,
258  const labelList& eFaces,
259  const surfaceLocation& start,
260  const label excludeEdgeI,
261  const label excludePointi,
262  const surfaceLocation& end,
263  const plane& cutPlane
264  );
265 
266 
267 public:
268 
269  // OBJ writing
270 
271  //- Write pointField to OBJ format file
272  static void writeOBJ
273  (
274  const fileName& fName,
275  const pointField& pts
276  );
277 
278  //- Write vertex subset to OBJ format file
279  static void writeOBJ
280  (
281  const triSurface& surf,
282  const fileName& fName,
283  const boolList& markedVerts
284  );
285 
286 
287  // Additional addressing
288 
289  //- Get all triangles using edge endpoint
290  static void getVertexTriangles
291  (
292  const triSurface& surf,
293  const label edgeI,
294  labelList& edgeTris
295  );
296 
297  //- Get all vertices (local numbering) connected to vertices of edge
299  (
300  const triSurface& surf,
301  const edge& e
302  );
303 
304  //- Get face connected to edge not facei
305  static label otherFace
306  (
307  const triSurface& surf,
308  const label facei,
309  const label edgeI
310  );
311 
312  //- Get the two edges on facei counterclockwise after edgeI
313  static void otherEdges
314  (
315  const triSurface& surf,
316  const label facei,
317  const label edgeI,
318  label& e1,
319  label& e2
320  );
321 
322  //- Get the two vertices (local numbering) on facei counterclockwise
323  // vertI
324  static void otherVertices
325  (
326  const triSurface& surf,
327  const label facei,
328  const label vertI,
329  label& vert1I,
330  label& vert2I
331  );
332 
333  //- Get edge opposite vertex (local numbering)
334  static label oppositeEdge
335  (
336  const triSurface& surf,
337  const label facei,
338  const label vertI
339  );
340 
341  //- Get vertex (local numbering) opposite edge
342  static label oppositeVertex
343  (
344  const triSurface& surf,
345  const label facei,
346  const label edgeI
347  );
348 
349  //- Returns edge label connecting v1, v2 (local numbering)
350  static label getEdge
351  (
352  const triSurface& surf,
353  const label vert1I,
354  const label vert2I
355  );
356 
357  //- Return index of triangle (or -1) using all three edges
358  static label getTriangle
359  (
360  const triSurface& surf,
361  const label e0I,
362  const label e1I,
363  const label e2I
364  );
365 
366  // Coarsening
367 
368  //- Create new triSurface by collapsing edges to edge mids.
370  (
371  const triSurface& surf,
372  const labelList& collapsableEdges
373  );
374 
375 
376  //- Face collapse status.
377  // anyEdge: any edge can be collapsed
378  // noEdge: no edge can be collapsed
379  // collapsed: already collapsed
380  // >0: edge label that can be collapsed
381  static const label ANYEDGE;
382  static const label NOEDGE;
383  static const label COLLAPSED;
384 
385  //- Create new triSurface by collapsing edges to specified
386  // positions. faceStatus allows
387  // explicit control over which faces need to be protected (see above).
388  // faceStatus gets updated to protect collapsing already collapsed
389  // faces.
391  (
392  const triSurface& surf,
393  const labelList& collapsableEdges,
394  const pointField& edgeMids,
395  labelList& faceStatus
396  );
397 
398 
399  // Refinement
400 
401  //- Refine edges by splitting to opposite vertex
402  static triSurface greenRefine
403  (
404  const triSurface& surf,
405  const labelList& refineEdges
406  );
407 
408  //- Refine face by splitting all edges. Neighbouring face is
409  // greenRefine'd.
411  (
412  const triSurface& surf,
413  const labelList& refineFaces
414  );
415 
416 
417  // Geometric
418 
419  //- Returns element in edgeIndices with minimum length
420  static label minEdge
421  (
422  const triSurface& surf,
423  const labelList& edgeIndices
424  );
425 
426  //- Returns element in edgeIndices with minimum length
427  static label maxEdge
428  (
429  const triSurface& surf,
430  const labelList& edgeIndices
431  );
432 
433  //- Merge points within distance
434  static triSurface mergePoints
435  (
436  const triSurface& surf,
437  const scalar mergeTol
438  );
439 
440  //- Triangle (unit) normal. If nearest point to triangle on edge use
441  // edge normal (calculated on the fly); if on vertex use vertex normal.
442  // Uses planarTol.
443  static vector surfaceNormal
444  (
445  const triSurface& surf,
446  const label nearestFacei,
447  const point& nearestPt
448  );
449 
450  //- On which side of surface
451  enum sideType
452  {
453  UNKNOWN, // cannot be determined (e.g. non-manifold)
454  INSIDE, // inside
455  OUTSIDE // outside
456  };
457 
458  //- If nearest point is on edgeI, determine on which side of surface
459  // sample is.
460  static sideType edgeSide
461  (
462  const triSurface& surf,
463  const point& sample,
464  const point& nearestPoint,
465  const label edgeI
466  );
467 
468  //- Given nearest point (to sample) on surface determines which side
469  // sample is. Uses either face normal, edge normal or point normal
470  // (non-trivial). Uses triangle::classify.
471  static sideType surfaceSide
472  (
473  const triSurface& surf,
474  const point& sample,
475  const label nearestFacei
476  );
477 
478  // Triangulation of faces
479 
480  //- Simple triangulation of (selected patches of) boundaryMesh. Needs
481  // polyMesh (or polyBoundaryMesh) since only at this level are the
482  // triangles on neighbouring patches connected.
483  // Return faceMap from triI to faceI
484  static triSurface triangulate
485  (
486  const polyBoundaryMesh& mBesh,
487  const labelHashSet& includePatches,
489  const bool verbose = false
490  );
491 
492 
493  static triSurface triangulate
494  (
495  const polyBoundaryMesh& bMesh,
496  const labelHashSet& includePatches,
497  const boundBox& bBox,
498  const bool verbose = false
499  );
500 
501 
502  //- Face-centre triangulation of (selected patches of) boundaryMesh.
503  // Needs
504  // polyMesh (or polyBoundaryMesh) since only at this level are the
505  // triangles on neighbouring patches connected.
507  (
508  const polyBoundaryMesh& mBesh,
509  const labelHashSet& includePatches,
510  const bool verbose = false
511  );
512 
513 
514  // Triangulation and interpolation
515 
516  //- Do unconstrained Delaunay of points. Returns triSurface with 3D
517  // points with z=0. All triangles in region 0.
518  static triSurface delaunay2D(const List<vector2D>&);
519 
520  //- Calculate linear interpolation weights for point (guaranteed to be
521  // inside triangle)
522  static void calcInterpolationWeights
523  (
524  const triPointRef& tri,
525  const point& p,
526  FixedList<scalar, 3>& weights
527  );
528 
529  // Calculate weighting factors from samplePts to triangle it is in.
530  // Uses linear search to find triangle.
531  // Vertices are:
532  // (a b c) : vertices of the triangle abc the point is in
533  // or if the point is outside all triangles:
534  // (a b -1) : the edge ab the point is nearest to.
535  // (a -1 -1) : the vertex a the point is nearest to
536  static void calcInterpolationWeights
537  (
538  const triSurface& s,
539  const pointField& samplePts,
540  List<FixedList<label, 3>>& allVerts,
541  List<FixedList<scalar, 3>>& allWeights
542  );
543 
544 
545  // Curvature
546 
547  //- Weighting for normals of faces attached to vertex
548  static scalar vertexNormalWeight
549  (
550  const triFace& f,
551  const label pI,
552  const vector& fN,
553  const UList<point>& points
554  );
555 
556  //- Weighted average of normals of attached faces
557  static tmp<vectorField> vertexNormals(const triSurface& surf);
558 
559  //- Local coordinate-system for each point normal
561  (
562  const triSurface& surf,
563  const vectorField& pointNormals
564  );
565 
566  //- Surface curvatures at the vertex points
568  (
569  const triSurface& surf,
570  const vectorField& pointNormals,
571  const triadField& pointTriads
572  );
573 
574  //- Surface curvatures at the vertex points
576  (
577  const triSurface& surf
578  );
579 
580  //- Write surface curvature at the vertex points and return the field
582  (
583  const Time& runTime,
584  const word& basename,
585  const triSurface& surf
586  );
587 
588 
589  // Closeness
590 
591  //- Check and write internal/external closeness fields
593  (
594  const Time& runTime,
595  const word& basename,
596  const triSurface& surf,
597  const scalar internalAngleTolerance = 45,
598  const scalar externalAngleTolerance = 10
599  );
600 
601 
602  // Feature Proximity
603 
604  //- Calculate feature proximity
606  (
607  const triSurface& surf,
608  const scalar searchDistance
609  );
610 
611  //- Check and write internal/external closeness fields
612  static void writeFeatureProximity
613  (
614  const Time& runTime,
615  const word& basename,
616  const triSurface& surf,
617  const bool writeVTK,
618  const scalar searchDistance
619  );
620 
621 
622  // Surface checking functionality
623 
624  //- Check single triangle for (topological) validity
625  static bool validTri
626  (
627  const triSurface&,
628  const label facei,
629  const bool verbose = true
630  );
631 
632  //- Check single triangle for (topological) validity
633  static bool validTri
634  (
635  const MeshedSurface<face>&,
636  const label facei,
637  const bool verbose = true
638  );
639 
640 
641  // Tracking
642 
643  //- Test point on plane of triangle to see if on edge or point or inside
645  (
646  const triSurface&,
647  const label triI,
648  const point& trianglePoint
649  );
650 
651  //- Track on surface to get closer to point.
652  // Possible situations:
653  // - 1. reached endpoint
654  // - 2. reached edge (normal situation)
655  // - 3. reached end of surface (edge on single face)
656  // Input:
657  // - starting position+triangle/edge/point (so has to be on surface!)
658  // - (optional) previous position+triangle/edge/point to prevent
659  // going back. Set index (of triangle/edge/point) to -1 if not
660  // used.
661  // - end position+triangle/edge/point (so has to be on surface!)
662  // - plane to follow. Has to go through end point!
663  // Returns:
664  // - true if end point reached (situation 1)
665  // - new position+triangle/edge/point
666  // Caller has to check for situation 3 by checking that triangle()
667  // is not set.
669  (
670  const triSurface&,
671  const surfaceLocation& start,
672  const surfaceLocation& end,
673  const plane& cutPlane
674  );
675 
676  //- Track from edge to edge across surface. Uses trackToEdge.
677  // Not really useful by itself, more example of how to use trackToEdge.
678  // endInfo should be location on surface.
679  // hitInfo should be initialised to starting location (on surface as
680  // well). Upon return is set to end location.
681  static void track
682  (
683  const triSurface&,
684  const surfaceLocation& endInfo,
685  const plane& cutPlane,
686  surfaceLocation& hitInfo
687  );
688 };
689 
690 
691 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
692 
693 } // End namespace Foam
694 
695 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
696 
697 #endif
698 
699 // ************************************************************************* //
runTime
engineTime & runTime
Definition: createEngineTime.H:13
Foam::triSurfaceTools::writeOBJ
static void writeOBJ(const fileName &fName, const pointField &pts)
Write pointField to OBJ format file.
Definition: triSurfaceTools.C:1202
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::Time
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:73
Foam::faceMap
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
Definition: blockMeshMergeTopological.C:94
Foam::triSurfaceTools::getEdge
static label getEdge(const triSurface &surf, const label vert1I, const label vert2I)
Returns edge label connecting v1, v2 (local numbering)
Definition: triSurfaceTools.C:1467
boolList.H
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::fileName
A class for handling file names.
Definition: fileName.H:73
s
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.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))
Definition: gmvOutputSpray.H:25
Foam::polyBoundaryMesh
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
Definition: polyBoundaryMesh.H:63
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::triSurfaceTools::surfaceNormal
static vector surfaceNormal(const triSurface &surf, const label nearestFacei, const point &nearestPt)
Triangle (unit) normal. If nearest point to triangle on edge use.
Definition: triSurfaceTools.C:1955
Foam::DynamicList
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:55
Foam::triSurfaceTools::triangulateFaceCentre
triSurface triangulateFaceCentre(const polyBoundaryMesh &mBesh, const labelHashSet &includePatches, const bool verbose=false)
Face-centre triangulation of (selected patches of) boundaryMesh.
Definition: triSurfaceTools.C:2376
Foam::triSurfaceTools::otherFace
static label otherFace(const triSurface &surf, const label facei, const label edgeI)
Get face connected to edge not facei.
Definition: triSurfaceTools.C:1321
Foam::triSurfaceTools::maxEdge
static label maxEdge(const triSurface &surf, const labelList &edgeIndices)
Returns element in edgeIndices with minimum length.
Definition: triSurfaceTools.C:1871
Foam::triSurfaceTools::delaunay2D
static triSurface delaunay2D(const List< vector2D > &)
Do unconstrained Delaunay of points. Returns triSurface with 3D.
Definition: triSurfaceTools.C:2472
Foam::triSurfaceTools::ANYEDGE
static const label ANYEDGE
Face collapse status.
Definition: triSurfaceTools.H:380
Foam::triSurfaceTools::collapseEdges
static triSurface collapseEdges(const triSurface &surf, const labelList &collapsableEdges)
Create new triSurface by collapsing edges to edge mids.
Definition: triSurfaceTools.C:1535
Foam::triSurfaceTools::mergePoints
static triSurface mergePoints(const triSurface &surf, const scalar mergeTol)
Merge points within distance.
Definition: triSurfaceTools.C:1898
Foam::triSurfaceTools::vertexTriads
static tmp< triadField > vertexTriads(const triSurface &surf, const vectorField &pointNormals)
Local coordinate-system for each point normal.
Definition: triSurfaceCurvature.C:110
Foam::surfaceLocation
Contains information about location on a triSurface.
Definition: surfaceLocation.H:71
Foam::edge
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:63
Foam::Map< label >
Foam::triSurfaceTools::writeCurvature
static tmp< scalarField > writeCurvature(const Time &runTime, const word &basename, const triSurface &surf)
Write surface curvature at the vertex points and return the field.
Definition: triSurfaceCurvature.C:349
Pair.H
Foam::triSurfaceTools::oppositeEdge
static label oppositeEdge(const triSurface &surf, const label facei, const label vertI)
Get edge opposite vertex (local numbering)
Definition: triSurfaceTools.C:1411
Foam::triSurfaceTools::triangulate
static triSurface triangulate(const polyBoundaryMesh &mBesh, const labelHashSet &includePatches, labelList &faceMap, const bool verbose=false)
Simple triangulation of (selected patches of) boundaryMesh. Needs.
Definition: triSurfaceTools.C:2192
triPointRef.H
Foam::HashSet< label, Hash< label > >
Foam::writeVTK
void writeVTK(OFstream &os, const Type &value)
Foam::plane
Geometric class that creates a 3D plane and can return the intersection point between a line and the ...
Definition: plane.H:89
Foam::triSurfaceTools::COLLAPSED
static const label COLLAPSED
Definition: triSurfaceTools.H:382
Foam::triSurfaceTools::otherVertices
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.
Definition: triSurfaceTools.C:1375
Foam::triangle
A triangle primitive used to calculate face normals and swept volumes.
Definition: triangle.H:59
Map.H
Foam::triSurfaceTools::getTriangle
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.
Definition: triSurfaceTools.C:1490
Foam::Field< vector >
Foam::triSurface
Triangulated surface description with patch information.
Definition: triSurface.H:76
Foam::triSurfaceTools::otherEdges
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.
Definition: triSurfaceTools.C:1346
Foam::triSurfaceTools::redGreenRefine
static triSurface redGreenRefine(const triSurface &surf, const labelList &refineFaces)
Refine face by splitting all edges. Neighbouring face is.
Definition: triSurfaceTools.C:1742
vector2D.H
Foam::triSurfaceTools::vertexNormalWeight
static scalar vertexNormalWeight(const triFace &f, const label pI, const vector &fN, const UList< point > &points)
Weighting for normals of faces attached to vertex.
Definition: triSurfaceCurvature.C:43
HashSet.H
Foam::triSurfaceTools::calcInterpolationWeights
static void calcInterpolationWeights(const triPointRef &tri, const point &p, FixedList< scalar, 3 > &weights)
Calculate linear interpolation weights for point (guaranteed to be.
Definition: triSurfaceTools.C:2538
newPointi
label newPointi
Definition: readKivaGrid.H:496
Foam::toPoint
PointFrompoint toPoint(const Foam::point &p)
Definition: pointConversion.H:82
Foam::triSurfaceTools::vertexNormals
static tmp< vectorField > vertexNormals(const triSurface &surf)
Weighted average of normals of attached faces.
Definition: triSurfaceCurvature.C:66
stdFoam::end
constexpr auto end(C &c) -> decltype(c.end())
Return iterator to the end of the container c.
Definition: stdFoam.H:121
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::triSurfaceTools::NOEDGE
static const label NOEDGE
Definition: triSurfaceTools.H:381
Foam::triSurfaceTools::UNKNOWN
Definition: triSurfaceTools.H:452
Foam::triSurfaceTools::surfaceSide
static sideType surfaceSide(const triSurface &surf, const point &sample, const label nearestFacei)
Given nearest point (to sample) on surface determines which side.
Definition: triSurfaceTools.C:2036
Foam::triSurfaceTools::getVertexTriangles
static void getVertexTriangles(const triSurface &surf, const label edgeI, labelList &edgeTris)
Get all triangles using edge endpoint.
Definition: triSurfaceTools.C:1249
Foam::triSurfaceTools::oppositeVertex
static label oppositeVertex(const triSurface &surf, const label facei, const label edgeI)
Get vertex (local numbering) opposite edge.
Definition: triSurfaceTools.C:1440
Foam::triSurfaceTools::edgeSide
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.
Definition: triSurfaceTools.C:1999
pointField.H
Foam::triSurfaceTools::INSIDE
Definition: triSurfaceTools.H:453
Foam::triSurfaceTools::trackToEdge
static surfaceLocation trackToEdge(const triSurface &, const surfaceLocation &start, const surfaceLocation &end, const plane &cutPlane)
Track on surface to get closer to point.
Definition: triSurfaceTools.C:2913
Foam::triSurfaceTools::classify
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.
Definition: triSurfaceTools.C:2869
Foam::triSurfaceTools::writeFeatureProximity
static void writeFeatureProximity(const Time &runTime, const word &basename, const triSurface &surf, const bool writeVTK, const scalar searchDistance)
Check and write internal/external closeness fields.
Foam::triSurfaceTools::writeCloseness
static Pair< tmp< scalarField > > writeCloseness(const Time &runTime, const word &basename, const triSurface &surf, const scalar internalAngleTolerance=45, const scalar externalAngleTolerance=10)
Check and write internal/external closeness fields.
Definition: triSurfaceCloseness.C:96
Foam::Pair
An ordered pair of two objects of type <T> with first() and second() elements.
Definition: Pair.H:54
Foam::triFace
A triangular face using a FixedList of labels corresponding to mesh vertices.
Definition: triFace.H:69
Foam::triSurfaceTools::minEdge
static label minEdge(const triSurface &surf, const labelList &edgeIndices)
Returns element in edgeIndices with minimum length.
Definition: triSurfaceTools.C:1844
f
labelList f(nPoints)
Foam::triSurfaceTools::getVertexVertices
static labelList getVertexVertices(const triSurface &surf, const edge &e)
Get all vertices (local numbering) connected to vertices of edge.
Definition: triSurfaceTools.C:1290
Foam::Vector< scalar >
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: BitOps.H:63
Foam::triSurfaceTools
A collection of tools for triSurface.
Definition: triSurfaceTools.H:86
vectorField.H
Foam::FixedList< scalar, 3 >
Foam::triSurfaceTools::OUTSIDE
Definition: triSurfaceTools.H:454
triadFieldFwd.H
Forward declarations of Field<T> triad specialisation.
Foam::UList
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: HashTable.H:103
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::boundBox
A bounding box defined in terms of min/max extrema points.
Definition: boundBox.H:63
Foam::triSurfaceTools::featureProximity
scalarField featureProximity(const triSurface &surf, const scalar searchDistance)
Calculate feature proximity.
DynamicList.H
FixedList.H
Foam::triSurfaceTools::curvatures
static tmp< scalarField > curvatures(const triSurface &surf, const vectorField &pointNormals, const triadField &pointTriads)
Surface curvatures at the vertex points.
Definition: triSurfaceCurvature.C:147
Foam::triSurfaceTools::validTri
static bool validTri(const triSurface &, const label facei, const bool verbose=true)
Check single triangle for (topological) validity.
Definition: triSurfaceTools.C:2696
Foam::MeshedSurface
A surface geometry mesh with zone information, not to be confused with the similarly named surfaceMes...
Definition: triSurfaceTools.H:80
Foam::triSurfaceTools::sideType
sideType
On which side of surface.
Definition: triSurfaceTools.H:450
Foam::triSurfaceTools::track
static void track(const triSurface &, const surfaceLocation &endInfo, const plane &cutPlane, surfaceLocation &hitInfo)
Track from edge to edge across surface. Uses trackToEdge.
Definition: triSurfaceTools.C:2991
surfaceLocation.H
Foam::PrimitivePatch
A list of faces which address into the list of points.
Definition: PrimitivePatch.H:79
sample
Minimal example by using system/controlDict.functions: