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 | Copyright (C) 2011-2015 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 "triPoints.H"
45 #include "NamedEnum.H"
46 
47 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
48 
49 namespace Foam
50 {
51 
52 /*---------------------------------------------------------------------------*\
53  Class faceAreaIntersect Declaration
54 \*---------------------------------------------------------------------------*/
55 
57 {
58 public:
59 
61  {
63  tmMesh
64  };
65 
67 
68 
69 private:
70 
71  // Private data
72 
73  //- Reference to the points of sideA
74  const pointField& pointsA_;
75 
76  //- Reference to the points of sideB
77  const pointField& pointsB_;
78 
79  //- Flag to reverse B faces
80  const bool reverseB_;
81 
82 
83  // Static data members
84 
85  static scalar tol;
86 
87 
88  // Private Member Functions
89 
90  //- Get triPoints from face
91  inline triPoints getTriPoints
92  (
93  const pointField& points,
94  const face& f,
95  const bool reverse
96  ) const;
97 
98  //- Set triPoints into tri list
99  inline void setTriPoints
100  (
101  const point& a,
102  const point& b,
103  const point& c,
104  label& count,
106  ) const;
107 
108  //- Decompose face into triangle fan
109  inline void triangleFan
110  (
111  const face& f,
112  DynamicList<face>& faces
113  ) const;
114 
115  //- Return point of intersection between plane and triangle edge
116  inline point planeIntersection
117  (
118  const FixedList<scalar, 3>& d,
119  const triPoints& t,
120  const label negI,
121  const label posI
122  ) const;
123 
124  //- Return triangle area
125  inline scalar triArea(const triPoints& t) const;
126 
127 
128  //- Slice triangle with plane and generate new cut sub-triangles
129  void triSliceWithPlane
130  (
131  const triPoints& tri,
132  const plane& p,
134  label& nTris,
135  const scalar len
136  );
137 
138  //- Return area of intersection of triangles src and tgt
139  scalar triangleIntersect
140  (
141  const triPoints& src,
142  const triPoints& tgt,
143  const vector& n
144  );
145 
146 
147 public:
148 
149  // Constructors
150 
151  //- Construct from components
153  (
154  const pointField& pointsA,
155  const pointField& pointsB,
156  const bool reverseB = false
157  );
158 
159 
160  // Public Member Functions
161 
162  //- Fraction of local length scale to use as intersection tolerance
163  inline static scalar& tolerance();
164 
165  //- Return area of intersection of faceA with faceB
166  scalar calc
167  (
168  const face& faceA,
169  const face& faceB,
170  const vector& n,
171  const triangulationMode& triMode
172  );
173 };
174 
175 
176 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
177 
178 } // End namespace Foam
179 
180 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
181 
182 #include "faceAreaIntersectI.H"
183 
184 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
185 
186 #endif
187 
188 // ************************************************************************* //
Foam::faceAreaIntersect::triangleFan
void triangleFan(const face &f, DynamicList< face > &faces) const
Decompose face into triangle fan.
Definition: faceAreaIntersectI.H:71
Foam::faceAreaIntersect::pointsA_
const pointField & pointsA_
Reference to the points of sideA.
Definition: faceAreaIntersect.H:73
p
p
Definition: pEqn.H:62
Foam::faceAreaIntersect::pointsB_
const pointField & pointsB_
Reference to the points of sideB.
Definition: faceAreaIntersect.H:76
Foam::faceAreaIntersect::triangulationModeNames_
static const NamedEnum< triangulationMode, 2 > triangulationModeNames_
Definition: faceAreaIntersect.H:65
Foam::faceAreaIntersect::triangulationMode
triangulationMode
Definition: faceAreaIntersect.H:59
Foam::DynamicList
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:56
NamedEnum.H
Foam::faceAreaIntersect::reverseB_
const bool reverseB_
Flag to reverse B faces.
Definition: faceAreaIntersect.H:79
Foam::faceAreaIntersect::setTriPoints
void setTriPoints(const point &a, const point &b, const point &c, label &count, FixedList< triPoints, 10 > &tris) const
Set triPoints into tri list.
Definition: faceAreaIntersectI.H:29
face.H
Foam::faceAreaIntersect::triangleIntersect
scalar triangleIntersect(const triPoints &src, const triPoints &tgt, const vector &n)
Return area of intersection of triangles src and tgt.
Definition: faceAreaIntersect.C:223
triPoints.H
Foam::faceAreaIntersect::triSliceWithPlane
void triSliceWithPlane(const triPoints &tri, const plane &p, FixedList< triPoints, 10 > &tris, label &nTris, const scalar len)
Slice triangle with plane and generate new cut sub-triangles.
Definition: faceAreaIntersect.C:48
Foam::faceAreaIntersect::planeIntersection
point planeIntersection(const FixedList< scalar, 3 > &d, const triPoints &t, const label negI, const label posI) const
Return point of intersection between plane and triangle edge.
Definition: faceAreaIntersectI.H:94
Foam::plane
Geometric class that creates a 2D plane and can return the intersection point between a line and the ...
Definition: plane.H:60
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::triPoints
Triangle storage. Null constructable (unfortunately triangle<point, point> is not)
Definition: triPoints.H:49
Foam::label
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
Foam::constant::physicoChemical::b
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:28
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
faceAreaIntersectI.H
plane.H
Foam::faceAreaIntersect::getTriPoints
triPoints getTriPoints(const pointField &points, const face &f, const bool reverse) const
Get triPoints from face.
Definition: faceAreaIntersectI.H:45
Foam::faceAreaIntersect
Face intersection class.
Definition: faceAreaIntersect.H:55
Foam::faceAreaIntersect::tmFan
@ tmFan
Definition: faceAreaIntersect.H:61
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
pointField.H
Foam::faceAreaIntersect::triArea
scalar triArea(const triPoints &t) const
Return triangle area.
Definition: faceAreaIntersectI.H:105
f
labelList f(nPoints)
Foam::Vector< scalar >
Foam::faceAreaIntersect::tolerance
static scalar & tolerance()
Fraction of local length scale to use as intersection tolerance.
Definition: faceAreaIntersectI.H:113
Foam::FixedList
A 1D vector of objects of type <T> with a fixed size <Size>.
Definition: FixedList.H:53
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::faceAreaIntersect::calc
scalar calc(const face &faceA, const face &faceB, const vector &n, const triangulationMode &triMode)
Return area of intersection of faceA with faceB.
Definition: faceAreaIntersect.C:329
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:75
Foam::faceAreaIntersect::faceAreaIntersect
faceAreaIntersect(const pointField &pointsA, const pointField &pointsB, const bool reverseB=false)
Construct from components.
Definition: faceAreaIntersect.C:314
FixedList.H
Foam::faceAreaIntersect::tol
static scalar tol
Definition: faceAreaIntersect.H:84
Foam::NamedEnum< triangulationMode, 2 >
Foam::faceAreaIntersect::tmMesh
@ tmMesh
Definition: faceAreaIntersect.H:62
Foam::reverse
void reverse(UList< T > &, const label n)
Definition: UListI.H:322