FreeFOAM The Cross-Platform CFD Toolkit
Hosted by SourceForge:
Get FreeFOAM at SourceForge.net.
            Fast, secure and Free Open Source software downloads

collapseBase.C

Go to the documentation of this file.
00001 /*---------------------------------------------------------------------------*\
00002   =========                 |
00003   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
00004    \\    /   O peration     |
00005     \\  /    A nd           | Copyright (C) 1991-2010 OpenCFD Ltd.
00006      \\/     M anipulation  |
00007 -------------------------------------------------------------------------------
00008 License
00009     This file is part of OpenFOAM.
00010 
00011     OpenFOAM is free software: you can redistribute it and/or modify it
00012     under the terms of the GNU General Public License as published by
00013     the Free Software Foundation, either version 3 of the License, or
00014     (at your option) any later version.
00015 
00016     OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
00017     ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
00018     FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
00019     for more details.
00020 
00021     You should have received a copy of the GNU General Public License
00022     along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.
00023 
00024 \*---------------------------------------------------------------------------*/
00025 
00026 #include "collapseBase.H"
00027 #include <meshTools/triSurfaceTools.H>
00028 #include <OpenFOAM/argList.H>
00029 #include <OpenFOAM/OFstream.H>
00030 #include <OpenFOAM/SubList.H>
00031 #include <OpenFOAM/labelPair.H>
00032 #include <meshTools/meshTools.H>
00033 #include <OpenFOAM/OSspecific.H>
00034 
00035 using namespace Foam;
00036 
00037 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00038 
00039 // Dump collapse region to .obj file
00040 static void writeRegionOBJ
00041 (
00042     const triSurface& surf,
00043     const label regionI,
00044     const labelList& collapseRegion,
00045     const labelList& outsideVerts
00046 )
00047 {
00048     fileName dir("regions");
00049 
00050     mkDir(dir);
00051     fileName regionName(dir / "region_" + name(regionI) + ".obj");
00052 
00053     Pout<< "Dumping region " << regionI << " to file " << regionName << endl;
00054 
00055     boolList include(surf.size(), false);
00056 
00057     forAll(collapseRegion, faceI)
00058     {
00059         if (collapseRegion[faceI] == regionI)
00060         {
00061             include[faceI] = true;
00062         }
00063     }
00064 
00065     labelList pointMap, faceMap;
00066 
00067     triSurface regionSurf(surf.subsetMesh(include, pointMap, faceMap));
00068 
00069     //Pout<< "Region " << regionI << " surface:" << nl;
00070     //regionSurf.writeStats(Pout);
00071 
00072     regionSurf.write(regionName);
00073 
00074 
00075     // Dump corresponding outside vertices.
00076     fileName pointsName(dir / "regionPoints_" + name(regionI) + ".obj");
00077 
00078     Pout<< "Dumping region " << regionI << " points to file " << pointsName
00079         << endl;
00080 
00081     OFstream str(pointsName);
00082 
00083     forAll(outsideVerts, i)
00084     {
00085         meshTools::writeOBJ(str, surf.localPoints()[outsideVerts[i]]);
00086     }
00087 }
00088 
00089 
00090 // Split triangle into multiple triangles because edge e being split
00091 // into multiple edges.
00092 static void splitTri
00093 (
00094     const labelledTri& f,
00095     const edge& e,
00096     const labelList& splitPoints,
00097     DynamicList<labelledTri>& tris
00098 )
00099 {
00100     label oldNTris = tris.size();
00101 
00102     label fp = findIndex(f, e[0]);
00103     label fp1 = (fp+1)%3;
00104     label fp2 = (fp1+1)%3;
00105 
00106     if (f[fp1] == e[1])
00107     {
00108         // Split triangle along fp to fp1
00109         tris.append(labelledTri(f[fp2], f[fp], splitPoints[0], f.region()));
00110 
00111         for (label i = 1; i < splitPoints.size(); i++)
00112         {
00113             tris.append
00114             (
00115                 labelledTri
00116                 (
00117                     f[fp2],
00118                     splitPoints[i-1],
00119                     splitPoints[i],
00120                     f.region()
00121                 )
00122             );
00123         }
00124 
00125         tris.append
00126         (
00127             labelledTri
00128             (
00129                 f[fp2],
00130                 splitPoints[splitPoints.size()-1],
00131                 f[fp1],
00132                 f.region()
00133             )
00134         );
00135     }
00136     else if (f[fp2] == e[1])
00137     {
00138         // Split triangle along fp2 to fp. (Reverse order of splitPoints)
00139 
00140         tris.append
00141         (
00142             labelledTri
00143             (
00144                 f[fp1],
00145                 f[fp2],
00146                 splitPoints[splitPoints.size()-1],
00147                 f.region()
00148             )
00149         );
00150 
00151         for (label i = splitPoints.size()-1; i > 0; --i)
00152         {
00153             tris.append
00154             (
00155                 labelledTri
00156                 (
00157                     f[fp1],
00158                     splitPoints[i],
00159                     splitPoints[i-1],
00160                     f.region()
00161                 )
00162             );
00163         }
00164 
00165         tris.append
00166         (
00167             labelledTri
00168             (
00169                 f[fp1],
00170                 splitPoints[0],
00171                 f[fp],
00172                 f.region()
00173             )
00174         );
00175     }
00176     else
00177     {
00178         FatalErrorIn("splitTri")
00179             << "Edge " << e << " not part of triangle " << f
00180             << " fp:" << fp
00181             << " fp1:" << fp1
00182             << " fp2:" << fp2
00183             << abort(FatalError);
00184     }
00185 
00186     Pout<< "Split face " << f << " along edge " << e
00187         << " into triangles:" << endl;
00188 
00189     for (label i = oldNTris; i < tris.size(); i++)
00190     {
00191         Pout<< "   " << tris[i] << nl;
00192     }
00193 }
00194 
00195 
00196 // Insert scalar into sortedVerts/sortedWeights so the weights are in
00197 // incrementing order.
00198 static bool insertSorted
00199 (   
00200     const label vertI,
00201     const scalar weight,
00202 
00203     labelList& sortedVerts,
00204     scalarField& sortedWeights
00205 )
00206 {
00207     if (findIndex(sortedVerts, vertI) != -1)
00208     {
00209         FatalErrorIn("insertSorted") << "Trying to insert vertex " << vertI
00210             << " which is already in list of sorted vertices "
00211             << sortedVerts << abort(FatalError);
00212     }
00213 
00214     if (weight <= 0 || weight >= 1)
00215     {
00216         FatalErrorIn("insertSorted") << "Trying to insert vertex " << vertI
00217             << " with illegal weight " << weight
00218             << " into list of sorted vertices "
00219             << sortedVerts << abort(FatalError);
00220     }
00221 
00222 
00223     label insertI = sortedVerts.size();
00224 
00225     forAll(sortedVerts, sortedI)
00226     {
00227         scalar w = sortedWeights[sortedI];
00228 
00229         if (mag(w - weight) < SMALL)
00230         {
00231             WarningIn("insertSorted")
00232                 << "Trying to insert weight " << weight << " which is close to"
00233                 << " existing weight " << w << " in " << sortedWeights
00234                 << endl;
00235         }
00236 
00237         if (w > weight)
00238         {
00239             // Start inserting before sortedI.
00240             insertI = sortedI;
00241             break;
00242         }
00243     }
00244 
00245     
00246     label sz = sortedWeights.size();
00247 
00248     sortedWeights.setSize(sz + 1);
00249     sortedVerts.setSize(sz + 1);
00250 
00251     // Leave everything up to (not including) insertI intact.
00252 
00253     // Make space by copying from insertI up.
00254     for (label i = sz-1; i >= insertI; --i)
00255     {
00256         sortedWeights[i+1] = sortedWeights[i];
00257         sortedVerts[i+1] = sortedVerts[i];
00258     }
00259     sortedWeights[insertI] = weight;
00260     sortedVerts[insertI] = vertI;
00261 
00262     return true;
00263 }
00264 
00265 
00266 // Mark all faces that are going to be collapsed.
00267 // faceToEdge: per face -1 or the base edge of the face.
00268 static void markCollapsedFaces
00269 (
00270     const triSurface& surf,
00271     const scalar minLen,
00272     labelList& faceToEdge
00273 )
00274 {
00275     faceToEdge.setSize(surf.size());
00276     faceToEdge = -1;
00277 
00278     const pointField& localPoints = surf.localPoints();
00279     const labelListList& edgeFaces = surf.edgeFaces();
00280 
00281     forAll(edgeFaces, edgeI)
00282     {
00283         const edge& e = surf.edges()[edgeI];
00284 
00285         const labelList& eFaces = surf.edgeFaces()[edgeI];
00286 
00287         forAll(eFaces, i)
00288         {
00289             label faceI = eFaces[i];
00290 
00291             // Check distance of vertex to edge.
00292             label opposite0 =
00293                 triSurfaceTools::oppositeVertex
00294                 (   
00295                     surf,
00296                     faceI,
00297                     edgeI
00298                 );
00299 
00300             pointHit pHit =
00301                 e.line(localPoints).nearestDist
00302                 (
00303                     localPoints[opposite0]
00304                 );
00305 
00306             if (pHit.hit() && pHit.distance() < minLen)
00307             {
00308                 // Remove faceI and split all other faces using this
00309                 // edge. This is done by 'replacing' the edgeI with the
00310                 // opposite0 vertex
00311                 Pout<< "Splitting face " << faceI << " since distance "
00312                     << pHit.distance()
00313                     << " from vertex " << opposite0
00314                     << " to edge " << edgeI
00315                     << "  points "
00316                     << localPoints[e[0]]
00317                     << localPoints[e[1]]
00318                     << " is too small" << endl;
00319 
00320                 // Mark face as being collapsed
00321                 if (faceToEdge[faceI] != -1)
00322                 {
00323                     FatalErrorIn("markCollapsedFaces")
00324                         << "Cannot collapse face " << faceI << " since "
00325                         << " is marked to be collapsed both to edge "
00326                         << faceToEdge[faceI] << " and " << edgeI
00327                         << abort(FatalError);
00328                 }
00329 
00330                 faceToEdge[faceI] = edgeI;
00331             }
00332         }
00333     }
00334 }
00335 
00336 
00337 // Recurse through collapsed faces marking all of them with regionI (in
00338 // collapseRegion)
00339 static void markRegion
00340 (
00341     const triSurface& surf,
00342     const labelList& faceToEdge,
00343     const label regionI,
00344     const label faceI,
00345     labelList& collapseRegion
00346 )
00347 {
00348     if (faceToEdge[faceI] == -1 || collapseRegion[faceI] != -1)
00349     {
00350         FatalErrorIn("markRegion")
00351             << "Problem : crossed into uncollapsed/regionized face"
00352             << abort(FatalError);
00353     }
00354 
00355     collapseRegion[faceI] = regionI;
00356 
00357     // Recurse across edges to collapsed neighbours
00358 
00359     const labelList& fEdges = surf.faceEdges()[faceI];
00360 
00361     forAll(fEdges, fEdgeI)
00362     {
00363         label edgeI = fEdges[fEdgeI];
00364 
00365         const labelList& eFaces = surf.edgeFaces()[edgeI];
00366 
00367         forAll(eFaces, i)
00368         {
00369             label nbrFaceI = eFaces[i];
00370 
00371             if (faceToEdge[nbrFaceI] != -1)
00372             {
00373                 if (collapseRegion[nbrFaceI] == -1)
00374                 {
00375                     markRegion
00376                     (
00377                         surf,
00378                         faceToEdge,
00379                         regionI,
00380                         nbrFaceI,
00381                         collapseRegion
00382                     );
00383                 }
00384                 else if (collapseRegion[nbrFaceI] != regionI)
00385                 {
00386                     FatalErrorIn("markRegion")
00387                         << "Edge:" << edgeI << " between face " << faceI
00388                         << " with region " << regionI
00389                         << " and face " << nbrFaceI
00390                         << " with region " << collapseRegion[nbrFaceI]
00391                         << endl;
00392                 }
00393             }
00394         }
00395     }
00396 }
00397 
00398 
00399 // Mark every face with region (in collapseRegion) (or -1).
00400 // Return number of regions.
00401 static label markRegions
00402 (
00403     const triSurface& surf,
00404     const labelList& faceToEdge,
00405     labelList& collapseRegion
00406 )
00407 {
00408     label regionI = 0;
00409 
00410     forAll(faceToEdge, faceI)
00411     {
00412         if (collapseRegion[faceI] == -1 && faceToEdge[faceI] != -1)
00413         {
00414             Pout<< "markRegions : Marking region:" << regionI
00415                 << " starting from face " << faceI << endl;
00416 
00417             // Collapsed face. Mark connected region with current region number
00418             markRegion(surf, faceToEdge, regionI++, faceI, collapseRegion);
00419         }
00420     }
00421     return regionI;
00422 }
00423 
00424 
00425 // Type of region.
00426 // -1  : edge inbetween uncollapsed faces.
00427 // -2  : edge inbetween collapsed faces
00428 // >=0 : edge inbetween uncollapsed and collapsed region. Returns region.
00429 static label edgeType
00430 (
00431     const triSurface& surf,
00432     const labelList& collapseRegion,
00433     const label edgeI
00434 )
00435 {
00436     const labelList& eFaces = surf.edgeFaces()[edgeI];
00437 
00438     // Detect if edge is inbetween collapseRegion and non-collapse face
00439     bool usesUncollapsed = false;
00440     label usesRegion = -1;
00441 
00442     forAll(eFaces, i)
00443     {
00444         label faceI = eFaces[i];
00445 
00446         label region = collapseRegion[faceI];
00447 
00448         if (region == -1)
00449         {
00450             usesUncollapsed = true;
00451         }
00452         else if (usesRegion == -1)
00453         {
00454             usesRegion = region;
00455         }
00456         else if (usesRegion != region)
00457         {
00458             FatalErrorIn("edgeRegion") << abort(FatalError);
00459         }
00460         else
00461         {
00462             // Equal regions.
00463         }
00464     }
00465 
00466     if (usesUncollapsed)
00467     {
00468         if (usesRegion == -1)
00469         {
00470             // uncollapsed faces only.
00471             return -1;
00472         }
00473         else
00474         {
00475             // between uncollapsed and collapsed.
00476             return usesRegion;
00477         }
00478     }
00479     else
00480     {
00481         if (usesRegion == -1)
00482         {
00483             FatalErrorIn("edgeRegion") << abort(FatalError); 
00484             return -2;
00485         }
00486         else
00487         {
00488             return -2;
00489         }
00490     }
00491 }
00492 
00493 
00494 // Get points on outside edge of region (= outside points)
00495 static labelListList getOutsideVerts
00496 (
00497     const triSurface& surf,
00498     const labelList& collapseRegion,
00499     const label nRegions
00500 )
00501 {
00502     const labelListList& edgeFaces = surf.edgeFaces();
00503 
00504     // Per region all the outside vertices.
00505     labelListList outsideVerts(nRegions);
00506 
00507     forAll(edgeFaces, edgeI)
00508     {
00509         // Detect if edge is inbetween collapseRegion and non-collapse face
00510         label regionI = edgeType(surf, collapseRegion, edgeI);
00511 
00512         if (regionI >= 0)
00513         {
00514             // Edge borders both uncollapsed face and collapsed face on region
00515             // usesRegion.
00516 
00517             const edge& e = surf.edges()[edgeI];
00518 
00519             labelList& regionVerts = outsideVerts[regionI];
00520 
00521             // Add both edge points to regionVerts.
00522             forAll(e, eI)
00523             {
00524                 label v = e[eI];
00525 
00526                 if (findIndex(regionVerts, v) == -1)
00527                 {
00528                     label sz = regionVerts.size();
00529                     regionVerts.setSize(sz+1);
00530                     regionVerts[sz] = v;
00531                 }
00532             }
00533         }
00534     }
00535 
00536     return outsideVerts;
00537 }
00538 
00539 
00540 // n^2 search for furthest removed point pair.
00541 static labelPair getSpanPoints
00542 (
00543     const triSurface& surf,
00544     const labelList& outsideVerts
00545 )
00546 {
00547     const pointField& localPoints = surf.localPoints();
00548 
00549     scalar maxDist = -GREAT;
00550     labelPair maxPair;
00551 
00552     forAll(outsideVerts, i)
00553     {
00554         label v0 = outsideVerts[i];
00555 
00556         for (label j = i+1; j < outsideVerts.size(); j++)
00557         {
00558             label v1 = outsideVerts[j];
00559 
00560             scalar d = mag(localPoints[v0] - localPoints[v1]);
00561 
00562             if (d > maxDist)
00563             {
00564                 maxDist = d;
00565                 maxPair[0] = v0;
00566                 maxPair[1] = v1;
00567             }
00568         }
00569     }
00570 
00571     return maxPair;
00572 }
00573 
00574 
00575 // Project all non-span points onto the span edge.
00576 static void projectNonSpanPoints
00577 (
00578     const triSurface& surf,
00579     const labelList& outsideVerts,
00580     const labelPair& spanPair,
00581     labelList& sortedVertices,
00582     scalarField& sortedWeights
00583 )
00584 {
00585     const point& p0 = surf.localPoints()[spanPair[0]];
00586     const point& p1 = surf.localPoints()[spanPair[1]];
00587 
00588     forAll(outsideVerts, i)
00589     {
00590         label v = outsideVerts[i];
00591 
00592         if (v != spanPair[0] && v != spanPair[1])
00593         {
00594             // Is a non-span point. Project onto spanning edge.
00595 
00596             pointHit pHit =
00597                 linePointRef(p0, p1).nearestDist
00598                 (
00599                     surf.localPoints()[v]
00600                 );
00601 
00602             if (!pHit.hit())
00603             {
00604                 FatalErrorIn("projectNonSpanPoints")
00605                     << abort(FatalError);
00606             }
00607 
00608             scalar w = mag(pHit.hitPoint() - p0) / mag(p1 - p0);
00609 
00610             insertSorted(v, w, sortedVertices, sortedWeights);
00611         }
00612     }
00613 }
00614 
00615 
00616 // Slice part of the orderVertices (and optionally reverse) for this edge.
00617 static void getSplitVerts
00618 (
00619     const triSurface& surf,
00620     const label regionI,
00621     const labelPair& spanPoints,
00622     const labelList& orderedVerts,
00623     const scalarField& orderedWeights,
00624     const label edgeI,
00625 
00626     labelList& splitVerts,
00627     scalarField& splitWeights
00628 )
00629 {
00630     const edge& e = surf.edges()[edgeI];
00631     const label sz = orderedVerts.size();
00632 
00633     if (e[0] == spanPoints[0])
00634     {
00635         // Edge in same order as spanPoints&orderedVerts. Keep order.
00636 
00637         if (e[1] == spanPoints[1])
00638         {
00639             // Copy all.
00640             splitVerts = orderedVerts;
00641             splitWeights = orderedWeights;
00642         }
00643         else
00644         {
00645             // Copy upto (but not including) e[1]
00646             label i1 = findIndex(orderedVerts, e[1]);
00647             splitVerts = SubList<label>(orderedVerts, i1, 0);
00648             splitWeights = SubList<scalar>(orderedWeights, i1, 0);
00649         }
00650     }
00651     else if (e[0] == spanPoints[1])
00652     {
00653         // Reverse.
00654 
00655         if (e[1] == spanPoints[0])
00656         {
00657             // Copy all.
00658             splitVerts = orderedVerts;
00659             reverse(splitVerts);
00660             splitWeights = orderedWeights;
00661             reverse(splitWeights);
00662         }
00663         else
00664         {
00665             // Copy downto (but not including) e[1]
00666 
00667             label i1 = findIndex(orderedVerts, e[1]);
00668             splitVerts = SubList<label>(orderedVerts, sz-(i1+1), i1+1);
00669             reverse(splitVerts);
00670             splitWeights = SubList<scalar>(orderedWeights, sz-(i1+1), i1+1);
00671             reverse(splitWeights);
00672         }
00673     }
00674     else if (e[1] == spanPoints[0])
00675     {
00676         // Reverse.
00677 
00678         // Copy upto (but not including) e[0]
00679 
00680         label i0 = findIndex(orderedVerts, e[0]);
00681         splitVerts = SubList<label>(orderedVerts, i0, 0);
00682         reverse(splitVerts);
00683         splitWeights = SubList<scalar>(orderedWeights, i0, 0);
00684         reverse(splitWeights);
00685     }
00686     else if (e[1] == spanPoints[1])
00687     {
00688         // Copy from (but not including) e[0] to end
00689 
00690         label i0 = findIndex(orderedVerts, e[0]);
00691         splitVerts = SubList<label>(orderedVerts, sz-(i0+1), i0+1);
00692         splitWeights = SubList<scalar>(orderedWeights, sz-(i0+1), i0+1);
00693     }
00694     else
00695     {
00696         label i0 = findIndex(orderedVerts, e[0]);
00697         label i1 = findIndex(orderedVerts, e[1]);
00698 
00699         if (i0 == -1 || i1 == -1)
00700         {
00701             FatalErrorIn("getSplitVerts")
00702                 << "Did not find edge in projected vertices." << nl
00703                 << "region:" << regionI << nl
00704                 << "spanPoints:" << spanPoints
00705                 << "  coords:" << surf.localPoints()[spanPoints[0]]
00706                 << surf.localPoints()[spanPoints[1]] << nl
00707                 << "edge:" << edgeI
00708                 << "  verts:" << e
00709                 << "  coords:" << surf.localPoints()[e[0]]
00710                 << surf.localPoints()[e[1]] << nl
00711                 << "orderedVerts:" << orderedVerts << nl
00712                 << abort(FatalError);
00713         }
00714 
00715         if (i0 < i1)
00716         {
00717             splitVerts = SubList<label>(orderedVerts, i1-i0-1, i0+1);
00718             splitWeights = SubList<scalar>(orderedWeights, i1-i0-1, i0+1);
00719         }
00720         else
00721         {
00722             splitVerts = SubList<label>(orderedVerts, i0-i1-1, i1+1);
00723             reverse(splitVerts);
00724             splitWeights = SubList<scalar>(orderedWeights, i0-i1-1, i1+1);
00725             reverse(splitWeights);
00726         }
00727     }
00728 }
00729 
00730 
00731 label collapseBase(triSurface& surf, const scalar minLen)
00732 {
00733     label nTotalSplit = 0;
00734 
00735     label iter = 0;
00736 
00737     while (true)
00738     {
00739         // Detect faces to collapse
00740         // ~~~~~~~~~~~~~~~~~~~~~~~~
00741 
00742         // -1 or edge the face is collapsed onto.
00743         labelList faceToEdge(surf.size(), -1);
00744 
00745         // Calculate faceToEdge (face collapses)
00746         markCollapsedFaces(surf, minLen, faceToEdge);
00747 
00748 
00749         // Find regions of connected collapsed faces
00750         // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
00751 
00752         // per face -1 or region 
00753         labelList collapseRegion(surf.size(), -1);
00754         
00755         label nRegions = markRegions(surf, faceToEdge, collapseRegion);
00756 
00757         Pout<< "Detected " << nRegions << " regions of faces to be collapsed"
00758             << nl << endl;
00759 
00760         // Pick up all vertices on outside of region
00761         labelListList outsideVerts
00762         (
00763             getOutsideVerts(surf, collapseRegion, nRegions)
00764         );
00765 
00766         // For all regions determine maximum distance between points
00767         List<labelPair> spanPoints(nRegions);
00768         labelListList orderedVertices(nRegions);
00769         List<scalarField> orderedWeights(nRegions);
00770 
00771         forAll(spanPoints, regionI)
00772         {
00773             spanPoints[regionI] = getSpanPoints(surf, outsideVerts[regionI]);
00774 
00775             Pout<< "For region " << regionI << " found extrema at points "
00776                 << surf.localPoints()[spanPoints[regionI][0]]
00777                 << surf.localPoints()[spanPoints[regionI][1]]
00778                 << endl;
00779 
00780             // Project all non-span points onto the span edge.
00781             projectNonSpanPoints
00782             (
00783                 surf,
00784                 outsideVerts[regionI],
00785                 spanPoints[regionI], 
00786                 orderedVertices[regionI],
00787                 orderedWeights[regionI]
00788             );
00789 
00790             Pout<< "For region:" << regionI
00791                 << " span:" << spanPoints[regionI]
00792                 << " orderedVerts:" << orderedVertices[regionI]
00793                 << " orderedWeights:" << orderedWeights[regionI]
00794                 << endl;
00795 
00796             writeRegionOBJ
00797             (
00798                 surf,
00799                 regionI,
00800                 collapseRegion,
00801                 outsideVerts[regionI]
00802             );
00803 
00804             Pout<< endl;
00805         }
00806 
00807 
00808 
00809         // Actually split the edges
00810         // ~~~~~~~~~~~~~~~~~~~~~~~~
00811 
00812 
00813         const List<labelledTri>& localFaces = surf.localFaces();
00814         const edgeList& edges = surf.edges();
00815 
00816         label nSplit = 0;
00817 
00818         // Storage for new triangles.
00819         DynamicList<labelledTri> newTris(surf.size());
00820 
00821         // Whether face has been dealt with (either copied/split or deleted)
00822         boolList faceHandled(surf.size(), false);
00823 
00824 
00825         forAll(edges, edgeI)
00826         {
00827             const edge& e = edges[edgeI];
00828 
00829             // Detect if edge is inbetween collapseRegion and non-collapse face
00830             label regionI = edgeType(surf, collapseRegion, edgeI);
00831 
00832             if (regionI == -2)
00833             {
00834                 // inbetween collapsed faces. nothing needs to be done.
00835             }
00836             else if (regionI == -1)
00837             {
00838                 // edge inbetween uncollapsed faces. Handle these later on.
00839             }
00840             else
00841             {
00842                 // some faces around edge are collapsed.
00843 
00844                 // Find additional set of points on edge to be used to split
00845                 // the remaining faces.
00846 
00847                 labelList splitVerts;
00848                 scalarField splitWeights;
00849                 getSplitVerts
00850                 (
00851                     surf,
00852                     regionI,
00853                     spanPoints[regionI],
00854                     orderedVertices[regionI],
00855                     orderedWeights[regionI],
00856                     edgeI,
00857 
00858                     splitVerts,
00859                     splitWeights
00860                 );
00861 
00862                 if (splitVerts.size())
00863                 {
00864                     // Split edge using splitVerts. All non-collapsed triangles 
00865                     // using edge will get split.
00866 
00867 
00868                     {
00869                         const pointField& localPoints = surf.localPoints();
00870                         Pout<< "edge " << edgeI << ' ' << e
00871                             << "  points "
00872                             << localPoints[e[0]] << ' ' << localPoints[e[1]]
00873                             << " split into edges with extra points:"
00874                             << endl;
00875                         forAll(splitVerts, i)
00876                         {
00877                             Pout<< "    " << splitVerts[i] << " weight "
00878                                 << splitWeights[i] << nl;
00879                         }
00880                     }
00881 
00882                     const labelList& eFaces = surf.edgeFaces()[edgeI];
00883 
00884                     forAll(eFaces, i)
00885                     {
00886                         label faceI = eFaces[i];
00887 
00888                         if (!faceHandled[faceI] && faceToEdge[faceI] == -1)
00889                         {
00890                             // Split face to use vertices.
00891                             splitTri
00892                             (
00893                                 localFaces[faceI],
00894                                 e,
00895                                 splitVerts,
00896                                 newTris
00897                             );
00898 
00899                             faceHandled[faceI] = true;
00900 
00901                             nSplit++;
00902                         }
00903                     }
00904                 }
00905             }
00906         }
00907 
00908         // Copy all unsplit faces
00909         forAll(faceHandled, faceI)
00910         {
00911             if (!faceHandled[faceI] && faceToEdge[faceI] == -1)
00912             {
00913                 newTris.append(localFaces[faceI]);
00914             }
00915         }
00916 
00917         Pout<< "collapseBase : splitting " << nSplit << " triangles"
00918             << endl;
00919 
00920         nTotalSplit += nSplit;
00921 
00922         if (nSplit == 0)
00923         {
00924             break;
00925         }
00926 
00927         // Pack the triangles
00928         newTris.shrink();
00929 
00930         Pout<< "Resetting surface from " << surf.size() << " to "
00931             << newTris.size() << " triangles" << endl;
00932         surf = triSurface(newTris, surf.patches(), surf.localPoints());
00933 
00934         {
00935             fileName fName("bla" + name(iter) + ".obj");
00936             Pout<< "Writing surf to " << fName << endl;
00937             surf.write(fName);
00938         }
00939 
00940         iter++;
00941     }
00942 
00943     // Remove any unused vertices
00944     surf = triSurface(surf.localFaces(), surf.patches(), surf.localPoints());
00945 
00946     return nTotalSplit;
00947 }
00948 
00949 
00950 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines