collapseBase.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 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "collapseBase.H"
29 #include "triSurfaceTools.H"
30 #include "argList.H"
31 #include "OFstream.H"
32 #include "SubList.H"
33 #include "labelPair.H"
34 #include "meshTools.H"
35 #include "OSspecific.H"
36 
37 using namespace Foam;
38 
39 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
40 
42 //static void writeRegionOBJ
43 //(
44 // const triSurface& surf,
45 // const label regionI,
46 // const labelList& collapseRegion,
47 // const labelList& outsideVerts
48 //)
49 //{
50 // fileName dir("regions");
51 //
52 // mkDir(dir);
53 // fileName regionName(dir / "region_" + name(regionI) + ".obj");
54 //
55 // Pout<< "Dumping region " << regionI << " to file " << regionName << endl;
56 //
57 // boolList include(surf.size(), false);
58 //
59 // forAll(collapseRegion, facei)
60 // {
61 // if (collapseRegion[facei] == regionI)
62 // {
63 // include[facei] = true;
64 // }
65 // }
66 //
67 // triSurface regionSurf(surf.subsetMesh(include));
68 //
69 // Pout<< "Region " << regionI << " surface:" << nl;
70 // regionSurf.writeStats(Pout);
71 //
72 // regionSurf.write(regionName);
73 //
74 //
75 // // Dump corresponding outside vertices.
76 // fileName pointsName(dir / "regionPoints_" + name(regionI) + ".obj");
77 //
78 // Pout<< "Dumping region " << regionI << " points to file " << pointsName
79 // << endl;
80 //
81 // OFstream str(pointsName);
82 //
83 // forAll(outsideVerts, i)
84 // {
85 // meshTools::writeOBJ(str, surf.localPoints()[outsideVerts[i]]);
86 // }
87 //}
88 
89 
90 // Split triangle into multiple triangles because edge e being split
91 // into multiple edges.
92 static void splitTri
93 (
94  const labelledTri& f,
95  const edge& e,
96  const labelList& splitPoints,
98 )
99 {
100  //label oldNTris = tris.size();
101 
102  label fp = f.find(e[0]);
103  label fp1 = f.fcIndex(fp);
104  label fp2 = f.fcIndex(fp1);
105 
106  if (f[fp1] == e[1])
107  {
108  // Split triangle along fp to fp1
109  tris.append(labelledTri(f[fp2], f[fp], splitPoints[0], f.region()));
110 
111  for (label i = 1; i < splitPoints.size(); i++)
112  {
113  tris.append
114  (
116  (
117  f[fp2],
118  splitPoints[i-1],
119  splitPoints[i],
120  f.region()
121  )
122  );
123  }
124 
125  tris.append
126  (
128  (
129  f[fp2],
130  splitPoints.last(),
131  f[fp1],
132  f.region()
133  )
134  );
135  }
136  else if (f[fp2] == e[1])
137  {
138  // Split triangle along fp2 to fp. (Reverse order of splitPoints)
139 
140  tris.append
141  (
143  (
144  f[fp1],
145  f[fp2],
146  splitPoints.last(),
147  f.region()
148  )
149  );
150 
151  for (label i = splitPoints.size()-1; i > 0; --i)
152  {
153  tris.append
154  (
156  (
157  f[fp1],
158  splitPoints[i],
159  splitPoints[i-1],
160  f.region()
161  )
162  );
163  }
164 
165  tris.append
166  (
168  (
169  f[fp1],
170  splitPoints[0],
171  f[fp],
172  f.region()
173  )
174  );
175  }
176  else
177  {
179  << "Edge " << e << " not part of triangle " << f
180  << " fp:" << fp
181  << " fp1:" << fp1
182  << " fp2:" << fp2
183  << abort(FatalError);
184  }
185 
186  //Pout<< "Split face " << f << " along edge " << e
187  // << " into triangles:" << endl;
188  //
189  //for (label i = oldNTris; i < tris.size(); i++)
190  //{
191  // Pout<< " " << tris[i] << nl;
192  //}
193 }
194 
195 
196 // Insert scalar into sortedVerts/sortedWeights so the weights are in
197 // incrementing order.
198 static bool insertSorted
199 (
200  const label vertI,
201  const scalar weight,
202 
203  labelList& sortedVerts,
204  scalarField& sortedWeights
205 )
206 {
207  if (sortedVerts.found(vertI))
208  {
210  << " which is already in list of sorted vertices "
211  << sortedVerts << abort(FatalError);
212  }
213 
214  if (weight <= 0 || weight >= 1)
215  {
217  << " with illegal weight " << weight
218  << " into list of sorted vertices "
219  << sortedVerts << abort(FatalError);
220  }
221 
222 
223  label insertI = sortedVerts.size();
224 
225  forAll(sortedVerts, sortedI)
226  {
227  scalar w = sortedWeights[sortedI];
228 
229  if (mag(w - weight) < SMALL)
230  {
232  << "Trying to insert weight " << weight << " which is close to"
233  << " existing weight " << w << " in " << sortedWeights
234  << endl;
235  }
236 
237  if (w > weight)
238  {
239  // Start inserting before sortedI.
240  insertI = sortedI;
241  break;
242  }
243  }
244 
245 
246  label sz = sortedWeights.size();
247 
248  sortedWeights.setSize(sz + 1);
249  sortedVerts.setSize(sz + 1);
250 
251  // Leave everything up to (not including) insertI intact.
252 
253  // Make space by copying from insertI up.
254  for (label i = sz-1; i >= insertI; --i)
255  {
256  sortedWeights[i+1] = sortedWeights[i];
257  sortedVerts[i+1] = sortedVerts[i];
258  }
259  sortedWeights[insertI] = weight;
260  sortedVerts[insertI] = vertI;
261 
262  return true;
263 }
264 
265 
266 // Is triangle candidate for collapse? Small height or small quality
267 bool isSliver
268 (
269  const triSurface& surf,
270  const scalar minLen,
271  const scalar minQuality,
272  const label facei,
273  const label edgeI
274 )
275 {
276  const pointField& localPoints = surf.localPoints();
277 
278  // Check
279  // - opposite vertex projects onto base edge
280  // - normal distance is small
281  // - or triangle quality is small
282 
283  label opposite0 =
285  (
286  surf,
287  facei,
288  edgeI
289  );
290 
291  const edge& e = surf.edges()[edgeI];
292  const labelledTri& f = surf[facei];
293 
294  pointHit pHit =
295  e.line(localPoints).nearestDist
296  (
297  localPoints[opposite0]
298  );
299 
300  if
301  (
302  pHit.hit()
303  && (
304  pHit.distance() < minLen
305  || f.tri(surf.points()).quality() < minQuality
306  )
307  )
308  {
309  // Remove facei and split all other faces using this
310  // edge. This is done by 'replacing' the edgeI with the
311  // opposite0 vertex
312  //Pout<< "Splitting face " << facei << " since distance "
313  // << pHit.distance()
314  // << " from vertex " << opposite0
315  // << " to edge " << edgeI
316  // << " points "
317  // << localPoints[e[0]]
318  // << localPoints[e[1]]
319  // << " is too small or triangle quality "
320  // << f.tri(surf.points()).quality()
321  // << " too small." << endl;
322 
323  return true;
324  }
325 
326  return false;
327 }
328 
329 
330 // Mark all faces that are going to be collapsed.
331 // faceToEdge: per face -1 or the base edge of the face.
332 static void markCollapsedFaces
333 (
334  const triSurface& surf,
335  const scalar minLen,
336  const scalar minQuality,
337  labelList& faceToEdge
338 )
339 {
340  faceToEdge.setSize(surf.size());
341  faceToEdge = -1;
342 
343  const labelListList& edgeFaces = surf.edgeFaces();
344 
345  forAll(edgeFaces, edgeI)
346  {
347  const labelList& eFaces = surf.edgeFaces()[edgeI];
348 
349  forAll(eFaces, i)
350  {
351  label facei = eFaces[i];
352 
353  bool isCandidate = isSliver(surf, minLen, minQuality, facei, edgeI);
354 
355  if (isCandidate)
356  {
357  // Mark face as being collapsed
358  if (faceToEdge[facei] != -1)
359  {
361  << "Cannot collapse face " << facei << " since "
362  << " is marked to be collapsed both to edge "
363  << faceToEdge[facei] << " and " << edgeI
364  << abort(FatalError);
365  }
366 
367  faceToEdge[facei] = edgeI;
368  }
369  }
370  }
371 }
372 
373 
374 // Recurse through collapsed faces marking all of them with regionI (in
375 // collapseRegion)
376 static void markRegion
377 (
378  const triSurface& surf,
379  const labelList& faceToEdge,
380  const label regionI,
381  const label facei,
382  labelList& collapseRegion
383 )
384 {
385  if (faceToEdge[facei] == -1 || collapseRegion[facei] != -1)
386  {
388  << "Problem : crossed into uncollapsed/regionized face"
389  << abort(FatalError);
390  }
391 
392  collapseRegion[facei] = regionI;
393 
394  // Recurse across edges to collapsed neighbours
395 
396  const labelList& fEdges = surf.faceEdges()[facei];
397 
398  forAll(fEdges, fEdgeI)
399  {
400  label edgeI = fEdges[fEdgeI];
401 
402  const labelList& eFaces = surf.edgeFaces()[edgeI];
403 
404  forAll(eFaces, i)
405  {
406  label nbrFacei = eFaces[i];
407 
408  if (faceToEdge[nbrFacei] != -1)
409  {
410  if (collapseRegion[nbrFacei] == -1)
411  {
412  markRegion
413  (
414  surf,
415  faceToEdge,
416  regionI,
417  nbrFacei,
418  collapseRegion
419  );
420  }
421  else if (collapseRegion[nbrFacei] != regionI)
422  {
424  << "Edge:" << edgeI << " between face " << facei
425  << " with region " << regionI
426  << " and face " << nbrFacei
427  << " with region " << collapseRegion[nbrFacei]
428  << endl;
429  }
430  }
431  }
432  }
433 }
434 
435 
436 // Mark every face with region (in collapseRegion) (or -1).
437 // Return number of regions.
438 static label markRegions
439 (
440  const triSurface& surf,
441  const labelList& faceToEdge,
442  labelList& collapseRegion
443 )
444 {
445  label regionI = 0;
446 
447  forAll(faceToEdge, facei)
448  {
449  if (collapseRegion[facei] == -1 && faceToEdge[facei] != -1)
450  {
451  //Pout<< "markRegions : Marking region:" << regionI
452  // << " starting from face " << facei << endl;
453 
454  // Collapsed face. Mark connected region with current region number
455  markRegion(surf, faceToEdge, regionI++, facei, collapseRegion);
456  }
457  }
458  return regionI;
459 }
460 
461 
462 // Type of region.
463 // -1 : edge inbetween uncollapsed faces.
464 // -2 : edge inbetween collapsed faces
465 // >=0 : edge inbetween uncollapsed and collapsed region. Returns region.
466 static label edgeType
467 (
468  const triSurface& surf,
469  const labelList& collapseRegion,
470  const label edgeI
471 )
472 {
473  const labelList& eFaces = surf.edgeFaces()[edgeI];
474 
475  // Detect if edge is inbetween collapseRegion and non-collapse face
476  bool usesUncollapsed = false;
477  label usesRegion = -1;
478 
479  forAll(eFaces, i)
480  {
481  label facei = eFaces[i];
482 
483  label region = collapseRegion[facei];
484 
485  if (region == -1)
486  {
487  usesUncollapsed = true;
488  }
489  else if (usesRegion == -1)
490  {
491  usesRegion = region;
492  }
493  else if (usesRegion != region)
494  {
496  }
497  else
498  {
499  // Equal regions.
500  }
501  }
502 
503  if (usesUncollapsed)
504  {
505  if (usesRegion == -1)
506  {
507  // uncollapsed faces only.
508  return -1;
509  }
510  else
511  {
512  // between uncollapsed and collapsed.
513  return usesRegion;
514  }
515  }
516  else
517  {
518  if (usesRegion == -1)
519  {
521  return -2;
522  }
523  else
524  {
525  return -2;
526  }
527  }
528 }
529 
530 
531 // Get points on outside edge of region (= outside points)
532 static labelListList getOutsideVerts
533 (
534  const triSurface& surf,
535  const labelList& collapseRegion,
536  const label nRegions
537 )
538 {
539  const labelListList& edgeFaces = surf.edgeFaces();
540 
541  // Per region all the outside vertices.
542  labelListList outsideVerts(nRegions);
543 
544  forAll(edgeFaces, edgeI)
545  {
546  // Detect if edge is inbetween collapseRegion and non-collapse face
547  label regionI = edgeType(surf, collapseRegion, edgeI);
548 
549  if (regionI >= 0)
550  {
551  // Edge borders both uncollapsed face and collapsed face on region
552  // usesRegion.
553 
554  const edge& e = surf.edges()[edgeI];
555 
556  labelList& regionVerts = outsideVerts[regionI];
557 
558  // Add both edge points to regionVerts.
559  forAll(e, eI)
560  {
561  label v = e[eI];
562 
563  if (!regionVerts.found(v))
564  {
565  label sz = regionVerts.size();
566  regionVerts.setSize(sz+1);
567  regionVerts[sz] = v;
568  }
569  }
570  }
571  }
572 
573  return outsideVerts;
574 }
575 
576 
577 // n^2 search for furthest removed point pair.
578 static labelPair getSpanPoints
579 (
580  const triSurface& surf,
581  const labelList& outsideVerts
582 )
583 {
584  const pointField& localPoints = surf.localPoints();
585 
586  scalar maxDist = -GREAT;
587  labelPair maxPair;
588 
589  forAll(outsideVerts, i)
590  {
591  label v0 = outsideVerts[i];
592 
593  for (label j = i+1; j < outsideVerts.size(); j++)
594  {
595  label v1 = outsideVerts[j];
596 
597  scalar d = mag(localPoints[v0] - localPoints[v1]);
598 
599  if (d > maxDist)
600  {
601  maxDist = d;
602  maxPair[0] = v0;
603  maxPair[1] = v1;
604  }
605  }
606  }
607 
608  return maxPair;
609 }
610 
611 
612 // Project all non-span points onto the span edge.
613 static void projectNonSpanPoints
614 (
615  const triSurface& surf,
616  const labelList& outsideVerts,
617  const labelPair& spanPair,
618  labelList& sortedVertices,
619  scalarField& sortedWeights
620 )
621 {
622  const point& p0 = surf.localPoints()[spanPair[0]];
623  const point& p1 = surf.localPoints()[spanPair[1]];
624 
625  forAll(outsideVerts, i)
626  {
627  label v = outsideVerts[i];
628 
629  if (v != spanPair[0] && v != spanPair[1])
630  {
631  // Is a non-span point. Project onto spanning edge.
632 
633  pointHit pHit =
634  linePointRef(p0, p1).nearestDist
635  (
636  surf.localPoints()[v]
637  );
638 
639  if (!pHit.hit())
640  {
642  << abort(FatalError);
643  }
644 
645  scalar w = mag(pHit.hitPoint() - p0) / mag(p1 - p0);
646 
647  insertSorted(v, w, sortedVertices, sortedWeights);
648  }
649  }
650 }
651 
652 
653 // Slice part of the orderVertices (and optionally reverse) for this edge.
654 static void getSplitVerts
655 (
656  const triSurface& surf,
657  const label regionI,
658  const labelPair& spanPoints,
659  const labelList& orderedVerts,
660  const scalarField& orderedWeights,
661  const label edgeI,
662 
663  labelList& splitVerts,
664  scalarField& splitWeights
665 )
666 {
667  const edge& e = surf.edges()[edgeI];
668  const label sz = orderedVerts.size();
669 
670  if (e[0] == spanPoints[0])
671  {
672  // Edge in same order as spanPoints&orderedVerts. Keep order.
673 
674  if (e[1] == spanPoints[1])
675  {
676  // Copy all.
677  splitVerts = orderedVerts;
678  splitWeights = orderedWeights;
679  }
680  else
681  {
682  // Copy upto (but not including) e[1]
683  label i1 = orderedVerts.find(e[1]);
684  splitVerts = SubList<label>(orderedVerts, i1, 0);
685  splitWeights = SubList<scalar>(orderedWeights, i1, 0);
686  }
687  }
688  else if (e[0] == spanPoints[1])
689  {
690  // Reverse.
691 
692  if (e[1] == spanPoints[0])
693  {
694  // Copy all.
695  splitVerts = orderedVerts;
696  reverse(splitVerts);
697  splitWeights = orderedWeights;
698  reverse(splitWeights);
699  }
700  else
701  {
702  // Copy downto (but not including) e[1]
703 
704  label i1 = orderedVerts.find(e[1]);
705  splitVerts = SubList<label>(orderedVerts, sz-(i1+1), i1+1);
706  reverse(splitVerts);
707  splitWeights = SubList<scalar>(orderedWeights, sz-(i1+1), i1+1);
708  reverse(splitWeights);
709  }
710  }
711  else if (e[1] == spanPoints[0])
712  {
713  // Reverse.
714 
715  // Copy upto (but not including) e[0]
716 
717  label i0 = orderedVerts.find(e[0]);
718  splitVerts = SubList<label>(orderedVerts, i0, 0);
719  reverse(splitVerts);
720  splitWeights = SubList<scalar>(orderedWeights, i0, 0);
721  reverse(splitWeights);
722  }
723  else if (e[1] == spanPoints[1])
724  {
725  // Copy from (but not including) e[0] to end
726 
727  label i0 = orderedVerts.find(e[0]);
728  splitVerts = SubList<label>(orderedVerts, sz-(i0+1), i0+1);
729  splitWeights = SubList<scalar>(orderedWeights, sz-(i0+1), i0+1);
730  }
731  else
732  {
733  label i0 = orderedVerts.find(e[0]);
734  label i1 = orderedVerts.find(e[1]);
735 
736  if (i0 == -1 || i1 == -1)
737  {
739  << "Did not find edge in projected vertices." << nl
740  << "region:" << regionI << nl
741  << "spanPoints:" << spanPoints
742  << " coords:" << surf.localPoints()[spanPoints[0]]
743  << surf.localPoints()[spanPoints[1]] << nl
744  << "edge:" << edgeI
745  << " verts:" << e
746  << " coords:" << surf.localPoints()[e[0]]
747  << surf.localPoints()[e[1]] << nl
748  << "orderedVerts:" << orderedVerts << nl
749  << abort(FatalError);
750  }
751 
752  if (i0 < i1)
753  {
754  splitVerts = SubList<label>(orderedVerts, i1-i0-1, i0+1);
755  splitWeights = SubList<scalar>(orderedWeights, i1-i0-1, i0+1);
756  }
757  else
758  {
759  splitVerts = SubList<label>(orderedVerts, i0-i1-1, i1+1);
760  reverse(splitVerts);
761  splitWeights = SubList<scalar>(orderedWeights, i0-i1-1, i1+1);
762  reverse(splitWeights);
763  }
764  }
765 }
766 
767 
768 label collapseBase
769 (
770  triSurface& surf,
771  const scalar minLen,
772  const scalar minQuality
773 )
774 {
775  label nTotalSplit = 0;
776 
777  label iter = 0;
778 
779  while (true)
780  {
781  // Detect faces to collapse
782  // ~~~~~~~~~~~~~~~~~~~~~~~~
783 
784  // -1 or edge the face is collapsed onto.
785  labelList faceToEdge(surf.size(), -1);
786 
787  // Calculate faceToEdge (face collapses)
788  markCollapsedFaces(surf, minLen, minQuality, faceToEdge);
789 
790 
791  // Find regions of connected collapsed faces
792  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
793 
794  // per face -1 or region
795  labelList collapseRegion(surf.size(), -1);
796 
797  label nRegions = markRegions(surf, faceToEdge, collapseRegion);
798 
799  //Pout<< "Detected " << nRegions << " regions of faces to be collapsed"
800  // << nl << endl;
801 
802  // Pick up all vertices on outside of region
803  labelListList outsideVerts
804  (
805  getOutsideVerts(surf, collapseRegion, nRegions)
806  );
807 
808  // For all regions determine maximum distance between points
809  List<labelPair> spanPoints(nRegions);
810  labelListList orderedVertices(nRegions);
811  List<scalarField> orderedWeights(nRegions);
812 
813  forAll(spanPoints, regionI)
814  {
815  spanPoints[regionI] = getSpanPoints(surf, outsideVerts[regionI]);
816 
817  //Pout<< "For region " << regionI << " found extrema at points "
818  // << surf.localPoints()[spanPoints[regionI][0]]
819  // << surf.localPoints()[spanPoints[regionI][1]]
820  // << endl;
821 
822  // Project all non-span points onto the span edge.
823  projectNonSpanPoints
824  (
825  surf,
826  outsideVerts[regionI],
827  spanPoints[regionI],
828  orderedVertices[regionI],
829  orderedWeights[regionI]
830  );
831 
832  //Pout<< "For region:" << regionI
833  // << " span:" << spanPoints[regionI]
834  // << " orderedVerts:" << orderedVertices[regionI]
835  // << " orderedWeights:" << orderedWeights[regionI]
836  // << endl;
837 
838  //writeRegionOBJ
839  //(
840  // surf,
841  // regionI,
842  // collapseRegion,
843  // outsideVerts[regionI]
844  //);
845 
846  //Pout<< endl;
847  }
848 
849 
850 
851  // Actually split the edges
852  // ~~~~~~~~~~~~~~~~~~~~~~~~
853 
854 
855  const List<labelledTri>& localFaces = surf.localFaces();
856  const edgeList& edges = surf.edges();
857 
858  label nSplit = 0;
859 
860  // Storage for new triangles.
861  DynamicList<labelledTri> newTris(surf.size());
862 
863  // Whether face has been dealt with (either copied/split or deleted)
864  boolList faceHandled(surf.size(), false);
865 
866 
867  forAll(edges, edgeI)
868  {
869  const edge& e = edges[edgeI];
870 
871  // Detect if edge is inbetween collapseRegion and non-collapse face
872  label regionI = edgeType(surf, collapseRegion, edgeI);
873 
874  if (regionI == -2)
875  {
876  // inbetween collapsed faces. nothing needs to be done.
877  }
878  else if (regionI == -1)
879  {
880  // edge inbetween uncollapsed faces. Handle these later on.
881  }
882  else
883  {
884  // some faces around edge are collapsed.
885 
886  // Find additional set of points on edge to be used to split
887  // the remaining faces.
888 
889  labelList splitVerts;
890  scalarField splitWeights;
891  getSplitVerts
892  (
893  surf,
894  regionI,
895  spanPoints[regionI],
896  orderedVertices[regionI],
897  orderedWeights[regionI],
898  edgeI,
899 
900  splitVerts,
901  splitWeights
902  );
903 
904  if (splitVerts.size())
905  {
906  // Split edge using splitVerts. All non-collapsed triangles
907  // using edge will get split.
908 
909  //{
910  // const pointField& localPoints = surf.localPoints();
911  // Pout<< "edge " << edgeI << ' ' << e
912  // << " points "
913  // << localPoints[e[0]] << ' ' << localPoints[e[1]]
914  // << " split into edges with extra points:"
915  // << endl;
916  // forAll(splitVerts, i)
917  // {
918  // Pout<< " " << splitVerts[i] << " weight "
919  // << splitWeights[i] << nl;
920  // }
921  //}
922 
923  const labelList& eFaces = surf.edgeFaces()[edgeI];
924 
925  forAll(eFaces, i)
926  {
927  label facei = eFaces[i];
928 
929  if (!faceHandled[facei] && faceToEdge[facei] == -1)
930  {
931  // Split face to use vertices.
932  splitTri
933  (
934  localFaces[facei],
935  e,
936  splitVerts,
937  newTris
938  );
939 
940  faceHandled[facei] = true;
941 
942  nSplit++;
943  }
944  }
945  }
946  }
947  }
948 
949  // Copy all unsplit faces
950  forAll(faceHandled, facei)
951  {
952  if (!faceHandled[facei] && faceToEdge[facei] == -1)
953  {
954  newTris.append(localFaces[facei]);
955  }
956  }
957 
958  Info<< "collapseBase : collapsing " << nSplit
959  << " triangles by splitting their base edge."
960  << endl;
961 
962  nTotalSplit += nSplit;
963 
964  if (nSplit == 0)
965  {
966  break;
967  }
968 
969  // Pack the triangles
970  newTris.shrink();
971 
972  //Pout<< "Resetting surface from " << surf.size() << " to "
973  // << newTris.size() << " triangles" << endl;
974  surf = triSurface(newTris, surf.patches(), surf.localPoints());
975 
976  //{
977  // fileName fName("bla" + name(iter) + ".obj");
978  // Pout<< "Writing surf to " << fName << endl;
979  // surf.write(fName);
980  //}
981 
982  iter++;
983  }
984 
985  // Remove any unused vertices
986  surf = triSurface(surf.localFaces(), surf.patches(), surf.localPoints());
987 
988  return nTotalSplit;
989 }
990 
991 
992 // ************************************************************************* //
Foam::PrimitivePatch::points
const Field< point_type > & points() const noexcept
Definition: PrimitivePatch.H:295
Foam::reverse
void reverse(UList< T > &list, const label n)
Definition: UListI.H:442
meshTools.H
OSspecific.H
Functions used by OpenFOAM that are specific to POSIX compliant operating systems and need to be repl...
SubList.H
Foam::PrimitivePatch::edgeFaces
const labelListList & edgeFaces() const
Definition: PrimitivePatch.C:255
Foam::PrimitivePatch::edges
const edgeList & edges() const
Definition: PrimitivePatch.C:176
Foam::PointHit
Describes the interaction of a face and a point. It carries the info of a successful hit and (if succ...
Definition: PointHit.H:47
Foam::DynamicList
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:51
Foam::SubList
A List obtained as a section of another List.
Definition: SubList.H:50
Foam::edge
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:59
Foam::PointHit::distance
scalar distance() const noexcept
Definition: PointHit.H:133
Foam::PointHit::hitPoint
const point_type & hitPoint() const
Definition: PointHit.H:139
Foam::endl
Ostream & endl(Ostream &os)
Definition: Ostream.H:381
Foam::PrimitivePatch::faceEdges
const labelListList & faceEdges() const
Definition: PrimitivePatch.C:268
Foam::triSurface::patches
const geometricSurfacePatchList & patches() const noexcept
Definition: triSurface.H:395
forAll
#define forAll(list, i)
Definition: stdFoam.H:349
OFstream.H
Foam::Field
Generic templated field type.
Definition: Field.H:59
Foam::triSurface
Triangulated surface description with patch information.
Definition: triSurface.H:72
Foam::Info
messageStream Info
Foam::List::setSize
void setSize(const label n)
Definition: List.H:218
Foam::DynamicList::append
DynamicList< T, SizeMin > & append(const T &val)
Definition: DynamicListI.H:504
argList.H
Foam::FatalError
error FatalError
Foam::PrimitivePatch::localFaces
const List< face_type > & localFaces() const
Definition: PrimitivePatch.C:310
Foam
Definition: atmBoundaryLayer.C:26
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:139
Foam::triSurfaceTools::oppositeVertex
static label oppositeVertex(const triSurface &surf, const label facei, const label edgeI)
Definition: triSurfaceTools.C:1433
Foam::PrimitivePatch::localPoints
const Field< point_type > & localPoints() const
Definition: PrimitivePatch.C:352
collapseBase.H
Routines collapse sliver triangles by splitting the base edge.
Foam::linePointRef
line< point, const point & > linePointRef
A line using referred points.
Definition: linePointRef.H:41
collapseBase
label collapseBase(triSurface &surf, const scalar minLen, const scalar minQuality)
FatalErrorInFunction
#define FatalErrorInFunction
Definition: error.H:465
Foam::nl
constexpr char nl
Definition: Ostream.H:424
Foam::Pair
An ordered pair of two objects of type <T> with first() and second() elements.
Definition: Pair.H:50
f
labelList f(nPoints)
Foam::Vector< scalar >
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: BitOps.H:58
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::constant::electromagnetic::e
const dimensionedScalar e
Definition: createFields.H:11
Foam::labelledTri
A triFace with additional (region) index.
Definition: labelledTri.H:53
p0
const volScalarField & p0
Definition: EEqn.H:36
Foam::PointHit::hit
bool hit() const noexcept
Definition: PointHit.H:115
triSurfaceTools.H
WarningInFunction
#define WarningInFunction
Definition: messageStream.H:353
labelPair.H