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 | Copyright (C) 2011-2015 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "face.H"
27 #include "triFace.H"
28 #include "triPointRef.H"
29 #include "mathematicalConstants.H"
30 #include "Swap.H"
31 #include "ConstCirculator.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 const char* const Foam::face::typeName = "face";
36 
37 
38 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
39 
41 Foam::face::calcEdges(const pointField& points) const
42 {
43  tmp<vectorField> tedges(new vectorField(size()));
44  vectorField& edges = tedges();
45 
46  forAll(*this, i)
47  {
48  label ni = fcIndex(i);
49 
50  point thisPt = points[operator[](i)];
51  point nextPt = points[operator[](ni)];
52 
53  vector vec(nextPt - thisPt);
54  vec /= Foam::mag(vec) + VSMALL;
55 
56  edges[i] = vec;
57  }
58 
59  return tedges;
60 }
61 
62 
63 Foam::scalar Foam::face::edgeCos
64 (
65  const vectorField& edges,
66  const label index
67 ) const
68 {
69  label leftEdgeI = left(index);
70  label rightEdgeI = right(index);
71 
72  // Note negate on left edge to get correct left-pointing edge.
73  return -(edges[leftEdgeI] & edges[rightEdgeI]);
74 }
75 
76 
77 Foam::label Foam::face::mostConcaveAngle
78 (
79  const pointField& points,
80  const vectorField& edges,
81  scalar& maxAngle
82 ) const
83 {
84  vector n(normal(points));
85 
86  label index = 0;
87  maxAngle = -GREAT;
88 
89  forAll(edges, i)
90  {
91  label leftEdgeI = left(i);
92  label rightEdgeI = right(i);
93 
94  vector edgeNormal = edges[rightEdgeI] ^ edges[leftEdgeI];
95 
96  scalar edgeCos = edges[leftEdgeI] & edges[rightEdgeI];
97  scalar edgeAngle = acos(max(-1.0, min(1.0, edgeCos)));
98 
99  scalar angle;
100 
101  if ((edgeNormal & n) > 0)
102  {
103  // Concave angle.
104  angle = constant::mathematical::pi + edgeAngle;
105  }
106  else
107  {
108  // Convex angle. Note '-' to take into account that rightEdge
109  // and leftEdge are head-to-tail connected.
110  angle = constant::mathematical::pi - edgeAngle;
111  }
112 
113  if (angle > maxAngle)
114  {
115  maxAngle = angle;
116  index = i;
117  }
118  }
119 
120  return index;
121 }
122 
123 
124 Foam::label Foam::face::split
125 (
126  const face::splitMode mode,
127  const pointField& points,
128  label& triI,
129  label& quadI,
130  faceList& triFaces,
131  faceList& quadFaces
132 ) const
133 {
134  label oldIndices = (triI + quadI);
135 
136  if (size() <= 2)
137  {
139  (
140  "face::split"
141  "(const face::splitMode, const pointField&, label&, label&"
142  ", faceList&, faceList&)"
143  )
144  << "Serious problem: asked to split a face with < 3 vertices"
145  << abort(FatalError);
146  }
147  if (size() == 3)
148  {
149  // Triangle. Just copy.
150  if (mode == COUNTTRIANGLE || mode == COUNTQUAD)
151  {
152  triI++;
153  }
154  else
155  {
156  triFaces[triI++] = *this;
157  }
158  }
159  else if (size() == 4)
160  {
161  if (mode == COUNTTRIANGLE)
162  {
163  triI += 2;
164  }
165  if (mode == COUNTQUAD)
166  {
167  quadI++;
168  }
169  else if (mode == SPLITTRIANGLE)
170  {
171  // Start at point with largest internal angle.
172  const vectorField edges(calcEdges(points));
173 
174  scalar minAngle;
175  label startIndex = mostConcaveAngle(points, edges, minAngle);
176 
177  label nextIndex = fcIndex(startIndex);
178  label splitIndex = fcIndex(nextIndex);
179 
180  // Create triangles
181  face triFace(3);
182  triFace[0] = operator[](startIndex);
183  triFace[1] = operator[](nextIndex);
184  triFace[2] = operator[](splitIndex);
185 
186  triFaces[triI++] = triFace;
187 
188  triFace[0] = operator[](splitIndex);
189  triFace[1] = operator[](fcIndex(splitIndex));
190  triFace[2] = operator[](startIndex);
191 
192  triFaces[triI++] = triFace;
193  }
194  else
195  {
196  quadFaces[quadI++] = *this;
197  }
198  }
199  else
200  {
201  // General case. Like quad: search for largest internal angle.
202 
203  const vectorField edges(calcEdges(points));
204 
205  scalar minAngle = 1;
206  label startIndex = mostConcaveAngle(points, edges, minAngle);
207 
208  scalar bisectAngle = minAngle/2;
209  vector rightEdge = edges[right(startIndex)];
210 
211  //
212  // Look for opposite point which as close as possible bisects angle
213  //
214 
215  // split candidate starts two points away.
216  label index = fcIndex(fcIndex(startIndex));
217 
218  label minIndex = index;
219  scalar minDiff = constant::mathematical::pi;
220 
221  for (label i = 0; i < size() - 3; i++)
222  {
223  vector splitEdge
224  (
225  points[operator[](index)]
226  - points[operator[](startIndex)]
227  );
228  splitEdge /= Foam::mag(splitEdge) + VSMALL;
229 
230  const scalar splitCos = splitEdge & rightEdge;
231  const scalar splitAngle = acos(max(-1.0, min(1.0, splitCos)));
232  const scalar angleDiff = fabs(splitAngle - bisectAngle);
233 
234  if (angleDiff < minDiff)
235  {
236  minDiff = angleDiff;
237  minIndex = index;
238  }
239 
240  // Go to next candidate
241  index = fcIndex(index);
242  }
243 
244 
245  // Split into two subshapes.
246  // face1: startIndex to minIndex
247  // face2: minIndex to startIndex
248 
249  // Get sizes of the two subshapes
250  label diff = 0;
251  if (minIndex > startIndex)
252  {
253  diff = minIndex - startIndex;
254  }
255  else
256  {
257  // Folded around
258  diff = minIndex + size() - startIndex;
259  }
260 
261  label nPoints1 = diff + 1;
262  label nPoints2 = size() - diff + 1;
263 
264  // Collect face1 points
265  face face1(nPoints1);
266 
267  index = startIndex;
268  for (label i = 0; i < nPoints1; i++)
269  {
270  face1[i] = operator[](index);
271  index = fcIndex(index);
272  }
273 
274  // Collect face2 points
275  face face2(nPoints2);
276 
277  index = minIndex;
278  for (label i = 0; i < nPoints2; i++)
279  {
280  face2[i] = operator[](index);
281  index = fcIndex(index);
282  }
283 
284  // Split faces
285  face1.split(mode, points, triI, quadI, triFaces, quadFaces);
286  face2.split(mode, points, triI, quadI, triFaces, quadFaces);
287  }
288 
289  return (triI + quadI - oldIndices);
290 }
291 
292 
293 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
294 
296 :
297  labelList(f)
298 {}
299 
300 
301 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
302 
303 int Foam::face::compare(const face& a, const face& b)
304 {
305  // Basic rule: we assume that the sequence of labels in each list
306  // will be circular in the same order (but not necessarily in the
307  // same direction or from the same starting point).
308 
309  // Trivial reject: faces are different size
310  label sizeA = a.size();
311  label sizeB = b.size();
312 
313  if (sizeA != sizeB || sizeA == 0)
314  {
315  return 0;
316  }
317  else if (sizeA == 1)
318  {
319  if (a[0] == b[0])
320  {
321  return 1;
322  }
323  else
324  {
325  return 0;
326  }
327  }
328 
329  ConstCirculator<face> aCirc(a);
330  ConstCirculator<face> bCirc(b);
331 
332  // Rotate face b until its element matches the starting element of face a.
333  do
334  {
335  if (aCirc() == bCirc())
336  {
337  // Set bCirc fulcrum to its iterator and increment the iterators
338  bCirc.setFulcrumToIterator();
339  ++aCirc;
340  ++bCirc;
341 
342  break;
343  }
344  } while (bCirc.circulate(CirculatorBase::CLOCKWISE));
345 
346  // If the circulator has stopped then faces a and b do not share a matching
347  // point. Doesn't work on matching, single element face.
348  if (!bCirc.circulate())
349  {
350  return 0;
351  }
352 
353  // Look forwards around the faces for a match
354  do
355  {
356  if (aCirc() != bCirc())
357  {
358  break;
359  }
360  }
361  while
362  (
365  );
366 
367  // If the circulator has stopped then faces a and b matched.
368  if (!aCirc.circulate())
369  {
370  return 1;
371  }
372  else
373  {
374  // Reset the circulators back to their fulcrum
375  aCirc.setIteratorToFulcrum();
376  bCirc.setIteratorToFulcrum();
377  ++aCirc;
378  --bCirc;
379  }
380 
381  // Look backwards around the faces for a match
382  do
383  {
384  if (aCirc() != bCirc())
385  {
386  break;
387  }
388  }
389  while
390  (
393  );
394 
395  // If the circulator has stopped then faces a and b matched.
396  if (!aCirc.circulate())
397  {
398  return -1;
399  }
400 
401  return 0;
402 }
403 
404 
405 bool Foam::face::sameVertices(const face& a, const face& b)
406 {
407  label sizeA = a.size();
408  label sizeB = b.size();
409 
410  // Trivial reject: faces are different size
411  if (sizeA != sizeB)
412  {
413  return false;
414  }
415  // Check faces with a single vertex
416  else if (sizeA == 1)
417  {
418  if (a[0] == b[0])
419  {
420  return true;
421  }
422  else
423  {
424  return false;
425  }
426  }
427 
428  forAll(a, i)
429  {
430  // Count occurrences of a[i] in a
431  label aOcc = 0;
432  forAll(a, j)
433  {
434  if (a[i] == a[j]) aOcc++;
435  }
436 
437  // Count occurrences of a[i] in b
438  label bOcc = 0;
439  forAll(b, j)
440  {
441  if (a[i] == b[j]) bOcc++;
442  }
443 
444  // Check if occurrences of a[i] in a and b are the same
445  if (aOcc != bOcc) return false;
446  }
447 
448  return true;
449 }
450 
451 
452 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
453 
455 {
456  if (size() > 1)
457  {
458  label ci = 0;
459  for (label i=1; i<size(); i++)
460  {
461  if (operator[](i) != operator[](ci))
462  {
463  operator[](++ci) = operator[](i);
464  }
465  }
466 
467  if (operator[](ci) != operator[](0))
468  {
469  ci++;
470  }
471 
472  setSize(ci);
473  }
474 
475  return size();
476 }
477 
478 
480 {
481  const label n = size();
482 
483  if (n > 2)
484  {
485  for (label i=1; i < (n+1)/2; ++i)
486  {
487  Swap(operator[](i), operator[](n-i));
488  }
489  }
490 }
491 
492 
494 {
495  // Calculate the centre by breaking the face into triangles and
496  // area-weighted averaging their centres
497 
498  const label nPoints = size();
499 
500  // If the face is a triangle, do a direct calculation
501  if (nPoints == 3)
502  {
503  return
504  (1.0/3.0)
505  *(
506  points[operator[](0)]
507  + points[operator[](1)]
508  + points[operator[](2)]
509  );
510  }
511 
512 
513  point centrePoint = point::zero;
514  for (label pI=0; pI<nPoints; ++pI)
515  {
516  centrePoint += points[operator[](pI)];
517  }
518  centrePoint /= nPoints;
519 
520  scalar sumA = 0;
521  vector sumAc = vector::zero;
522 
523  for (label pI=0; pI<nPoints; ++pI)
524  {
525  const point& nextPoint = points[operator[]((pI + 1) % nPoints)];
526 
527  // Calculate 3*triangle centre
528  const vector ttc
529  (
530  points[operator[](pI)]
531  + nextPoint
532  + centrePoint
533  );
534 
535  // Calculate 2*triangle area
536  const scalar ta = Foam::mag
537  (
538  (points[operator[](pI)] - centrePoint)
539  ^ (nextPoint - centrePoint)
540  );
541 
542  sumA += ta;
543  sumAc += ta*ttc;
544  }
545 
546  if (sumA > VSMALL)
547  {
548  return sumAc/(3.0*sumA);
549  }
550  else
551  {
552  return centrePoint;
553  }
554 }
555 
556 
558 {
559  const label nPoints = size();
560 
561  // Calculate the normal by summing the face triangle normals.
562  // Changed to deal with small concavity by using a central decomposition
563  //
564 
565  // If the face is a triangle, do a direct calculation to avoid round-off
566  // error-related problems
567  //
568  if (nPoints == 3)
569  {
570  return triPointRef
571  (
572  p[operator[](0)],
573  p[operator[](1)],
574  p[operator[](2)]
575  ).normal();
576  }
577 
578  label pI;
579 
580  point centrePoint = vector::zero;
581  for (pI = 0; pI < nPoints; ++pI)
582  {
583  centrePoint += p[operator[](pI)];
584  }
585  centrePoint /= nPoints;
586 
588 
589  point nextPoint = centrePoint;
590 
591  for (pI = 0; pI < nPoints; ++pI)
592  {
593  if (pI < nPoints - 1)
594  {
595  nextPoint = p[operator[](pI + 1)];
596  }
597  else
598  {
599  nextPoint = p[operator[](0)];
600  }
601 
602  // Note: for best accuracy, centre point always comes last
603  //
604  n += triPointRef
605  (
606  p[operator[](pI)],
607  nextPoint,
608  centrePoint
609  ).normal();
610  }
611 
612  return n;
613 }
614 
615 
617 {
618  // Reverse the label list and return
619  // The starting points of the original and reverse face are identical.
620 
621  const labelList& f = *this;
622  labelList newList(size());
623 
624  newList[0] = f[0];
625 
626  for (label pointI = 1; pointI < newList.size(); pointI++)
627  {
628  newList[pointI] = f[size() - pointI];
629  }
630 
631  return face(xferMove(newList));
632 }
633 
634 
636 {
637  const labelList& f = *this;
638 
639  forAll(f, localIdx)
640  {
641  if (f[localIdx] == globalIndex)
642  {
643  return localIdx;
644  }
645  }
646 
647  return -1;
648 }
649 
650 
651 Foam::scalar Foam::face::sweptVol
652 (
653  const pointField& oldPoints,
654  const pointField& newPoints
655 ) const
656 {
657  // This Optimization causes a small discrepancy between the swept-volume of
658  // opposite faces of complex cells with triangular faces opposing polygons.
659  // It could be used without problem for tetrahedral cells
660  // if (size() == 3)
661  // {
662  // return
663  // (
664  // triPointRef
665  // (
666  // oldPoints[operator[](0)],
667  // oldPoints[operator[](1)],
668  // oldPoints[operator[](2)]
669  // ).sweptVol
670  // (
671  // triPointRef
672  // (
673  // newPoints[operator[](0)],
674  // newPoints[operator[](1)],
675  // newPoints[operator[](2)]
676  // )
677  // )
678  // );
679  // }
680 
681  scalar sv = 0;
682 
683  // Calculate the swept volume by breaking the face into triangles and
684  // summing their swept volumes.
685  // Changed to deal with small concavity by using a central decomposition
686 
687  point centreOldPoint = centre(oldPoints);
688  point centreNewPoint = centre(newPoints);
689 
690  label nPoints = size();
691 
692  for (label pi=0; pi<nPoints-1; ++pi)
693  {
694  // Note: for best accuracy, centre point always comes last
695  sv += triPointRef
696  (
697  centreOldPoint,
698  oldPoints[operator[](pi)],
699  oldPoints[operator[](pi + 1)]
700  ).sweptVol
701  (
703  (
704  centreNewPoint,
705  newPoints[operator[](pi)],
706  newPoints[operator[](pi + 1)]
707  )
708  );
709  }
710 
711  sv += triPointRef
712  (
713  centreOldPoint,
714  oldPoints[operator[](nPoints-1)],
715  oldPoints[operator[](0)]
716  ).sweptVol
717  (
719  (
720  centreNewPoint,
721  newPoints[operator[](nPoints-1)],
722  newPoints[operator[](0)]
723  )
724  );
725 
726  return sv;
727 }
728 
729 
731 (
732  const pointField& p,
733  const point& refPt,
734  scalar density
735 ) const
736 {
737  // If the face is a triangle, do a direct calculation
738  if (size() == 3)
739  {
740  return triPointRef
741  (
742  p[operator[](0)],
743  p[operator[](1)],
744  p[operator[](2)]
745  ).inertia(refPt, density);
746  }
747 
748  const point ctr = centre(p);
749 
750  tensor J = tensor::zero;
751 
752  forAll(*this, i)
753  {
754  J += triPointRef
755  (
756  p[operator[](i)],
757  p[operator[](fcIndex(i))],
758  ctr
759  ).inertia(refPt, density);
760  }
761 
762  return J;
763 }
764 
765 
767 {
768  const labelList& points = *this;
769 
770  edgeList e(points.size());
771 
772  for (label pointI = 0; pointI < points.size() - 1; ++pointI)
773  {
774  e[pointI] = edge(points[pointI], points[pointI + 1]);
775  }
776 
777  // Add last edge
778  e.last() = edge(points.last(), points[0]);
779 
780  return e;
781 }
782 
783 
785 {
786  forAll(*this, i)
787  {
788  if (operator[](i) == e.start())
789  {
790  if (operator[](rcIndex(i)) == e.end())
791  {
792  // Reverse direction
793  return -1;
794  }
795  else if (operator[](fcIndex(i)) == e.end())
796  {
797  // Forward direction
798  return 1;
799  }
800 
801  // No match
802  return 0;
803  }
804  else if (operator[](i) == e.end())
805  {
806  if (operator[](rcIndex(i)) == e.start())
807  {
808  // Forward direction
809  return 1;
810  }
811  else if (operator[](fcIndex(i)) == e.start())
812  {
813  // Reverse direction
814  return -1;
815  }
816 
817  // No match
818  return 0;
819  }
820  }
821 
822  // Not found
823  return 0;
824 }
825 
826 
827 // Number of triangles directly known from number of vertices
829 {
830  return nTriangles();
831 }
832 
833 
835 (
836  const pointField& points,
837  label& triI,
838  faceList& triFaces
839 ) const
840 {
841  label quadI = 0;
842  faceList quadFaces;
843 
844  return split(SPLITTRIANGLE, points, triI, quadI, triFaces, quadFaces);
845 }
846 
847 
849 (
850  const pointField& points,
851  label& triI,
852  label& quadI
853 ) const
854 {
855  faceList triFaces;
856  faceList quadFaces;
857 
858  return split(COUNTQUAD, points, triI, quadI, triFaces, quadFaces);
859 }
860 
861 
863 (
864  const pointField& points,
865  label& triI,
866  label& quadI,
867  faceList& triFaces,
868  faceList& quadFaces
869 ) const
870 {
871  return split(SPLITQUAD, points, triI, quadI, triFaces, quadFaces);
872 }
873 
874 
876 {
877  const edgeList& eds = f.edges();
878 
879  label longestEdgeI = -1;
880  scalar longestEdgeLength = -SMALL;
881 
882  forAll(eds, edI)
883  {
884  scalar edgeLength = eds[edI].mag(pts);
885 
886  if (edgeLength > longestEdgeLength)
887  {
888  longestEdgeI = edI;
889  longestEdgeLength = edgeLength;
890  }
891  }
892 
893  return longestEdgeI;
894 }
895 
896 
897 // ************************************************************************* //
label which(const label globalIndex) const
Navigation through face vertices.
Definition: face.C:635
label trianglesQuads(const pointField &points, label &triI, label &quadI, faceList &triFaces, faceList &quadFaces) const
Split into triangles and quads.
Definition: face.C:863
vector point
Point is a vector.
Definition: point.H:41
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:28
dimensioned< scalar > mag(const dimensioned< Type > &)
face reverseFace() const
Return face with reverse direction.
Definition: face.C:616
labelList f(nPoints)
Walks over a container as if it were circular. The container must have the following members defined:...
static const Tensor zero
Definition: Tensor.H:80
T & last()
Return the last element of the list.
Definition: UListI.H:131
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Definition: globalIndex.H:63
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:59
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:76
label fcIndex(const label i) const
Return the forward circular index, i.e. the next index.
Definition: UListI.H:58
face triFace(3)
edgeList edges() const
Return edges in face point ordering,.
Definition: face.C:766
List< face > faceList
Definition: faceListFwd.H:43
Xfer< T > xferMove(T &)
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:58
label collapse()
Collapse face by removing duplicate point labels.
Definition: face.C:454
point centre(const pointField &) const
Centre point of face.
Definition: face.C:493
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
A triangle primitive used to calculate face normals and swept volumes.
Definition: triangle.H:57
triangle< point, const point & > triPointRef
Definition: triPointRef.H:44
label n
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
scalar sweptVol(const pointField &oldPoints, const pointField &newPoints) const
Return the volume swept out by the face when its points move.
Definition: face.C:652
T & operator[](const label)
Return element of UList.
Definition: UListI.H:163
static bool sameVertices(const face &, const face &)
Return true if the faces have the same vertices.
Definition: face.C:405
const double e
Elementary charge.
Definition: doubleFloat.H:78
void setSize(const label)
Reset size of List.
Definition: List.C:318
tensor inertia(const pointField &, const point &refPt=vector::zero, scalar density=1.0) const
Return the inertia tensor, with optional reference.
Definition: face.C:731
void setIteratorToFulcrum()
Set the iterator to the current position of the fulcrum.
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
Field< vector > vectorField
Specialisation of Field<T> for vector.
#define forAll(list, i)
Definition: UList.H:421
Swap its arguments.
dimensionedScalar acos(const dimensionedScalar &ds)
void setFulcrumToIterator()
Set the fulcrum to the current position of the iterator.
errorManip< error > abort(error &err)
Definition: errorManip.H:131
face()
Construct null.
Definition: faceI.H:44
label nTrianglesQuads(const pointField &points, label &nTris, label &nQuads) const
Number of triangles and quads after splitting.
Definition: face.C:849
void Swap(T &a, T &b)
Definition: Swap.H:43
label longestEdge(const face &f, const pointField &pts)
Find the longest edge on a face. Face point labels index into pts.
Definition: face.C:875
#define FatalErrorIn(functionName)
Report an error message using Foam::FatalError.
Definition: error.H:314
void flip()
Flip the face in-place.
Definition: face.C:479
int edgeDirection(const edge &) const
Return the edge direction on the face.
Definition: face.C:784
vector normal(const pointField &) const
Vector normal; magnitude is equal to area of face.
Definition: face.C:557
error FatalError
label end() const
Return end vertex label.
Definition: edgeI.H:92
static const Vector zero
Definition: Vector.H:80
label size() const
Return the number of elements in the UList.
Definition: ListI.H:83
bool circulate(const CirculatorBase::direction dir=NONE)
Circulate around the list in the given direction.
A triangular face using a FixedList of labels corresponding to mesh vertices.
Definition: triFace.H:68
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
label rcIndex(const label i) const
Return the reverse circular index, i.e. the previous index.
Definition: UListI.H:65
label start() const
Return start vertex label.
Definition: edgeI.H:81
static const char *const typeName
Definition: face.H:143
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
label nPoints
label nTriangles() const
Number of triangles after splitting.
Definition: faceI.H:130
label triangles(const pointField &points, label &triI, faceList &triFaces) const
Split into triangles using existing points.
Definition: face.C:835
A class for managing temporary objects.
Definition: PtrList.H:118
static int compare(const face &, const face &)
Compare faces.
Definition: face.C:303