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 -------------------------------------------------------------------------------
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::faceAreaIntersect
29 
30 Description
31  Face intersection class
32  - calculates intersection area by sub-dividing face into triangles
33  and cutting
34 
35 SourceFiles
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 "triPoints.H"
48 #include "Enum.H"
49 #include "searchableSurface.H"
50 
51 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
52 
53 namespace Foam
54 {
55 
56 /*---------------------------------------------------------------------------*\
57  Class faceAreaIntersect Declaration
58 \*---------------------------------------------------------------------------*/
59 
61 {
62 public:
63 
65  {
67  tmMesh
68  };
69 
71 
72 
73 private:
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  //- Return triangle area
135  inline scalar triArea(const triPoints& t) const;
136 
137  //- Return triangle centre
138  inline vector triCentroid(const triPoints& t) const;
139 
140  //- Slice triangle with plane and generate new cut sub-triangles
141  void triSliceWithPlane
142  (
143  const triPoints& tri,
144  const plane& pln,
146  label& nTris,
147  const scalar len
148  ) const;
149 
150  //- Return area of intersection of triangles src and tgt
151  void triangleIntersect
152  (
153  const triPoints& src,
154  const point& tgt0,
155  const point& tgt1,
156  const point& tgt2,
157  const vector& n,
158  scalar& area,
159  vector& centroid
160  ) const;
161 
162 
163 public:
164 
165  // Constructors
166 
167  //- Construct from components
169  (
170  const pointField& pointsA,
171  const pointField& pointsB,
172  const DynamicList<face>& trisA,
173  const DynamicList<face>& trisB,
174  const bool reverseB = false,
175  const bool cacheTriangulation = false
176  );
177 
178 
179  // Public Member Functions
180 
181  //- Fraction of local length scale to use as intersection tolerance
182  inline static scalar& tolerance();
183 
184  //- Triangulate a face using the given triangulation mode
185  static void triangulate
186  (
187  const face& f,
188  const pointField& points,
189  const triangulationMode& triMode,
190  faceList& faceTris
191  );
192 
193  //- Const access to the cacheTriangulation flag
194  inline bool cacheTriangulation() const;
195 
196  //- Const access to the triangulation
197  inline const DynamicList<triPoints> triangles() const;
198 
199  //- Non-const access to the triangulation
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
222  (
223  const face& faceA,
224  const face& faceB,
225  const vector& n,
226  const scalar threshold
227  ) const;
228 };
229 
230 
231 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
232 
233 } // End namespace Foam
234 
235 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
236 
237 #include "faceAreaIntersectI.H"
238 
239 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
240 
241 #endif
242 
243 // ************************************************************************* //
Foam::reverse
void reverse(UList< T > &list, const label n)
Definition: UListI.H:449
Foam::Enum< triangulationMode >
searchableSurface.H
Foam::faceAreaIntersect::triangulationMode
triangulationMode
Definition: faceAreaIntersect.H:63
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::faceAreaIntersect::triangulate
static void triangulate(const face &f, const pointField &points, const triangulationMode &triMode, faceList &faceTris)
Triangulate a face using the given triangulation mode.
Definition: faceAreaIntersect.C:408
face.H
triPoints.H
Foam::faceAreaIntersect::calc
void calc(const face &faceA, const face &faceB, const vector &n, scalar &area, vector &centroid) const
Definition: faceAreaIntersect.C:451
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::faceAreaIntersect::triangulationModeNames_
static const Enum< triangulationMode > triangulationModeNames_
Definition: faceAreaIntersect.H:69
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::triPoints
Triangle storage. Null constructable (unfortunately triangle<point, point> is not)
Definition: triPoints.H:52
Foam::faceAreaIntersect::faceAreaIntersect
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.
Definition: faceAreaIntersect.C:386
Foam::constant::physicoChemical::b
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
Foam::Field< vector >
faceAreaIntersectI.H
plane.H
Foam::faceAreaIntersect
Face intersection class.
Definition: faceAreaIntersect.H:59
Foam::faceAreaIntersect::tmFan
Definition: faceAreaIntersect.H:65
Foam::faceAreaIntersect::overlaps
bool overlaps(const face &faceA, const face &faceB, const vector &n, const scalar threshold) const
Return area of intersection of faceA with faceB.
Definition: faceAreaIntersect.C:512
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
pointField.H
Foam::faceAreaIntersect::cacheTriangulation
bool cacheTriangulation() const
Const access to the cacheTriangulation flag.
Definition: faceAreaIntersectI.H:133
f
labelList f(nPoints)
Foam::BitOps::count
unsigned int count(const UList< bool > &bools, const bool val=true)
Count number of 'true' entries.
Definition: BitOps.H:77
Foam::Vector< scalar >
Foam::List< face >
Foam::faceAreaIntersect::tolerance
static scalar & tolerance()
Fraction of local length scale to use as intersection tolerance.
Definition: faceAreaIntersectI.H:127
Foam::FixedList
A 1D vector of objects of type <T> with a fixed length <N>.
Definition: HashTable.H:104
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::fieldTypes::area
const wordList area
Standard area field types (scalar, vector, tensor, etc)
Foam::constant::universal::c
const dimensionedScalar c
Speed of light in a vacuum.
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:72
FixedList.H
Foam::faceAreaIntersect::triangles
const DynamicList< triPoints > triangles() const
Const access to the triangulation.
Definition: faceAreaIntersectI.H:140
Foam::faceAreaIntersect::triangleFan
static void triangleFan(const face &f, DynamicList< face > &faces)
Decompose face into triangle fan.
Definition: faceAreaIntersectI.H:32
Foam::faceAreaIntersect::tmMesh
Definition: faceAreaIntersect.H:66
Enum.H