triFace.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-2016 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  //- Magnitude of face area
123  inline scalar mag(const pointField&) const;
124 
125  //- Vector normal; magnitude is equal to area of face
126  inline vector normal(const pointField&) const;
127 
128  //- Number of triangles after splitting
129  inline label nTriangles() const;
130 
131  //- Return face with reverse direction
132  // The starting points of the original and reverse face are identical.
133  inline triFace reverseFace() const;
134 
135  //- Return swept-volume
136  inline scalar sweptVol
137  (
138  const pointField& oldPoints,
139  const pointField& newPoints
140  ) const;
141 
142  //- Return the inertia tensor, with optional reference
143  // point and density specification
144  inline tensor inertia
145  (
146  const pointField&,
147  const point& refPt = vector::zero,
148  scalar density = 1.0
149  ) const;
150 
151  //- Return point intersection with a ray starting at p,
152  // with direction q.
153  inline pointHit ray
154  (
155  const point& p,
156  const vector& q,
157  const pointField& points,
160  ) const;
161 
162  //- Fast intersection with a ray.
163  inline pointHit intersection
164  (
165  const point& p,
166  const vector& q,
167  const pointField& points,
168  const intersection::algorithm alg,
169  const scalar tol = 0.0
170  ) const;
171 
172  inline pointHit intersection
173  (
174  const point& p,
175  const vector& q,
176  const point& ctr,
177  const pointField& points,
178  const intersection::algorithm alg,
179  const scalar tol = 0.0
180  ) const;
181 
182  //- Return nearest point to face
183  inline pointHit nearestPoint
184  (
185  const point& p,
186  const pointField& points
187  ) const;
188 
189 
190  //- Return nearest point to face and classify it:
191  // + near point (nearType=POINT, nearLabel=0, 1, 2)
192  // + near edge (nearType=EDGE, nearLabel=0, 1, 2)
193  // Note: edges are counted from starting vertex so
194  // e.g. edge n is from f[n] to f[0], where the face has n + 1
195  // points
197  (
198  const point& p,
199  const pointField& points,
200  label& nearType,
201  label& nearLabel
202  ) const;
203 
204  //- Return number of edges
205  inline label nEdges() const;
206 
207  //- Return edges in face point ordering,
208  // i.e. edges()[0] is edge between [0] and [1]
209  inline edgeList edges() const;
210 
211  //- Return n-th face edge
212  inline edge faceEdge(const label n) const;
213 
214  //- Return the edge direction on the face
215  // Returns:
216  // - +1: forward (counter-clockwise) on the face
217  // - -1: reverse (clockwise) on the face
218  // - 0: edge not found on the face
219  inline int edgeDirection(const edge&) const;
220 
221  //- Compare triFaces
222  // Returns:
223  // - 0: different
224  // - +1: identical
225  // - -1: same face, but different orientation
226  static inline int compare(const triFace&, const triFace&);
227 
228  // Friend Operators
229 
230  inline friend bool operator==(const triFace&, const triFace&);
231  inline friend bool operator!=(const triFace&, const triFace&);
232 };
233 
234 
235 //- Hash specialization for hashing triFace - a commutative hash value.
236 // Hash incrementally.
237 template<>
238 inline unsigned Hash<triFace>::operator()(const triFace& t, unsigned seed) const
239 {
240  // Fortunately we don't need this very often
241  const uLabel t0(t[0]);
242  const uLabel t1(t[1]);
243  const uLabel t2(t[2]);
244 
245  const uLabel val = (t0*t1*t2 + t0+t1+t2);
246 
247  return Hash<uLabel>()(val, seed);
248 }
249 
250 
251 //- Hash specialization for hashing triFace - a commutative hash value.
252 template<>
253 inline unsigned Hash<triFace>::operator()(const triFace& t) const
254 {
255  return Hash<triFace>::operator()(t, 0);
256 }
257 
258 
259 template<>
260 inline bool contiguous<triFace>() {return true;}
261 
262 
263 //- Hash specialization to offset faces in ListListOps::combineOffset
264 template<>
265 class offsetOp<triFace>
266 {
267 
268 public:
269 
270  inline triFace operator()
271  (
272  const triFace& x,
273  const label offset
274  ) const
275  {
276  triFace result(x);
277 
278  forAll(x, xI)
279  {
280  result[xI] = x[xI] + offset;
281  }
282  return result;
283  }
284 };
285 
286 
287 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
288 
289 } // End namespace Foam
290 
291 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
292 
293 #include "triFaceI.H"
294 
295 #ifdef NoRepository
296  #include "triFaceTemplates.C"
297 #endif
298 
299 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
300 
301 #endif
302 
303 // ************************************************************************* //
pointHit nearestPointClassify(const point &p, const pointField &points, label &nearType, label &nearLabel) const
Return nearest point to face and classify it:
Definition: triFaceI.H:301
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:428
A triangle primitive used to calculate face normals and swept volumes.
Definition: triangle.H:58
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
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:197
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
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:318
unsigned operator()(const PrimitiveType &p, unsigned seed) const
Definition: Hash.H:61
int edgeDirection(const edge &) const
Return the edge direction on the face.
Definition: triFaceI.H:345
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:264
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:58
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
scalar mag(const pointField &) const
Magnitude of face area.
Definition: triFaceI.H:174
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:61
label nEdges() const
Return number of edges.
Definition: triFaceI.H:312
scalar sweptVol(const pointField &oldPoints, const pointField &newPoints) const
Return swept-volume.
Definition: triFaceI.H:205
triFace()
Construct null.
Definition: triFaceI.H:64
vector normal(const pointField &) const
Vector normal; magnitude is equal to area of face.
Definition: triFaceI.H:181
edge faceEdge(const label n) const
Return n-th face edge.
Definition: triFaceI.H:335
pointHit nearestPoint(const point &p, const pointField &points) const
Return nearest point to face.
Definition: triFaceI.H:291
const dimensionedScalar c
Speed of light in a vacuum.
Hash function class for primitives. All non-primitives used to hash entries on hash tables likely nee...
Definition: Hash.H:54
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
pointHit ray(const point &p, const vector &q, const pointField &points, const intersection::algorithm=intersection::FULL_RAY, const intersection::direction dir=intersection::VECTOR) const
Return point intersection with a ray starting at p,.
Definition: triFaceI.H:250
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:238
bool operator!=(const particle &, const particle &)
Definition: particle.C:1106
bool contiguous< triFace >()
Definition: triFace.H:259
label nTriangles() const
Number of triangles after splitting.
Definition: triFaceI.H:191
face triFaceFace() const
Return triangle as a face.
Definition: triFaceI.H:140
Namespace for OpenFOAM.