triFaceI.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-2018 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 \*---------------------------------------------------------------------------*/
25 
26 #include "IOstreams.H"
27 #include "face.H"
28 #include "triPointRef.H"
29 #include "Swap.H"
30 
31 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
32 
33 inline int Foam::triFace::compare(const triFace& a, const triFace& b)
34 {
35  if
36  (
37  (a[0] == b[0] && a[1] == b[1] && a[2] == b[2])
38  || (a[0] == b[1] && a[1] == b[2] && a[2] == b[0])
39  || (a[0] == b[2] && a[1] == b[0] && a[2] == b[1])
40  )
41  {
42  // identical
43  return 1;
44  }
45  else if
46  (
47  (a[0] == b[2] && a[1] == b[1] && a[2] == b[0])
48  || (a[0] == b[1] && a[1] == b[0] && a[2] == b[2])
49  || (a[0] == b[0] && a[1] == b[2] && a[2] == b[1])
50  )
51  {
52  // same face, but reversed orientation
53  return -1;
54  }
55  else
56  {
57  return 0;
58  }
59 }
60 
61 
62 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
63 
65 {}
66 
67 
69 (
70  const label a,
71  const label b,
72  const label c
73 )
74 {
75  operator[](0) = a;
76  operator[](1) = b;
77  operator[](2) = c;
78 }
79 
80 
82 :
83  FixedList<label, 3>(lst)
84 {}
85 
86 
88 :
89  FixedList<label, 3>(is)
90 {}
91 
92 
93 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
94 
96 {
97  // we cannot resize a FixedList, so mark duplicates with '-1'
98  // (the lower vertex is retained)
99  // catch any '-1' (eg, if called twice)
100 
101  label n = 3;
102  if (operator[](0) == operator[](1) || operator[](1) == -1)
103  {
104  operator[](1) = -1;
105  n--;
106  }
107  else if (operator[](1) == operator[](2) || operator[](2) == -1)
108  {
109  operator[](2) = -1;
110  n--;
111  }
112  if (operator[](0) == operator[](2))
113  {
114  operator[](2) = -1;
115  n--;
116  }
117 
118  return n;
119 }
120 
121 
122 inline void Foam::triFace::flip()
123 {
124  Swap(operator[](1), operator[](2));
125 }
126 
127 
129 {
130  pointField p(3);
131 
132  p[0] = points[operator[](0)];
133  p[1] = points[operator[](1)];
134  p[2] = points[operator[](2)];
135 
136  return p;
137 }
138 
139 
141 {
142  Foam::face f(3);
143 
144  f[0] = operator[](0);
145  f[1] = operator[](1);
146  f[2] = operator[](2);
147 
148  return f;
149 }
150 
151 
153 {
154  return triPointRef
155  (
156  points[operator[](0)],
157  points[operator[](1)],
158  points[operator[](2)]
159  );
160 }
161 
162 
163 inline Foam::point Foam::triFace::centre(const pointField& points) const
164 {
165  return (1.0/3.0)*
166  (
167  points[operator[](0)]
168  + points[operator[](1)]
169  + points[operator[](2)]
170  );
171 }
172 
173 
174 inline Foam::vector Foam::triFace::area(const pointField& points) const
175 {
176  return 0.5*
177  (
178  (points[operator[](1)] - points[operator[](0)])
179  ^(points[operator[](2)] - points[operator[](0)])
180  );
181 }
182 
183 
184 inline Foam::scalar Foam::triFace::mag(const pointField& points) const
185 {
186  return ::Foam::mag(area(points));
187 }
188 
189 
190 inline Foam::vector Foam::triFace::normal(const pointField& points) const
191 {
192  const vector a = area(points);
193  const scalar maga = Foam::mag(a);
194  return maga > 0 ? a/maga : Zero;
195 }
196 
197 
199 {
200  return 1;
201 }
202 
203 
205 {
206  // The starting points of the original and reverse face are identical.
207  return triFace(operator[](0), operator[](2), operator[](1));
208 }
209 
210 
211 inline Foam::scalar Foam::triFace::sweptVol
212 (
213  const pointField& opts,
214  const pointField& npts
215 ) const
216 {
217  return (1.0/6.0)*
218  (
219  (
220  (npts[operator[](0)] - opts[operator[](0)])
221  & (
222  (opts[operator[](1)] - opts[operator[](0)])
223  ^ (opts[operator[](2)] - opts[operator[](0)])
224  )
225  )
226  + (
227  (npts[operator[](1)] - opts[operator[](1)])
228  & (
229  (opts[operator[](2)] - opts[operator[](1)])
230  ^ (npts[operator[](0)] - opts[operator[](1)])
231  )
232  )
233  + (
234  (opts[operator[](2)] - npts[operator[](2)])
235  & (
236  (npts[operator[](1)] - npts[operator[](2)])
237  ^ (npts[operator[](0)] - npts[operator[](2)])
238  )
239  )
240  );
241 }
242 
243 
245 (
246  const pointField& points,
247  const point& refPt,
248  scalar density
249 ) const
250 {
251  // a triangle, do a direct calculation
252  return this->tri(points).inertia(refPt, density);
253 }
254 
255 
257 (
258  const point& p,
259  const vector& q,
260  const pointField& points,
261  const intersection::algorithm alg,
262  const intersection::direction dir
263 ) const
264 {
265  return this->tri(points).ray(p, q, alg, dir);
266 }
267 
268 
269 
271 (
272  const point& p,
273  const vector& q,
274  const pointField& points,
275  const intersection::algorithm alg,
276  const scalar tol
277 ) const
278 {
279  return this->tri(points).intersection(p, q, alg, tol);
280 }
281 
282 
284 (
285  const point& p,
286  const vector& q,
287  const point& ctr,
288  const pointField& points,
289  const intersection::algorithm alg,
290  const scalar tol
291 ) const
292 {
293  return intersection(p, q, points, alg, tol);
294 }
295 
296 
298 (
299  const point& p,
300  const pointField& points
301 ) const
302 {
303  return this->tri(points).nearestPoint(p);
304 }
305 
306 
308 (
309  const point& p,
310  const pointField& points,
311  label& nearType,
312  label& nearLabel
313 ) const
314 {
315  return this->tri(points).nearestPointClassify(p, nearType, nearLabel);
316 }
317 
318 
320 {
321  return 3;
322 }
323 
324 
326 {
327  edgeList e(3);
328 
329  e[0].start() = operator[](0);
330  e[0].end() = operator[](1);
331 
332  e[1].start() = operator[](1);
333  e[1].end() = operator[](2);
334 
335  e[2].start() = operator[](2);
336  e[2].end() = operator[](0);
337 
338  return e;
339 }
340 
341 
343 {
344  return edge(operator[](n), operator[](fcIndex(n)));
345 }
346 
347 
348 // return
349 // - +1: forward (counter-clockwise) on the face
350 // - -1: reverse (clockwise) on the face
351 // - 0: edge not found on the face
352 inline int Foam::triFace::edgeDirection(const edge& e) const
353 {
354  if
355  (
356  (operator[](0) == e.start() && operator[](1) == e.end())
357  || (operator[](1) == e.start() && operator[](2) == e.end())
358  || (operator[](2) == e.start() && operator[](0) == e.end())
359  )
360  {
361  return 1;
362  }
363  else if
364  (
365  (operator[](0) == e.end() && operator[](1) == e.start())
366  || (operator[](1) == e.end() && operator[](2) == e.start())
367  || (operator[](2) == e.end() && operator[](0) == e.start())
368  )
369  {
370  return -1;
371  }
372  else
373  {
374  return 0;
375  }
376 }
377 
378 
379 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
380 
381 inline bool Foam::operator==(const triFace& a, const triFace& b)
382 {
383  return triFace::compare(a,b) != 0;
384 }
385 
386 
387 inline bool Foam::operator!=(const triFace& a, const triFace& b)
388 {
389  return triFace::compare(a,b) == 0;
390 }
391 
392 
393 // ************************************************************************* //
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
A triangle primitive used to calculate face areas and swept volumes.
Definition: triangle.H:57
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
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
pointHit nearestPointClassify(const point &p, label &nearType, label &nearLabel) const
Find the nearest point to p on the triangle and classify it:
Definition: triangleI.H:518
edgeList edges() const
Return edges in face point ordering,.
Definition: triFaceI.H:325
vector area(const pointField &) const
Return vector area.
Definition: triFaceI.H:174
iterator end()
Return an iterator to end traversing the UList.
Definition: UListI.H:224
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
Useful combination of include files which define Sin, Sout and Serr and the use of IO streams general...
pointHit nearestPoint(const point &p) const
Return nearest point to p on triangle.
Definition: triangleI.H:667
void Swap(T &a, T &b)
Definition: Swap.H:43
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 > &)
label & operator[](const label)
Return element of FixedList.
Definition: FixedListI.H:247
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
static const zero Zero
Definition: zero.H:97
triPointRef tri(const pointField &) const
Return the triangle.
Definition: triFaceI.H:152
tensor inertia(PointRef refPt=vector::zero, scalar density=1.0) const
Return the inertia tensor, with optional reference.
Definition: triangleI.H:205
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 fcIndex(const label i) const
Return the forward circular index, i.e. the next index.
Definition: FixedListI.H:125
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
labelList f(nPoints)
triFace()
Construct null.
Definition: triFaceI.H:64
pointHit ray(const point &p, const vector &q, const intersection::algorithm=intersection::algorithm::fullRay, const intersection::direction dir=intersection::direction::vector) const
Return point intersection with a ray.
Definition: triangleI.H:310
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
label end() const
Return end vertex label.
Definition: edgeI.H:92
pointHit nearestPoint(const point &p, const pointField &points) const
Return nearest point to face.
Definition: triFaceI.H:298
Swap its arguments.
triangle< point, const point & > triPointRef
Definition: triPointRef.H:44
void flip()
Flip the face in-place.
Definition: triFaceI.H:122
dimensioned< scalar > mag(const dimensioned< Type > &)
label n
const doubleScalar e
Elementary charge.
Definition: doubleScalar.H:105
volScalarField & p
pointHit intersection(const point &p, const vector &q, const intersection::algorithm alg, const scalar tol=0.0) const
Fast intersection with a ray.
Definition: triangleI.H:431
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:1205
label nTriangles() const
Number of triangles after splitting.
Definition: triFaceI.H:198
label start() const
Return start vertex label.
Definition: edgeI.H:81
face triFaceFace() const
Return triangle as a face.
Definition: triFaceI.H:140