All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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-2019 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 
46 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
47 
48 namespace Foam
49 {
50 
51 /*---------------------------------------------------------------------------*\
52  Class faceAreaIntersect Declaration
53 \*---------------------------------------------------------------------------*/
54 
56 {
57 public:
58 
60 
62  {
64  tmMesh
65  };
66 
68 
69 
70 private:
71 
72  // Private Data
73 
74  //- Reference to the points of sideA
75  const pointField& pointsA_;
76 
77  //- Reference to the points of sideB
78  const pointField& pointsB_;
79 
80  //- Flag to reverse B faces
81  const bool reverseB_;
82 
83 
84  // Static Data Members
85 
86  static scalar tol;
87 
88 
89  // Private Member Functions
90 
91  //- Get triPoints from face
92  inline triPoints getTriPoints
93  (
94  const pointField& points,
95  const face& f,
96  const bool reverse
97  ) const;
98 
99  //- Set triPoints into tri list
100  inline void setTriPoints
101  (
102  const point& a,
103  const point& b,
104  const point& c,
105  label& count,
107  ) const;
108 
109  //- Return point of intersection between plane and triangle edge
110  inline point planeIntersection
111  (
112  const FixedList<scalar, 3>& d,
113  const triPoints& t,
114  const label negI,
115  const label posI
116  ) const;
117 
118  //- Return triangle area
119  inline scalar triArea(const triPoints& t) const;
120 
121 
122  //- Slice triangle with plane and generate new cut sub-triangles
123  void triSliceWithPlane
124  (
125  const triPoints& tri,
126  const plane& p,
128  label& nTris,
129  const scalar len
130  );
131 
132  //- Return area of intersection of triangles src and tgt
133  scalar triangleIntersect
134  (
135  const triPoints& src,
136  const triPoints& tgt,
137  const vector& n
138  );
139 
140 
141 public:
142 
143  // Constructors
144 
145  //- Construct from components
147  (
148  const pointField& pointsA,
149  const pointField& pointsB,
150  const bool reverseB = false
151  );
152 
153 
154  // Public Member Functions
155 
156  //- Fraction of local length scale to use as intersection tolerance
157  inline static scalar& tolerance();
158 
159  //- Triangulate a face using the given triangulation mode
160  static void triangulate
161  (
162  const face& f,
163  const pointField& points,
164  const triangulationMode& triMode,
165  faceList& faceTris
166  );
167 
168  //- Return area of intersection of faceA with faceB
169  scalar calc
170  (
171  const face& faceA,
172  const face& faceB,
173  const vector& n,
174  const triangulationMode& triMode
175  );
176 };
177 
178 
179 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
180 
181 } // End namespace Foam
182 
183 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
184 
185 #include "faceAreaIntersectI.H"
186 
187 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
188 
189 #endif
190 
191 // ************************************************************************* //
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
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_
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
const dimensionedScalar & b
Wien displacement law constant: default SI units: [m K].
Definition: createFields.H:27
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:322
Face intersection class.
labelList f(nPoints)
FixedList< point, 3 > triPoints
label n
static void triangulate(const face &f, const pointField &points, const triangulationMode &triMode, faceList &faceTris)
Triangulate a face using the given triangulation mode.
volScalarField & p
Namespace for OpenFOAM.