triFaceI.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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2013 OpenFOAM Foundation
9  Copyright (C) 2017-2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "IOstreams.H"
30 #include "face.H"
31 #include "triPointRef.H"
32 #include "Swap.H"
33 
34 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
35 
36 inline int Foam::triFace::compare(const triFace& a, const triFace& b)
37 {
38  if
39  (
40  (a[0] == b[0] && a[1] == b[1] && a[2] == b[2])
41  || (a[0] == b[1] && a[1] == b[2] && a[2] == b[0])
42  || (a[0] == b[2] && a[1] == b[0] && a[2] == b[1])
43  )
44  {
45  // identical
46  return 1;
47  }
48  else if
49  (
50  (a[0] == b[2] && a[1] == b[1] && a[2] == b[0])
51  || (a[0] == b[1] && a[1] == b[0] && a[2] == b[2])
52  || (a[0] == b[0] && a[1] == b[2] && a[2] == b[1])
53  )
54  {
55  // same face, but reversed orientation
56  return -1;
57  }
58 
59  return 0;
60 }
61 
62 
63 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
64 
66 :
67  FixedList<label, 3>(-1)
68 {}
69 
70 
72 (
73  const label a,
74  const label b,
75  const label c
76 )
77 {
78  operator[](0) = a;
79  operator[](1) = b;
80  operator[](2) = c;
81 }
82 
83 
84 inline Foam::triFace::triFace(std::initializer_list<label> list)
85 :
86  FixedList<label, 3>(list)
87 {}
88 
89 
91 :
92  FixedList<label, 3>(list)
93 {}
94 
95 
97 (
98  const labelUList& list,
99  const FixedList<label, 3>& indices
100 )
101 :
102  FixedList<label, 3>(list, indices)
103 {}
104 
105 
107 :
108  FixedList<label, 3>(is)
109 {}
110 
111 
112 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
113 
114 inline Foam::label Foam::triFace::collapse()
115 {
116  // Cannot resize FixedList, so mark duplicates with '-1'
117  // (the lower vertex is retained)
118  // catch any '-1' (eg, if called multiple times)
119 
120  label n = 3;
121  if (operator[](0) == operator[](1) || operator[](1) == -1)
122  {
123  operator[](1) = -1;
124  n--;
125  }
126  else if (operator[](1) == operator[](2) || operator[](2) == -1)
127  {
128  operator[](2) = -1;
129  n--;
130  }
131  if (operator[](0) == operator[](2))
132  {
133  operator[](2) = -1;
134  n--;
135  }
136 
137  return n;
138 }
139 
140 
141 inline void Foam::triFace::flip()
142 {
143  Swap(operator[](1), operator[](2));
144 }
145 
146 
148 {
149  pointField p(3);
150 
151  p[0] = points[operator[](0)];
152  p[1] = points[operator[](1)];
153  p[2] = points[operator[](2)];
154 
155  return p;
156 }
157 
158 
160 {
161  return Foam::face(*this);
162 }
163 
164 
166 {
167  return triPointRef
168  (
169  points[operator[](0)],
170  points[operator[](1)],
171  points[operator[](2)]
172  );
173 }
174 
175 
177 {
178  return (1.0/3.0)*
179  (
180  points[operator[](0)]
181  + points[operator[](1)]
182  + points[operator[](2)]
183  );
184 }
185 
186 
188 {
189  // As per triPointRef(...).areaNormal()
190  return 0.5*
191  (
192  (points[operator[](1)] - points[operator[](0)])
193  ^(points[operator[](2)] - points[operator[](0)])
194  );
195 }
196 
197 
199 {
200  const vector n(areaNormal(points));
201  const scalar s(Foam::mag(n));
202  return s < ROOTVSMALL ? Zero : n/s;
203 }
204 
205 
206 inline Foam::scalar Foam::triFace::mag(const UList<point>& points) const
207 {
208  return ::Foam::mag(areaNormal(points));
209 }
210 
211 
212 inline Foam::label Foam::triFace::nTriangles() const
213 {
214  return 1;
215 }
216 
217 
219 {
220  // The starting points of the original and reverse face are identical.
221  return triFace(operator[](0), operator[](2), operator[](1));
222 }
223 
224 
225 inline bool Foam::triFace::found(const label pointLabel) const
226 {
227  return FixedList<label, 3>::found(pointLabel);
228 }
229 
230 
231 inline Foam::label Foam::triFace::which(const label pointLabel) const
232 {
233  return FixedList<label, 3>::find(pointLabel);
234 }
235 
236 
237 inline Foam::scalar Foam::triFace::sweptVol
238 (
239  const UList<point>& opts,
240  const UList<point>& npts
241 ) const
242 {
243  return (1.0/6.0)*
244  (
245  (
246  (npts[operator[](0)] - opts[operator[](0)])
247  & (
248  (opts[operator[](1)] - opts[operator[](0)])
249  ^ (opts[operator[](2)] - opts[operator[](0)])
250  )
251  )
252  + (
253  (npts[operator[](1)] - opts[operator[](1)])
254  & (
255  (opts[operator[](2)] - opts[operator[](1)])
256  ^ (npts[operator[](0)] - opts[operator[](1)])
257  )
258  )
259  + (
260  (opts[operator[](2)] - npts[operator[](2)])
261  & (
262  (npts[operator[](1)] - npts[operator[](2)])
263  ^ (npts[operator[](0)] - npts[operator[](2)])
264  )
265  )
266  );
267 }
268 
269 
271 (
272  const UList<point>& points,
273  const point& refPt,
274  scalar density
275 ) const
276 {
277  // a triangle, do a direct calculation
278  return this->tri(points).inertia(refPt, density);
279 }
280 
281 
283 (
284  const point& p,
285  const vector& q,
286  const UList<point>& points,
287  const intersection::algorithm alg,
288  const intersection::direction dir
289 ) const
290 {
291  return this->tri(points).ray(p, q, alg, dir);
292 }
293 
294 
295 
297 (
298  const point& p,
299  const vector& q,
300  const UList<point>& points,
301  const intersection::algorithm alg,
302  const scalar tol
303 ) const
304 {
305  return this->tri(points).intersection(p, q, alg, tol);
306 }
307 
308 
310 (
311  const point& p,
312  const vector& q,
313  const point& ctr,
314  const UList<point>& points,
315  const intersection::algorithm alg,
316  const scalar tol
317 ) const
318 {
319  return intersection(p, q, points, alg, tol);
320 }
321 
322 
324 (
325  const point& p,
326  const UList<point>& points
327 ) const
328 {
329  return this->tri(points).nearestPoint(p);
330 }
331 
332 
334 (
335  const point& p,
336  const UList<point>& points,
337  label& nearType,
338  label& nearLabel
339 ) const
340 {
341  return this->tri(points).nearestPointClassify(p, nearType, nearLabel);
342 }
343 
344 
345 inline int Foam::triFace::sign
346 (
347  const point& p,
348  const UList<point>& points,
349  const scalar tol
350 ) const
351 {
352  return this->tri(points).sign(p, tol);
353 }
354 
355 
356 inline Foam::label Foam::triFace::nEdges() const
357 {
358  return 3;
359 }
360 
361 
363 {
364  edgeList e(3);
365 
366  e[0].first() = operator[](0);
367  e[0].second() = operator[](1);
368 
369  e[1].first() = operator[](1);
370  e[1].second() = operator[](2);
371 
372  e[2].first() = operator[](2);
373  e[2].second() = operator[](0);
374 
375  return e;
376 }
377 
378 
379 inline Foam::edge Foam::triFace::faceEdge(const label n) const
380 {
381  return edge(operator[](n), operator[](fcIndex(n)));
382 }
383 
384 
385 // return
386 // - +1: forward (counter-clockwise) on the face
387 // - -1: reverse (clockwise) on the face
388 // - 0: edge not found on the face
389 inline int Foam::triFace::edgeDirection(const edge& e) const
390 {
391  if
392  (
393  (operator[](0) == e.first() && operator[](1) == e.second())
394  || (operator[](1) == e.first() && operator[](2) == e.second())
395  || (operator[](2) == e.first() && operator[](0) == e.second())
396  )
397  {
398  return 1;
399  }
400  else if
401  (
402  (operator[](0) == e.second() && operator[](1) == e.first())
403  || (operator[](1) == e.second() && operator[](2) == e.first())
404  || (operator[](2) == e.second() && operator[](0) == e.first())
405  )
406  {
407  return -1;
408  }
409 
410  return 0;
411 }
412 
413 
414 // * * * * * * * * * * * * * * * Global Operators * * * * * * * * * * * * * //
415 
416 inline bool Foam::operator==(const triFace& a, const triFace& b)
417 {
418  return triFace::compare(a,b) != 0;
419 }
420 
421 
422 inline bool Foam::operator!=(const triFace& a, const triFace& b)
423 {
424  return triFace::compare(a,b) == 0;
425 }
426 
427 
428 // ************************************************************************* //
Foam::triFace::inertia
tensor inertia(const UList< point > &points, const point &refPt=vector::zero, scalar density=1.0) const
Return the inertia tensor, with optional reference.
Definition: triFaceI.H:271
Foam::intersection::direction
direction
Definition: intersection.H:66
Foam::Tensor< scalar >
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::triFace::reverseFace
triFace reverseFace() const
Return face with reverse direction.
Definition: triFaceI.H:218
IOstreams.H
Useful combination of include files which define Sin, Sout and Serr and the use of IO streams general...
s
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Definition: gmvOutputSpray.H:25
Foam::PointHit
Describes the interaction of a face and a point. It carries the info of a successful hit and (if succ...
Definition: PointHit.H:53
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::triFace::areaNormal
vector areaNormal(const UList< point > &points) const
The area normal - with magnitude equal to area of face.
Definition: triFaceI.H:187
Foam::edge
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:63
Foam::triFace::found
bool found(const label pointLabel) const
Return true if the point label is found in face.
Definition: triFaceI.H:225
face.H
Foam::triFace::triFaceFace
face triFaceFace() const
Return triangle as a face.
Definition: triFaceI.H:159
Foam::FixedList::find
label find(const T &val, label pos=0) const
Find index of the first occurrence of the value.
Definition: FixedList.C:35
Foam::triFace::unitNormal
vector unitNormal(const UList< point > &points) const
The unit normal.
Definition: triFaceI.H:198
triPointRef.H
Foam::Swap
void Swap(DynamicList< T, SizeMin1 > &a, DynamicList< T, SizeMin2 > &b)
Definition: DynamicListI.H:913
Foam::triFace::tri
triPointRef tri(const UList< point > &points) const
Return the triangle.
Definition: triFaceI.H:165
Foam::triFace::flip
void flip()
Flip the face in-place.
Definition: triFaceI.H:141
Foam::triFace::collapse
label collapse()
'Collapse' face by marking duplicate point labels.
Definition: triFaceI.H:114
Swap.H
Swap arguments as per std::swap, but in Foam namespace.
Foam::operator!=
bool operator!=(const eddy &a, const eddy &b)
Definition: eddy.H:235
Foam::triangle
A triangle primitive used to calculate face normals and swept volumes.
Definition: triangle.H:61
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::triFace::which
label which(const label pointLabel) const
Find local index on face for the point label.
Definition: triFaceI.H:231
Foam::triFace::points
pointField points(const UList< point > &points) const
Return the points corresponding to this face.
Definition: triFaceI.H:147
Foam::constant::physicoChemical::b
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
Foam::Field< vector >
Foam::operator==
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
Foam::triFace::sign
int sign(const point &p, const UList< point > &points, const scalar tol=SMALL) const
The sign for which side of the face plane the point is on.
Definition: triFaceI.H:346
Foam::Istream
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:61
Foam::triFace::edges
edgeList edges() const
Return edges in face point ordering,.
Definition: triFaceI.H:362
Foam::intersection::algorithm
algorithm
Definition: intersection.H:72
Foam::triFace::faceEdge
edge faceEdge(const label n) const
Return n-th face edge.
Definition: triFaceI.H:379
Foam::triFace::sweptVol
scalar sweptVol(const UList< point > &opts, const UList< point > &npts) const
Return swept-volume from old-points to new-points.
Definition: triFaceI.H:238
Foam::triFace::centre
point centre(const UList< point > &points) const
Return centre (centroid)
Definition: triFaceI.H:176
Foam::triFace::nearestPoint
pointHit nearestPoint(const point &p, const UList< point > &points) const
Return nearest point to face.
Definition: triFaceI.H:324
Foam::triFace::nearestPointClassify
pointHit nearestPointClassify(const point &p, const UList< point > &points, label &nearType, label &nearLabel) const
Return nearest point to face and classify it:
Definition: triFaceI.H:334
Foam::triFace
A triangular face using a FixedList of labels corresponding to mesh vertices.
Definition: triFace.H:69
Foam::triFace::ray
pointHit ray(const point &p, const vector &q, const UList< point > &points, const intersection::algorithm=intersection::FULL_RAY, const intersection::direction dir=intersection::VECTOR) const
Return point intersection with a ray starting at p, in direction q.
Definition: triFaceI.H:283
Foam::Vector< scalar >
Foam::triFace::mag
scalar mag(const UList< point > &points) const
Magnitude of face area.
Definition: triFaceI.H:206
Foam::List< edge >
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::FixedList
A 1D vector of objects of type <T> with a fixed length <N>.
Definition: HashTable.H:104
Foam::triFace::intersection
pointHit intersection(const point &p, const vector &q, const UList< point > &points, const intersection::algorithm alg, const scalar tol=0.0) const
Fast intersection with a ray.
Definition: triFaceI.H:297
Foam::UList< label >
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::intersection
Foam::intersection.
Definition: intersection.H:52
Foam::constant::universal::c
const dimensionedScalar c
Speed of light in a vacuum.
Foam::triFace::edgeDirection
int edgeDirection(const edge &e) const
Return the edge direction on the face.
Definition: triFaceI.H:389
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:72
Foam::FixedList::found
bool found(const T &val, label pos=0) const
True if the value if found in the list.
Definition: FixedListI.H:313
Foam::triFace::triFace
triFace()
Default construct, with invalid point labels (-1)
Definition: triFaceI.H:65
Foam::triFace::compare
static int compare(const triFace &a, const triFace &b)
Compare triFaces.
Definition: triFaceI.H:36
Foam::triFace::nTriangles
label nTriangles() const
Number of triangles after splitting.
Definition: triFaceI.H:212
Foam::triFace::nEdges
label nEdges() const
Return number of edges.
Definition: triFaceI.H:356
triFace
face triFace(3)
Foam::triPointRef
triangle< point, const point & > triPointRef
A triangle using referred points.
Definition: triangle.H:73