All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
triFace.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::triFace
26 
27 Description
28  A triangular face using a FixedList of labels corresponding to mesh
29  vertices.
30 
31 See also
32  Foam::face, Foam::triangle
33 
34 SourceFiles
35  triFaceI.H
36  triFaceTemplates.C
37 
38 \*---------------------------------------------------------------------------*/
39 
40 #ifndef triFace_H
41 #define triFace_H
42 
43 #include "FixedList.H"
44 #include "edgeList.H"
45 #include "pointHit.H"
46 #include "intersection.H"
47 #include "pointField.H"
48 #include "triPointRef.H"
49 #include "ListListOps.H"
50 
51 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
52 
53 namespace Foam
54 {
55 
56 // Forward declaration of friend functions and operators
57 
58 class face;
59 class triFace;
60 
61 inline bool operator==(const triFace&, const triFace&);
62 inline bool operator!=(const triFace&, const triFace&);
63 
64 
65 /*---------------------------------------------------------------------------*\
66  Class triFace Declaration
67 \*---------------------------------------------------------------------------*/
68 
69 class triFace
70 :
71  public FixedList<label, 3>
72 {
73 
74 public:
75 
76  // Constructors
77 
78  //- Construct null
79  inline triFace();
80 
81  //- Construct from three point labels
82  inline triFace
83  (
84  const label a,
85  const label b,
86  const label c
87  );
88 
89  //- Construct from a list of labels
90  explicit inline triFace(const labelUList&);
91 
92  //- Construct from Istream
93  inline triFace(Istream&);
94 
95 
96  // Member Functions
97 
98  //- Collapse face by removing duplicate point labels
99  // return the collapsed size, set collapsed point labels to -1
100  inline label collapse();
101 
102  //- Flip the face in-place.
103  // The starting points of the original and flipped face are identical.
104  inline void flip();
105 
106  //- Return the points corresponding to this face
107  inline pointField points(const pointField&) const;
108 
109  //- Return triangle as a face
110  inline face triFaceFace() const;
111 
112  //- Return the triangle
113  inline triPointRef tri(const pointField&) const;
114 
115  //- Return centre (centroid)
116  inline point centre(const pointField&) const;
117 
118  //- Calculate average value at centroid of face
119  template<class Type>
120  Type average(const pointField&, const Field<Type>&) const;
121 
122  //- Return vector area
123  inline vector area(const pointField&) const;
124 
125  //- Return scalar magnitude
126  inline scalar mag(const pointField&) const;
127 
128  //- Return unit normal
129  inline vector normal(const pointField&) const;
130 
131  //- Number of triangles after splitting
132  inline label nTriangles() const;
133 
134  //- Return face with reverse direction
135  // The starting points of the original and reverse face are identical.
136  inline triFace reverseFace() const;
137 
138  //- Return swept-volume
139  inline scalar sweptVol
140  (
141  const pointField& oldPoints,
142  const pointField& newPoints
143  ) const;
144 
145  //- Return the inertia tensor, with optional reference
146  // point and density specification
147  inline tensor inertia
148  (
149  const pointField&,
150  const point& refPt = vector::zero,
151  scalar density = 1.0
152  ) const;
153 
154  //- Return point intersection with a ray starting at p,
155  // with direction q.
156  inline pointHit ray
157  (
158  const point& p,
159  const vector& q,
160  const pointField& points,
163  ) const;
164 
165  //- Fast intersection with a ray.
166  inline pointHit intersection
167  (
168  const point& p,
169  const vector& q,
170  const pointField& points,
171  const intersection::algorithm alg,
172  const scalar tol = 0.0
173  ) const;
174 
175  inline pointHit intersection
176  (
177  const point& p,
178  const vector& q,
179  const point& ctr,
180  const pointField& points,
181  const intersection::algorithm alg,
182  const scalar tol = 0.0
183  ) const;
184 
185  //- Return nearest point to face
186  inline pointHit nearestPoint
187  (
188  const point& p,
189  const pointField& points
190  ) const;
191 
192 
193  //- Return nearest point to face and classify it:
194  // + near point (nearType=POINT, nearLabel=0, 1, 2)
195  // + near edge (nearType=EDGE, nearLabel=0, 1, 2)
196  // Note: edges are counted from starting vertex so
197  // e.g. edge n is from f[n] to f[0], where the face has n + 1
198  // points
200  (
201  const point& p,
202  const pointField& points,
203  label& nearType,
204  label& nearLabel
205  ) const;
206 
207  //- Return number of edges
208  inline label nEdges() const;
209 
210  //- Return edges in face point ordering,
211  // i.e. edges()[0] is edge between [0] and [1]
212  inline edgeList edges() const;
213 
214  //- Return n-th face edge
215  inline edge faceEdge(const label n) const;
216 
217  //- Return the edge direction on the face
218  // Returns:
219  // - +1: forward (counter-clockwise) on the face
220  // - -1: reverse (clockwise) on the face
221  // - 0: edge not found on the face
222  inline int edgeDirection(const edge&) const;
223 
224  //- Compare triFaces
225  // Returns:
226  // - 0: different
227  // - +1: identical
228  // - -1: same face, but different orientation
229  static inline int compare(const triFace&, const triFace&);
230 
231  // Friend Operators
232 
233  inline friend bool operator==(const triFace&, const triFace&);
234  inline friend bool operator!=(const triFace&, const triFace&);
235 };
236 
237 
238 //- Hash specialisation for hashing triFace - a commutative hash value.
239 // Hash incrementally.
240 template<>
241 inline unsigned Hash<triFace>::operator()(const triFace& t, unsigned seed) const
242 {
243  // Fortunately we don't need this very often
244  const uLabel t0(t[0]);
245  const uLabel t1(t[1]);
246  const uLabel t2(t[2]);
247 
248  const uLabel val = (t0*t1*t2 + t0+t1+t2);
249 
250  return Hash<uLabel>()(val, seed);
251 }
252 
253 
254 //- Hash specialisation for hashing triFace - a commutative hash value.
255 template<>
256 inline unsigned Hash<triFace>::operator()(const triFace& t) const
257 {
258  return Hash<triFace>::operator()(t, 0);
259 }
260 
261 
262 template<>
263 inline bool contiguous<triFace>() {return true;}
264 
265 
266 //- Hash specialisation to offset faces in ListListOps::combineOffset
267 template<>
268 class offsetOp<triFace>
269 {
270 
271 public:
272 
273  inline triFace operator()
274  (
275  const triFace& x,
276  const label offset
277  ) const
278  {
279  triFace result(x);
280 
281  forAll(x, xI)
282  {
283  result[xI] = x[xI] + offset;
284  }
285  return result;
286  }
287 };
288 
289 
290 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
291 
292 } // End namespace Foam
293 
294 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
295 
296 #include "triFaceI.H"
297 
298 #ifdef NoRepository
299  #include "triFaceTemplates.C"
300 #endif
301 
302 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
303 
304 #endif
305 
306 // ************************************************************************* //
pointHit nearestPointClassify(const point &p, const pointField &points, label &nearType, label &nearLabel) const
Return nearest point to face and classify it:
Definition: triFaceI.H:308
label collapse()
Collapse face by removing duplicate point labels.
Definition: triFaceI.H:95
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
A triangle primitive used to calculate face areas and swept volumes.
Definition: triangle.H:57
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
Type average(const pointField &, const Field< Type > &) const
Calculate average value at centroid of face.
A 1D vector of objects of type <T> with a fixed size <Size>.
Definition: FixedList.H:54
triFace reverseFace() const
Return face with reverse direction.
Definition: triFaceI.H:204
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
pointHit ray(const point &p, const vector &q, const pointField &points, const intersection::algorithm=intersection::algorithm::fullRay, const intersection::direction dir=intersection::direction::vector) const
Return point intersection with a ray starting at p,.
Definition: triFaceI.H:257
const dimensionedScalar b
Wien displacement law constant: default SI units: [m K].
Definition: createFields.H:27
point centre(const pointField &) const
Return centre (centroid)
Definition: triFaceI.H:163
friend bool operator==(const triFace &, const triFace &)
friend Ostream & operator(Ostream &, const FixedList< label, Size > &)
Write FixedList to Ostream.
edgeList edges() const
Return edges in face point ordering,.
Definition: triFaceI.H:325
unsigned operator()(const PrimitiveType &p, unsigned seed) const
Definition: Hash.H:59
vector area(const pointField &) const
Return vector area.
Definition: triFaceI.H:174
const dimensionedScalar c
Speed of light in a vacuum.
int edgeDirection(const edge &) const
Return the edge direction on the face.
Definition: triFaceI.H:352
pointHit intersection(const point &p, const vector &q, const pointField &points, const intersection::algorithm alg, const scalar tol=0.0) const
Fast intersection with a ray.
Definition: triFaceI.H:271
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:58
scalar mag(const pointField &) const
Return scalar magnitude.
Definition: triFaceI.H:184
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
A triangular face using a FixedList of labels corresponding to mesh vertices.
Definition: triFace.H:68
static int compare(const triFace &, const triFace &)
Compare triFaces.
Definition: triFaceI.H:33
pointField points(const pointField &) const
Return the points corresponding to this face.
Definition: triFaceI.H:128
face triFace(3)
friend bool operator!=(const triFace &, const triFace &)
triPointRef tri(const pointField &) const
Return the triangle.
Definition: triFaceI.H:152
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: HashTable.H:60
label nEdges() const
Return number of edges.
Definition: triFaceI.H:319
scalar sweptVol(const pointField &oldPoints, const pointField &newPoints) const
Return swept-volume.
Definition: triFaceI.H:212
triFace()
Construct null.
Definition: triFaceI.H:64
vector normal(const pointField &) const
Return unit normal.
Definition: triFaceI.H:190
edge faceEdge(const label n) const
Return n-th face edge.
Definition: triFaceI.H:342
pointHit nearestPoint(const point &p, const pointField &points) const
Return nearest point to face.
Definition: triFaceI.H:298
Hash function class for primitives. All non-primitives used to hash entries on hash tables likely nee...
Definition: Hash.H:52
void flip()
Flip the face in-place.
Definition: triFaceI.H:122
uintWM_LABEL_SIZE_t uLabel
A uLabel is an uint32_t or uint64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: uLabel.H:59
label n
volScalarField & p
tensor inertia(const pointField &, const point &refPt=vector::zero, scalar density=1.0) const
Return the inertia tensor, with optional reference.
Definition: triFaceI.H:245
bool operator!=(const particle &, const particle &)
Definition: particle.C:1282
bool contiguous< triFace >()
Definition: triFace.H:262
label nTriangles() const
Number of triangles after splitting.
Definition: triFaceI.H:198
face triFaceFace() const
Return triangle as a face.
Definition: triFaceI.H:140
Namespace for OpenFOAM.