Loading...
Searching...
No Matches
face.C
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-2016 OpenFOAM Foundation
9 Copyright (C) 2021-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 "face.H"
30#include "triFace.H"
31#include "triangle.H"
33#include "Circulator.H"
34
35// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36
37const char* const Foam::face::typeName = "face";
38
39
40// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
41
43Foam::face::calcEdgeVectors(const UList<point>& points) const
44{
45 auto tedgeVecs = tmp<vectorField>::New(size());
46 auto& edgeVecs = tedgeVecs.ref();
47
48 forAll(edgeVecs, i)
49 {
50 edgeVecs[i] = vector(points[nextLabel(i)] - points[thisLabel(i)]);
51 edgeVecs[i].normalise();
52 }
53
54 return tedgeVecs;
55}
56
57
58Foam::label Foam::face::mostConcaveAngle
59(
60 const UList<point>& points,
61 const vectorField& edges,
62 scalar& maxAngle
63) const
64{
65 vector n(areaNormal(points));
66
67 label index = 0;
68 maxAngle = -GREAT;
69
70 forAll(edges, i)
71 {
72 const vector& leftEdge = edges[rcIndex(i)];
73 const vector& rightEdge = edges[i];
74
75 vector edgeNormal = (rightEdge ^ leftEdge);
76
77 // NOTE: is -ve angle since left edge pointing in other direction
78 scalar edgeCos = (leftEdge & rightEdge);
79 scalar edgeAngle = acos(clamp(edgeCos, -1, 1));
80
81 scalar angle;
82
83 if ((edgeNormal & n) > 0)
84 {
85 // Concave angle.
86 angle = constant::mathematical::pi + edgeAngle;
87 }
88 else
89 {
90 // Convex angle. Note '-' to take into account that rightEdge
91 // and leftEdge are head-to-tail connected.
92 angle = constant::mathematical::pi - edgeAngle;
93 }
94
95 if (maxAngle < angle)
96 {
97 maxAngle = angle;
98 index = i;
99 }
100 }
101
102 return index;
103}
104
105
106Foam::label Foam::face::split
107(
108 const face::splitMode mode,
109 const UList<point>& points,
110 label& triI,
111 label& quadI,
112 faceList& triFaces,
113 faceList& quadFaces
114) const
115{
116 const label oldIndices = (triI + quadI);
117
118 if (size() < 3)
119 {
121 << "Serious problem: asked to split a face with < 3 vertices"
122 << abort(FatalError);
123 }
124 else if (size() == 3)
125 {
126 // Triangle. Just copy.
127 if (mode == COUNTTRIANGLE || mode == COUNTQUAD)
128 {
129 triI++;
130 }
131 else
132 {
133 triFaces[triI++] = *this;
134 }
135 }
136 else if (size() == 4)
137 {
138 if (mode == COUNTTRIANGLE)
139 {
140 triI += 2;
141 }
142 if (mode == COUNTQUAD)
143 {
144 quadI++;
145 }
146 else if (mode == SPLITTRIANGLE)
147 {
148 // Start at point with largest internal angle.
149 const vectorField edges(calcEdgeVectors(points));
150
151 scalar minAngle;
152 label startIndex = mostConcaveAngle(points, edges, minAngle);
153
154 label nextIndex = fcIndex(startIndex);
155 label splitIndex = fcIndex(nextIndex);
156
157 // Create triangles
158 face triFace(3);
159 triFace[0] = operator[](startIndex);
160 triFace[1] = operator[](nextIndex);
161 triFace[2] = operator[](splitIndex);
162
163 triFaces[triI++] = triFace;
164
165 triFace[0] = operator[](splitIndex);
166 triFace[1] = operator[](fcIndex(splitIndex));
167 triFace[2] = operator[](startIndex);
168
169 triFaces[triI++] = triFace;
170 }
171 else
172 {
173 quadFaces[quadI++] = *this;
174 }
175 }
176 else
177 {
178 // General case. Like quad: search for largest internal angle.
179
180 const vectorField edges(calcEdgeVectors(points));
181
182 scalar minAngle = 1;
183 label startIndex = mostConcaveAngle(points, edges, minAngle);
184
185 scalar bisectAngle = minAngle/2;
186 const vector& rightEdge = edges[startIndex];
187
188 //
189 // Look for opposite point which as close as possible bisects angle
190 //
191
192 // split candidate starts two points away.
193 label index = fcIndex(fcIndex(startIndex));
194
195 label minIndex = index;
196 scalar minDiff = constant::mathematical::pi;
197
198 for (label i = 0; i < size() - 3; i++)
199 {
200 vector splitEdge
201 (
202 points[operator[](index)]
203 - points[operator[](startIndex)]
204 );
205 splitEdge.normalise();
206
207 const scalar splitCos = splitEdge & rightEdge;
208 const scalar splitAngle = acos(clamp(splitCos, -1, 1));
209 const scalar angleDiff = fabs(splitAngle - bisectAngle);
210
211 if (angleDiff < minDiff)
212 {
213 minDiff = angleDiff;
214 minIndex = index;
215 }
216
217 // Go to next candidate
218 index = fcIndex(index);
219 }
220
221
222 // Split into two subshapes.
223 // face1: startIndex to minIndex
224 // face2: minIndex to startIndex
225
226 // Get sizes of the two subshapes
227 label diff = 0;
228 if (minIndex > startIndex)
229 {
230 diff = minIndex - startIndex;
231 }
232 else
233 {
234 // Folded around
235 diff = minIndex + size() - startIndex;
236 }
237
238 label nPoints1 = diff + 1;
239 label nPoints2 = size() - diff + 1;
240
241 // Collect face1 points
242 face face1(nPoints1);
243
244 index = startIndex;
245 for (label i = 0; i < nPoints1; i++)
246 {
247 face1[i] = operator[](index);
248 index = fcIndex(index);
249 }
250
251 // Collect face2 points
252 face face2(nPoints2);
253
254 index = minIndex;
255 for (label i = 0; i < nPoints2; i++)
256 {
257 face2[i] = operator[](index);
258 index = fcIndex(index);
259 }
260
261 // Split faces
262 face1.split(mode, points, triI, quadI, triFaces, quadFaces);
263 face2.split(mode, points, triI, quadI, triFaces, quadFaces);
264 }
266 return (triI + quadI - oldIndices);
267}
268
269
270// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
271
274 labelList(f)
275{}
276
277
278// * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
279
280int Foam::face::compare(const face& a, const face& b)
281{
282 // Basic rule: we assume that the sequence of labels in each list
283 // will be circular in the same order (but not necessarily in the
284 // same direction or from the same starting point).
285
286 const label sizeA = a.size();
287 const label sizeB = b.size();
288
289 if (sizeA != sizeB)
290 {
291 // Trivial reject: faces have different sizes
292 return 0;
293 }
294 else if (sizeA == 0)
295 {
296 // Both faces with zero vertices. Always identical
297 return 1;
298 }
299 else if (sizeA == 1)
300 {
301 // Both faces with a single vertex. Simple check
302 return (a[0] == b[0] ? 1 : 0);
303 }
304
305 ConstCirculator<face> aCirc(a);
307
308 // Rotate face b until its element matches the starting element of face a.
309 do
310 {
311 if (aCirc() == bCirc())
312 {
313 // Set bCirc fulcrum to its iterator and increment the iterators
314 bCirc.setFulcrumToIterator();
315 ++aCirc;
316 ++bCirc;
317
318 break;
319 }
320 } while (bCirc.circulate(CirculatorBase::CLOCKWISE));
321
322 // If the circulator has stopped then faces a and b do not share a matching
323 // point. Doesn't work on matching, single element face.
324 if (!bCirc.circulate())
325 {
326 return 0;
327 }
328
329 // Look forwards around the faces for a match
330 do
331 {
332 if (*aCirc != *bCirc)
333 {
334 break;
335 }
336 }
337 while
338 (
339 aCirc.circulate(CirculatorBase::CLOCKWISE),
340 bCirc.circulate(CirculatorBase::CLOCKWISE)
341 );
342
343 // If the circulator has stopped then faces a and b matched.
344 if (!aCirc.circulate())
345 {
346 return 1;
347 }
348 else
349 {
350 // Reset the circulators back to their fulcrum
351 aCirc.setIteratorToFulcrum();
352 bCirc.setIteratorToFulcrum();
353 ++aCirc;
354 --bCirc;
355 }
356
357 // Look backwards around the faces for a match
358 do
359 {
360 if (*aCirc != *bCirc)
361 {
362 break;
363 }
364 }
365 while
366 (
367 aCirc.circulate(CirculatorBase::CLOCKWISE),
368 bCirc.circulate(CirculatorBase::ANTICLOCKWISE)
369 );
370
371 // If the circulator has stopped then faces a and b matched.
372 if (!aCirc.circulate())
373 {
374 return -1;
375 }
376
377 return 0;
378}
379
380
381bool Foam::face::sameVertices(const face& a, const face& b)
382{
383 const label sizeA = a.size();
384 const label sizeB = b.size();
385
386 // Trivial reject: faces are different size
387 if (sizeA != sizeB)
388 {
389 return false;
390 }
391 // Trivial: face with a single vertex
392 else if (sizeA == 1)
393 {
394 return (a[0] == b[0]);
395 }
396
397 forAll(a, i)
398 {
399 // Count occurrences of a[i] in a
400 label aOcc = 0;
401 forAll(a, j)
402 {
403 if (a[i] == a[j]) aOcc++;
404 }
405
406 // Count occurrences of a[i] in b
407 label bOcc = 0;
408 forAll(b, j)
409 {
410 if (a[i] == b[j]) bOcc++;
411 }
412
413 // Check if occurrences of a[i] in a and b are the same
414 if (aOcc != bOcc) return false;
415 }
416
417 return true;
418}
419
420
421unsigned Foam::face::symmhash_code(const UList<label>& f, unsigned seed)
422{
424
425 label len = f.size();
426
427 if (!len)
428 {
429 // Trivial: zero-sized
430 return 0;
431 }
432 else if (len == 1)
433 {
434 // Trivial: single vertex
435 return op(f[0], seed);
436 }
437
438 // Find location of the min vertex
439 label pivot = 0;
440 for (label i = 1; i < len; ++i)
441 {
442 if (f[pivot] > f[i])
443 {
444 pivot = i;
445 }
446 }
447
448 // Use next lowest value for deciding direction to circulate
449 if (f.fcValue(pivot) < f.rcValue(pivot))
450 {
451 // Forward circulate
452 while (len--)
453 {
454 seed = op(f[pivot], seed);
455 pivot = f.fcIndex(pivot);
456 }
457 }
458 else
459 {
460 // Reverse circulate
461 while (len--)
462 {
463 seed = op(f[pivot], seed);
464 pivot = f.rcIndex(pivot);
465 }
466 }
468 return seed;
469}
470
471
472// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
473
474Foam::label Foam::face::collapse()
475{
476 if (size() > 1)
477 {
478 label ci = 0;
479 for (label i=1; i<size(); i++)
480 {
481 if (operator[](i) != operator[](ci))
482 {
483 operator[](++ci) = operator[](i);
484 }
485 }
486
487 if (operator[](ci) != operator[](0))
488 {
489 ++ci;
490 }
491
493 }
494
495 return size();
496}
497
498
499void Foam::face::flip()
500{
501 const label n = size();
502
503 if (n > 2)
504 {
505 for (label i=1; i < (n+1)/2; ++i)
507 std::swap(operator[](i), operator[](n-i));
508 }
509 }
510}
511
512
514{
515 // Calculate the centre by breaking the face into triangles and
516 // area-weighted averaging their centres
517
518 const label nPoints = size();
519
520 // If the face is a triangle, do a direct calculation
521 if (nPoints == 3)
522 {
524 (
525 points[operator[](0)],
526 points[operator[](1)],
527 points[operator[](2)]
528 );
529 }
530 // else if (nPoints == 0)
531 // {
532 // // Should not happen
533 // return Zero;
534 // }
535
536
537 // Note: use double precision to avoid overflows when summing
538
539 // Estimated centre by averaging the face points
540 const point centrePoint(this->average(points));
541
542 doubleScalar sumA(0);
543 doubleVector sumAc(Zero);
544
545 for (label pI = 0; pI < nPoints; ++pI)
546 {
547 const label nextPti(pI == nPoints-1 ? 0 : pI+1);
548 const auto& thisPoint = points[operator[](pI)];
549 const auto& nextPoint = points[operator[](nextPti)];
550
551 // Calculate 3*triangle centre
552 const vector ttc(thisPoint + nextPoint + centrePoint);
553
554 // Calculate 2*triangle area
555 const scalar ta = Foam::mag
556 (
557 (thisPoint - centrePoint)
558 ^ (nextPoint - centrePoint)
559 );
560
561 sumA += ta;
562 sumAc += ta*ttc;
563 }
564
565 if (sumA > VSMALL)
566 {
567 return sumAc/(3.0*sumA);
568 }
569 else
570 {
571 return centrePoint;
572 }
573}
574
575
576Foam::vector Foam::face::areaNormal(const UList<point>& points) const
577{
578 const label nPoints = size();
579
580 // Calculate the area normal by summing the face triangle area normals.
581 // Changed to deal with small concavity by using a central decomposition
582
583 // If the face is a triangle, do a direct calculation
584 if (nPoints == 3)
585 {
587 (
588 points[operator[](0)],
589 points[operator[](1)],
590 points[operator[](2)]
591 );
592 }
593 // else if (nPoints == 0)
594 // {
595 // // Should not happen
596 // return Zero;
597 // }
598
599 // Estimated centre by averaging the face points
600 const point centrePoint(this->average(points));
601
602 vector n(Zero);
603
604 for (label pI = 0; pI < nPoints; ++pI)
605 {
606 const label nextPti(pI == nPoints-1 ? 0 : pI+1);
607 const auto& thisPoint = points[operator[](pI)];
608 const auto& nextPoint = points[operator[](nextPti)];
609
610 // Note: for best accuracy, centre point always comes last
611
612 n += triPointRef::areaNormal(thisPoint, nextPoint, centrePoint);
613 }
614
615 return n;
616}
617
618
620{
621 // Reverse the label list and return
622 // The starting points of the original and reverse face are identical.
623
624 const labelUList& origFace = *this;
625 const label len = origFace.size();
626
627 face newFace(len);
628 if (len)
629 {
630 newFace[0] = origFace[0];
631 for (label i=1; i < len; ++i)
632 {
633 newFace[i] = origFace[len - i];
635 }
636
637 return newFace;
638}
639
640
641Foam::scalar Foam::face::sweptVol
642(
643 const UList<point>& oldPoints,
644 const UList<point>& newPoints
645) const
646{
647 // This Optimization causes a small discrepancy between the swept-volume of
648 // opposite faces of complex cells with triangular faces opposing polygons.
649 // It could be used without problem for tetrahedral cells
650 // if (size() == 3)
651 // {
652 // return
653 // (
654 // triPointRef
655 // (
656 // oldPoints[operator[](0)],
657 // oldPoints[operator[](1)],
658 // oldPoints[operator[](2)]
659 // ).sweptVol
660 // (
661 // triPointRef
662 // (
663 // newPoints[operator[](0)],
664 // newPoints[operator[](1)],
665 // newPoints[operator[](2)]
666 // )
667 // )
668 // );
669 // }
670
671 scalar sv = 0;
672
673 // Calculate the swept volume by breaking the face into triangles and
674 // summing their swept volumes.
675 // Changed to deal with small concavity by using a central decomposition
676
677 const point centreOldPoint = centre(oldPoints);
678 const point centreNewPoint = centre(newPoints);
679
680 const label nPoints = size();
681
682 for (label pi=0; pi<nPoints-1; ++pi)
683 {
684 // Note: for best accuracy, centre point always comes last
685 sv += triPointRef
686 (
687 centreOldPoint,
688 oldPoints[operator[](pi)],
689 oldPoints[operator[](pi + 1)]
690 ).sweptVol
691 (
693 (
694 centreNewPoint,
695 newPoints[operator[](pi)],
696 newPoints[operator[](pi + 1)]
697 )
698 );
699 }
700
701 sv += triPointRef
702 (
703 centreOldPoint,
704 oldPoints[operator[](nPoints-1)],
705 oldPoints[operator[](0)]
706 ).sweptVol
707 (
709 (
710 centreNewPoint,
711 newPoints[operator[](nPoints-1)],
712 newPoints[operator[](0)]
714 );
715
716 return sv;
717}
718
719
721(
722 const UList<point>& p,
723 const point& refPt,
724 scalar density
725) const
726{
727 // If the face is a triangle, do a direct calculation
728 if (size() == 3)
729 {
730 return triPointRef
731 (
732 p[operator[](0)],
733 p[operator[](1)],
734 p[operator[](2)]
735 ).inertia(refPt, density);
736 }
737
738 const point ctr = centre(p);
739
740 tensor J = Zero;
741
742 forAll(*this, i)
743 {
744 J += triPointRef
745 (
746 p[operator[](i)],
747 p[operator[](fcIndex(i))],
748 ctr
749 ).inertia(refPt, density);
750 }
751
752 return J;
753}
754
755
757{
758 const labelList& verts = *this;
759 const label nVerts = verts.size();
760
761 edgeList theEdges(nVerts);
762
763 // Last edge closes the polygon
764 theEdges.back().first() = verts.back();
765 theEdges.back().second() = verts.front();
766
767 for (label verti = 0; verti < nVerts - 1; ++verti)
768 {
769 theEdges[verti].first() = verts[verti];
770 theEdges[verti].second() = verts[verti + 1];
771 }
772
773 return theEdges;
774}
775
776
778{
779 const labelList& verts = *this;
780 const label nVerts = verts.size();
781
782 edgeList theEdges(nVerts);
783
784 // First edge closes the polygon
785 theEdges.front().first() = verts.front();
786 theEdges.front().second() = verts.back();
787
788 for (label verti = 1; verti < nVerts; ++verti)
789 {
790 theEdges[verti].first() = verts[nVerts - verti];
791 theEdges[verti].second() = verts[nVerts - verti - 1];
792 }
793
794 return theEdges;
795}
796
797
798Foam::label Foam::face::find(const Foam::edge& e) const
799{
800 const label idx = labelList::find(e.first());
801
802 // NB: the point index is simultaneously the edge index.
803 // ie, face point 0 starts edge 0, point 1 starts edge 1, ...
804
805 if (idx != -1)
806 {
807 if (e.second() == nextLabel(idx))
808 {
809 // Current edge matches forward
810 return idx;
811 }
812 if (e.second() == prevLabel(idx))
813 {
814 // Previous edge matches reverse
815 return rcIndex(idx);
817 }
818
819 return -1; // Not found
820}
821
822
823int Foam::face::edgeDirection(const Foam::edge& e) const
824{
825 const label idx = labelList::find(e.first());
826
827 if (idx != -1)
828 {
829 if (e.second() == nextLabel(idx))
830 {
831 // Current edge matches forward. Could encode: (idx+1)
832 return 1;
833 }
834 if (e.second() == prevLabel(idx))
835 {
836 // Previous edge matches reverse. Could encode: -(rcIndex(idx)+1)
837 return -1;
839 }
840
841 return 0; // Not found
842}
843
845Foam::label Foam::face::nTriangles(const UList<point>&) const
846{
847 return nTriangles();
848}
849
850
851Foam::label Foam::face::triangles
852(
853 const UList<point>& points,
854 label& triI,
855 faceList& triFaces
856) const
857{
858 label quadI = 0;
859 faceList quadFaces;
860
861 return split(SPLITTRIANGLE, points, triI, quadI, triFaces, quadFaces);
862}
863
864
866(
867 const UList<point>& points,
868 label& triI,
869 label& quadI
870) const
871{
872 faceList triFaces;
873 faceList quadFaces;
874
875 return split(COUNTQUAD, points, triI, quadI, triFaces, quadFaces);
876}
877
878
880(
881 const UList<point>& points,
882 label& triI,
883 label& quadI,
884 faceList& triFaces,
885 faceList& quadFaces
886) const
887{
888 return split(SPLITQUAD, points, triI, quadI, triFaces, quadFaces);
889}
890
891
892Foam::label Foam::face::longestEdge(const UList<point>& pts) const
893{
894 const labelList& verts = *this;
895 const label nVerts = verts.size();
896
897 // Last edge closes the polygon. Use it to initialize loop
898 label longest = nVerts - 1;
899 scalar longestLen = pts[verts.front()].distSqr(pts[verts.back()]);
900
901 // Examine other edges
902 for (label edgei = 0; edgei < nVerts - 1; ++edgei)
903 {
904 scalar edgeLen = pts[verts[edgei]].distSqr(pts[verts[edgei+1]]);
905
906 if (longestLen < edgeLen)
907 {
908 longest = edgei;
909 longestLen = edgeLen;
910 }
911 }
912
913 return longest;
914}
915
916
917// ************************************************************************* //
constexpr scalar pi(M_PI)
label n
void setFulcrumToIterator()
Set the fulcrum to the current position of the iterator.
bool circulate(const CirculatorBase::direction dir=CirculatorBase::NONE)
Circulate around the list in the given direction.
void setIteratorToFulcrum()
Set the iterator to the current position of the fulcrum.
Like Foam::Circulator, but with const-access iterators.
Definition Circulator.H:413
const T & first() const noexcept
Access the first element.
Definition Pair.H:137
const T & second() const noexcept
Access the second element.
Definition Pair.H:147
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
const T & rcValue(const label i) const
Return reverse circular value (ie, previous value in the list).
Definition UListI.H:127
T & back()
Access last element of the list, position [size()-1].
Definition UListI.H:253
const T & fcValue(const label i) const
Return forward circular value (ie, next value in the list).
Definition UListI.H:113
label rcIndex(const label i) const noexcept
T & front()
Access first element of the list, position [0].
Definition UListI.H:239
label size() const noexcept
Definition UList.H:706
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
label & operator[](const label i)
label find(const label &val) const
label fcIndex(const label i) const noexcept
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
label trianglesQuads(const UList< point > &points, label &triI, label &quadI, faceList &triFaces, faceList &quadFaces) const
Split into triangles and quads.
Definition face.C:873
void flip()
Flip the face in-place.
Definition face.C:492
static unsigned symmhash_code(const UList< label > &f, unsigned seed=0)
The symmetric hash value for face.
Definition face.C:414
scalar sweptVol(const UList< point > &oldPoints, const UList< point > &newPoints) const
Return the volume swept out by the face when its points move.
Definition face.C:635
label triangles(const UList< point > &points, label &triI, faceList &triFaces) const
Split into triangles using existing points.
Definition face.C:845
label nTrianglesQuads(const UList< point > &points, label &nTris, label &nQuads) const
Number of triangles and quads after splitting.
Definition face.C:859
static int compare(const face &a, const face &b)
Compare faces.
Definition face.C:273
tensor inertia(const UList< point > &p, const point &refPt=vector::zero, scalar density=1.0) const
Return the inertia tensor, with optional reference point and density specification.
Definition face.C:714
label find(const Foam::edge &e) const
Find the edge within the face.
Definition face.C:791
int edgeDirection(const Foam::edge &e) const
Test the edge direction on the face.
Definition face.C:816
label longestEdge(const UList< point > &pts) const
Find the longest edge on a face.
Definition face.C:885
point centre(const UList< point > &points) const
Centre point of face.
Definition face.C:506
pointField points(const UList< point > &pts) const
Return the points corresponding to this face.
Definition faceI.H:80
label nextLabel(const label i) const
Next vertex on face.
Definition faceI.H:208
label thisLabel(const label i) const
The vertex on face - identical to operator[], but with naming similar to nextLabel(),...
Definition faceI.H:202
edgeList rcEdges() const
Return list of edges in reverse walk order.
Definition face.C:770
label collapse()
Collapse face by removing duplicate vertex labels.
Definition face.C:467
face reverseFace() const
Return face with reverse direction.
Definition face.C:612
static bool sameVertices(const face &a, const face &b)
True if the faces have all the same vertices.
Definition face.C:374
edgeList edges() const
Return list of edges in forward walk order.
Definition face.C:749
vector areaNormal(const UList< point > &p) const
The area normal - with magnitude equal to area of face.
Definition face.C:569
label nTriangles() const noexcept
Number of triangles after splitting.
Definition faceI.H:220
constexpr face() noexcept=default
Default construct.
label prevLabel(const label i) const
Previous vertex on face.
Definition faceI.H:214
static const char *const typeName
Definition face.H:137
Tensor of scalars, i.e. Tensor<scalar>.
A class for managing temporary objects.
Definition tmp.H:75
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition tmp.H:215
A triangular face using a FixedList of labels corresponding to mesh vertices.
Definition triFace.H:68
scalar sweptVol(const triangle &t) const
Return swept-volume.
Definition triangleI.H:393
tensor inertia(PointRef refPt=Zero, scalar density=1.0) const
Return the inertia tensor, with optional reference.
Definition triangleI.H:412
volScalarField & p
static bool split(const std::string &line, std::string &key, std::string &val)
Definition cpuInfo.C:32
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
const pointField & points
label nPoints
constexpr scalar pi(M_PI)
List< edge > edgeList
List of edge.
Definition edgeList.H:32
List< label > labelList
A List of labels.
Definition List.H:62
dimensionSet clamp(const dimensionSet &a, const dimensionSet &range)
mode_t mode(const fileName &name, const bool followLink=true)
Return the file mode, normally following symbolic links.
Definition POSIX.C:775
List< face > faceList
List of faces.
Definition faceListFwd.H:41
Tensor< scalar > tensor
Definition symmTensor.H:57
double doubleScalar
A typedef for double.
Definition scalarFwd.H:48
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
scalar diff(const triad &A, const triad &B)
Return a quantity of the difference between two triads.
Definition triad.C:373
Vector< double > doubleVector
Definition vector.H:54
static scalar angleDiff(scalar angle1, scalar angle2)
errorManip< error > abort(error &err)
Definition errorManip.H:139
Field< vector > vectorField
Specialisation of Field<T> for vector.
vector point
Point is a vector.
Definition point.H:37
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
triangle< point, const point & > triPointRef
A triangle using referred points.
Definition triangleFwd.H:46
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
Vector< scalar > vector
Definition vector.H:57
UList< label > labelUList
A UList of labels.
Definition UList.H:75
dimensioned< Type > average(const DimensionedField< Type, GeoMesh > &f1, const label comm)
dimensionedScalar acos(const dimensionedScalar &ds)
points setSize(newPointi)
face triFace(3)
labelList f(nPoints)
const pointField & pts
volScalarField & b
volScalarField & e
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
Hash function class. The default definition is for primitives. Non-primitives used to hash entries on...
Definition Hash.H:48