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-2018 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 
85 :
86  FixedList<label, 3>(lst)
87 {}
88 
89 
90 inline Foam::triFace::triFace(std::initializer_list<label> lst)
91 :
92  FixedList<label, 3>(lst)
93 {}
94 
95 
97 :
98  FixedList<label, 3>(is)
99 {}
100 
101 
102 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
103 
105 {
106  // Cannot resize FixedList, so mark duplicates with '-1'
107  // (the lower vertex is retained)
108  // catch any '-1' (eg, if called multiple times)
109 
110  label n = 3;
111  if (operator[](0) == operator[](1) || operator[](1) == -1)
112  {
113  operator[](1) = -1;
114  n--;
115  }
116  else if (operator[](1) == operator[](2) || operator[](2) == -1)
117  {
118  operator[](2) = -1;
119  n--;
120  }
121  if (operator[](0) == operator[](2))
122  {
123  operator[](2) = -1;
124  n--;
125  }
126 
127  return n;
128 }
129 
130 
131 inline void Foam::triFace::flip()
132 {
133  Swap(operator[](1), operator[](2));
134 }
135 
136 
138 {
139  pointField p(3);
140 
141  p[0] = points[operator[](0)];
142  p[1] = points[operator[](1)];
143  p[2] = points[operator[](2)];
144 
145  return p;
146 }
147 
148 
150 {
151  return Foam::face(*this);
152 }
153 
154 
156 {
157  return triPointRef
158  (
159  points[operator[](0)],
160  points[operator[](1)],
161  points[operator[](2)]
162  );
163 }
164 
165 
167 {
168  return (1.0/3.0)*
169  (
170  points[operator[](0)]
171  + points[operator[](1)]
172  + points[operator[](2)]
173  );
174 }
175 
176 
178 {
179  // As per triPointRef(...).areaNormal()
180  return 0.5*
181  (
182  (points[operator[](1)] - points[operator[](0)])
183  ^(points[operator[](2)] - points[operator[](0)])
184  );
185 }
186 
187 
189 {
190  const vector n(areaNormal(points));
191  const scalar s(Foam::mag(n));
192  return s < ROOTVSMALL ? Zero : n/s;
193 }
194 
195 
196 inline Foam::scalar Foam::triFace::mag(const UList<point>& points) const
197 {
198  return ::Foam::mag(areaNormal(points));
199 }
200 
201 
203 {
204  return 1;
205 }
206 
207 
209 {
210  // The starting points of the original and reverse face are identical.
211  return triFace(operator[](0), operator[](2), operator[](1));
212 }
213 
214 
215 inline bool Foam::triFace::found(const label pointLabel) const
216 {
217  return FixedList<label, 3>::found(pointLabel);
218 }
219 
220 
221 inline Foam::label Foam::triFace::which(const label pointLabel) const
222 {
223  return FixedList<label, 3>::find(pointLabel);
224 }
225 
226 
227 inline Foam::scalar Foam::triFace::sweptVol
228 (
229  const UList<point>& opts,
230  const UList<point>& npts
231 ) const
232 {
233  return (1.0/6.0)*
234  (
235  (
236  (npts[operator[](0)] - opts[operator[](0)])
237  & (
238  (opts[operator[](1)] - opts[operator[](0)])
239  ^ (opts[operator[](2)] - opts[operator[](0)])
240  )
241  )
242  + (
243  (npts[operator[](1)] - opts[operator[](1)])
244  & (
245  (opts[operator[](2)] - opts[operator[](1)])
246  ^ (npts[operator[](0)] - opts[operator[](1)])
247  )
248  )
249  + (
250  (opts[operator[](2)] - npts[operator[](2)])
251  & (
252  (npts[operator[](1)] - npts[operator[](2)])
253  ^ (npts[operator[](0)] - npts[operator[](2)])
254  )
255  )
256  );
257 }
258 
259 
261 (
262  const UList<point>& points,
263  const point& refPt,
264  scalar density
265 ) const
266 {
267  // a triangle, do a direct calculation
268  return this->tri(points).inertia(refPt, density);
269 }
270 
271 
273 (
274  const point& p,
275  const vector& q,
276  const UList<point>& points,
277  const intersection::algorithm alg,
278  const intersection::direction dir
279 ) const
280 {
281  return this->tri(points).ray(p, q, alg, dir);
282 }
283 
284 
285 
287 (
288  const point& p,
289  const vector& q,
290  const UList<point>& points,
291  const intersection::algorithm alg,
292  const scalar tol
293 ) const
294 {
295  return this->tri(points).intersection(p, q, alg, tol);
296 }
297 
298 
300 (
301  const point& p,
302  const vector& q,
303  const point& ctr,
304  const UList<point>& points,
305  const intersection::algorithm alg,
306  const scalar tol
307 ) const
308 {
309  return intersection(p, q, points, alg, tol);
310 }
311 
312 
314 (
315  const point& p,
316  const UList<point>& points
317 ) const
318 {
319  return this->tri(points).nearestPoint(p);
320 }
321 
322 
324 (
325  const point& p,
326  const UList<point>& points,
327  label& nearType,
328  label& nearLabel
329 ) const
330 {
331  return this->tri(points).nearestPointClassify(p, nearType, nearLabel);
332 }
333 
334 
335 inline int Foam::triFace::sign
336 (
337  const point& p,
338  const UList<point>& points,
339  const scalar tol
340 ) const
341 {
342  return this->tri(points).sign(p, tol);
343 }
344 
345 
347 {
348  return 3;
349 }
350 
351 
353 {
354  edgeList e(3);
355 
356  e[0].first() = operator[](0);
357  e[0].second() = operator[](1);
358 
359  e[1].first() = operator[](1);
360  e[1].second() = operator[](2);
361 
362  e[2].first() = operator[](2);
363  e[2].second() = operator[](0);
364 
365  return e;
366 }
367 
368 
370 {
371  return edge(operator[](n), operator[](fcIndex(n)));
372 }
373 
374 
375 // return
376 // - +1: forward (counter-clockwise) on the face
377 // - -1: reverse (clockwise) on the face
378 // - 0: edge not found on the face
379 inline int Foam::triFace::edgeDirection(const edge& e) const
380 {
381  if
382  (
383  (operator[](0) == e.first() && operator[](1) == e.second())
384  || (operator[](1) == e.first() && operator[](2) == e.second())
385  || (operator[](2) == e.first() && operator[](0) == e.second())
386  )
387  {
388  return 1;
389  }
390  else if
391  (
392  (operator[](0) == e.second() && operator[](1) == e.first())
393  || (operator[](1) == e.second() && operator[](2) == e.first())
394  || (operator[](2) == e.second() && operator[](0) == e.first())
395  )
396  {
397  return -1;
398  }
399 
400  return 0;
401 }
402 
403 
404 // * * * * * * * * * * * * * * * Global Operators * * * * * * * * * * * * * //
405 
406 inline bool Foam::operator==(const triFace& a, const triFace& b)
407 {
408  return triFace::compare(a,b) != 0;
409 }
410 
411 
412 inline bool Foam::operator!=(const triFace& a, const triFace& b)
413 {
414  return triFace::compare(a,b) == 0;
415 }
416 
417 
418 // ************************************************************************* //
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:261
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:208
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< point >
Foam::Zero
static constexpr const zero Zero
Global zero.
Definition: zero.H:128
Foam::triFace::areaNormal
vector areaNormal(const UList< point > &points) const
The area normal - with magnitude equal to area of face.
Definition: triFaceI.H:177
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:215
face.H
Foam::triFace::triFaceFace
face triFaceFace() const
Return triangle as a face.
Definition: triFaceI.H:149
Foam::triFace::unitNormal
vector unitNormal(const UList< point > &points) const
The unit normal.
Definition: triFaceI.H:188
triPointRef.H
Foam::Swap
void Swap(DynamicList< T, SizeMin1 > &a, DynamicList< T, SizeMin2 > &b)
Definition: DynamicListI.H:909
Foam::triFace::tri
triPointRef tri(const UList< point > &points) const
Return the triangle.
Definition: triFaceI.H:155
Foam::triFace::flip
void flip()
Flip the face in-place.
Definition: triFaceI.H:131
Foam::triFace::collapse
label collapse()
'Collapse' face by marking duplicate point labels.
Definition: triFaceI.H:104
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::FixedList::found
bool found(const T &val, const label start=0) const
True if the value if found in the list. Linear search.
Definition: FixedListI.H:305
Foam::triangle
A triangle primitive used to calculate face normals and swept volumes.
Definition: triangle.H:62
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:221
Foam::triFace::points
pointField points(const UList< point > &points) const
Return the points corresponding to this face.
Definition: triFaceI.H:137
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:62
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:336
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:352
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:369
Foam::FixedList::find
label find(const T &val, const label start=0) const
Find index of the first occurence of the value.
Definition: FixedList.C:36
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:228
Foam::triFace::centre
point centre(const UList< point > &points) const
Return centre (centroid)
Definition: triFaceI.H:166
Foam::triFace::nearestPoint
pointHit nearestPoint(const point &p, const UList< point > &points) const
Return nearest point to face.
Definition: triFaceI.H:314
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:324
Foam::triFace
A triangular face using a FixedList of labels corresponding to mesh vertices.
Definition: triFace.H:70
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:273
Foam::Vector< scalar >
Foam::triFace::mag
scalar mag(const UList< point > &points) const
Magnitude of face area.
Definition: triFaceI.H:196
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:287
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:379
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:74
Foam::triFace::triFace
triFace()
Construct null 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:202
Foam::triFace::nEdges
label nEdges() const
Return number of edges.
Definition: triFaceI.H:346
triFace
face triFace(3)
Foam::triPointRef
triangle< point, const point & > triPointRef
Definition: triangle.H:78