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

cellCuts.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 "cellCuts.H"
00027 #include <OpenFOAM/polyMesh.H>
00028 #include <OpenFOAM/Time.H>
00029 #include <OpenFOAM/ListOps.H>
00030 #include <dynamicMesh/cellLooper.H>
00031 #include <dynamicMesh/refineCell.H>
00032 #include <meshTools/meshTools.H>
00033 #include <dynamicMesh/geomCellLooper.H>
00034 #include <OpenFOAM/OFstream.H>
00035 #include <OpenFOAM/plane.H>
00036 
00037 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
00038 
00039 defineTypeNameAndDebug(Foam::cellCuts, 0);
00040 
00041 
00042 // * * * * * * * * * * * * * Private Static Functions  * * * * * * * * * * * //
00043 
00044 // Find val in first nElems elements of list.
00045 Foam::label Foam::cellCuts::findPartIndex
00046 (
00047     const labelList& elems,
00048     const label nElems,
00049     const label val
00050 )
00051 {
00052     for (label i = 0; i < nElems; i++)
00053     {
00054         if (elems[i] == val)
00055         {
00056             return i;
00057         }
00058     }
00059     return -1;
00060 }
00061 
00062 
00063 Foam::boolList Foam::cellCuts::expand
00064 (
00065     const label size,
00066     const labelList& labels
00067 )
00068 {
00069     boolList result(size, false);
00070 
00071     forAll(labels, labelI)
00072     {
00073         result[labels[labelI]] = true;
00074     }
00075     return result;
00076 }
00077 
00078 
00079 Foam::scalarField Foam::cellCuts::expand
00080 (
00081     const label size,
00082     const labelList& labels,
00083     const scalarField& weights
00084 )
00085 {
00086     scalarField result(size, -GREAT);
00087 
00088     forAll(labels, labelI)
00089     {
00090         result[labels[labelI]] = weights[labelI];
00091     }
00092     return result;
00093 }
00094 
00095 
00096 // Find first point in lst not in map.
00097 Foam::label Foam::cellCuts::firstUnique
00098 (
00099     const labelList& lst,
00100     const Map<label>& map
00101 )
00102 {
00103     forAll(lst, i)
00104     {
00105         if (!map.found(lst[i]))
00106         {
00107             return i;
00108         }
00109     }
00110     return -1;
00111 }
00112 
00113 
00114 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
00115 
00116 // Write cell and raw cuts on any of the elements
00117 void Foam::cellCuts::writeUncutOBJ
00118 (
00119     const fileName& dir,
00120     const label cellI
00121 ) const
00122 {
00123     //- Cell edges
00124     OFstream cutsStream(dir / "cell_" + name(cellI) + ".obj");
00125 
00126     Pout<< "Writing cell for time " <<  mesh().time().timeName()
00127         << " to " << cutsStream.name() << nl;
00128 
00129     meshTools::writeOBJ
00130     (
00131         cutsStream,
00132         mesh().cells(),
00133         mesh().faces(),
00134         mesh().points(),
00135         labelList(1, cellI)
00136     );
00137 
00138     //- Loop cutting cell in two
00139     OFstream cutStream(dir / "cellCuts_" + name(cellI) + ".obj");
00140 
00141     Pout<< "Writing raw cuts on cell for time " <<  mesh().time().timeName()
00142         << " to " << cutStream.name() << nl;
00143 
00144     const labelList& cPoints = mesh().cellPoints()[cellI];
00145 
00146     forAll(cPoints, i)
00147     {
00148         label pointI = cPoints[i];
00149         if (pointIsCut_[pointI])
00150         {
00151             meshTools::writeOBJ(cutStream, mesh().points()[pointI]);
00152         }
00153     }
00154 
00155     const pointField& pts = mesh().points();
00156 
00157     const labelList& cEdges = mesh().cellEdges()[cellI];
00158 
00159     forAll(cEdges, i)
00160     {
00161         label edgeI = cEdges[i];
00162 
00163         if (edgeIsCut_[edgeI])
00164         {
00165             const edge& e = mesh().edges()[edgeI];
00166 
00167             const scalar w = edgeWeight_[edgeI];
00168 
00169             meshTools::writeOBJ(cutStream, w*pts[e[1]] + (1-w)*pts[e[0]]);
00170         }
00171     }
00172 }
00173 
00174 
00175 void Foam::cellCuts::writeOBJ
00176 (
00177     const fileName& dir,
00178     const label cellI,
00179     const pointField& loopPoints,
00180     const labelList& anchors
00181 ) const
00182 {
00183     //- Cell edges
00184     OFstream cutsStream(dir / "cell_" + name(cellI) + ".obj");
00185 
00186     Pout<< "Writing cell for time " <<  mesh().time().timeName()
00187         << " to " << cutsStream.name() << nl;
00188 
00189     meshTools::writeOBJ
00190     (
00191         cutsStream,
00192         mesh().cells(),
00193         mesh().faces(),
00194         mesh().points(),
00195         labelList(1, cellI)
00196     );
00197 
00198 
00199     //- Loop cutting cell in two
00200     OFstream loopStream(dir / "cellLoop_" + name(cellI) + ".obj");
00201 
00202     Pout<< "Writing loop for time " <<  mesh().time().timeName()
00203         << " to " << loopStream.name() << nl;
00204 
00205     label vertI = 0;
00206 
00207     writeOBJ(loopStream, loopPoints, vertI);
00208 
00209 
00210     //- Anchors for cell
00211     OFstream anchorStream(dir / "anchors_" + name(cellI) + ".obj");
00212 
00213     Pout<< "Writing anchors for time " <<  mesh().time().timeName()
00214         << " to " << anchorStream.name() << endl;
00215 
00216     forAll(anchors, i)
00217     {
00218         meshTools::writeOBJ(anchorStream, mesh().points()[anchors[i]]);
00219     }
00220 }
00221 
00222 
00223 // Find face on cell using the two edges.
00224 Foam::label Foam::cellCuts::edgeEdgeToFace
00225 (
00226     const label cellI,
00227     const label edgeA,
00228     const label edgeB
00229 ) const
00230 {
00231     const labelList& cFaces = mesh().cells()[cellI];
00232 
00233     forAll(cFaces, cFaceI)
00234     {
00235         label faceI = cFaces[cFaceI];
00236 
00237         const labelList& fEdges = mesh().faceEdges()[faceI];
00238 
00239         if
00240         (
00241             findIndex(fEdges, edgeA) != -1
00242          && findIndex(fEdges, edgeB) != -1
00243         )
00244         {
00245            return faceI;
00246         }
00247     }
00248 
00249     // Coming here means the loop is illegal since the two edges
00250     // are not shared by a face. We just mark loop as invalid and continue.
00251 
00252     WarningIn
00253     (
00254         "Foam::cellCuts::edgeEdgeToFace"
00255         "(const label cellI, const label edgeA,"
00256         "const label edgeB) const"
00257     )   << "cellCuts : Cannot find face on cell "
00258         << cellI << " that has both edges " << edgeA << ' ' << edgeB << endl
00259         << "faces : " << cFaces << endl
00260         << "edgeA : " << mesh().edges()[edgeA] << endl
00261         << "edgeB : " << mesh().edges()[edgeB] << endl
00262         << "Marking the loop across this cell as invalid" << endl;
00263 
00264     return -1;
00265 }
00266 
00267 
00268 // Find face on cell using an edge and a vertex.
00269 Foam::label Foam::cellCuts::edgeVertexToFace
00270 (
00271     const label cellI,
00272     const label edgeI,
00273     const label vertI
00274 ) const
00275 {
00276     const labelList& cFaces = mesh().cells()[cellI];
00277 
00278     forAll(cFaces, cFaceI)
00279     {
00280         label faceI = cFaces[cFaceI];
00281 
00282         const face& f = mesh().faces()[faceI];
00283 
00284         const labelList& fEdges = mesh().faceEdges()[faceI];
00285 
00286         if
00287         (
00288             findIndex(fEdges, edgeI) != -1
00289          && findIndex(f, vertI) != -1
00290         )
00291         {
00292            return faceI;
00293         }
00294     }
00295 
00296     WarningIn
00297     (
00298         "Foam::cellCuts::edgeVertexToFace"
00299         "(const label cellI, const label edgeI, "
00300         "const label vertI) const"
00301     )   << "cellCuts : Cannot find face on cell "
00302         << cellI << " that has both edge " << edgeI << " and vertex "
00303         << vertI << endl
00304         << "faces : " << cFaces << endl
00305         << "edge : " << mesh().edges()[edgeI] << endl
00306         << "Marking the loop across this cell as invalid" << endl;
00307 
00308     return -1;
00309 }
00310 
00311 
00312 // Find face using two vertices (guaranteed not to be along edge)
00313 Foam::label Foam::cellCuts::vertexVertexToFace
00314 (
00315     const label cellI,
00316     const label vertA,
00317     const label vertB
00318 ) const
00319 {
00320     const labelList& cFaces = mesh().cells()[cellI];
00321 
00322     forAll(cFaces, cFaceI)
00323     {
00324         label faceI = cFaces[cFaceI];
00325 
00326         const face& f = mesh().faces()[faceI];
00327 
00328         if
00329         (
00330             findIndex(f, vertA) != -1
00331          && findIndex(f, vertB) != -1
00332         )
00333         {
00334            return faceI;
00335         }
00336     }
00337 
00338     WarningIn("Foam::cellCuts::vertexVertexToFace")
00339         << "cellCuts : Cannot find face on cell "
00340         << cellI << " that has vertex " << vertA << " and vertex "
00341         << vertB << endl
00342         << "faces : " << cFaces << endl
00343         << "Marking the loop across this cell as invalid" << endl;
00344 
00345     return -1;
00346 }
00347 
00348 
00349 void Foam::cellCuts::calcFaceCuts() const
00350 {
00351     if (faceCutsPtr_)
00352     {
00353         FatalErrorIn("cellCuts::calcFaceCuts()")
00354             << "faceCuts already calculated" << abort(FatalError);
00355     }
00356 
00357     const faceList& faces = mesh().faces();
00358 
00359 
00360     faceCutsPtr_ = new labelListList(mesh().nFaces());
00361     labelListList& faceCuts = *faceCutsPtr_;
00362 
00363     for (label faceI = 0; faceI < mesh().nFaces(); faceI++)
00364     {
00365         const face& f = faces[faceI];
00366 
00367         // Big enough storage (possibly all points and all edges cut). Shrink
00368         // later on.
00369         labelList& cuts = faceCuts[faceI];
00370 
00371         cuts.setSize(2*f.size());
00372 
00373         label cutI = 0;
00374 
00375         // Do point-edge-point walk over face and collect all cuts.
00376         // Problem is that we want to start from one of the endpoints of a
00377         // string of connected cuts; we don't want to start somewhere in the
00378         // middle.
00379 
00380         // Pass1: find first point cut not preceeded by a cut.
00381         label startFp = -1;
00382 
00383         forAll(f, fp)
00384         {
00385             label v0 = f[fp];
00386 
00387             if (pointIsCut_[v0])
00388             {
00389                 label vMin1 = f[f.rcIndex(fp)];
00390 
00391                 if
00392                 (
00393                     !pointIsCut_[vMin1]
00394                  && !edgeIsCut_[findEdge(faceI, v0, vMin1)]
00395                 )
00396                 {
00397                     cuts[cutI++] = vertToEVert(v0);
00398                     startFp = f.fcIndex(fp);
00399                     break;
00400                 }
00401             }
00402         }
00403 
00404         // Pass2: first edge cut not preceeded by point cut
00405         if (startFp == -1)
00406         {
00407             forAll(f, fp)
00408             {
00409                 label fp1 = f.fcIndex(fp);
00410 
00411                 label v0 = f[fp];
00412                 label v1 = f[fp1];
00413 
00414                 label edgeI = findEdge(faceI, v0, v1);
00415 
00416                 if (edgeIsCut_[edgeI] && !pointIsCut_[v0])
00417                 {
00418                     cuts[cutI++] = edgeToEVert(edgeI);
00419                     startFp = fp1;
00420                     break;
00421                 }
00422             }
00423         }
00424 
00425         // Pass3: nothing found so far. Either face is not cut at all or all
00426         // vertices are cut. Start from 0.
00427         if (startFp == -1)
00428         {
00429             startFp = 0;
00430         }
00431 
00432         // Store all cuts starting from startFp;
00433         label fp = startFp;
00434 
00435         bool allVerticesCut = true;
00436 
00437         forAll(f, i)
00438         {
00439             label fp1 = f.fcIndex(fp);
00440 
00441             // Get the three items: current vertex, next vertex and edge
00442             // inbetween
00443             label v0 = f[fp];
00444             label v1 = f[fp1];
00445             label edgeI = findEdge(faceI, v0, v1);
00446 
00447             if (pointIsCut_[v0])
00448             {
00449                 cuts[cutI++] = vertToEVert(v0);
00450             }
00451             else
00452             {
00453                 // For check if all vertices have been cut (= illegal)
00454                 allVerticesCut = false;
00455             }
00456 
00457             if (edgeIsCut_[edgeI])
00458             {
00459                 cuts[cutI++] = edgeToEVert(edgeI);
00460             }
00461 
00462             fp = fp1;
00463         }
00464 
00465 
00466         if (allVerticesCut)
00467         {
00468             WarningIn("Foam::cellCuts::calcFaceCuts() const")
00469                 << "Face " << faceI << " vertices " << f
00470                 << " has all its vertices cut. Not cutting face." << endl;
00471 
00472             cutI = 0;
00473         }
00474 
00475         // Remove duplicate starting point
00476         if (cutI > 1 && cuts[cutI-1] == cuts[0])
00477         {
00478             cutI--;
00479         }
00480         cuts.setSize(cutI);
00481     }
00482 }
00483 
00484 
00485 // Find edge on face using two vertices
00486 Foam::label Foam::cellCuts::findEdge
00487 (
00488     const label faceI,
00489     const label v0,
00490     const label v1
00491 ) const
00492 {
00493     const edgeList& edges = mesh().edges();
00494 
00495     const labelList& fEdges = mesh().faceEdges()[faceI];
00496 
00497     forAll(fEdges, i)
00498     {
00499         const edge& e = edges[fEdges[i]];
00500 
00501         if
00502         (
00503             (e[0] == v0 && e[1] == v1)
00504          || (e[0] == v1 && e[1] == v0)
00505         )
00506         {
00507             return fEdges[i];
00508         }
00509     }
00510 
00511     return -1;
00512 }
00513 
00514 
00515 // Check if there is a face on the cell on which all cuts are.
00516 Foam::label Foam::cellCuts::loopFace
00517 (
00518     const label cellI,
00519     const labelList& loop
00520 ) const
00521 {
00522     const cell& cFaces = mesh().cells()[cellI];
00523 
00524     forAll(cFaces, cFaceI)
00525     {
00526         label faceI = cFaces[cFaceI];
00527 
00528         const labelList& fEdges = mesh().faceEdges()[faceI];
00529         const face& f = mesh().faces()[faceI];
00530 
00531         bool allOnFace = true;
00532 
00533         forAll(loop, i)
00534         {
00535             label cut = loop[i];
00536 
00537             if (isEdge(cut))
00538             {
00539                 if (findIndex(fEdges, getEdge(cut)) == -1)
00540                 {
00541                     // Edge not on face. Skip face.
00542                     allOnFace = false;
00543                     break;
00544                 }
00545             }
00546             else
00547             {
00548                 if (findIndex(f, getVertex(cut)) == -1)
00549                 {
00550                     // Vertex not on face. Skip face.
00551                     allOnFace = false;
00552                     break;
00553                 }
00554             }
00555         }
00556 
00557         if (allOnFace)
00558         {
00559             // Found face where all elements of loop are on the face.
00560             return faceI;
00561         }
00562     }
00563     return -1;
00564 }
00565 
00566 
00567 // From point go into connected face
00568 bool Foam::cellCuts::walkPoint
00569 (
00570     const label cellI,
00571     const label startCut,
00572 
00573     const label exclude0,
00574     const label exclude1,
00575 
00576     const label otherCut,
00577 
00578     label& nVisited,
00579     labelList& visited
00580 ) const
00581 {
00582     label vertI = getVertex(otherCut);
00583 
00584     const labelList& pFaces = mesh().pointFaces()[vertI];
00585 
00586     forAll(pFaces, pFaceI)
00587     {
00588         label otherFaceI = pFaces[pFaceI];
00589 
00590         if
00591         (
00592             otherFaceI != exclude0
00593          && otherFaceI != exclude1
00594          && meshTools::faceOnCell(mesh(), cellI, otherFaceI)
00595         )
00596         {
00597             label oldNVisited = nVisited;
00598 
00599             bool foundLoop =
00600                 walkCell
00601                 (
00602                     cellI,
00603                     startCut,
00604                     otherFaceI,
00605                     otherCut,
00606                     nVisited,
00607                     visited
00608                 );
00609 
00610             if (foundLoop)
00611             {
00612                 return true;
00613             }
00614 
00615             // No success. Restore state and continue
00616             nVisited = oldNVisited;
00617         }
00618     }
00619     return false;
00620 }
00621 
00622 
00623 // Cross cut (which is edge on faceI) onto next face
00624 bool Foam::cellCuts::crossEdge
00625 (
00626     const label cellI,
00627     const label startCut,
00628     const label faceI,
00629     const label otherCut,
00630 
00631     label& nVisited,
00632     labelList& visited
00633 ) const
00634 {
00635     // Cross edge to other face
00636     label edgeI = getEdge(otherCut);
00637 
00638     label otherFaceI = meshTools::otherFace(mesh(), cellI, faceI, edgeI);
00639 
00640     // Store old state
00641     label oldNVisited = nVisited;
00642 
00643     // Recurse to otherCut
00644     bool foundLoop =
00645         walkCell
00646         (
00647             cellI,
00648             startCut,
00649             otherFaceI,
00650             otherCut,
00651             nVisited,
00652             visited
00653         );
00654 
00655     if (foundLoop)
00656     {
00657         return true;
00658     }
00659     else
00660     {
00661         // No success. Restore state (i.e. backtrack)
00662         nVisited = oldNVisited;
00663 
00664         return false;
00665     }
00666 }
00667 
00668 
00669 bool Foam::cellCuts::addCut
00670 (
00671     const label cellI,
00672     const label cut,
00673     label& nVisited,
00674     labelList& visited
00675 ) const
00676 {
00677     // Check for duplicate cuts.
00678     if (findPartIndex(visited, nVisited, cut) != -1)
00679     {
00680         // Truncate (copy of) visited for ease of printing.
00681         labelList truncVisited(visited);
00682         truncVisited.setSize(nVisited);
00683 
00684         Pout<< "For cell " << cellI << " : trying to add duplicate cut " << cut;
00685         labelList cuts(1, cut);
00686         writeCuts(Pout, cuts, loopWeights(cuts));
00687 
00688         Pout<< " to path:";
00689         writeCuts(Pout, truncVisited, loopWeights(truncVisited));
00690         Pout<< endl;
00691 
00692         return false;
00693     }
00694     else
00695     {
00696         visited[nVisited++] = cut;
00697 
00698         return true;
00699     }
00700 }
00701 
00702 
00703 // Walk across faceI, storing cuts as you go. Returns last two cuts visisted.
00704 // Returns true if valid walk.
00705 bool Foam::cellCuts::walkFace
00706 (
00707     const label cellI,
00708     const label startCut,
00709     const label faceI,
00710     const label cut,
00711 
00712     label& lastCut,
00713     label& beforeLastCut,
00714     label& nVisited,
00715     labelList& visited
00716 ) const
00717 {
00718     const labelList& fCuts = faceCuts()[faceI];
00719 
00720     if (fCuts.size() < 2)
00721     {
00722         return false;
00723     }
00724 
00725     // Easy case : two cuts.
00726     if (fCuts.size() == 2)
00727     {
00728         if (fCuts[0] == cut)
00729         {
00730             if (!addCut(cellI, cut, nVisited, visited))
00731             {
00732                 return false;
00733             }
00734 
00735             beforeLastCut = cut;
00736             lastCut = fCuts[1];
00737 
00738             return true;
00739         }
00740         else
00741         {
00742             if (!addCut(cellI, cut, nVisited, visited))
00743             {
00744                 return false;
00745             }
00746 
00747             beforeLastCut = cut;
00748             lastCut = fCuts[0];
00749 
00750             return true;
00751         }
00752     }
00753 
00754     // Harder case: more than 2 cuts on face.
00755     // Should be start or end of string of cuts. Store all but last cut.
00756     if (fCuts[0] == cut)
00757     {
00758         // Walk forward
00759         for (label i = 0; i < fCuts.size()-1; i++)
00760         {
00761             if (!addCut(cellI, fCuts[i], nVisited, visited))
00762             {
00763                 return false;
00764             }
00765         }
00766         beforeLastCut = fCuts[fCuts.size()-2];
00767         lastCut = fCuts[fCuts.size()-1];
00768     }
00769     else if (fCuts[fCuts.size()-1] == cut)
00770     {
00771         for (label i = fCuts.size()-1; i >= 1; --i)
00772         {
00773             if (!addCut(cellI, fCuts[i], nVisited, visited))
00774             {
00775                 return false;
00776             }
00777         }
00778         beforeLastCut = fCuts[1];
00779         lastCut = fCuts[0];
00780     }
00781     else
00782     {
00783         WarningIn("Foam::cellCuts::walkFace")
00784             << "In middle of cut. cell:" << cellI << " face:" << faceI
00785             << " cuts:" << fCuts << " current cut:" << cut << endl;
00786 
00787         return false;
00788     }
00789 
00790     return true;
00791 }
00792 
00793 
00794 
00795 // Walk across cuts (cut edges or cut vertices) of cell. Stops when hit cut
00796 // already visited. Returns true when loop of 3 or more vertices found.
00797 bool Foam::cellCuts::walkCell
00798 (
00799     const label cellI,
00800     const label startCut,
00801     const label faceI,
00802     const label cut,
00803 
00804     label& nVisited,
00805     labelList& visited
00806 ) const
00807 {
00808     // Walk across face, storing cuts. Return last two cuts visited.
00809     label lastCut = -1;
00810     label beforeLastCut = -1;
00811 
00812 
00813     if (debug & 2)
00814     {
00815         Pout<< "For cell:" << cellI << " walked across face " << faceI
00816             << " from cut ";
00817         labelList cuts(1, cut);
00818         writeCuts(Pout, cuts, loopWeights(cuts));
00819         Pout<< endl;
00820     }
00821 
00822     bool validWalk = walkFace
00823     (
00824         cellI,
00825         startCut,
00826         faceI,
00827         cut,
00828 
00829         lastCut,
00830         beforeLastCut,
00831         nVisited,
00832         visited
00833     );
00834 
00835     if (!validWalk)
00836     {
00837         return false;
00838     }
00839 
00840     if (debug & 2)
00841     {
00842         Pout<< "    to last cut ";
00843         labelList cuts(1, lastCut);
00844         writeCuts(Pout, cuts, loopWeights(cuts));
00845         Pout<< endl;
00846     }
00847 
00848     // Check if starting point reached.
00849     if (lastCut == startCut)
00850     {
00851         if (nVisited >= 3)
00852         {
00853             if (debug & 2)
00854             {
00855                 // Truncate (copy of) visited for ease of printing.
00856                 labelList truncVisited(visited);
00857                 truncVisited.setSize(nVisited);
00858 
00859                 Pout<< "For cell " << cellI << " : found closed path:";
00860                 writeCuts(Pout, truncVisited, loopWeights(truncVisited));
00861                 Pout<< " closed by " << lastCut << endl;
00862             }
00863 
00864             return true;
00865         }
00866         else
00867         {
00868             // Loop but too short. Probably folded back on itself.
00869             return false;
00870         }
00871     }
00872 
00873 
00874     // Check last cut and one before that to determine type.
00875 
00876     // From beforeLastCut to lastCut:
00877     //  - from edge to edge
00878     //      (- always not along existing edge)
00879     //  - from edge to vertex
00880     //      - not along existing edge
00881     //      - along existing edge: not allowed?
00882     //  - from vertex to edge
00883     //      - not along existing edge
00884     //      - along existing edge. Not allowed. See above.
00885     //  - from vertex to vertex
00886     //      - not along existing edge
00887     //      - along existing edge
00888 
00889     if (isEdge(beforeLastCut))
00890     {
00891         if (isEdge(lastCut))
00892         {
00893             // beforeLastCut=edge, lastCut=edge.
00894 
00895             // Cross lastCut (=edge) into face not faceI
00896             return crossEdge
00897             (
00898                 cellI,
00899                 startCut,
00900                 faceI,
00901                 lastCut,
00902                 nVisited,
00903                 visited
00904             );
00905         }
00906         else
00907         {
00908             // beforeLastCut=edge, lastCut=vertex.
00909 
00910             // Go from lastCut to all connected faces (except faceI)
00911             return walkPoint
00912             (
00913                 cellI,
00914                 startCut,
00915                 faceI,          // exclude0: don't cross back on current face
00916                 -1,             // exclude1
00917                 lastCut,
00918                 nVisited,
00919                 visited
00920             );
00921         }
00922     }
00923     else
00924     {
00925         if (isEdge(lastCut))
00926         {
00927             // beforeLastCut=vertex, lastCut=edge.
00928             return crossEdge
00929             (
00930                 cellI,
00931                 startCut,
00932                 faceI,
00933                 lastCut,
00934                 nVisited,
00935                 visited
00936             );
00937         }
00938         else
00939         {
00940             // beforeLastCut=vertex, lastCut=vertex. Check if along existing
00941             // edge.
00942             label edgeI =
00943                 findEdge
00944                 (
00945                     faceI,
00946                     getVertex(beforeLastCut),
00947                     getVertex(lastCut)
00948                 );
00949 
00950             if (edgeI != -1)
00951             {
00952                 // Cut along existing edge. So is in fact on two faces.
00953                 // Get faces on both sides of the edge to make
00954                 // sure we dont fold back on to those.
00955 
00956                 label f0, f1;
00957                 meshTools::getEdgeFaces(mesh(), cellI, edgeI, f0, f1);
00958 
00959                 return walkPoint
00960                 (
00961                     cellI,
00962                     startCut,
00963                     f0,
00964                     f1,
00965                     lastCut,
00966                     nVisited,
00967                     visited
00968                 );
00969 
00970             }
00971             else
00972             {
00973                 // Cross cut across face.
00974                 return walkPoint
00975                 (
00976                     cellI,
00977                     startCut,
00978                     faceI,      // face to exclude
00979                     -1,         // face to exclude
00980                     lastCut,
00981                     nVisited,
00982                     visited
00983                 );
00984             }
00985         }
00986     }
00987 }
00988 
00989 
00990 // Determine for every cut cell the loop (= face) it is cut by. Done by starting
00991 // from a cut edge or cut vertex and walking across faces, from cut to cut,
00992 // until starting cut hit.
00993 // If multiple loops are possible across a cell circumference takes the first
00994 // one found.
00995 void Foam::cellCuts::calcCellLoops(const labelList& cutCells)
00996 {
00997     // Calculate cuts per face.
00998     const labelListList& allFaceCuts = faceCuts();
00999 
01000     // Per cell the number of faces with valid cuts. Is used as quick
01001     // rejection to see if cell can be cut.
01002     labelList nCutFaces(mesh().nCells(), 0);
01003 
01004     forAll(allFaceCuts, faceI)
01005     {
01006         const labelList& fCuts = allFaceCuts[faceI];
01007 
01008         if (fCuts.size() == mesh().faces()[faceI].size())
01009         {
01010             // Too many cuts on face. WalkCell would get very upset so disable.
01011             nCutFaces[mesh().faceOwner()[faceI]] = labelMin;
01012 
01013             if (mesh().isInternalFace(faceI))
01014             {
01015                 nCutFaces[mesh().faceNeighbour()[faceI]] = labelMin;
01016             }
01017         }
01018         else if (fCuts.size() >= 2)
01019         {
01020             // Could be valid cut. Update count for owner and neighbour.
01021             nCutFaces[mesh().faceOwner()[faceI]]++;
01022 
01023             if (mesh().isInternalFace(faceI))
01024             {
01025                 nCutFaces[mesh().faceNeighbour()[faceI]]++;
01026             }
01027         }
01028     }
01029 
01030 
01031     // Stack of visited cuts (nVisited used as stack pointer)
01032     // Size big enough.
01033     labelList visited(mesh().nPoints());
01034 
01035     forAll(cutCells, i)
01036     {
01037         label cellI = cutCells[i];
01038 
01039         bool validLoop = false;
01040 
01041         // Quick rejection: has enough faces that are cut?
01042         if (nCutFaces[cellI] >= 3)
01043         {
01044             const labelList& cFaces = mesh().cells()[cellI];
01045 
01046             if (debug & 2)
01047             {
01048                 Pout<< "cell:" << cellI << " cut faces:" << endl;
01049                 forAll(cFaces, i)
01050                 {
01051                     label faceI = cFaces[i];
01052                     const labelList& fCuts = allFaceCuts[faceI];
01053 
01054                     Pout<< "    face:" << faceI << " cuts:";
01055                     writeCuts(Pout, fCuts, loopWeights(fCuts));
01056                     Pout<< endl;
01057                 }
01058             }
01059 
01060             label nVisited = 0;
01061 
01062             // Determine the first cut face to start walking from.
01063             forAll(cFaces, cFaceI)
01064             {
01065                 label faceI = cFaces[cFaceI];
01066 
01067                 const labelList& fCuts = allFaceCuts[faceI];
01068 
01069                 // Take first or last cut of multiple on face.
01070                 // Note that in calcFaceCuts
01071                 // we have already made sure this is the start or end of a
01072                 // string of cuts.
01073                 if (fCuts.size() >= 2)
01074                 {
01075                     // Try walking from start of fCuts.
01076                     nVisited = 0;
01077 
01078                     if (debug & 2)
01079                     {
01080                         Pout<< "cell:" << cellI
01081                             << " start walk at face:" << faceI
01082                             << " cut:";
01083                         labelList cuts(1, fCuts[0]);
01084                         writeCuts(Pout, cuts, loopWeights(cuts));
01085                         Pout<< endl;
01086                     }
01087 
01088                     validLoop =
01089                         walkCell
01090                         (
01091                             cellI,
01092                             fCuts[0],
01093                             faceI,
01094                             fCuts[0],
01095 
01096                             nVisited,
01097                             visited
01098                         );
01099 
01100                     if (validLoop)
01101                     {
01102                         break;
01103                     }
01104 
01105                     // No need to try and walk from end of cuts since
01106                     // all paths already tried by walkCell.
01107                 }
01108             }
01109 
01110             if (validLoop)
01111             {
01112                 // Copy nVisited elements out of visited (since visited is
01113                 // never truncated for efficiency reasons)
01114 
01115                 labelList& loop = cellLoops_[cellI];
01116 
01117                 loop.setSize(nVisited);
01118 
01119                 for (label i = 0; i < nVisited; i++)
01120                 {
01121                     loop[i] = visited[i];
01122                 }
01123             }
01124             else
01125             {
01126                 // Invalid loop. Leave cellLoops_[cellI] zero size which
01127                 // flags this.
01128                 Pout<< "calcCellLoops(const labelList&) : did not find valid"
01129                     << " loop for cell " << cellI << endl;
01130                 // Dump cell and cuts on cell.
01131                 writeUncutOBJ(".", cellI);
01132 
01133                 cellLoops_[cellI].setSize(0);
01134             }
01135         }
01136         else
01137         {
01138             //Pout<< "calcCellLoops(const labelList&) : did not find valid"
01139             //    << " loop for cell " << cellI << " since not enough cut faces"
01140             //    << endl;
01141             cellLoops_[cellI].setSize(0);
01142         }
01143     }
01144 }
01145 
01146 
01147 // Walk unset edges of single cell from starting point and marks visited
01148 // edges and vertices with status.
01149 void Foam::cellCuts::walkEdges
01150 (
01151     const label cellI,
01152     const label pointI,
01153     const label status,
01154 
01155     Map<label>& edgeStatus,
01156     Map<label>& pointStatus
01157 ) const
01158 {
01159     if (pointStatus.insert(pointI, status))
01160     {
01161         // First visit to pointI
01162 
01163         const labelList& pEdges = mesh().pointEdges()[pointI];
01164 
01165         forAll(pEdges, pEdgeI)
01166         {
01167             label edgeI = pEdges[pEdgeI];
01168 
01169             if
01170             (
01171                 meshTools::edgeOnCell(mesh(), cellI, edgeI)
01172              && edgeStatus.insert(edgeI, status)
01173             )
01174             {
01175                 // First visit to edgeI so recurse.
01176 
01177                 label v2 = mesh().edges()[edgeI].otherVertex(pointI);
01178 
01179                 walkEdges(cellI, v2, status, edgeStatus, pointStatus);
01180             }
01181         }
01182     }
01183 }
01184 
01185 
01186 // Invert anchor point selection.
01187 Foam::labelList Foam::cellCuts::nonAnchorPoints
01188 (
01189     const labelList& cellPoints,
01190     const labelList& anchorPoints,
01191     const labelList& loop
01192 ) const
01193 {
01194     labelList newElems(cellPoints.size());
01195     label newElemI = 0;
01196 
01197     forAll(cellPoints, i)
01198     {
01199         label pointI = cellPoints[i];
01200 
01201         if
01202         (
01203             findIndex(anchorPoints, pointI) == -1
01204          && findIndex(loop, vertToEVert(pointI)) == -1
01205         )
01206         {
01207             newElems[newElemI++] = pointI;
01208         }
01209     }
01210 
01211     newElems.setSize(newElemI);
01212 
01213     return newElems;
01214 }
01215 
01216 
01217 //- Check anchor points on 'outside' of loop
01218 bool Foam::cellCuts::loopAnchorConsistent
01219 (
01220     const label cellI,
01221     const pointField& loopPts,
01222     const labelList& anchorPoints
01223 ) const
01224 {
01225     // Create identity face for ease of calculation of normal etc.
01226     face f(identity(loopPts.size()));
01227 
01228     vector normal = f.normal(loopPts);
01229     point ctr = f.centre(loopPts);
01230 
01231 
01232     // Get average position of anchor points.
01233     vector avg(vector::zero);
01234 
01235     forAll(anchorPoints, ptI)
01236     {
01237         avg += mesh().points()[anchorPoints[ptI]];
01238     }
01239     avg /= anchorPoints.size();
01240 
01241 
01242     if (((avg - ctr) & normal) > 0)
01243     {
01244         return true;
01245     }
01246     else
01247     {
01248         return false;
01249     }
01250 }
01251 
01252 
01253 // Determines set of anchor points given a loop. The loop should split the
01254 // cell into (one or) two sets of vertices. The set of vertices that is
01255 // on the 'normal' side of the loop is the anchor set.
01256 // Returns true if valid set, false otherwise.
01257 bool Foam::cellCuts::calcAnchors
01258 (
01259     const label cellI,
01260     const labelList& loop,
01261     const pointField& loopPts,
01262 
01263     labelList& anchorPoints
01264 ) const
01265 {
01266     const edgeList& edges = mesh().edges();
01267 
01268     const labelList& cPoints = mesh().cellPoints()[cellI];
01269     const labelList& cEdges = mesh().cellEdges()[cellI];
01270     const cell& cFaces = mesh().cells()[cellI];
01271 
01272     // Points on loop
01273 
01274     // Status of point:
01275     // 0 - on loop
01276     // 1 - point set 1
01277     // 2 - point set 2
01278     Map<label> pointStatus(2*cPoints.size());
01279     Map<label> edgeStatus(2*cEdges.size());
01280 
01281     // Mark loop vertices
01282     forAll(loop, i)
01283     {
01284         label cut = loop[i];
01285 
01286         if (isEdge(cut))
01287         {
01288             edgeStatus.insert(getEdge(cut), 0);
01289         }
01290         else
01291         {
01292             pointStatus.insert(getVertex(cut), 0);
01293         }
01294     }
01295     // Since edges between two cut vertices have not been marked, mark them
01296     // explicitly
01297     forAll(cEdges, i)
01298     {
01299         label edgeI = cEdges[i];
01300         const edge& e = edges[edgeI];
01301 
01302         if (pointStatus.found(e[0]) && pointStatus.found(e[1]))
01303         {
01304             edgeStatus.insert(edgeI, 0);
01305         }
01306     }
01307 
01308 
01309     // Find uncut starting vertex
01310     label uncutIndex = firstUnique(cPoints, pointStatus);
01311 
01312     if (uncutIndex == -1)
01313     {
01314         WarningIn("Foam::cellCuts::calcAnchors")
01315             << "Invalid loop " << loop << " for cell " << cellI << endl
01316             << "Can not find point on cell which is not cut by loop."
01317             << endl;
01318 
01319         writeOBJ(".", cellI, loopPts, labelList(0));
01320 
01321         return false;
01322     }
01323 
01324     // Walk unset vertices and edges and mark with 1 in pointStatus, edgeStatus
01325     walkEdges(cellI, cPoints[uncutIndex], 1, edgeStatus, pointStatus);
01326 
01327     // Find new uncut starting vertex
01328     uncutIndex = firstUnique(cPoints, pointStatus);
01329 
01330     if (uncutIndex == -1)
01331     {
01332         // All vertices either in loop or in anchor. So split is along single
01333         // face.
01334         WarningIn("Foam::cellCuts::calcAnchors")
01335             << "Invalid loop " << loop << " for cell " << cellI << endl
01336             << "All vertices of cell are either in loop or in anchor set"
01337             << endl;
01338 
01339         writeOBJ(".", cellI, loopPts, labelList(0));
01340 
01341         return false;
01342     }
01343 
01344     // Walk unset vertices and edges and mark with 2. These are the
01345     // pointset 2.
01346     walkEdges(cellI, cPoints[uncutIndex], 2, edgeStatus, pointStatus);
01347 
01348     // Collect both sets in lists.
01349     DynamicList<label> connectedPoints(cPoints.size());
01350     DynamicList<label> otherPoints(cPoints.size());
01351 
01352     forAllConstIter(Map<label>, pointStatus, iter)
01353     {
01354         if (iter() == 1)
01355         {
01356             connectedPoints.append(iter.key());
01357         }
01358         else if (iter() == 2)
01359         {
01360             otherPoints.append(iter.key());
01361         }
01362     }
01363     connectedPoints.shrink();
01364     otherPoints.shrink();
01365 
01366     // Check that all points have been used.
01367     uncutIndex = firstUnique(cPoints, pointStatus);
01368 
01369     if (uncutIndex != -1)
01370     {
01371         WarningIn("Foam::cellCuts::calcAnchors")
01372             << "Invalid loop " << loop << " for cell " << cellI
01373             << " since it splits the cell into more than two cells" << endl;
01374 
01375         writeOBJ(".", cellI, loopPts, connectedPoints);
01376 
01377         return false;
01378     }
01379 
01380 
01381     // Check that both parts (connectedPoints, otherPoints) have enough faces.
01382     labelHashSet connectedFaces(2*cFaces.size());
01383     labelHashSet otherFaces(2*cFaces.size());
01384 
01385     forAllConstIter(Map<label>, pointStatus, iter)
01386     {
01387         label pointI = iter.key();
01388 
01389         const labelList& pFaces = mesh().pointFaces()[pointI];
01390 
01391         if (iter() == 1)
01392         {
01393             forAll(pFaces, pFaceI)
01394             {
01395                 if (meshTools::faceOnCell(mesh(), cellI, pFaces[pFaceI]))
01396                 {
01397                     connectedFaces.insert(pFaces[pFaceI]);
01398                 }
01399             }
01400         }
01401         else if (iter() == 2)
01402         {
01403             forAll(pFaces, pFaceI)
01404             {
01405                 if (meshTools::faceOnCell(mesh(), cellI, pFaces[pFaceI]))
01406                 {
01407                     otherFaces.insert(pFaces[pFaceI]);
01408                 }
01409             }
01410         }
01411     }
01412 
01413     if (connectedFaces.size() < 3)
01414     {
01415         WarningIn("Foam::cellCuts::calcAnchors")
01416             << "Invalid loop " << loop << " for cell " << cellI
01417             << " since would have too few faces on one side." << nl
01418             << "All faces:" << cFaces << endl;
01419 
01420         writeOBJ(".", cellI, loopPts, connectedPoints);
01421 
01422         return false;
01423     }
01424 
01425     if (otherFaces.size() < 3)
01426     {
01427         WarningIn("Foam::cellCuts::calcAnchors")
01428             << "Invalid loop " << loop << " for cell " << cellI
01429             << " since would have too few faces on one side." << nl
01430             << "All faces:" << cFaces << endl;
01431 
01432         writeOBJ(".", cellI, loopPts, otherPoints);
01433 
01434         return false;
01435     }
01436 
01437 
01438     // Check that faces are split into two regions and not more.
01439     // When walking across the face the only transition of pointStatus is
01440     // from set1 to loop to set2 or back. Not allowed is from set1 to loop to
01441     // set1.
01442     {
01443         forAll(cFaces, i)
01444         {
01445             label faceI = cFaces[i];
01446 
01447             const face& f = mesh().faces()[faceI];
01448 
01449             bool hasSet1 = false;
01450             bool hasSet2 = false;
01451 
01452             label prevStat = pointStatus[f[0]];
01453 
01454             forAll(f, fp)
01455             {
01456                 label v0 = f[fp];
01457                 label pStat = pointStatus[v0];
01458 
01459                 if (pStat == prevStat)
01460                 {
01461                 }
01462                 else if (pStat == 0)
01463                 {
01464                     // Loop.
01465                 }
01466                 else if (pStat == 1)
01467                 {
01468                     if (hasSet1)
01469                     {
01470                         // Second occurence of set1.
01471                         WarningIn("Foam::cellCuts::calcAnchors")
01472                             << "Invalid loop " << loop << " for cell " << cellI
01473                             << " since face " << f << " would be split into"
01474                             << " more than two faces" << endl;
01475 
01476                         writeOBJ(".", cellI, loopPts, otherPoints);
01477 
01478                         return false;
01479                     }
01480 
01481                     hasSet1 = true;
01482                 }
01483                 else if (pStat == 2)
01484                 {
01485                     if (hasSet2)
01486                     {
01487                         // Second occurence of set1.
01488                         WarningIn("Foam::cellCuts::calcAnchors")
01489                             << "Invalid loop " << loop << " for cell " << cellI
01490                             << " since face " << f << " would be split into"
01491                             << " more than two faces" << endl;
01492 
01493                         writeOBJ(".", cellI, loopPts, otherPoints);
01494 
01495                         return false;
01496                     }
01497 
01498                     hasSet2 = true;
01499                 }
01500                 else
01501                 {
01502                     FatalErrorIn("Foam::cellCuts::calcAnchors")
01503                         << abort(FatalError);
01504                 }
01505 
01506                 prevStat = pStat;
01507 
01508 
01509                 label v1 = f.nextLabel(fp);
01510                 label edgeI = findEdge(faceI, v0, v1);
01511 
01512                 label eStat = edgeStatus[edgeI];
01513 
01514                 if (eStat == prevStat)
01515                 {
01516                 }
01517                 else if (eStat == 0)
01518                 {
01519                     // Loop.
01520                 }
01521                 else if (eStat == 1)
01522                 {
01523                     if (hasSet1)
01524                     {
01525                         // Second occurence of set1.
01526                         WarningIn("Foam::cellCuts::calcAnchors")
01527                             << "Invalid loop " << loop << " for cell " << cellI
01528                             << " since face " << f << " would be split into"
01529                             << " more than two faces" << endl;
01530 
01531                         writeOBJ(".", cellI, loopPts, otherPoints);
01532 
01533                         return false;
01534                     }
01535 
01536                     hasSet1 = true;
01537                 }
01538                 else if (eStat == 2)
01539                 {
01540                     if (hasSet2)
01541                     {
01542                         // Second occurence of set1.
01543                         WarningIn("Foam::cellCuts::calcAnchors")
01544                             << "Invalid loop " << loop << " for cell " << cellI
01545                             << " since face " << f << " would be split into"
01546                             << " more than two faces" << endl;
01547 
01548                         writeOBJ(".", cellI, loopPts, otherPoints);
01549 
01550                         return false;
01551                     }
01552 
01553                     hasSet2 = true;
01554                 }
01555                 prevStat = eStat;
01556             }
01557         }
01558     }
01559 
01560 
01561 
01562 
01563     // Check which one of point sets to use.
01564     bool loopOk = loopAnchorConsistent(cellI, loopPts, connectedPoints);
01565 
01566     //if (debug)
01567     {
01568         // Additional check: are non-anchor points on other side?
01569         bool otherLoopOk = loopAnchorConsistent(cellI, loopPts, otherPoints);
01570 
01571         if (loopOk == otherLoopOk)
01572         {
01573             // Both sets of points are supposedly on the same side as the
01574             // loop normal. Oops.
01575 
01576             WarningIn("Foam::cellCuts::calcAnchors")
01577                 << " For cell:" << cellI
01578                 << " achorpoints and nonanchorpoints are geometrically"
01579                 << " on same side!" << endl
01580                 << "cellPoints:" << cPoints << endl
01581                 << "loop:" << loop << endl
01582                 << "anchors:" << connectedPoints << endl
01583                 << "otherPoints:" << otherPoints << endl;
01584 
01585             writeOBJ(".", cellI, loopPts, connectedPoints);
01586         }
01587     }
01588 
01589     if (loopOk)
01590     {
01591         // connectedPoints on 'outside' of loop so these are anchor points
01592         anchorPoints.transfer(connectedPoints);
01593         connectedPoints.clear();
01594     }
01595     else
01596     {
01597         anchorPoints.transfer(otherPoints);
01598         otherPoints.clear();
01599     }
01600     return true;
01601 }
01602 
01603 
01604 Foam::pointField Foam::cellCuts::loopPoints
01605 (
01606     const labelList& loop,
01607     const scalarField& loopWeights
01608 ) const
01609 {
01610     pointField loopPts(loop.size());
01611 
01612     forAll(loop, fp)
01613     {
01614         loopPts[fp] = coord(loop[fp], loopWeights[fp]);
01615     }
01616     return loopPts;
01617 }
01618 
01619 
01620 // Returns weights of loop. Inverse of loopPoints.
01621 Foam::scalarField Foam::cellCuts::loopWeights(const labelList& loop) const
01622 {
01623     scalarField weights(loop.size());
01624 
01625     forAll(loop, fp)
01626     {
01627         label cut = loop[fp];
01628 
01629         if (isEdge(cut))
01630         {
01631             label edgeI = getEdge(cut);
01632 
01633             weights[fp] = edgeWeight_[edgeI];
01634         }
01635         else
01636         {
01637             weights[fp] = -GREAT;
01638         }
01639     }
01640     return weights;
01641 }
01642 
01643 
01644 // Check if cut edges in loop are compatible with ones in edgeIsCut_
01645 bool Foam::cellCuts::validEdgeLoop
01646 (
01647     const labelList& loop,
01648     const scalarField& loopWeights
01649 ) const
01650 {
01651     forAll(loop, fp)
01652     {
01653         label cut = loop[fp];
01654 
01655         if (isEdge(cut))
01656         {
01657             label edgeI = getEdge(cut);
01658 
01659             // Check: cut compatible only if can be snapped to existing one.
01660             if (edgeIsCut_[edgeI])
01661             {
01662                 scalar edgeLen = mesh().edges()[edgeI].mag(mesh().points());
01663 
01664                 if
01665                 (
01666                     mag(loopWeights[fp] - edgeWeight_[edgeI])
01667                   > geomCellLooper::snapTol()*edgeLen
01668                 )
01669                 {
01670                     // Positions differ too much->would create two cuts.
01671                     return false;
01672                 }
01673             }
01674         }
01675     }
01676     return true;
01677 }
01678 
01679 
01680 // Counts cuts on face. Includes cuts through vertices and through edges.
01681 // Assumes that if edge is cut both in edgeIsCut and in loop that the position
01682 // of the cut is the same.
01683 Foam::label Foam::cellCuts::countFaceCuts
01684 (
01685     const label faceI,
01686     const labelList& loop
01687 ) const
01688 {
01689     label nCuts = 0;
01690 
01691     // Count cut vertices
01692     const face& f = mesh().faces()[faceI];
01693 
01694     forAll(f, fp)
01695     {
01696         label vertI = f[fp];
01697 
01698         // Vertex already cut or mentioned in current loop.
01699         if
01700         (
01701             pointIsCut_[vertI]
01702          || (findIndex(loop, vertToEVert(vertI)) != -1)
01703         )
01704         {
01705             nCuts++;
01706         }
01707     }
01708 
01709     // Count cut edges.
01710     const labelList& fEdges = mesh().faceEdges()[faceI];
01711 
01712     forAll(fEdges, fEdgeI)
01713     {
01714         label edgeI = fEdges[fEdgeI];
01715 
01716         // Edge already cut or mentioned in current loop.
01717         if
01718         (
01719             edgeIsCut_[edgeI]
01720          || (findIndex(loop, edgeToEVert(edgeI)) != -1)
01721         )
01722         {
01723             nCuts++;
01724         }
01725     }
01726 
01727     return nCuts;
01728 }
01729 
01730 
01731 // Determine compatibility of loop with existing cut pattern. Does not use
01732 // cut-addressing (faceCuts_, cutCuts_)
01733 bool Foam::cellCuts::conservativeValidLoop
01734 (
01735     const label cellI,
01736     const labelList& loop
01737 ) const
01738 {
01739 
01740     if (loop.size() < 2)
01741     {
01742         return false;
01743     }
01744 
01745     forAll(loop, cutI)
01746     {
01747         if (isEdge(loop[cutI]))
01748         {
01749             label edgeI = getEdge(loop[cutI]);
01750 
01751             if (edgeIsCut_[edgeI])
01752             {
01753                 // edge compatibility already checked.
01754             }
01755             else
01756             {
01757                 // Quick rejection: vertices of edge itself cannot be cut.
01758                 const edge& e = mesh().edges()[edgeI];
01759 
01760                 if (pointIsCut_[e.start()] || pointIsCut_[e.end()])
01761                 {
01762                     return false;
01763                 }
01764 
01765 
01766                 // Check faces using this edge
01767                 const labelList& eFaces = mesh().edgeFaces()[edgeI];
01768 
01769                 forAll(eFaces, eFaceI)
01770                 {
01771                     label nCuts = countFaceCuts(eFaces[eFaceI], loop);
01772 
01773                     if (nCuts > 2)
01774                     {
01775                         return false;
01776                     }
01777                 }
01778             }
01779         }
01780         else
01781         {
01782             // Vertex cut
01783 
01784             label vertI = getVertex(loop[cutI]);
01785 
01786             if (!pointIsCut_[vertI])
01787             {
01788                 // New cut through vertex.
01789 
01790                 // Check edges using vertex.
01791                 const labelList& pEdges = mesh().pointEdges()[vertI];
01792 
01793                 forAll(pEdges, pEdgeI)
01794                 {
01795                     label edgeI = pEdges[pEdgeI];
01796 
01797                     if (edgeIsCut_[edgeI])
01798                     {
01799                         return false;
01800                     }
01801                 }
01802 
01803                 // Check faces using vertex.
01804                 const labelList& pFaces = mesh().pointFaces()[vertI];
01805 
01806                 forAll(pFaces, pFaceI)
01807                 {
01808                     label nCuts = countFaceCuts(pFaces[pFaceI], loop);
01809 
01810                     if (nCuts > 2)
01811                     {
01812                         return false;
01813                     }
01814                 }
01815             }
01816         }
01817     }
01818 
01819     // All ok.
01820     return true;
01821 }
01822 
01823 
01824 // Determine compatibility of loop with existing cut pattern. Does not use
01825 // derived cut-addressing (faceCuts), only pointIsCut, edgeIsCut.
01826 // Adds any cross-cuts found to newFaceSplitCut and sets cell points on
01827 // one side of the loop in anchorPoints.
01828 bool Foam::cellCuts::validLoop
01829 (
01830     const label cellI,
01831     const labelList& loop,
01832     const scalarField& loopWeights,
01833 
01834     Map<edge>& newFaceSplitCut,
01835     labelList& anchorPoints
01836 ) const
01837 {
01838     if (loop.size() < 2)
01839     {
01840         return false;
01841     }
01842 
01843     if (debug & 4)
01844     {
01845         // Allow as fallback the 'old' loop checking where only a single
01846         // cut per face is allowed.
01847         if (!conservativeValidLoop(cellI, loop))
01848         {
01849             return  false;
01850         }
01851     }
01852 
01853     forAll(loop, fp)
01854     {
01855         label cut = loop[fp];
01856         label nextCut = loop[(fp+1) % loop.size()];
01857 
01858         // Label (if any) of face cut (so cut not along existing edge)
01859         label meshFaceI = -1;
01860 
01861         if (isEdge(cut))
01862         {
01863             label edgeI = getEdge(cut);
01864 
01865             // Look one cut ahead to find if it is along existing edge.
01866 
01867             if (isEdge(nextCut))
01868             {
01869                 // From edge to edge -> cross cut
01870                 label nextEdgeI = getEdge(nextCut);
01871 
01872                 // Find face and mark as to be split.
01873                 meshFaceI = edgeEdgeToFace(cellI, edgeI, nextEdgeI);
01874 
01875                 if (meshFaceI == -1)
01876                 {
01877                     // Can't find face using both cut edges.
01878                     return false;
01879                 }
01880             }
01881             else
01882             {
01883                 // From edge to vertex -> cross cut only if no existing edge.
01884 
01885                 label nextVertI = getVertex(nextCut);
01886                 const edge& e = mesh().edges()[edgeI];
01887 
01888                 if (e.start() != nextVertI && e.end() != nextVertI)
01889                 {
01890                     // New edge. Find face and mark as to be split.
01891                     meshFaceI = edgeVertexToFace(cellI, edgeI, nextVertI);
01892 
01893                     if (meshFaceI == -1)
01894                     {
01895                         // Can't find face. Ilegal.
01896                         return false;
01897                     }
01898                 }
01899             }
01900         }
01901         else
01902         {
01903             // Cut is vertex.
01904             label vertI = getVertex(cut);
01905 
01906             if (isEdge(nextCut))
01907             {
01908                 // From vertex to edge -> cross cut only if no existing edge.
01909                 label nextEdgeI = getEdge(nextCut);
01910 
01911                 const edge& nextE = mesh().edges()[nextEdgeI];
01912 
01913                 if (nextE.start() != vertI && nextE.end() != vertI)
01914                 {
01915                     // Should be cross cut. Find face.
01916                     meshFaceI = edgeVertexToFace(cellI, nextEdgeI, vertI);
01917 
01918                     if (meshFaceI == -1)
01919                     {
01920                         return false;
01921                     }
01922                 }
01923             }
01924             else
01925             {
01926                 // From vertex to vertex -> cross cut only if no existing edge.
01927                 label nextVertI = getVertex(nextCut);
01928 
01929                 if (meshTools::findEdge(mesh(), vertI, nextVertI) == -1)
01930                 {
01931                     // New edge. Find face.
01932                     meshFaceI = vertexVertexToFace(cellI, vertI, nextVertI);
01933 
01934                     if (meshFaceI == -1)
01935                     {
01936                         return false;
01937                     }
01938                 }
01939             }
01940         }
01941 
01942         if (meshFaceI != -1)
01943         {
01944             // meshFaceI is cut across along cut-nextCut (so not along existing
01945             // edge). Check if this is compatible with existing pattern.
01946             edge cutEdge(cut, nextCut);
01947 
01948             Map<edge>::const_iterator iter = faceSplitCut_.find(meshFaceI);
01949 
01950             if (iter == faceSplitCut_.end())
01951             {
01952                 // Face not yet cut so insert.
01953                 newFaceSplitCut.insert(meshFaceI, cutEdge);
01954             }
01955             else
01956             {
01957                 // Face already cut. Ok if same edge.
01958                 if (iter() != cutEdge)
01959                 {
01960                     return false;
01961                 }
01962             }
01963         }
01964     }
01965 
01966     // Is there a face on which all cuts are?
01967     label faceContainingLoop = loopFace(cellI, loop);
01968 
01969     if (faceContainingLoop != -1)
01970     {
01971         WarningIn("Foam::cellCuts::validLoop")
01972             << "Found loop on cell " << cellI << " with all points"
01973             << " on face " << faceContainingLoop << endl;
01974 
01975         //writeOBJ(".", cellI, loopPoints(loop, loopWeights), labelList(0));
01976 
01977         return false;
01978     }
01979 
01980     // Calculate anchor points
01981     // Final success is determined by whether anchor points can be determined.
01982     return calcAnchors
01983     (
01984         cellI,
01985         loop,
01986         loopPoints(loop, loopWeights),
01987         anchorPoints
01988     );
01989 }
01990 
01991 
01992 // Update basic cut information (pointIsCut, edgeIsCut) from cellLoops.
01993 // Assumes cellLoops_ and edgeWeight_ already set and consistent.
01994 // Does not use any other information.
01995 void Foam::cellCuts::setFromCellLoops()
01996 {
01997     // 'Uncut' edges/vertices that are not used in loops.
01998     pointIsCut_ = false;
01999 
02000     edgeIsCut_ = false;
02001 
02002     faceSplitCut_.clear();
02003 
02004     forAll(cellLoops_, cellI)
02005     {
02006         const labelList& loop = cellLoops_[cellI];
02007 
02008         if (loop.size())
02009         {
02010             // Storage for cross-face cuts
02011             Map<edge> faceSplitCuts(loop.size());
02012             // Storage for points on one side of cell.
02013             labelList anchorPoints;
02014 
02015             if
02016             (
02017                !validLoop
02018                 (
02019                     cellI,
02020                     loop,
02021                     loopWeights(loop),
02022                     faceSplitCuts,
02023                     anchorPoints
02024                 )
02025             )
02026             {
02027                 //writeOBJ(".", cellI, loopPoints(cellI), anchorPoints);
02028 
02029                 //FatalErrorIn("cellCuts::setFromCellLoops()")
02030                 WarningIn("cellCuts::setFromCellLoops")
02031                     << "Illegal loop " << loop
02032                     << " when recreating cut-addressing"
02033                     << " from existing cellLoops for cell " << cellI
02034                     << endl;
02035                     //<< abort(FatalError);
02036 
02037                 cellLoops_[cellI].setSize(0);
02038                 cellAnchorPoints_[cellI].setSize(0);
02039             }
02040             else
02041             {
02042                 // Copy anchor points.
02043                 cellAnchorPoints_[cellI].transfer(anchorPoints);
02044 
02045                 // Copy faceSplitCuts into overall faceSplit info.
02046                 forAllConstIter(Map<edge>, faceSplitCuts, iter)
02047                 {
02048                     faceSplitCut_.insert(iter.key(), iter());
02049                 }
02050 
02051                 // Update edgeIsCut, pointIsCut information
02052                 forAll(loop, cutI)
02053                 {
02054                     label cut = loop[cutI];
02055 
02056                     if (isEdge(cut))
02057                     {
02058                         edgeIsCut_[getEdge(cut)] = true;
02059                     }
02060                     else
02061                     {
02062                         pointIsCut_[getVertex(cut)] = true;
02063                     }
02064                 }
02065             }
02066         }
02067     }
02068 
02069     // Reset edge weights
02070     forAll(edgeIsCut_, edgeI)
02071     {
02072         if (!edgeIsCut_[edgeI])
02073         {
02074             // Weight not used. Set to illegal value to make any use fall over.
02075             edgeWeight_[edgeI] = -GREAT;
02076         }
02077     }
02078 }
02079 
02080 
02081 // Upate basic cut information from single cellLoop. Returns true if loop
02082 // was valid.
02083 bool Foam::cellCuts::setFromCellLoop
02084 (
02085     const label cellI,
02086     const labelList& loop,
02087     const scalarField& loopWeights
02088 )
02089 {
02090     // Dump loop for debugging.
02091     if (debug)
02092     {
02093         OFstream str("last_cell.obj");
02094 
02095         str<< "# edges of cell " << cellI << nl;
02096 
02097         meshTools::writeOBJ
02098         (
02099             str,
02100             mesh().cells(),
02101             mesh().faces(),
02102             mesh().points(),
02103             labelList(1, cellI)
02104         );
02105 
02106 
02107         OFstream loopStr("last_loop.obj");
02108 
02109         loopStr<< "# looppoints for cell " << cellI << nl;
02110 
02111         pointField pointsOfLoop = loopPoints(loop, loopWeights);
02112 
02113         forAll(pointsOfLoop, i)
02114         {
02115             meshTools::writeOBJ(loopStr, pointsOfLoop[i]);
02116         }
02117 
02118         str << 'l';
02119 
02120         forAll(pointsOfLoop, i)
02121         {
02122             str << ' ' << i + 1;
02123         }
02124         str << ' ' << 1 << nl;
02125     }
02126 
02127     bool okLoop = false;
02128 
02129     if (validEdgeLoop(loop, loopWeights))
02130     {
02131         // Storage for cuts across face
02132         Map<edge> faceSplitCuts(loop.size());
02133         // Storage for points on one side of cell
02134         labelList anchorPoints;
02135 
02136         okLoop =
02137             validLoop(cellI, loop, loopWeights, faceSplitCuts, anchorPoints);
02138 
02139         if (okLoop)
02140         {
02141             // Valid loop. Copy cellLoops and anchorPoints
02142             cellLoops_[cellI] = loop;
02143             cellAnchorPoints_[cellI].transfer(anchorPoints);
02144 
02145             // Copy split cuts
02146             forAllConstIter(Map<edge>, faceSplitCuts, iter)
02147             {
02148                 faceSplitCut_.insert(iter.key(), iter());
02149             }
02150 
02151 
02152             // Update edgeIsCut, pointIsCut information
02153             forAll(loop, cutI)
02154             {
02155                 label cut = loop[cutI];
02156 
02157                 if (isEdge(cut))
02158                 {
02159                     label edgeI = getEdge(cut);
02160 
02161                     edgeIsCut_[edgeI] = true;
02162 
02163                     edgeWeight_[edgeI] = loopWeights[cutI];
02164                 }
02165                 else
02166                 {
02167                     label vertI = getVertex(cut);
02168 
02169                     pointIsCut_[vertI] = true;
02170                 }
02171             }
02172         }
02173     }
02174 
02175     return okLoop;
02176 }
02177 
02178 
02179 // Update basic cut information from cellLoops. Checks for consistency with
02180 // existing cut pattern.
02181 void Foam::cellCuts::setFromCellLoops
02182 (
02183     const labelList& cellLabels,
02184     const labelListList& cellLoops,
02185     const List<scalarField>& cellLoopWeights
02186 )
02187 {
02188     // 'Uncut' edges/vertices that are not used in loops.
02189     pointIsCut_ = false;
02190 
02191     edgeIsCut_ = false;
02192 
02193     forAll(cellLabels, cellLabelI)
02194     {
02195         label cellI = cellLabels[cellLabelI];
02196 
02197         const labelList& loop = cellLoops[cellLabelI];
02198 
02199         if (loop.size())
02200         {
02201             const scalarField& loopWeights = cellLoopWeights[cellLabelI];
02202 
02203             if (setFromCellLoop(cellI, loop, loopWeights))
02204             {
02205                 // Valid loop. Call above will have upated all already.
02206             }
02207             else
02208             {
02209                 // Clear cellLoops
02210                 cellLoops_[cellI].setSize(0);
02211             }
02212         }
02213     }
02214 }
02215 
02216 
02217 // Cut cells and update basic cut information from cellLoops. Checks each loop
02218 // for consistency with existing cut pattern.
02219 void Foam::cellCuts::setFromCellCutter
02220 (
02221     const cellLooper& cellCutter,
02222     const List<refineCell>& refCells
02223 )
02224 {
02225     // 'Uncut' edges/vertices that are not used in loops.
02226     pointIsCut_ = false;
02227 
02228     edgeIsCut_ = false;
02229 
02230     // storage for loop of cuts (cut vertices and/or cut edges)
02231     labelList cellLoop;
02232     scalarField cellLoopWeights;
02233 
02234     // For debugging purposes
02235     DynamicList<label> invalidCutCells(2);
02236     DynamicList<labelList> invalidCutLoops(2);
02237     DynamicList<scalarField> invalidCutLoopWeights(2);
02238 
02239 
02240     forAll(refCells, refCellI)
02241     {
02242         const refineCell& refCell = refCells[refCellI];
02243 
02244         label cellI = refCell.cellNo();
02245 
02246         const vector& refDir = refCell.direction();
02247 
02248 
02249         // Cut cell. Determines cellLoop and cellLoopWeights
02250         bool goodCut =
02251             cellCutter.cut
02252             (
02253                 refDir,
02254                 cellI,
02255 
02256                 pointIsCut_,
02257                 edgeIsCut_,
02258                 edgeWeight_,
02259 
02260                 cellLoop,
02261                 cellLoopWeights
02262             );
02263 
02264         // Check whether edge refinement is on a per face basis compatible with
02265         // current pattern.
02266         if (goodCut)
02267         {
02268             if (setFromCellLoop(cellI, cellLoop, cellLoopWeights))
02269             {
02270                 // Valid loop. Will have updated all info already.
02271             }
02272             else
02273             {
02274                 cellLoops_[cellI].setSize(0);
02275 
02276                 // Discarded by validLoop
02277                 if (debug)
02278                 {
02279                     invalidCutCells.append(cellI);
02280                     invalidCutLoops.append(cellLoop);
02281                     invalidCutLoopWeights.append(cellLoopWeights);
02282                 }
02283             }
02284         }
02285         else
02286         {
02287             // Clear cellLoops
02288             cellLoops_[cellI].setSize(0);
02289         }
02290     }
02291 
02292     if (debug && invalidCutCells.size())
02293     {
02294         invalidCutCells.shrink();
02295         invalidCutLoops.shrink();
02296         invalidCutLoopWeights.shrink();
02297 
02298         fileName cutsFile("invalidLoopCells.obj");
02299 
02300         Pout<< "cellCuts : writing inValidLoops cells to " << cutsFile << endl;
02301 
02302         OFstream cutsStream(cutsFile);
02303 
02304         meshTools::writeOBJ
02305         (
02306             cutsStream,
02307             mesh().cells(),
02308             mesh().faces(),
02309             mesh().points(),
02310             invalidCutCells
02311         );
02312 
02313         fileName loopsFile("invalidLoops.obj");
02314 
02315         Pout<< "cellCuts : writing inValidLoops loops to " << loopsFile << endl;
02316 
02317         OFstream loopsStream(loopsFile);
02318 
02319         label vertI = 0;
02320 
02321         forAll(invalidCutLoops, i)
02322         {
02323             writeOBJ
02324             (
02325                 loopsStream,
02326                 loopPoints(invalidCutLoops[i], invalidCutLoopWeights[i]),
02327                 vertI
02328             );
02329         }
02330     }
02331 }
02332 
02333 
02334 // Same as one before but cut plane prescribed (instead of just normal)
02335 void Foam::cellCuts::setFromCellCutter
02336 (
02337     const cellLooper& cellCutter,
02338     const labelList& cellLabels,
02339     const List<plane>& cellCutPlanes
02340 )
02341 {
02342     // 'Uncut' edges/vertices that are not used in loops.
02343     pointIsCut_ = false;
02344 
02345     edgeIsCut_ = false;
02346 
02347     // storage for loop of cuts (cut vertices and/or cut edges)
02348     labelList cellLoop;
02349     scalarField cellLoopWeights;
02350 
02351     // For debugging purposes
02352     DynamicList<label> invalidCutCells(2);
02353     DynamicList<labelList> invalidCutLoops(2);
02354     DynamicList<scalarField> invalidCutLoopWeights(2);
02355 
02356 
02357     forAll(cellLabels, i)
02358     {
02359         label cellI = cellLabels[i];
02360 
02361         // Cut cell. Determines cellLoop and cellLoopWeights
02362         bool goodCut =
02363             cellCutter.cut
02364             (
02365                 cellCutPlanes[i],
02366                 cellI,
02367 
02368                 pointIsCut_,
02369                 edgeIsCut_,
02370                 edgeWeight_,
02371 
02372                 cellLoop,
02373                 cellLoopWeights
02374             );
02375 
02376         // Check whether edge refinement is on a per face basis compatible with
02377         // current pattern.
02378         if (goodCut)
02379         {
02380             if (setFromCellLoop(cellI, cellLoop, cellLoopWeights))
02381             {
02382                 // Valid loop. Will have updated all info already.
02383             }
02384             else
02385             {
02386                 cellLoops_[cellI].setSize(0);
02387 
02388                 // Discarded by validLoop
02389                 if (debug)
02390                 {
02391                     invalidCutCells.append(cellI);
02392                     invalidCutLoops.append(cellLoop);
02393                     invalidCutLoopWeights.append(cellLoopWeights);
02394                 }
02395             }
02396         }
02397         else
02398         {
02399             // Clear cellLoops
02400             cellLoops_[cellI].setSize(0);
02401         }
02402     }
02403 
02404     if (debug && invalidCutCells.size())
02405     {
02406         invalidCutCells.shrink();
02407         invalidCutLoops.shrink();
02408         invalidCutLoopWeights.shrink();
02409 
02410         fileName cutsFile("invalidLoopCells.obj");
02411 
02412         Pout<< "cellCuts : writing inValidLoops cells to " << cutsFile << endl;
02413 
02414         OFstream cutsStream(cutsFile);
02415 
02416         meshTools::writeOBJ
02417         (
02418             cutsStream,
02419             mesh().cells(),
02420             mesh().faces(),
02421             mesh().points(),
02422             invalidCutCells
02423         );
02424 
02425         fileName loopsFile("invalidLoops.obj");
02426 
02427         Pout<< "cellCuts : writing inValidLoops loops to " << loopsFile << endl;
02428 
02429         OFstream loopsStream(loopsFile);
02430 
02431         label vertI = 0;
02432 
02433         forAll(invalidCutLoops, i)
02434         {
02435             writeOBJ
02436             (
02437                 loopsStream,
02438                 loopPoints(invalidCutLoops[i], invalidCutLoopWeights[i]),
02439                 vertI
02440             );
02441         }
02442     }
02443 }
02444 
02445 
02446 // Set orientation of loops
02447 void Foam::cellCuts::orientPlanesAndLoops()
02448 {
02449     // Determine anchorPoints if not yet done by validLoop.
02450     forAll(cellLoops_, cellI)
02451     {
02452         const labelList& loop = cellLoops_[cellI];
02453 
02454         if (loop.size() && cellAnchorPoints_[cellI].empty())
02455         {
02456             // Leave anchor points empty if illegal loop.
02457             calcAnchors
02458             (
02459                 cellI,
02460                 loop,
02461                 loopPoints(cellI),
02462                 cellAnchorPoints_[cellI]
02463             );
02464         }
02465     }
02466 
02467     if (debug & 2)
02468     {
02469         Pout<< "cellAnchorPoints:" << endl;
02470     }
02471     forAll(cellAnchorPoints_, cellI)
02472     {
02473         if (cellLoops_[cellI].size())
02474         {
02475             if (cellAnchorPoints_[cellI].empty())
02476             {
02477                 FatalErrorIn("orientPlanesAndLoops()")
02478                     << "No anchor points for cut cell " << cellI << endl
02479                     << "cellLoop:" << cellLoops_[cellI] << abort(FatalError);
02480             }
02481 
02482             if (debug & 2)
02483             {
02484                 Pout<< "    cell:" << cellI << " anchored at "
02485                     << cellAnchorPoints_[cellI] << endl;
02486             }
02487         }
02488     }
02489 
02490     // Calculate number of valid cellLoops
02491     nLoops_ = 0;
02492 
02493     forAll(cellLoops_, cellI)
02494     {
02495         if (cellLoops_[cellI].size())
02496         {
02497             nLoops_++;
02498         }
02499     }
02500 }
02501 
02502 
02503 // Do all: calculate addressing, calculate loops splitting cells
02504 void Foam::cellCuts::calcLoopsAndAddressing(const labelList& cutCells)
02505 {
02506     // Sanity check on weights
02507     forAll(edgeIsCut_, edgeI)
02508     {
02509         if (edgeIsCut_[edgeI])
02510         {
02511             scalar weight = edgeWeight_[edgeI];
02512 
02513             if (weight < 0 || weight > 1)
02514             {
02515                 FatalErrorIn
02516                 (
02517                     "cellCuts::calcLoopsAndAddressing(const labelList&)"
02518                 )   << "Weight out of range [0,1]. Edge " << edgeI
02519                     << " verts:" << mesh().edges()[edgeI]
02520                     << " weight:" << weight << abort(FatalError);
02521             }
02522         }
02523         else
02524         {
02525             // Weight not used. Set to illegal value to make any use fall over.
02526             edgeWeight_[edgeI] = -GREAT;
02527         }
02528     }
02529 
02530 
02531     // Calculate faces that split cells in two
02532     calcCellLoops(cutCells);
02533 
02534     if (debug & 2)
02535     {
02536         Pout<< "-- cellLoops --" << endl;
02537         forAll(cellLoops_, cellI)
02538         {
02539             const labelList& loop = cellLoops_[cellI];
02540 
02541             if (loop.size())
02542             {
02543                 Pout<< "cell:" << cellI << "  ";
02544                 writeCuts(Pout, loop, loopWeights(loop));
02545                 Pout<< endl;
02546             }
02547         }
02548     }
02549 
02550     // Redo basic cut information (pointIsCut, edgeIsCut, faceSplitCut)
02551     // using cellLoop only.
02552     setFromCellLoops();
02553 }
02554 
02555 
02556 void Foam::cellCuts::check() const
02557 {
02558     // Check weights for unsnapped values
02559     forAll(edgeIsCut_, edgeI)
02560     {
02561         if (edgeIsCut_[edgeI])
02562         {
02563             if
02564             (
02565                 edgeWeight_[edgeI] < geomCellLooper::snapTol()
02566              || edgeWeight_[edgeI] > (1 - geomCellLooper::snapTol())
02567             )
02568             {
02569                 // Should have been snapped.
02570                 //FatalErrorIn("cellCuts::check()")
02571                 WarningIn("cellCuts::check()")
02572                     << "edge:" << edgeI << " vertices:"
02573                     << mesh().edges()[edgeI]
02574                     << " weight:" << edgeWeight_[edgeI] << " should have been"
02575                     << " snapped to one of its endpoints"
02576                     //<< abort(FatalError);
02577                     << endl;
02578             }
02579         }
02580         else
02581         {
02582             if (edgeWeight_[edgeI] > - 1)
02583             {
02584                 FatalErrorIn("cellCuts::check()")
02585                     << "edge:" << edgeI << " vertices:"
02586                     << mesh().edges()[edgeI]
02587                     << " weight:" << edgeWeight_[edgeI] << " is not cut but"
02588                     << " its weight is not marked invalid"
02589                     << abort(FatalError);
02590             }
02591         }
02592     }
02593 
02594     // Check that all elements of cellloop are registered
02595     forAll(cellLoops_, cellI)
02596     {
02597         const labelList& loop = cellLoops_[cellI];
02598 
02599         forAll(loop, i)
02600         {
02601             label cut = loop[i];
02602 
02603             if
02604             (
02605                 (isEdge(cut) && !edgeIsCut_[getEdge(cut)])
02606              || (!isEdge(cut) && !pointIsCut_[getVertex(cut)])
02607             )
02608             {
02609                 labelList cuts(1, cut);
02610                 writeCuts(Pout, cuts, loopWeights(cuts));
02611 
02612                 FatalErrorIn("cellCuts::check()")
02613                     << "cell:" << cellI << " loop:"
02614                     << loop
02615                     << " cut:" << cut << " is not marked as cut"
02616                     << abort(FatalError);
02617             }
02618         }
02619     }
02620 
02621     // Check that no elements of cell loop are anchor point.
02622     forAll(cellLoops_, cellI)
02623     {
02624         const labelList& anchors = cellAnchorPoints_[cellI];
02625 
02626         const labelList& loop = cellLoops_[cellI];
02627 
02628         if (loop.size() && anchors.empty())
02629         {
02630             FatalErrorIn("cellCuts::check()")
02631                 << "cell:" << cellI << " loop:" << loop
02632                 << " has no anchor points"
02633                 << abort(FatalError);
02634         }
02635 
02636 
02637         forAll(loop, i)
02638         {
02639             label cut = loop[i];
02640 
02641             if
02642             (
02643                 !isEdge(cut)
02644              && findIndex(anchors, getVertex(cut)) != -1
02645             )
02646             {
02647                 FatalErrorIn("cellCuts::check()")
02648                     << "cell:" << cellI << " loop:" << loop
02649                     << " anchor points:" << anchors
02650                     << " anchor:" << getVertex(cut) << " is part of loop"
02651                     << abort(FatalError);
02652             }
02653         }
02654     }
02655 
02656 
02657     // Check that cut faces have a neighbour that is cut.
02658     forAllConstIter(Map<edge>, faceSplitCut_, iter)
02659     {
02660         label faceI = iter.key();
02661 
02662         if (mesh().isInternalFace(faceI))
02663         {
02664             label own = mesh().faceOwner()[faceI];
02665             label nei = mesh().faceNeighbour()[faceI];
02666 
02667             if (cellLoops_[own].empty() && cellLoops_[nei].empty())
02668             {
02669                 FatalErrorIn("cellCuts::check()")
02670                     << "Internal face:" << faceI << " cut by " << iter()
02671                     << " has owner:" << own
02672                     << " and neighbour:" << nei
02673                     << " that are both uncut"
02674                     << abort(FatalError);
02675             }
02676         }
02677         else
02678         {
02679             label own = mesh().faceOwner()[faceI];
02680 
02681             if (cellLoops_[own].empty())
02682             {
02683                 FatalErrorIn("cellCuts::check()")
02684                     << "Boundary face:" << faceI << " cut by " << iter()
02685                     << " has owner:" << own
02686                     << " that is uncut"
02687                     << abort(FatalError);
02688             }
02689         }
02690     }
02691 }
02692 
02693 
02694 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
02695 
02696 // Construct from cells to cut and pattern of cuts
02697 Foam::cellCuts::cellCuts
02698 (
02699     const polyMesh& mesh,
02700     const labelList& cutCells,
02701     const labelList& meshVerts,
02702     const labelList& meshEdges,
02703     const scalarField& meshEdgeWeights
02704 )
02705 :
02706     edgeVertex(mesh),
02707     pointIsCut_(expand(mesh.nPoints(), meshVerts)),
02708     edgeIsCut_(expand(mesh.nEdges(), meshEdges)),
02709     edgeWeight_(expand(mesh.nEdges(), meshEdges, meshEdgeWeights)),
02710     faceCutsPtr_(NULL),
02711     faceSplitCut_(cutCells.size()),
02712     cellLoops_(mesh.nCells()),
02713     nLoops_(-1),
02714     cellAnchorPoints_(mesh.nCells())
02715 {
02716     if (debug)
02717     {
02718         Pout<< "cellCuts : constructor from cut verts and edges" << endl;
02719     }
02720 
02721     calcLoopsAndAddressing(cutCells);
02722 
02723     // Calculate planes and flip cellLoops if nessecary
02724     orientPlanesAndLoops();
02725 
02726     if (debug)
02727     {
02728         check();
02729     }
02730 
02731     clearOut();
02732 
02733     if (debug)
02734     {
02735         Pout<< "cellCuts : leaving constructor from cut verts and edges"
02736             << endl;
02737     }
02738 }
02739 
02740 
02741 // Construct from pattern of cuts. Finds out itself which cells are cut.
02742 // (can go wrong if e.g. all neighbours of cell are refined)
02743 Foam::cellCuts::cellCuts
02744 (
02745     const polyMesh& mesh,
02746     const labelList& meshVerts,
02747     const labelList& meshEdges,
02748     const scalarField& meshEdgeWeights
02749 )
02750 :
02751     edgeVertex(mesh),
02752     pointIsCut_(expand(mesh.nPoints(), meshVerts)),
02753     edgeIsCut_(expand(mesh.nEdges(), meshEdges)),
02754     edgeWeight_(expand(mesh.nEdges(), meshEdges, meshEdgeWeights)),
02755     faceCutsPtr_(NULL),
02756     faceSplitCut_(mesh.nFaces()/10 + 1),
02757     cellLoops_(mesh.nCells()),
02758     nLoops_(-1),
02759     cellAnchorPoints_(mesh.nCells())
02760 {
02761     if (debug)
02762     {
02763         Pout<< "cellCuts : constructor from cellLoops" << endl;
02764     }
02765 
02766     calcLoopsAndAddressing(identity(mesh.nCells()));
02767 
02768     // Calculate planes and flip cellLoops if nessecary
02769     orientPlanesAndLoops();
02770 
02771     if (debug)
02772     {
02773         check();
02774     }
02775 
02776     clearOut();
02777 
02778     if (debug)
02779     {
02780         Pout<< "cellCuts : leaving constructor from cellLoops" << endl;
02781     }
02782 }
02783 
02784 
02785 // Construct from complete cellLoops. Assumes correct cut pattern.
02786 // Only constructs cut-cut addressing and cellAnchorPoints
02787 Foam::cellCuts::cellCuts
02788 (
02789     const polyMesh& mesh,
02790     const labelList& cellLabels,
02791     const labelListList& cellLoops,
02792     const List<scalarField>& cellEdgeWeights
02793 )
02794 :
02795     edgeVertex(mesh),
02796     pointIsCut_(mesh.nPoints(), false),
02797     edgeIsCut_(mesh.nEdges(), false),
02798     edgeWeight_(mesh.nEdges(), -GREAT),
02799     faceCutsPtr_(NULL),
02800     faceSplitCut_(cellLabels.size()),
02801     cellLoops_(mesh.nCells()),
02802     nLoops_(-1),
02803     cellAnchorPoints_(mesh.nCells())
02804 {
02805     if (debug)
02806     {
02807         Pout<< "cellCuts : constructor from cellLoops" << endl;
02808     }
02809 
02810     // Update pointIsCut, edgeIsCut, faceSplitCut from cell loops.
02811     // Makes sure cuts are consistent
02812     setFromCellLoops(cellLabels, cellLoops, cellEdgeWeights);
02813 
02814     // Calculate planes and flip cellLoops if nessecary
02815     orientPlanesAndLoops();
02816 
02817     if (debug)
02818     {
02819         check();
02820     }
02821 
02822     clearOut();
02823 
02824     if (debug)
02825     {
02826         Pout<< "cellCuts : leaving constructor from cellLoops" << endl;
02827     }
02828 }
02829 
02830 
02831 // Construct from list of cells to cut and cell cutter.
02832 Foam::cellCuts::cellCuts
02833 (
02834     const polyMesh& mesh,
02835     const cellLooper& cellCutter,
02836     const List<refineCell>& refCells
02837 )
02838 :
02839     edgeVertex(mesh),
02840     pointIsCut_(mesh.nPoints(), false),
02841     edgeIsCut_(mesh.nEdges(), false),
02842     edgeWeight_(mesh.nEdges(), -GREAT),
02843     faceCutsPtr_(NULL),
02844     faceSplitCut_(refCells.size()),
02845     cellLoops_(mesh.nCells()),
02846     nLoops_(-1),
02847     cellAnchorPoints_(mesh.nCells())
02848 {
02849     if (debug)
02850     {
02851         Pout<< "cellCuts : constructor from cellCutter" << endl;
02852     }
02853 
02854     // Update pointIsCut, edgeIsCut, faceSplitCut from cell loops.
02855     // Makes sure cuts are consistent
02856     setFromCellCutter(cellCutter, refCells);
02857 
02858     // Calculate planes and flip cellLoops if nessecary
02859     orientPlanesAndLoops();
02860 
02861     if (debug)
02862     {
02863         check();
02864     }
02865 
02866     clearOut();
02867 
02868     if (debug)
02869     {
02870         Pout<< "cellCuts : leaving constructor from cellCutter" << endl;
02871     }
02872 }
02873 
02874 
02875 // Construct from list of cells to cut, plane to cut with and cell cutter.
02876 Foam::cellCuts::cellCuts
02877 (
02878     const polyMesh& mesh,
02879     const cellLooper& cellCutter,
02880     const labelList& cellLabels,
02881     const List<plane>& cutPlanes
02882 )
02883 :
02884     edgeVertex(mesh),
02885     pointIsCut_(mesh.nPoints(), false),
02886     edgeIsCut_(mesh.nEdges(), false),
02887     edgeWeight_(mesh.nEdges(), -GREAT),
02888     faceCutsPtr_(NULL),
02889     faceSplitCut_(cellLabels.size()),
02890     cellLoops_(mesh.nCells()),
02891     nLoops_(-1),
02892     cellAnchorPoints_(mesh.nCells())
02893 {
02894     if (debug)
02895     {
02896         Pout<< "cellCuts : constructor from cellCutter with prescribed plane"
02897             << endl;
02898     }
02899 
02900     // Update pointIsCut, edgeIsCut, faceSplitCut from cell loops.
02901     // Makes sure cuts are consistent
02902     setFromCellCutter(cellCutter, cellLabels, cutPlanes);
02903 
02904     // Calculate planes and flip cellLoops if nessecary
02905     orientPlanesAndLoops();
02906 
02907     if (debug)
02908     {
02909         check();
02910     }
02911 
02912     clearOut();
02913 
02914     if (debug)
02915     {
02916         Pout<< "cellCuts : leaving constructor from cellCutter with prescribed"
02917             << " plane" << endl;
02918     }
02919 }
02920 
02921 
02922 // Construct from components
02923 Foam::cellCuts::cellCuts
02924 (
02925     const polyMesh& mesh,
02926     const boolList& pointIsCut,
02927     const boolList& edgeIsCut,
02928     const scalarField& edgeWeight,
02929     const Map<edge>& faceSplitCut,
02930     const labelListList& cellLoops,
02931     const label nLoops,
02932     const labelListList& cellAnchorPoints
02933 )
02934 :
02935     edgeVertex(mesh),
02936     pointIsCut_(pointIsCut),
02937     edgeIsCut_(edgeIsCut),
02938     edgeWeight_(edgeWeight),
02939     faceCutsPtr_(NULL),
02940     faceSplitCut_(faceSplitCut),
02941     cellLoops_(cellLoops),
02942     nLoops_(nLoops),
02943     cellAnchorPoints_(cellAnchorPoints)
02944 {
02945     if (debug)
02946     {
02947         Pout<< "cellCuts : constructor from components" << endl;
02948         Pout<< "cellCuts : leaving constructor from components" << endl;
02949     }
02950 }
02951 
02952 
02953 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
02954 
02955 Foam::cellCuts::~cellCuts()
02956 {
02957     clearOut();
02958 }
02959 
02960 
02961 void Foam::cellCuts::clearOut()
02962 {
02963     deleteDemandDrivenData(faceCutsPtr_);
02964 }
02965 
02966 
02967 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
02968 
02969 Foam::pointField Foam::cellCuts::loopPoints(const label cellI) const
02970 {
02971     const labelList& loop = cellLoops_[cellI];
02972 
02973     pointField loopPts(loop.size());
02974 
02975     forAll(loop, fp)
02976     {
02977         label cut = loop[fp];
02978 
02979         if (isEdge(cut))
02980         {
02981             label edgeI = getEdge(cut);
02982 
02983             loopPts[fp] = coord(cut, edgeWeight_[edgeI]);
02984         }
02985         else
02986         {
02987             loopPts[fp] = coord(cut, -GREAT);
02988         }
02989     }
02990     return loopPts;
02991 }
02992 
02993 
02994 // Flip loop for cell
02995 void Foam::cellCuts::flip(const label cellI)
02996 {
02997     labelList& loop = cellLoops_[cellI];
02998 
02999     reverse(loop);
03000 
03001     // Reverse anchor point set.
03002     cellAnchorPoints_[cellI] =
03003         nonAnchorPoints
03004         (
03005             mesh().cellPoints()[cellI],
03006             cellAnchorPoints_[cellI],
03007             loop
03008         );
03009 }
03010 
03011 
03012 // Flip loop only
03013 void Foam::cellCuts::flipLoopOnly(const label cellI)
03014 {
03015     labelList& loop = cellLoops_[cellI];
03016 
03017     reverse(loop);
03018 }
03019 
03020 
03021 void Foam::cellCuts::writeOBJ
03022 (
03023     Ostream& os,
03024     const pointField& loopPts,
03025     label& vertI
03026 ) const
03027 {
03028     label startVertI = vertI;
03029 
03030     forAll(loopPts, fp)
03031     {
03032         const point& pt = loopPts[fp];
03033 
03034         os << "v " << pt.x() << ' ' << pt.y() << ' ' << pt.z() << endl;
03035 
03036         vertI++;
03037     }
03038 
03039     os << 'f';
03040     forAll(loopPts, fp)
03041     {
03042         os  << ' ' << startVertI + fp + 1;
03043     }
03044     os<< endl;
03045 }
03046 
03047 
03048 void Foam::cellCuts::writeOBJ(Ostream& os) const
03049 {
03050     label vertI = 0;
03051 
03052     forAll(cellLoops_, cellI)
03053     {
03054         writeOBJ(os, loopPoints(cellI), vertI);
03055     }
03056 }
03057 
03058 
03059 void Foam::cellCuts::writeCellOBJ(const fileName& dir, const label cellI) const
03060 {
03061     const labelList& anchors = cellAnchorPoints_[cellI];
03062 
03063     writeOBJ(dir, cellI, loopPoints(cellI), anchors);
03064 }
03065 
03066 
03067 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines