Loading...
Searching...
No Matches
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-2025 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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
32// * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
33
34inline int Foam::triFace::compare(const triFace& a, const triFace& b)
35{
36 if
37 (
38 (a[0] == b[0] && a[1] == b[1] && a[2] == b[2])
39 || (a[0] == b[1] && a[1] == b[2] && a[2] == b[0])
40 || (a[0] == b[2] && a[1] == b[0] && a[2] == b[1])
41 )
42 {
43 // identical
44 return 1;
45 }
46 else if
47 (
48 (a[0] == b[2] && a[1] == b[1] && a[2] == b[0])
49 || (a[0] == b[1] && a[1] == b[0] && a[2] == b[2])
50 || (a[0] == b[0] && a[1] == b[2] && a[2] == b[1])
51 )
52 {
53 // same face, but reversed orientation
54 return -1;
55 }
57 return 0;
58}
59
60
61// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
64:
65 FixedList<label, 3>(-1)
66{}
67
68
70(
71 const label p0,
72 const label p1,
73 const label p2
74) noexcept
76 a() = p0;
77 b() = p1;
78 c() = p2;
79}
80
82inline Foam::triFace::triFace(std::initializer_list<label> list)
83:
84 FixedList<label, 3>(list)
85{}
86
88inline Foam::triFace::triFace(const labelUList& list)
89:
90 FixedList<label, 3>(list)
91{}
92
93
95(
96 const labelUList& list,
97 const FixedList<label, 3>& indices
98)
99:
100 FixedList<label, 3>(list, indices)
101{}
102
103
106 FixedList<label, 3>(is)
107{}
108
109
110// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
111
112inline bool Foam::triFace::good() const noexcept
113{
114 return
115 (
116 a() >= 0 && a() != b()
117 && b() >= 0 && b() != c()
118 && c() >= 0 && c() != a()
119 );
120}
121
122
123inline Foam::label Foam::triFace::collapse()
124{
125 // Cannot resize FixedList, so mark duplicates with '-1'
126 // (the lower vertex is retained)
127 // catch any '-1' (eg, if called multiple times)
128
129 label n = 3;
130 if (operator[](0) == operator[](1) || operator[](1) == -1)
131 {
132 operator[](1) = -1;
133 n--;
134 }
135 else if (operator[](1) == operator[](2) || operator[](2) == -1)
136 {
137 operator[](2) = -1;
138 n--;
139 }
140 if (operator[](0) == operator[](2))
141 {
142 operator[](2) = -1;
143 n--;
144 }
145
146 return n;
147}
148
150inline void Foam::triFace::flip()
151{
152 std::swap(get<1>(), get<2>());
153}
154
155
157{
158 pointField p(3);
159
160 p[0] = pts[a()];
161 p[1] = pts[b()];
162 p[2] = pts[c()];
163
164 return p;
165}
166
169{
170 return Foam::face(*this);
171}
172
175{
176 return triPointRef(points[a()], points[b()], points[c()]);
177}
178
181{
182 return triPointRef::centre(points[a()], points[b()], points[c()]);
183}
184
189}
190
195}
196
198inline Foam::scalar Foam::triFace::mag(const UList<point>& points) const
199{
200 return ::Foam::mag(areaNormal(points));
201}
202
203
204inline Foam::scalar Foam::triFace::magSqr(const UList<point>& points) const
205{
206 return ::Foam::magSqr(areaNormal(points));
207}
208
209
214}
215
217inline Foam::label Foam::triFace::nTriangles() const noexcept
218{
219 return 1;
220}
221
222
224{
225 // The starting points of the original and reverse face are identical.
226 return triFace(a(), c(), b());
227}
228
230inline Foam::label Foam::triFace::which(const label vertex) const
231{
232 return FixedList<label, 3>::find(vertex);
233}
234
236inline Foam::label Foam::triFace::thisLabel(const label i) const
237{
238 return operator[](i);
239}
240
242inline Foam::label Foam::triFace::nextLabel(const label i) const
243{
244 return operator[]((i == 2 ? 0 : i+1));
245}
246
248inline Foam::label Foam::triFace::prevLabel(const label i) const
249{
250 return operator[]((i ? i-1 : 2));
251}
252
253
254inline Foam::scalar Foam::triFace::sweptVol
255(
256 const UList<point>& opts,
257 const UList<point>& npts
258) const
259{
260 return (1.0/6.0)*
261 (
262 (
263 (npts[operator[](0)] - opts[operator[](0)])
264 & (
265 (opts[operator[](1)] - opts[operator[](0)])
266 ^ (opts[operator[](2)] - opts[operator[](0)])
267 )
268 )
269 + (
270 (npts[operator[](1)] - opts[operator[](1)])
271 & (
272 (opts[operator[](2)] - opts[operator[](1)])
273 ^ (npts[operator[](0)] - opts[operator[](1)])
274 )
275 )
276 + (
277 (opts[operator[](2)] - npts[operator[](2)])
278 & (
279 (npts[operator[](1)] - npts[operator[](2)])
280 ^ (npts[operator[](0)] - npts[operator[](2)])
281 )
282 )
283 );
284}
285
286
288(
289 const UList<point>& points,
290 const point& refPt,
291 scalar density
292) const
293{
294 // a triangle, do a direct calculation
295 return this->tri(points).inertia(refPt, density);
296}
297
298
300(
301 const point& p,
302 const vector& q,
303 const UList<point>& points,
304 const intersection::algorithm alg,
306) const
307{
308 return this->tri(points).ray(p, q, alg, dir);
309}
310
311
312
314(
315 const point& p,
316 const vector& q,
317 const UList<point>& points,
318 const intersection::algorithm alg,
319 const scalar tol
320) const
321{
322 return this->tri(points).intersection(p, q, alg, tol);
323}
324
325
327(
328 const point& p,
329 const vector& q,
330 const point& ctr,
331 const UList<point>& points,
332 const intersection::algorithm alg,
333 const scalar tol
334) const
335{
336 return intersection(p, q, points, alg, tol);
337}
338
339
341(
342 const point& p,
344) const
345{
346 return this->tri(points).nearestPoint(p);
347}
348
349
351(
352 const point& p,
353 const UList<point>& points,
354 label& nearType,
355 label& nearLabel
356) const
357{
358 return this->tri(points).nearestPointClassify(p, nearType, nearLabel);
359}
360
361
362inline int Foam::triFace::sign
363(
364 const point& p,
365 const UList<point>& points,
366 const scalar tol
367) const
368{
369 return this->tri(points).sign(p, tol);
370}
371
373inline Foam::label Foam::triFace::nEdges() const noexcept
374{
375 return 3;
376}
377
379inline Foam::edge Foam::triFace::edge(const label edgei) const
380{
381 return Foam::edge(thisLabel(edgei), nextLabel(edgei));
382}
383
384
386(
387 const label edgei,
389) const
390{
391 return vector(pts[nextLabel(edgei)] - pts[thisLabel(edgei)]);
392}
393
394
395inline Foam::edge Foam::triFace::rcEdge(const label edgei) const
396{
397 // Edge 0 (forward and reverse) always starts at [0]
398 // for consistency with face flipping
399 const label pointi = edgei ? (3 - edgei) : 0;
400 return Foam::edge(thisLabel(pointi), prevLabel(pointi));
401}
402
403
405(
406 const label edgei,
407 const UList<point>& pts
408) const
409{
410 // Edge 0 (forward and reverse) always starts at [0]
411 // for consistency with face flipping
412 const label pointi = edgei ? (3 - edgei) : 0;
413 return vector(pts[prevLabel(pointi)] - pts[thisLabel(pointi)]);
414}
415
416
418{
419 edgeList theEdges(3);
420
421 theEdges[0].first() = a();
422 theEdges[0].second() = b();
423
424 theEdges[1].first() = b();
425 theEdges[1].second() = c();
426
427 theEdges[2].first() = c();
428 theEdges[2].second() = a();
429
430 return theEdges;
431}
432
433
435{
436 edgeList theEdges(3);
437
438 theEdges[0].first() = a();
439 theEdges[0].second() = c();
440
441 theEdges[1].first() = c();
442 theEdges[1].second() = b();
443
444 theEdges[2].first() = b();
445 theEdges[2].second() = a();
446
447 return theEdges;
448}
449
450
451inline int Foam::triFace::edgeDirection(const Foam::edge& e) const
452{
453 if (e.first() == a())
454 {
455 if (e.second() == b()) return 1; // Forward edge 0 (encoded 1)
456 if (e.second() == c()) return -1; // Reverse edge 2 (encoded -3)
457 }
458 if (e.first() == b())
459 {
460 if (e.second() == c()) return 1; // Forward edge 1 (encoded 2)
461 if (e.second() == a()) return -1; // Reverse edge 0 (encoded -1)
462 }
463 if (e.first() == c())
464 {
465 if (e.second() == a()) return 1; // Forward edge 2 (encoded 3)
466 if (e.second() == b()) return -1; // Reverse edge 1 (encoded -2)
467 }
468
469 return 0; // Not found
470}
471
472
473inline Foam::label Foam::triFace::find(const Foam::edge& e) const
474{
475 if (e.first() == a())
476 {
477 if (e.second() == b()) return 0; // Forward edge 0
478 if (e.second() == c()) return 2; // Reverse edge 2
479 }
480 if (e.first() == b())
481 {
482 if (e.second() == c()) return 1; // Forward edge 1
483 if (e.second() == a()) return 0; // Reverse edge 0
484 }
485 if (e.first() == c())
486 {
487 if (e.second() == a()) return 2; // Forward edge 2
488 if (e.second() == b()) return 1; // Reverse edge 1
489 }
490
491 return -1; // Not found
492}
493
494
495inline bool Foam::triFace::contains(const Foam::edge& e) const
496{
497 // or (find(e) >= 0)
498 return (edgeDirection(e) != 0);
499}
500
501
502// * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
503
504inline void Foam::triFace::operator+=(const label vertexOffset)
505{
506 if (vertexOffset)
507 {
508 (*this)[0] += vertexOffset;
509 (*this)[1] += vertexOffset;
510 (*this)[2] += vertexOffset;
511 }
512}
513
514
515// * * * * * * * * * * * * * * * Global Operators * * * * * * * * * * * * * //
517inline bool Foam::operator==(const triFace& a, const triFace& b)
518{
519 return triFace::compare(a,b) != 0;
520}
521
522
523inline bool Foam::operator!=(const triFace& a, const triFace& b)
524{
525 return triFace::compare(a,b) == 0;
526}
527
528
529// ************************************************************************* //
Useful combination of include files which define Sin, Sout and Serr and the use of IO streams general...
label n
label & get() noexcept
Definition FixedListI.H:152
label & operator[](const label i)
Definition FixedListI.H:399
label find(const T &val) const
Find index of the first occurrence of the value.
Definition FixedList.C:42
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition Istream.H:60
An ordered pair of two objects of type <T> with first() and second() elements.
Definition Pair.H:66
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition UList.H:89
T & first()
Access first element of the list, position [0].
Definition UList.H:957
An edge is a list of two vertex labels. This can correspond to a directed graph edge or an edge on a ...
Definition edge.H:62
A face is a list of labels corresponding to mesh vertices.
Definition face.H:71
A triangular face using a FixedList of labels corresponding to mesh vertices.
Definition triFace.H:68
scalar sweptVol(const UList< point > &opts, const UList< point > &npts) const
Return swept-volume from old-points to new-points.
Definition triFaceI.H:248
vector unitNormal(const UList< point > &points) const
The unit normal.
Definition triFaceI.H:185
Foam::edge edge(const label edgei) const
Return i-th face edge (forward walk order).
Definition triFaceI.H:372
void flip()
Flip the face in-place.
Definition triFaceI.H:143
label b() const noexcept
The second vertex.
Definition triFace.H:133
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:293
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:344
pointHit nearestPoint(const point &p, const UList< point > &points) const
Return nearest point to face.
Definition triFaceI.H:334
scalar magSqr(const UList< point > &points) const
Magnitude squared of face area.
Definition triFaceI.H:197
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:307
bool contains(const Foam::edge &e) const
True if face contains(edge).
Definition triFaceI.H:488
label nEdges() const noexcept
Return number of edges == 3.
Definition triFaceI.H:366
triFace reverseFace() const
Return face with reverse direction.
Definition triFaceI.H:216
scalar mag(const UList< point > &points) const
Magnitude of face area.
Definition triFaceI.H:191
bool good() const noexcept
True if vertices are unique and non-negative.
Definition triFaceI.H:105
label find(const Foam::edge &e) const
Find the edge within the face.
Definition triFaceI.H:466
Pair< point > box(const UList< point > &points) const
The enclosing (bounding) box for the face.
Definition triFaceI.H:204
int edgeDirection(const Foam::edge &e) const
Test the edge direction on the face.
Definition triFaceI.H:444
void operator+=(const label vertexOffset)
Increment (offset) vertices by given amount.
Definition triFaceI.H:497
Foam::edge rcEdge(const label edgei) const
Return i-th face edge in reverse walk order.
Definition triFaceI.H:388
triFace()
Default construct, with invalid vertex labels (-1).
Definition triFaceI.H:56
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:281
label c() const noexcept
The third vertex.
Definition triFace.H:138
point centre(const UList< point > &points) const
Return centre (centroid).
Definition triFaceI.H:173
pointField points(const UList< point > &pts) const
Return the points corresponding to this face.
Definition triFaceI.H:149
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:356
static int compare(const triFace &a, const triFace &b)
Compare triFaces.
Definition triFaceI.H:27
label which(const label vertex) const
Find local vertex on face for the vertex label, same as find().
Definition triFaceI.H:223
label a() const noexcept
The first vertex.
Definition triFace.H:128
label nextLabel(const label i) const
Next vertex on face.
Definition triFaceI.H:235
label thisLabel(const label i) const
The vertex on face - identical to operator[], but with naming similar to nextLabel(),...
Definition triFaceI.H:229
edgeList rcEdges() const
Return list of edges in reverse walk order.
Definition triFaceI.H:427
label collapse()
'Collapse' face by marking duplicate vertex labels.
Definition triFaceI.H:116
face triFaceFace() const
Return triangle as a face.
Definition triFaceI.H:161
edgeList edges() const
Return list of edges in forward walk order.
Definition triFaceI.H:410
vector areaNormal(const UList< point > &points) const
The area normal - with magnitude equal to area of face.
Definition triFaceI.H:179
label nTriangles() const noexcept
Number of triangles after splitting == 1.
Definition triFaceI.H:210
label prevLabel(const label i) const
Previous vertex on face.
Definition triFaceI.H:241
triPointRef tri(const UList< point > &points) const
Return the triangle.
Definition triFaceI.H:167
volScalarField & p
const volScalarField & p0
Definition EEqn.H:36
const pointField & points
const dimensionedScalar c
Speed of light in a vacuum.
List< edge > edgeList
List of edge.
Definition edgeList.H:32
bool operator!=(const eddy &a, const eddy &b)
Definition eddy.H:297
PointHit< point > pointHit
A PointHit with a 3D point.
Definition pointHit.H:51
Tensor< scalar > tensor
Definition symmTensor.H:57
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
vector point
Point is a vector.
Definition point.H:37
const direction noexcept
Definition scalarImpl.H:265
triangle< point, const point & > triPointRef
A triangle using referred points.
Definition triangleFwd.H:46
vectorField pointField
pointField is a vectorField.
Vector< scalar > vector
Definition vector.H:57
UList< label > labelUList
A UList of labels.
Definition UList.H:75
face triFace(3)
const pointField & pts
volScalarField & b
volScalarField & e