Loading...
Searching...
No Matches
faceAreaIntersect.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-2017 OpenFOAM Foundation
9 Copyright (C) 2017-2018 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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
27Class
28 Foam::faceAreaIntersect
29
30Description
31 Face intersection class
32 - calculates intersection area by sub-dividing face into triangles
33 and cutting
34
35SourceFiles
36 faceAreaIntersect.C
37
38\*---------------------------------------------------------------------------*/
39
40#ifndef faceAreaIntersect_H
41#define faceAreaIntersect_H
42
43#include "pointField.H"
44#include "FixedList.H"
45#include "plane.H"
46#include "face.H"
47#include "triangle.H"
48#include "Enum.H"
49#include "searchableSurface.H"
50
51// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
52
53namespace Foam
54{
56/*---------------------------------------------------------------------------*\
57 Class faceAreaIntersect Declaration
58\*---------------------------------------------------------------------------*/
62public:
63
65 {
66 tmFan,
67 tmMesh
68 };
69
71
72
73private:
74
75 // Private data
76
77 //- Reference to the points of face A
78 const pointField& pointsA_;
79
80 //- Reference to the points of face B
81 const pointField& pointsB_;
82
83 //- Triangle decomposition of face A
84 const DynamicList<face>& trisA_;
85
86 //- Triangle decomposition of face B
87 const DynamicList<face>& trisB_;
88
89 //- Flag to reverse B faces
90 const bool reverseB_;
91
92 //- Flag to cache the final triangulation
93 bool cacheTriangulation_;
94
95 //- Final triangulation
96 mutable DynamicList<triPoints> triangles_;
97
98
99 // Static data members
100
101 //- Tolerance
102 static scalar tol;
103
104
105 // Private Member Functions
106
107 //- Get triPoints from face
108 inline triPoints getTriPoints
109 (
110 const pointField& points,
111 const face& f,
112 const bool reverse
113 ) const;
114
115 //- Set triPoints into tri list
116 inline void setTriPoints
117 (
118 const point& a,
119 const point& b,
120 const point& c,
121 label& count,
123 ) const;
124
125 //- Return point of intersection between plane and triangle edge
126 inline point planeIntersection
127 (
128 const FixedList<scalar, 3>& d,
129 const triPoints& t,
130 const label negI,
131 const label posI
132 ) const;
133
134 //- Slice triangle with plane and generate new cut sub-triangles
135 void triSliceWithPlane
136 (
137 const triPoints& tri,
138 const plane& pln,
140 label& nTris,
141 const scalar len
142 ) const;
143
144 //- Return area of intersection of triangles src and tgt
145 void triangleIntersect
146 (
147 const triPoints& src,
148 const point& tgt0,
149 const point& tgt1,
150 const point& tgt2,
151 const vector& n,
152 scalar& area,
153 vector& centroid
154 ) const;
155
156
157public:
158
159 // Constructors
160
161 //- Construct from components
163 (
164 const pointField& pointsA,
165 const pointField& pointsB,
166 const DynamicList<face>& trisA,
167 const DynamicList<face>& trisB,
168 const bool reverseB = false,
169 const bool cacheTriangulation = false
170 );
171
172
173 // Public Member Functions
174
175 //- Fraction of local length scale to use as intersection tolerance
176 inline static scalar& tolerance();
177
178 //- Triangulate a face using the given triangulation mode
179 static void triangulate
180 (
181 const face& f,
182 const pointField& points,
183 const triangulationMode& triMode,
184 faceList& faceTris
185 );
186
187 //- Const access to the cacheTriangulation flag
188 inline bool cacheTriangulation() const;
189
190 //- Const access to the triangulation
192 {
193 return triangles_;
194 }
195
196 //- Non-const access to the triangulation
197 DynamicList<triPoints>& triangles() noexcept
198 {
199 return triangles_;
200 }
201
202 //- Decompose face into triangle fan
203 static inline void triangleFan
204 (
205 const face& f,
206 DynamicList<face>& faces
207 );
208
209 //- Return area of intersection of faceA with faceB and effective
210 //- face centre
211 void calc
212 (
213 const face& faceA,
214 const face& faceB,
215 const vector& n,
216 scalar& area,
217 vector& centroid
218 ) const;
219
220 //- Return area of intersection of faceA with faceB
221 bool overlaps
223 const face& faceA,
224 const face& faceB,
225 const vector& n,
226 const scalar threshold
227 ) const;
228};
229
231// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
232
233} // End namespace Foam
234
235// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
236
237#include "faceAreaIntersectI.H"
238
239// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
240
241#endif
242
243// ************************************************************************* //
label n
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition DynamicList.H:68
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition Enum.H:57
A 1D vector of objects of type <T> with a fixed length <N>.
Definition FixedList.H:73
void calc(const face &faceA, const face &faceB, const vector &n, scalar &area, vector &centroid) const
Return area of intersection of faceA with faceB and effective face centre.
faceAreaIntersect(const pointField &pointsA, const pointField &pointsB, const DynamicList< face > &trisA, const DynamicList< face > &trisB, const bool reverseB=false, const bool cacheTriangulation=false)
Construct from components.
bool overlaps(const face &faceA, const face &faceB, const vector &n, const scalar threshold) const
Return area of intersection of faceA with faceB.
static void triangleFan(const face &f, DynamicList< face > &faces)
Decompose face into triangle fan.
static const Enum< triangulationMode > triangulationModeNames_
bool cacheTriangulation() const
Const access to the cacheTriangulation flag.
static void triangulate(const face &f, const pointField &points, const triangulationMode &triMode, faceList &faceTris)
Triangulate a face using the given triangulation mode.
static scalar & tolerance()
Fraction of local length scale to use as intersection tolerance.
DynamicList< triPoints > & triangles() noexcept
Non-const access to the triangulation.
const DynamicList< triPoints > & triangles() const noexcept
Const access to the triangulation.
A face is a list of labels corresponding to mesh vertices.
Definition face.H:71
Geometric class that creates a 3D plane and can return the intersection point between a line and the ...
Definition plane.H:91
Triangle point storage. Default constructable (triangle is not).
Definition triangle.H:77
const pointField & points
Namespace for OpenFOAM.
void reverse(UList< T > &list, const label n)
Reverse the first n elements of the list.
Definition UListI.H:539
List< face > faceList
List of faces.
Definition faceListFwd.H:41
vector point
Point is a vector.
Definition point.H:37
const direction noexcept
Definition scalarImpl.H:265
vectorField pointField
pointField is a vectorField.
Vector< scalar > vector
Definition vector.H:57
labelList f(nPoints)
volScalarField & b