faceAreaIntersect.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-2021 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::faceAreaIntersect
26 
27 Description
28  Face intersection class
29  - calculates intersection area by sub-dividing face into triangles
30  and cutting
31 
32 SourceFiles
33  faceAreaIntersect.C
34 
35 \*---------------------------------------------------------------------------*/
36 
37 #ifndef faceAreaIntersect_H
38 #define faceAreaIntersect_H
39 
40 #include "pointField.H"
41 #include "FixedList.H"
42 #include "plane.H"
43 #include "face.H"
44 #include "NamedEnum.H"
45 #include "triFaceList.H"
46 
47 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
48 
49 namespace Foam
50 {
51 
52 /*---------------------------------------------------------------------------*\
53  Class faceAreaIntersect Declaration
54 \*---------------------------------------------------------------------------*/
55 
57 {
58 public:
59 
61 
63  {
65  tmMesh
66  };
67 
69 
70 
71 private:
72 
73  // Private Data
74 
75  //- Reference to the points of sideA
76  const pointField& pointsA_;
77 
78  //- Reference to the points of sideB
79  const pointField& pointsB_;
80 
81  //- Flag to reverse B faces
82  const bool reverseB_;
83 
84 
85  // Static Data Members
86 
87  static scalar tol;
88 
89 
90  // Private Member Functions
91 
92  //- Get triPoints from face
93  inline triPoints getTriPoints
94  (
95  const pointField& points,
96  const triFace& f,
97  const bool reverse
98  ) const;
99 
100  //- Set triPoints into tri list
101  inline void setTriPoints
102  (
103  const point& a,
104  const point& b,
105  const point& c,
106  label& count,
108  ) const;
109 
110  //- Return point of intersection between plane and triangle edge
111  inline point planeIntersection
112  (
113  const FixedList<scalar, 3>& d,
114  const triPoints& t,
115  const label negI,
116  const label posI
117  ) const;
118 
119  //- Return triangle area
120  inline scalar triArea(const triPoints& t) const;
121 
122 
123  //- Slice triangle with plane and generate new cut sub-triangles
124  void triSliceWithPlane
125  (
126  const triPoints& tri,
127  const plane& p,
129  label& nTris,
130  const scalar len
131  );
132 
133  //- Return area of intersection of triangles src and tgt
134  scalar triangleIntersect
135  (
136  const triPoints& src,
137  const triPoints& tgt,
138  const vector& n
139  );
140 
141 
142 public:
143 
144  // Constructors
145 
146  //- Construct from components
148  (
149  const pointField& pointsA,
150  const pointField& pointsB,
151  const bool reverseB = false
152  );
153 
154 
155  // Public Member Functions
156 
157  //- Fraction of local length scale to use as intersection tolerance
158  inline static scalar& tolerance();
159 
160  //- Triangulate a face using the given triangulation mode
161  static void triangulate
162  (
163  const face& f,
164  const pointField& points,
165  const triangulationMode& triMode,
166  triFaceList& faceTris
167  );
168 
169  //- Return area of intersection of faceA with faceB
170  scalar calc
171  (
172  const face& faceA,
173  const face& faceB,
174  const vector& n,
175  const triangulationMode& triMode
176  );
177 };
178 
179 
180 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
181 
182 } // End namespace Foam
183 
184 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
185 
186 #include "faceAreaIntersectI.H"
187 
188 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
189 
190 #endif
191 
192 // ************************************************************************* //
static scalar & tolerance()
Fraction of local length scale to use as intersection tolerance.
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
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
A 1D vector of objects of type <T> with a fixed size <Size>.
Definition: FixedList.H:54
const dimensionedScalar b
Wien displacement law constant: default SI units: [m K].
Definition: createFields.H:27
label count(const ListType &l, typename ListType::const_reference x)
Count the number of occurrences of a value in a list.
Geometric class that creates a 2D plane and can return the intersection point between a line and the ...
Definition: plane.H:60
static const NamedEnum< triangulationMode, 2 > triangulationModeNames_
static void triangulate(const face &f, const pointField &points, const triangulationMode &triMode, triFaceList &faceTris)
Triangulate a face using the given triangulation mode.
const dimensionedScalar c
Speed of light in a vacuum.
faceAreaIntersect(const pointField &pointsA, const pointField &pointsB, const bool reverseB=false)
Construct from components.
const pointField & points
A triangular face using a FixedList of labels corresponding to mesh vertices.
Definition: triFace.H:68
scalar calc(const face &faceA, const face &faceB, const vector &n, const triangulationMode &triMode)
Return area of intersection of faceA with faceB.
void reverse(UList< T > &, const label n)
Definition: UListI.H:334
Face intersection class.
labelList f(nPoints)
FixedList< point, 3 > triPoints
label n
volScalarField & p
Namespace for OpenFOAM.