polygonTriangulate.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2021-2025 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 "polygonTriangulate.H"
27 #include "tensor2D.H"
28 
29 // * * * * * * * * * * * Private Static Member Functions * * * * * * * * * * //
30 
31 template<class PointField>
32 Foam::scalar Foam::polygonTriangulate::area
33 (
34  const triFace& triPoints,
35  const PointField& points,
36  const vector& normal
37 )
38 {
39  const point& o = points[triPoints[0]];
40  const vector ao = points[triPoints[1]] - o;
41  const vector ab = points[triPoints[2]] - o;
42 
43  return (ao ^ ab) & normal;
44 }
45 
46 
47 template<class PointField>
48 Foam::scalar Foam::polygonTriangulate::quality
49 (
50  const triFace& triPoints,
51  const PointField& points,
52  const vector& normal
53 )
54 {
55  const scalar An = area(triPoints, points, normal);
56 
57  scalar PSqr = 0;
58  forAll(triPoints, triEdgei)
59  {
60  const point& p0 = points[triPoints[triEdgei]];
61  const point& p1 = points[triPoints[(triEdgei + 1) % 3]];
62  PSqr += magSqr(p0 - p1);
63  }
64 
65  static const scalar equilateralAnByPSqr = sqrt(scalar(3))/12;
66 
67  return PSqr != 0 ? An/PSqr/equilateralAnByPSqr : 0;
68 }
69 
70 
71 template<class PointField>
72 bool Foam::polygonTriangulate::intersection
73 (
74  const edge& edgePointsA,
75  const edge& edgePointsB,
76  const PointField& points,
77  const vector& normal
78 )
79 {
80  const point& pointA0 = points[edgePointsA[0]];
81  const point& pointA1 = points[edgePointsA[1]];
82  const point& pointB0 = points[edgePointsB[0]];
83  const point& pointB1 = points[edgePointsB[1]];
84 
85  // Bounding sphere-based rejection for problematic co-linear cases
86  const point centreA = (pointA0 + pointA1)/2;
87  const point centreB = (pointB0 + pointB1)/2;
88  const scalar radiusSqrA = magSqr((pointA0 - pointA1)/2);
89  const scalar radiusSqrB = magSqr((pointB0 - pointB1)/2);
90 
91  if (magSqr(centreA - centreB) > radiusSqrA + radiusSqrB) return false;
92 
93  // Solve for the edge-local coordinates of the intersection. If these are
94  // both in the range 0 -> 1 then the edges intersect.
95  const vector tau0 = perpendicular(normal);
96  const vector tau1 = normal ^ tau0;
97 
98  const point2D point2DA0(tau0 & pointA0, tau1 & pointA0);
99  const point2D point2DA1(tau0 & pointA1, tau1 & pointA1);
100  const point2D point2DB0(tau0 & pointB0, tau1 & pointB0);
101  const point2D point2DB1(tau0 & pointB1, tau1 & pointB1);
102 
103  const tensor2D M =
104  tensor2D(point2DA0 - point2DA1, point2DB1 - point2DB0).T();
105 
106  const scalar detM = det(M);
107  const tensor2D detMInvM = cof(M);
108 
109  const vector2D magDetMTu = sign(detM)*detMInvM & (point2DA0 - point2DB0);
110 
111  return 0 <= cmptMin(magDetMTu) && cmptMax(magDetMTu) <= mag(detM);
112 }
113 
114 
115 template<class PointField>
116 Foam::label Foam::polygonTriangulate::nIntersections
117 (
118  const edge& edgePoints,
119  const PointField& points,
120  const vector& normal
121 )
122 {
123  label n = 0;
124 
125  forAll(points, pi)
126  {
127  const edge otherEdgePoints(pi, points.fcIndex(pi));
128 
129  if (edgePoints.commonVertex(otherEdgePoints) == -1)
130  {
131  n += intersection(edgePoints, otherEdgePoints, points, normal);
132  }
133  }
134 
135  return n;
136 }
137 
138 
139 template<class PointField>
140 Foam::scalar Foam::polygonTriangulate::angle
141 (
142  const label pointi,
143  const PointField& points,
144  const vector& normal
145 )
146 {
147  const point& o = points[pointi];
148  const vector oa = points[points.rcIndex(pointi)] - o;
149  const vector ob = points[points.fcIndex(pointi)] - o;
150 
151  const vector oaNegNNOa = oa - normal*(normal & oa);
152  const vector obNegNNOb = ob - normal*(normal & ob);
153 
154  return
155  atan2(normal & (oaNegNNOa ^ obNegNNOb), - oaNegNNOa & obNegNNOb)
157 }
158 
159 
160 template<class PointField>
161 bool Foam::polygonTriangulate::ear
162 (
163  const label pointi,
164  const PointField& points,
165  const vector& normal
166 )
167 {
168  const point& o = points[pointi];
169  const vector oa = points[points.rcIndex(pointi)] - o;
170  const vector ob = points[points.fcIndex(pointi)] - o;
171 
172  const tensor A = tensor(ob, oa, normal).T();
173  const tensor T(A.y() ^ A.z(), A.z() ^ A.x(), A.x() ^ A.y());
174  const scalar detA = det(A);
175 
176  for
177  (
178  label pointj = points.fcIndex(points.fcIndex(pointi));
179  pointj != points.rcIndex(pointi);
180  pointj = points.fcIndex(pointj)
181  )
182  {
183  const vector detAY = (points[pointj] - o) & T;
184 
185  if (detAY.x() > 0 && detAY.y() > 0 && detAY.x() + detAY.y() < detA)
186  {
187  return false;
188  }
189  }
190 
191  return true;
192 }
193 
194 
195 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
196 
198 (
200  const label n,
201  const scalar error
202 )
203 {
204  // Get random points on a unit disk
206  forAll(points, pointi)
207  {
208  const scalar theta =
210  const scalar r = sqrt(rndGen.sample01<scalar>());
211  points[pointi] = point(r*cos(theta), r*sin(theta), 0);
212  }
213 
214  // Initialise intersected polygon
215  labelList pointis(identityMap(n));
216 
217  // Reorder until non-intersected
218  bool incomplete = true;
219  while (incomplete)
220  {
221  incomplete = false;
222 
223  for (label piA0 = 0; piA0 < n; ++ piA0)
224  {
225  const label piA1 = points.fcIndex(piA0);
226 
227  for (label piB0 = piA0 + 2; piB0 < (piA0 == 0 ? n - 1 : n); ++ piB0)
228  {
229  const label piB1 = points.fcIndex(piB0);
230 
231  if
232  (
234  (
235  edge(pointis[piA0], pointis[piA1]),
236  edge(pointis[piB0], pointis[piB1]),
237  points,
238  vector(0, 0, 1)
239  )
240  )
241  {
242  SubList<label> subOrder(pointis, piB0 + 1 - piA1, piA1);
243  inplaceReverseList(subOrder);
244  incomplete = true;
245  }
246 
247  if (incomplete) break;
248  }
249 
250  if (incomplete) break;
251  }
252  }
253 
254  // Add noise and potential self-intersection
255  forAll(points, pointi)
256  {
257  const scalar theta =
259  const scalar r = sqrt(rndGen.sample01<scalar>());
260  points[pointi] =
261  (1 - error)*points[pointi]
262  + error*point(r*cos(theta), r*sin(theta), rndGen.sample01<scalar>());
263  }
264 
265  // Reorder and return
266  return List<point>(points, pointis);
267 }
268 
269 
270 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
271 
272 template<class PointField>
273 void Foam::polygonTriangulate::optimiseTriangulation
274 (
275  const label trii,
276  const PointField& points,
277  const vector& normal,
278  UList<triFace>& triPoints,
279  UList<FixedList<label, 3>>& triEdges,
280  UList<labelPair>& edgeTris
281 )
282 {
283  const triFace& t = triPoints[trii];
284  const scalar q = quality(t, points, normal);
285 
286  forAll(triEdges[trii], triEdgei)
287  {
288  const label edgei = triEdges[trii][triEdgei];
289 
290  const label otherTrii = edgeTris[edgei][edgeTris[edgei][0] == trii];
291 
292  if (otherTrii == -1) continue;
293 
294  const label otherTriEdgei = findIndex(triEdges[otherTrii], edgei);
295 
296  const triFace& otherT = triPoints[otherTrii];
297  const scalar otherQ = quality(otherT, points, normal);
298 
299  const label piA = triPoints[trii][(triEdgei + 1) % 3];
300  const label piB = triPoints[trii][(triEdgei + 2) % 3];
301  const label piC = triPoints[otherTrii][(otherTriEdgei + 1) % 3];
302  const label piD = triPoints[otherTrii][(otherTriEdgei + 2) % 3];
303 
304  const triFace flipT0(piA, piB, piD);
305  const triFace flipT1(piC, piD, piB);
306 
307  if
308  (
309  triPointsSet_.found(flipT0)
310  || triPointsSet_.found(flipT1)
311  ) continue;
312 
313  const scalar flipQ0 = quality(flipT0, points, normal);
314  const scalar flipQ1 = quality(flipT1, points, normal);
315 
316  if (min(q, otherQ) < min(flipQ0, flipQ1))
317  {
318  triPointsSet_.set(t);
319  triPointsSet_.set(otherT);
320 
321  const label eiA = triEdges[trii][(triEdgei + 1) % 3];
322  const label eiB = triEdges[trii][(triEdgei + 2) % 3];
323  const label eiC = triEdges[otherTrii][(otherTriEdgei + 1) % 3];
324  const label eiD = triEdges[otherTrii][(otherTriEdgei + 2) % 3];
325 
326  triPoints[trii] = {piA, piB, piD};
327  triEdges[trii] = {eiA, edgei, eiD};
328  edgeTris[eiD][edgeTris[eiD][0] != otherTrii] = trii;
329 
330  triPoints[otherTrii] = {piC, piD, piB};
331  triEdges[otherTrii] = {eiC, edgei, eiB};
332  edgeTris[eiB][edgeTris[eiB][0] != trii] = otherTrii;
333 
334  optimiseTriangulation
335  (
336  trii,
337  points,
338  normal,
339  triPoints,
340  triEdges,
341  edgeTris
342  );
343  optimiseTriangulation
344  (
345  otherTrii,
346  points,
347  normal,
348  triPoints,
349  triEdges,
350  edgeTris
351  );
352 
353  break;
354  }
355  }
356 }
357 
358 
359 void Foam::polygonTriangulate::simpleTriangulate
360 (
361  const UList<point>& allPoints,
362  const vector& normal,
363  UList<triFace>& triPoints,
364  UList<FixedList<label, 3>>& triEdges,
365  UList<labelPair>& edgeTris,
366  const bool optimal
367 )
368 {
369  // Get the polygon size
370  const label n = allPoints.size();
371 
372  // Clear the triangulation
373  triPoints = triFace(-1, -1, -1);
374  triEdges = FixedList<label, 3>({-1, -1, -1});
375  edgeTris = labelPair(-1, -1);
376 
377  // Initialise workspace
378  pointis_.setSize(n);
379  edges_.setSize(n);
380  angle_.setSize(n);
381  ear_.setSize(n);
382  forAll(pointis_, i)
383  {
384  pointis_[i] = i;
385  edges_[i] = i;
386  angle_[i] = angle(i, allPoints, normal);
387  ear_[i] = ear(i, allPoints, normal);
388  }
389 
390  // Create the indirect point list for the current un-triangulated polygon
391  UIndirectList<point> points(allPoints, pointis_);
392 
393  // Generate triangles
394  forAll(triPoints, trii)
395  {
396  // Find the ear with the smallest angle
397  scalar minEarAngle = vGreat;
398  label minEarAnglei = -1;
399  for (label i = 0; i < pointis_.size(); ++ i)
400  {
401  if (ear_[i] && minEarAngle > angle_[i])
402  {
403  minEarAngle = angle_[i];
404  minEarAnglei = i;
405  }
406  }
407 
408  // If nothing is an ear, then this face is degenerate. Default to
409  // the smallest angle, ignoring the ear test.
410  if (minEarAnglei == -1)
411  {
412  minEarAnglei = findMin(angle_);
413  }
414 
415  // Get the adjacent points
416  label minEarAnglei0 = pointis_.rcIndex(minEarAnglei);
417  label minEarAnglei1 = pointis_.fcIndex(minEarAnglei);
418 
419  // Add a triangle
420  triPoints[trii] = triFace
421  (
422  pointis_[minEarAnglei0],
423  pointis_[minEarAnglei],
424  pointis_[minEarAnglei1]
425  );
426  triEdges[trii] = FixedList<label, 3>
427  ({
428  edges_[minEarAnglei0],
429  edges_[minEarAnglei],
430  trii == triPoints.size() - 1 ? edges_[minEarAnglei1] : n + trii
431  });
432  forAll(triEdges[trii], triEdgei)
433  {
434  const label edgei = triEdges[trii][triEdgei];
435  edgeTris[edgei][edgeTris[edgei][0] != -1] = trii;
436  }
437 
438  // Remove the cut off point from the polygon
439  edges_[minEarAnglei0] = triEdges[trii][2];
440  for (label i = minEarAnglei; i < pointis_.size() - 1; ++ i)
441  {
442  pointis_[i] = pointis_[i + 1];
443  edges_[i] = edges_[i + 1];
444  angle_[i] = angle_[i + 1];
445  ear_[i] = ear_[i + 1];
446  }
447  pointis_.resize(pointis_.size() - 1);
448  edges_.resize(edges_.size() - 1);
449  angle_.resize(angle_.size() - 1);
450  ear_.resize(ear_.size() - 1);
451 
452  // Update connected points
453  if (minEarAnglei0 > minEarAnglei) -- minEarAnglei0;
454  angle_[minEarAnglei0] = angle(minEarAnglei0, points, normal);
455  ear_[minEarAnglei0] = ear(minEarAnglei0, points, normal);
456  if (minEarAnglei1 > minEarAnglei) -- minEarAnglei1;
457  angle_[minEarAnglei1] = angle(minEarAnglei1, points, normal);
458  ear_[minEarAnglei1] = ear(minEarAnglei1, points, normal);
459 
460  // Optimise the quality
461  if (optimal)
462  {
463  triPointsSet_.clear();
464  optimiseTriangulation
465  (
466  trii,
467  allPoints,
468  normal,
469  triPoints,
470  triEdges,
471  edgeTris
472  );
473  }
474  }
475 }
476 
477 
478 void Foam::polygonTriangulate::partitionTriangulate
479 (
480  const UList<point>& points,
481  const vector& normal,
482  const label spanPointi,
483  const label spanEdgei,
484  UList<triFace>& triPoints,
485  UList<FixedList<label, 3>>& triEdges,
486  UList<labelPair>& edgeTris,
487  const bool simple,
488  const bool optimal
489 )
490 {
491  // Get the polygon size
492  const label n = points.size();
493 
494  // Get the points of the spanning triangle
495  const label piA = (spanEdgei + 1) % n;
496  const label piP = spanPointi;
497  const label piB = spanEdgei;
498 
499  // Determine the sizes of the polygons connected to either side of the
500  // spanning triangle
501  const label nA = (piP - piA + n) % n + 1;
502  const label nB = (piB - piP + n) % n + 1;
503 
504  // Get the edges of the spanning triangle
505  const label eiA = nA >= 3 ? n : piA;
506  const label eiB = nB >= 3 ? n + (nA >= 3) : piP;
507  const label eiE = piB;
508 
509  // Add the spanning triangle
510  triPoints[0] = triFace(piA, piP, piB);
511  triEdges[0] = FixedList<label, 3>({eiA, eiB, eiE});
512 
513  // Rotate the point list so that it can be passed to sub-triangulations
514  // without indirect addressing
515  inplaceRotateList(const_cast<UList<point>&>(points), - piA);
516 
517  // Triangulate the sub-polygon connected to edge A
518  if (nA >= 3)
519  {
520  // Subset the polygon and triangulate
521  SubList<triFace> triPointsA(triPoints, nA - 2, 1);
522  SubList<FixedList<label, 3>> triEdgesA(triEdges, nA - 2, 1);
523  SubList<labelPair> edgeTrisA(edgeTris, 2*nA - 3);
524  triangulate
525  (
526  SubList<point>(points, nA),
527  normal,
528  triPointsA,
529  triEdgesA,
530  edgeTrisA,
531  simple,
532  optimal
533  );
534 
535  // Map the point labels back to the full polygon
536  forAll(triPointsA, triAi)
537  {
538  forAll(triPointsA[triAi], triPointAi)
539  {
540  label& pointi = triPointsA[triAi][triPointAi];
541  pointi = (pointi + piA) % n;
542  }
543  }
544 
545  // Map the edge labels back to the full polygon
546  forAll(triEdgesA, triAi)
547  {
548  forAll(triEdgesA[triAi], triEdgeAi)
549  {
550  label& edgei = triEdgesA[triAi][triEdgeAi];
551  edgei =
552  edgei >= nA ? max(eiA, eiB) + edgei - nA + 1
553  : edgei == nA - 1 ? eiA
554  : (piA + edgei) % n;
555  }
556  }
557 
558  // Map the tri labels back to the full polygon
559  forAll(edgeTrisA, edgeAi)
560  {
561  forAll(edgeTrisA[edgeAi], edgeTriAi)
562  {
563  label& trii = edgeTrisA[edgeAi][edgeTriAi];
564  trii = trii == -1 ? -1 : trii + 1;
565  }
566  }
567  }
568 
569  // Triangulate the sub-polygon connected to edge B
570  if (nB >= 3)
571  {
572  // Subset the polygon and triangulate
573  SubList<triFace> triPointsB(triPoints, nB - 2, nA - 1);
574  SubList<FixedList<label, 3>> triEdgesB(triEdges, nB - 2, nA - 1);
575  SubList<labelPair> edgeTrisB(edgeTris, 2*nB - 3, 2*nA - 3);
576  triangulate
577  (
578  SubList<point>(points, nB, nA - 1),
579  normal,
580  triPointsB,
581  triEdgesB,
582  edgeTrisB,
583  simple,
584  optimal
585  );
586 
587  // Map the point labels back to the full polygon
588  forAll(triPointsB, triBi)
589  {
590  forAll(triPointsB[triBi], triPointBi)
591  {
592  label& pointi = triPointsB[triBi][triPointBi];
593  pointi = (pointi + piA + nA - 1) % n;
594  }
595  }
596 
597  // Map the edge labels back to the full polygon
598  forAll(triEdgesB, triBi)
599  {
600  forAll(triEdgesB[triBi], triEdgeBi)
601  {
602  label& edgei = triEdgesB[triBi][triEdgeBi];
603  edgei =
604  edgei >= nB ? max(eiA, eiB) + nA - 3 + edgei - nB + 1
605  : edgei == nB - 1 ? eiB
606  : (piA + nA - 1 + edgei) % n;
607  }
608  }
609 
610  // Map the tri labels back to the full polygon
611  forAll(edgeTrisB, edgeBi)
612  {
613  forAll(edgeTrisB[edgeBi], edgeTriBi)
614  {
615  label& trii = edgeTrisB[edgeBi][edgeTriBi];
616  trii = trii == -1 ? -1 : trii + nA - 1;
617  }
618  }
619  }
620 
621  // Rotate the point list back
622  inplaceRotateList(const_cast<UList<point>&>(points), piA);
623 
624  // Reorder the edge-tris
625  {
626  // Swap the A-internal-edges and B-outer-edges to get all outer
627  // edges and all internal edges in contiguous blocks
628  SubList<labelPair> l(edgeTris, nA + nB - 3, nA - 1);
629  inplaceRotateList(l, - nA + 2);
630  }
631  {
632  // Rotate the outer edges right by one so that the intersection
633  // edge (currently at the end) moves adjacent to the outer edges
634  SubList<labelPair> l(edgeTris, nA + nB - 3, nA + nB - 2);
635  inplaceRotateList(l, 1);
636  }
637  {
638  // Rotate the outer edges so that they are in sequence with the
639  // original point ordering
640  SubList<labelPair> l(edgeTris, n);
641  inplaceRotateList(l, - n + piA);
642  }
643  if (nA >= 3 && nB >= 3)
644  {
645  // Move edge-B leftwards adjacent to edge-A as this is created
646  // before any of the other internal edges
647  const label fromi = n + nA - 2, toi = n + 1;
648  const labelPair temp = edgeTris[fromi];
649  for (label i = fromi; i > toi; -- i)
650  {
651  edgeTris[i] = edgeTris[i-1];
652  }
653  edgeTris[toi] = temp;
654  }
655 
656  // Add associations to the intersection triangle
657  edgeTris[eiA][edgeTris[eiA][0] != -1] = 0;
658  edgeTris[eiB][edgeTris[eiB][0] != -1] = 0;
659  edgeTris[eiE] = labelPair(0, -1);
660 }
661 
662 
663 void Foam::polygonTriangulate::complexTriangulate
664 (
665  const UList<point>& points,
666  const vector& normal,
667  UList<triFace>& triPoints,
668  UList<FixedList<label, 3>>& triEdges,
669  UList<labelPair>& edgeTris,
670  const bool optimal
671 )
672 {
673  // Get the polygon size
674  const label n = points.size();
675 
676  // Detect self intersections. When one is found, remove one of the
677  // intersected edges by adding a spanning triangle and recurse. Pick the
678  // spanning triangle that results in the smallest negative area.
679  for (label piA0 = 0; piA0 < n; ++ piA0)
680  {
681  const label piA1 = points.fcIndex(piA0);
682 
683  for (label piB0 = piA0 + 2; piB0 < (piA0 == 0 ? n - 1 : n); ++ piB0)
684  {
685  const label piB1 = points.fcIndex(piB0);
686 
687  if
688  (
689  intersection
690  (
691  edge(piA0, piA1),
692  edge(piB0, piB1),
693  points,
694  normal
695  )
696  )
697  {
698  // Get a bunch of unique spanning triangles that remove one of
699  // the intersected edges
700  FixedList<labelPair, 8> spanPointAndEdgeis
701  ({
702  labelPair(piA0, piB0),
703  labelPair(piA1, piB0),
704  labelPair(piB0, piA0),
705  labelPair(piB1, piA0),
706  labelPair(points.rcIndex(piA0), piA0),
707  labelPair(points.fcIndex(piA1), piA0),
708  labelPair(points.rcIndex(piB0), piB0),
709  labelPair(points.fcIndex(piB1), piB0)
710  });
711  forAll(spanPointAndEdgeis, spani)
712  {
713  for (label spanj = spani + 1; spanj < 8; ++ spanj)
714  {
715  if
716  (
717  spanPointAndEdgeis[spani]
718  == spanPointAndEdgeis[spanj]
719  )
720  {
721  spanPointAndEdgeis[spanj] = labelPair(-1, -1);
722  }
723  }
724  }
725 
726  // Find the spanning triangle that results in the smallest
727  // negative triangulation area
728  scalar spanSumNegA = - vGreat;
729  label spanPointi = -1, spanEdgei = -1;
730  forAll(spanPointAndEdgeis, spani)
731  {
732  if (spanPointAndEdgeis[spani] == labelPair(-1, -1))
733  continue;
734 
735  const label piP = spanPointAndEdgeis[spani].first();
736  const label ei = spanPointAndEdgeis[spani].second();
737  const label piA0 = points.rcIndex(ei);
738  const label piA = ei;
739  const label piB = points.fcIndex(ei);
740  const label piB1 = points.fcIndex(piB);
741 
742  // Only consider this spanning triangle if it decreases
743  // the number of intersections
744  const label netNIntersections =
745  - nIntersections(edge(piA, piB), points, normal)
746  + (piP == piA0 ? -1 : +1)
747  *nIntersections(edge(piP, piA), points, normal)
748  + (piP == piB1 ? -1 : +1)
749  *nIntersections(edge(piP, piB), points, normal);
750  if (netNIntersections >= 0) continue;
751 
752  // Triangulate
753  partitionTriangulate
754  (
755  points,
756  normal,
757  piP,
758  ei,
759  triPoints,
760  triEdges,
761  edgeTris,
762  false,
763  false
764  );
765 
766  // Compute the negative area. If it is smaller than the
767  // previous value then update the spanning triangle. If it
768  // is zero then break and use this triangle.
769  scalar sumNegA = 0;
770  forAll(triPoints, trii)
771  {
772  const scalar a = area(triPoints[trii], points, normal);
773  sumNegA += negPart(a);
774  }
775  if (sumNegA > spanSumNegA)
776  {
777  spanSumNegA = sumNegA;
778  spanPointi = piP;
779  spanEdgei = ei;
780  }
781  if (spanSumNegA == 0)
782  {
783  break;
784  }
785  }
786 
787  // If a suitable spanning triangle was found then partition the
788  // triangulation
789  if (spanPointi != -1)
790  {
791  partitionTriangulate
792  (
793  points,
794  normal,
795  spanPointi,
796  spanEdgei,
797  triPoints,
798  triEdges,
799  edgeTris,
800  false,
801  optimal
802  );
803 
804  return;
805  }
806  }
807  }
808  }
809 
810  // If there are no intersections then do simple polygon triangulation
811  simpleTriangulate(points, normal, triPoints, triEdges, edgeTris, optimal);
812 }
813 
814 
815 void Foam::polygonTriangulate::triangulate
816 (
817  const UList<point>& points,
818  const vector& normal,
819  UList<triFace>& triPoints,
820  UList<FixedList<label, 3>>& triEdges,
821  UList<labelPair>& edgeTris,
822  const bool simple,
823  const bool optimal
824 )
825 {
826  if (simple)
827  {
828  simpleTriangulate
829  (
830  points,
831  normal,
832  triPoints,
833  triEdges,
834  edgeTris,
835  optimal
836  );
837  }
838  else
839  {
840  complexTriangulate
841  (
842  points,
843  normal,
844  triPoints,
845  triEdges,
846  edgeTris,
847  optimal
848  );
849  }
850 }
851 
852 
853 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
854 
856 :
857  pointis_(),
858  edges_(),
859  angle_(),
860  ear_(),
861 
862  triPointsSet_(),
863 
864  points_(),
865  triPoints_(),
866  triEdges_(),
867  edgeTris_()
868 {}
869 
870 
871 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
872 
874 {}
875 
876 
877 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
878 
879 const Foam::UList<Foam::triFace>& Foam::polygonTriangulate::triangulate
880 (
881  const UList<point>& points,
882  const vector& normal,
883  const bool simple,
884  const bool optimal
885 )
886 {
887  // Get the polygon size
888  const label n = points.size();
889 
890  // Resize workspace
891  triPoints_.resize(n - 2);
892  triEdges_.resize(n - 2);
893  edgeTris_.resize(2*n - 3);
894 
895  // Run
896  triangulate
897  (
898  points,
899  normal,
900  triPoints_,
901  triEdges_,
902  edgeTris_,
903  simple,
904  optimal
905  );
906 
907  return triPoints_;
908 }
909 
910 
911 const Foam::UList<Foam::triFace>& Foam::polygonTriangulate::triangulate
912 (
913  const UList<point>& points,
914  const bool simple,
915  const bool optimal
916 )
917 {
918  return
919  triangulate
920  (
921  points,
923  simple,
924  optimal
925  );
926 }
927 
928 
929 // ************************************************************************* //
#define M(I)
label n
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: List.H:91
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
A List obtained as a section of another List.
Definition: SubList.H:56
Templated 2D tensor derived from VectorSpace adding construction from 4 components,...
Definition: Tensor2D.H:59
Tensor< Cmpt > T() const
Return transpose.
Definition: TensorI.H:331
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:74
label rcIndex(const label i) const
Return the reverse circular index, i.e. the previous index.
Definition: UListI.H:65
label fcIndex(const label i) const
Return the forward circular index, i.e. the next index.
Definition: UListI.H:58
const Cmpt & y() const
Definition: VectorI.H:81
const Cmpt & x() const
Definition: VectorI.H:75
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:61
label commonVertex(const edge &a) const
Return common vertex.
Definition: edgeI.H:122
Class to handle errors and exceptions in a simple, consistent stream-based manner.
Definition: error.H:71
static vector area(const PointField &ps)
Return vector area given face points.
Foam::intersection.
Definition: intersection.H:50
static List< point > randomPolygon(randomGenerator &rndGen, const label n, const scalar error)
Generate a random polygon for testing purposes.
polygonTriangulate()
Construct null.
const UList< triFace > & triPoints() const
Get the triangles' points.
Random number generator.
Type sample01()
Return a type with components uniformly distributed between.
A triangular face using a FixedList of labels corresponding to mesh vertices.
Definition: triFace.H:71
simpleControl simple(mesh)
const pointField & points
static const coefficient A("A", dimPressure, 611.21)
Pair< label > labelPair
Label pair.
Definition: labelPair.H:48
dimensionedScalar sign(const dimensionedScalar &ds)
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
Tensor< scalar > tensor
Tensor of scalars.
Definition: tensor.H:51
void inplaceReverseList(ListType &list)
Inplace reversal of a list using Swap.
dimensionedScalar sin(const dimensionedScalar &ds)
void T(LagrangianPatchField< Type > &f, const LagrangianPatchField< Type > &f1)
void cmptMin(FieldField< Field, typename FieldField< Field, Type >::cmptType > &cf, const FieldField< Field, Type > &f)
Tensor2D< scalar > tensor2D
Tensor2D or scalars.
Definition: tensor2D.H:49
void inplaceRotateList(ListType< DataType > &list, label n)
Inplace reversal of a list using the Reversal Block Swapping algorithm.
vector point
Point is a vector.
Definition: point.H:41
void det(LagrangianPatchField< scalar > &f, const LagrangianPatchField< tensor > &f1)
void mag(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
dimensionedScalar atan2(const dimensionedScalar &x, const dimensionedScalar &y)
dimensionedSymmTensor cof(const dimensionedSymmTensor &dt)
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
dimensionedScalar negPart(const dimensionedScalar &ds)
dimensionSet normalised(const dimensionSet &)
Definition: dimensionSet.C:507
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
void cmptMax(FieldField< Field, typename FieldField< Field, Type >::cmptType > &cf, const FieldField< Field, Type > &f)
label findMin(const ListType &, const label start=0)
Find index of min element (and less than given element).
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurrence of given element and return index,.
void magSqr(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
labelList identityMap(const label len)
Create identity map (map[i] == i) of given length.
Definition: ListOps.C:104
dimensionSet perpendicular(const dimensionSet &)
Definition: dimensionSet.C:513
void sqrt(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
dimensionedScalar cos(const dimensionedScalar &ds)
face triFace(3)
randomGenerator rndGen(653213)