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