1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 2011 OpenFOAM Foundation
7 -------------------------------------------------------------------------------
9 This file is part of OpenFOAM.
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.
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
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/>.
24 \*---------------------------------------------------------------------------*/
26 #include "collapseBase.H"
27 #include "triSurfaceTools.H"
31 #include "labelPair.H"
32 #include "meshTools.H"
33 #include "OSspecific.H"
37 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
39 // Dump collapse region to .obj file
40 static void writeRegionOBJ
42 const triSurface& surf,
44 const labelList& collapseRegion,
45 const labelList& outsideVerts
48 fileName dir("regions");
51 fileName regionName(dir / "region_" + name(regionI) + ".obj");
53 Pout<< "Dumping region " << regionI << " to file " << regionName << endl;
55 boolList include(surf.size(), false);
57 forAll(collapseRegion, faceI)
59 if (collapseRegion[faceI] == regionI)
61 include[faceI] = true;
65 labelList pointMap, faceMap;
67 triSurface regionSurf(surf.subsetMesh(include, pointMap, faceMap));
69 //Pout<< "Region " << regionI << " surface:" << nl;
70 //regionSurf.writeStats(Pout);
72 regionSurf.write(regionName);
75 // Dump corresponding outside vertices.
76 fileName pointsName(dir / "regionPoints_" + name(regionI) + ".obj");
78 Pout<< "Dumping region " << regionI << " points to file " << pointsName
81 OFstream str(pointsName);
83 forAll(outsideVerts, i)
85 meshTools::writeOBJ(str, surf.localPoints()[outsideVerts[i]]);
90 // Split triangle into multiple triangles because edge e being split
91 // into multiple edges.
96 const labelList& splitPoints,
97 DynamicList<labelledTri>& tris
100 label oldNTris = tris.size();
102 label fp = findIndex(f, e[0]);
103 label fp1 = f.fcIndex(fp);
104 label fp2 = f.fcIndex(fp1);
108 // Split triangle along fp to fp1
109 tris.append(labelledTri(f[fp2], f[fp], splitPoints[0], f.region()));
111 for (label i = 1; i < splitPoints.size(); i++)
136 else if (f[fp2] == e[1])
138 // Split triangle along fp2 to fp. (Reverse order of splitPoints)
151 for (label i = splitPoints.size()-1; i > 0; --i)
178 FatalErrorIn("splitTri")
179 << "Edge " << e << " not part of triangle " << f
183 << abort(FatalError);
186 Pout<< "Split face " << f << " along edge " << e
187 << " into triangles:" << endl;
189 for (label i = oldNTris; i < tris.size(); i++)
191 Pout<< " " << tris[i] << nl;
196 // Insert scalar into sortedVerts/sortedWeights so the weights are in
197 // incrementing order.
198 static bool insertSorted
203 labelList& sortedVerts,
204 scalarField& sortedWeights
207 if (findIndex(sortedVerts, vertI) != -1)
209 FatalErrorIn("insertSorted") << "Trying to insert vertex " << vertI
210 << " which is already in list of sorted vertices "
211 << sortedVerts << abort(FatalError);
214 if (weight <= 0 || weight >= 1)
216 FatalErrorIn("insertSorted") << "Trying to insert vertex " << vertI
217 << " with illegal weight " << weight
218 << " into list of sorted vertices "
219 << sortedVerts << abort(FatalError);
223 label insertI = sortedVerts.size();
225 forAll(sortedVerts, sortedI)
227 scalar w = sortedWeights[sortedI];
229 if (mag(w - weight) < SMALL)
231 WarningIn("insertSorted")
232 << "Trying to insert weight " << weight << " which is close to"
233 << " existing weight " << w << " in " << sortedWeights
239 // Start inserting before sortedI.
246 label sz = sortedWeights.size();
248 sortedWeights.setSize(sz + 1);
249 sortedVerts.setSize(sz + 1);
251 // Leave everything up to (not including) insertI intact.
253 // Make space by copying from insertI up.
254 for (label i = sz-1; i >= insertI; --i)
256 sortedWeights[i+1] = sortedWeights[i];
257 sortedVerts[i+1] = sortedVerts[i];
259 sortedWeights[insertI] = weight;
260 sortedVerts[insertI] = vertI;
266 // Mark all faces that are going to be collapsed.
267 // faceToEdge: per face -1 or the base edge of the face.
268 static void markCollapsedFaces
270 const triSurface& surf,
272 labelList& faceToEdge
275 faceToEdge.setSize(surf.size());
278 const pointField& localPoints = surf.localPoints();
279 const labelListList& edgeFaces = surf.edgeFaces();
281 forAll(edgeFaces, edgeI)
283 const edge& e = surf.edges()[edgeI];
285 const labelList& eFaces = surf.edgeFaces()[edgeI];
289 label faceI = eFaces[i];
291 // Check distance of vertex to edge.
293 triSurfaceTools::oppositeVertex
301 e.line(localPoints).nearestDist
303 localPoints[opposite0]
306 if (pHit.hit() && pHit.distance() < minLen)
308 // Remove faceI and split all other faces using this
309 // edge. This is done by 'replacing' the edgeI with the
311 Pout<< "Splitting face " << faceI << " since distance "
313 << " from vertex " << opposite0
314 << " to edge " << edgeI
318 << " is too small" << endl;
320 // Mark face as being collapsed
321 if (faceToEdge[faceI] != -1)
323 FatalErrorIn("markCollapsedFaces")
324 << "Cannot collapse face " << faceI << " since "
325 << " is marked to be collapsed both to edge "
326 << faceToEdge[faceI] << " and " << edgeI
327 << abort(FatalError);
330 faceToEdge[faceI] = edgeI;
337 // Recurse through collapsed faces marking all of them with regionI (in
339 static void markRegion
341 const triSurface& surf,
342 const labelList& faceToEdge,
345 labelList& collapseRegion
348 if (faceToEdge[faceI] == -1 || collapseRegion[faceI] != -1)
350 FatalErrorIn("markRegion")
351 << "Problem : crossed into uncollapsed/regionized face"
352 << abort(FatalError);
355 collapseRegion[faceI] = regionI;
357 // Recurse across edges to collapsed neighbours
359 const labelList& fEdges = surf.faceEdges()[faceI];
361 forAll(fEdges, fEdgeI)
363 label edgeI = fEdges[fEdgeI];
365 const labelList& eFaces = surf.edgeFaces()[edgeI];
369 label nbrFaceI = eFaces[i];
371 if (faceToEdge[nbrFaceI] != -1)
373 if (collapseRegion[nbrFaceI] == -1)
384 else if (collapseRegion[nbrFaceI] != regionI)
386 FatalErrorIn("markRegion")
387 << "Edge:" << edgeI << " between face " << faceI
388 << " with region " << regionI
389 << " and face " << nbrFaceI
390 << " with region " << collapseRegion[nbrFaceI]
399 // Mark every face with region (in collapseRegion) (or -1).
400 // Return number of regions.
401 static label markRegions
403 const triSurface& surf,
404 const labelList& faceToEdge,
405 labelList& collapseRegion
410 forAll(faceToEdge, faceI)
412 if (collapseRegion[faceI] == -1 && faceToEdge[faceI] != -1)
414 Pout<< "markRegions : Marking region:" << regionI
415 << " starting from face " << faceI << endl;
417 // Collapsed face. Mark connected region with current region number
418 markRegion(surf, faceToEdge, regionI++, faceI, collapseRegion);
426 // -1 : edge inbetween uncollapsed faces.
427 // -2 : edge inbetween collapsed faces
428 // >=0 : edge inbetween uncollapsed and collapsed region. Returns region.
429 static label edgeType
431 const triSurface& surf,
432 const labelList& collapseRegion,
436 const labelList& eFaces = surf.edgeFaces()[edgeI];
438 // Detect if edge is inbetween collapseRegion and non-collapse face
439 bool usesUncollapsed = false;
440 label usesRegion = -1;
444 label faceI = eFaces[i];
446 label region = collapseRegion[faceI];
450 usesUncollapsed = true;
452 else if (usesRegion == -1)
456 else if (usesRegion != region)
458 FatalErrorIn("edgeRegion") << abort(FatalError);
468 if (usesRegion == -1)
470 // uncollapsed faces only.
475 // between uncollapsed and collapsed.
481 if (usesRegion == -1)
483 FatalErrorIn("edgeRegion") << abort(FatalError);
494 // Get points on outside edge of region (= outside points)
495 static labelListList getOutsideVerts
497 const triSurface& surf,
498 const labelList& collapseRegion,
502 const labelListList& edgeFaces = surf.edgeFaces();
504 // Per region all the outside vertices.
505 labelListList outsideVerts(nRegions);
507 forAll(edgeFaces, edgeI)
509 // Detect if edge is inbetween collapseRegion and non-collapse face
510 label regionI = edgeType(surf, collapseRegion, edgeI);
514 // Edge borders both uncollapsed face and collapsed face on region
517 const edge& e = surf.edges()[edgeI];
519 labelList& regionVerts = outsideVerts[regionI];
521 // Add both edge points to regionVerts.
526 if (findIndex(regionVerts, v) == -1)
528 label sz = regionVerts.size();
529 regionVerts.setSize(sz+1);
540 // n^2 search for furthest removed point pair.
541 static labelPair getSpanPoints
543 const triSurface& surf,
544 const labelList& outsideVerts
547 const pointField& localPoints = surf.localPoints();
549 scalar maxDist = -GREAT;
552 forAll(outsideVerts, i)
554 label v0 = outsideVerts[i];
556 for (label j = i+1; j < outsideVerts.size(); j++)
558 label v1 = outsideVerts[j];
560 scalar d = mag(localPoints[v0] - localPoints[v1]);
575 // Project all non-span points onto the span edge.
576 static void projectNonSpanPoints
578 const triSurface& surf,
579 const labelList& outsideVerts,
580 const labelPair& spanPair,
581 labelList& sortedVertices,
582 scalarField& sortedWeights
585 const point& p0 = surf.localPoints()[spanPair[0]];
586 const point& p1 = surf.localPoints()[spanPair[1]];
588 forAll(outsideVerts, i)
590 label v = outsideVerts[i];
592 if (v != spanPair[0] && v != spanPair[1])
594 // Is a non-span point. Project onto spanning edge.
597 linePointRef(p0, p1).nearestDist
599 surf.localPoints()[v]
604 FatalErrorIn("projectNonSpanPoints")
605 << abort(FatalError);
608 scalar w = mag(pHit.hitPoint() - p0) / mag(p1 - p0);
610 insertSorted(v, w, sortedVertices, sortedWeights);
616 // Slice part of the orderVertices (and optionally reverse) for this edge.
617 static void getSplitVerts
619 const triSurface& surf,
621 const labelPair& spanPoints,
622 const labelList& orderedVerts,
623 const scalarField& orderedWeights,
626 labelList& splitVerts,
627 scalarField& splitWeights
630 const edge& e = surf.edges()[edgeI];
631 const label sz = orderedVerts.size();
633 if (e[0] == spanPoints[0])
635 // Edge in same order as spanPoints&orderedVerts. Keep order.
637 if (e[1] == spanPoints[1])
640 splitVerts = orderedVerts;
641 splitWeights = orderedWeights;
645 // Copy upto (but not including) e[1]
646 label i1 = findIndex(orderedVerts, e[1]);
647 splitVerts = SubList<label>(orderedVerts, i1, 0);
648 splitWeights = SubList<scalar>(orderedWeights, i1, 0);
651 else if (e[0] == spanPoints[1])
655 if (e[1] == spanPoints[0])
658 splitVerts = orderedVerts;
660 splitWeights = orderedWeights;
661 reverse(splitWeights);
665 // Copy downto (but not including) e[1]
667 label i1 = findIndex(orderedVerts, e[1]);
668 splitVerts = SubList<label>(orderedVerts, sz-(i1+1), i1+1);
670 splitWeights = SubList<scalar>(orderedWeights, sz-(i1+1), i1+1);
671 reverse(splitWeights);
674 else if (e[1] == spanPoints[0])
678 // Copy upto (but not including) e[0]
680 label i0 = findIndex(orderedVerts, e[0]);
681 splitVerts = SubList<label>(orderedVerts, i0, 0);
683 splitWeights = SubList<scalar>(orderedWeights, i0, 0);
684 reverse(splitWeights);
686 else if (e[1] == spanPoints[1])
688 // Copy from (but not including) e[0] to end
690 label i0 = findIndex(orderedVerts, e[0]);
691 splitVerts = SubList<label>(orderedVerts, sz-(i0+1), i0+1);
692 splitWeights = SubList<scalar>(orderedWeights, sz-(i0+1), i0+1);
696 label i0 = findIndex(orderedVerts, e[0]);
697 label i1 = findIndex(orderedVerts, e[1]);
699 if (i0 == -1 || i1 == -1)
701 FatalErrorIn("getSplitVerts")
702 << "Did not find edge in projected vertices." << nl
703 << "region:" << regionI << nl
704 << "spanPoints:" << spanPoints
705 << " coords:" << surf.localPoints()[spanPoints[0]]
706 << surf.localPoints()[spanPoints[1]] << nl
709 << " coords:" << surf.localPoints()[e[0]]
710 << surf.localPoints()[e[1]] << nl
711 << "orderedVerts:" << orderedVerts << nl
712 << abort(FatalError);
717 splitVerts = SubList<label>(orderedVerts, i1-i0-1, i0+1);
718 splitWeights = SubList<scalar>(orderedWeights, i1-i0-1, i0+1);
722 splitVerts = SubList<label>(orderedVerts, i0-i1-1, i1+1);
724 splitWeights = SubList<scalar>(orderedWeights, i0-i1-1, i1+1);
725 reverse(splitWeights);
731 label collapseBase(triSurface& surf, const scalar minLen)
733 label nTotalSplit = 0;
739 // Detect faces to collapse
740 // ~~~~~~~~~~~~~~~~~~~~~~~~
742 // -1 or edge the face is collapsed onto.
743 labelList faceToEdge(surf.size(), -1);
745 // Calculate faceToEdge (face collapses)
746 markCollapsedFaces(surf, minLen, faceToEdge);
749 // Find regions of connected collapsed faces
750 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
752 // per face -1 or region
753 labelList collapseRegion(surf.size(), -1);
755 label nRegions = markRegions(surf, faceToEdge, collapseRegion);
757 Pout<< "Detected " << nRegions << " regions of faces to be collapsed"
760 // Pick up all vertices on outside of region
761 labelListList outsideVerts
763 getOutsideVerts(surf, collapseRegion, nRegions)
766 // For all regions determine maximum distance between points
767 List<labelPair> spanPoints(nRegions);
768 labelListList orderedVertices(nRegions);
769 List<scalarField> orderedWeights(nRegions);
771 forAll(spanPoints, regionI)
773 spanPoints[regionI] = getSpanPoints(surf, outsideVerts[regionI]);
775 Pout<< "For region " << regionI << " found extrema at points "
776 << surf.localPoints()[spanPoints[regionI][0]]
777 << surf.localPoints()[spanPoints[regionI][1]]
780 // Project all non-span points onto the span edge.
784 outsideVerts[regionI],
786 orderedVertices[regionI],
787 orderedWeights[regionI]
790 Pout<< "For region:" << regionI
791 << " span:" << spanPoints[regionI]
792 << " orderedVerts:" << orderedVertices[regionI]
793 << " orderedWeights:" << orderedWeights[regionI]
801 outsideVerts[regionI]
809 // Actually split the edges
810 // ~~~~~~~~~~~~~~~~~~~~~~~~
813 const List<labelledTri>& localFaces = surf.localFaces();
814 const edgeList& edges = surf.edges();
818 // Storage for new triangles.
819 DynamicList<labelledTri> newTris(surf.size());
821 // Whether face has been dealt with (either copied/split or deleted)
822 boolList faceHandled(surf.size(), false);
827 const edge& e = edges[edgeI];
829 // Detect if edge is inbetween collapseRegion and non-collapse face
830 label regionI = edgeType(surf, collapseRegion, edgeI);
834 // inbetween collapsed faces. nothing needs to be done.
836 else if (regionI == -1)
838 // edge inbetween uncollapsed faces. Handle these later on.
842 // some faces around edge are collapsed.
844 // Find additional set of points on edge to be used to split
845 // the remaining faces.
847 labelList splitVerts;
848 scalarField splitWeights;
854 orderedVertices[regionI],
855 orderedWeights[regionI],
862 if (splitVerts.size())
864 // Split edge using splitVerts. All non-collapsed triangles
865 // using edge will get split.
869 const pointField& localPoints = surf.localPoints();
870 Pout<< "edge " << edgeI << ' ' << e
872 << localPoints[e[0]] << ' ' << localPoints[e[1]]
873 << " split into edges with extra points:"
875 forAll(splitVerts, i)
877 Pout<< " " << splitVerts[i] << " weight "
878 << splitWeights[i] << nl;
882 const labelList& eFaces = surf.edgeFaces()[edgeI];
886 label faceI = eFaces[i];
888 if (!faceHandled[faceI] && faceToEdge[faceI] == -1)
890 // Split face to use vertices.
899 faceHandled[faceI] = true;
908 // Copy all unsplit faces
909 forAll(faceHandled, faceI)
911 if (!faceHandled[faceI] && faceToEdge[faceI] == -1)
913 newTris.append(localFaces[faceI]);
917 Pout<< "collapseBase : splitting " << nSplit << " triangles"
920 nTotalSplit += nSplit;
927 // Pack the triangles
930 Pout<< "Resetting surface from " << surf.size() << " to "
931 << newTris.size() << " triangles" << endl;
932 surf = triSurface(newTris, surf.patches(), surf.localPoints());
935 fileName fName("bla" + name(iter) + ".obj");
936 Pout<< "Writing surf to " << fName << endl;
943 // Remove any unused vertices
944 surf = triSurface(surf.localFaces(), surf.patches(), surf.localPoints());
950 // ************************************************************************* //