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

meshCutAndRemove.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 "meshCutAndRemove.H"
00027 #include <OpenFOAM/polyMesh.H>
00028 #include <dynamicMesh/polyTopoChange.H>
00029 #include <dynamicMesh/polyAddFace.H>
00030 #include <dynamicMesh/polyAddPoint.H>
00031 #include <dynamicMesh/polyRemovePoint.H>
00032 #include <dynamicMesh/polyRemoveFace.H>
00033 #include <dynamicMesh/polyModifyFace.H>
00034 #include <dynamicMesh/cellCuts.H>
00035 #include <OpenFOAM/mapPolyMesh.H>
00036 #include <meshTools/meshTools.H>
00037 
00038 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
00039 
00040 defineTypeNameAndDebug(Foam::meshCutAndRemove, 0);
00041 
00042 // * * * * * * * * * * * * * Private Static Functions  * * * * * * * * * * * //
00043 
00044 // Returns -1 or index in elems1 of first shared element.
00045 Foam::label Foam::meshCutAndRemove::firstCommon
00046 (
00047     const labelList& elems1,
00048     const labelList& elems2
00049 )
00050 {
00051     forAll(elems1, elemI)
00052     {
00053         label index1 = findIndex(elems2, elems1[elemI]);
00054 
00055         if (index1 != -1)
00056         {
00057             return index1;
00058         }
00059     }
00060     return -1;
00061 }
00062 
00063 
00064 // Check if twoCuts at two consecutive position in cuts.
00065 bool Foam::meshCutAndRemove::isIn
00066 (
00067     const edge& twoCuts,
00068     const labelList& cuts
00069 )
00070 {
00071     label index = findIndex(cuts, twoCuts[0]);
00072 
00073     if (index == -1)
00074     {
00075         return false;
00076     }
00077 
00078     return
00079     (
00080         cuts[cuts.fcIndex(index)] == twoCuts[1]
00081      || cuts[cuts.rcIndex(index)] == twoCuts[1]
00082     );
00083 }
00084 
00085 
00086 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
00087 
00088 // Returns the cell in cellLabels that is cut. Or -1.
00089 Foam::label Foam::meshCutAndRemove::findCutCell
00090 (
00091     const cellCuts& cuts,
00092     const labelList& cellLabels
00093 ) const
00094 {
00095     forAll(cellLabels, labelI)
00096     {
00097         label cellI = cellLabels[labelI];
00098 
00099         if (cuts.cellLoops()[cellI].size())
00100         {
00101             return cellI;
00102         }
00103     }
00104     return -1;
00105 }
00106 
00107 
00108 //- Returns first pointI in pointLabels that uses an internal
00109 //  face. Used to find point to inflate cell/face from (has to be
00110 //  connected to internal face). Returns -1 (so inflate from nothing) if
00111 //  none found.
00112 Foam::label Foam::meshCutAndRemove::findInternalFacePoint
00113 (
00114     const labelList& pointLabels
00115 ) const
00116 {
00117     forAll(pointLabels, labelI)
00118     {
00119         label pointI = pointLabels[labelI];
00120 
00121         const labelList& pFaces = mesh().pointFaces()[pointI];
00122 
00123         forAll(pFaces, pFaceI)
00124         {
00125             label faceI = pFaces[pFaceI];
00126 
00127             if (mesh().isInternalFace(faceI))
00128             {
00129                 return pointI;
00130             }
00131         }
00132     }
00133 
00134     if (pointLabels.empty())
00135     {
00136         FatalErrorIn("meshCutAndRemove::findInternalFacePoint(const labelList&)")
00137             << "Empty pointLabels" << abort(FatalError);
00138     }
00139 
00140     return -1;
00141 }
00142 
00143 
00144 // Find point on face that is part of original mesh and that is point connected
00145 // to the patch
00146 Foam::label Foam::meshCutAndRemove::findPatchFacePoint
00147 (
00148     const face& f,
00149     const label exposedPatchI
00150 ) const
00151 {
00152     const labelListList& pointFaces = mesh().pointFaces();
00153     const polyBoundaryMesh& patches = mesh().boundaryMesh();
00154 
00155     forAll(f, fp)
00156     {
00157         label pointI = f[fp];
00158 
00159         if (pointI < mesh().nPoints())
00160         {
00161             const labelList& pFaces = pointFaces[pointI];
00162 
00163             forAll(pFaces, i)
00164             {
00165                 if (patches.whichPatch(pFaces[i]) == exposedPatchI)
00166                 {
00167                     return pointI;
00168                 }
00169             }
00170         }
00171     }
00172     return -1;
00173 }
00174 
00175 
00176 // Get new owner and neighbour of face. Checks anchor points to see if
00177 // cells have been removed.
00178 void Foam::meshCutAndRemove::faceCells
00179 (
00180     const cellCuts& cuts,
00181     const label exposedPatchI,
00182     const label faceI,
00183     label& own,
00184     label& nei,
00185     label& patchID
00186 ) const
00187 {
00188     const labelListList& anchorPts = cuts.cellAnchorPoints();
00189     const labelListList& cellLoops = cuts.cellLoops();
00190 
00191     const face& f = mesh().faces()[faceI];
00192 
00193     own = mesh().faceOwner()[faceI];
00194 
00195     if (cellLoops[own].size() && firstCommon(f, anchorPts[own]) == -1)
00196     {
00197         // owner has been split and this is the removed part.
00198         own = -1;
00199     }
00200 
00201     nei = -1;
00202 
00203     if (mesh().isInternalFace(faceI))
00204     {
00205         nei = mesh().faceNeighbour()[faceI];
00206 
00207         if (cellLoops[nei].size() && firstCommon(f, anchorPts[nei]) == -1)
00208         {
00209             nei = -1;
00210         }
00211     }
00212 
00213     patchID = mesh().boundaryMesh().whichPatch(faceI);
00214 
00215     if (patchID == -1 && (own == -1 || nei == -1))
00216     {
00217         // Face was internal but becomes external
00218         patchID = exposedPatchI;
00219     }
00220 }
00221 
00222 
00223 void Foam::meshCutAndRemove::getZoneInfo
00224 (
00225     const label faceI,
00226     label& zoneID,
00227     bool& zoneFlip
00228 ) const
00229 {
00230     zoneID = mesh().faceZones().whichZone(faceI);
00231 
00232     zoneFlip = false;
00233 
00234     if (zoneID >= 0)
00235     {
00236         const faceZone& fZone = mesh().faceZones()[zoneID];
00237 
00238         zoneFlip = fZone.flipMap()[fZone.whichFace(faceI)];
00239     }
00240 }
00241 
00242 
00243 // Adds a face from point.
00244 void Foam::meshCutAndRemove::addFace
00245 (
00246     polyTopoChange& meshMod,
00247     const label faceI,
00248     const label masterPointI,
00249     const face& newFace,
00250     const label own,
00251     const label nei,
00252     const label patchID
00253 )
00254 {
00255     label zoneID;
00256     bool zoneFlip;
00257 
00258     getZoneInfo(faceI, zoneID, zoneFlip);
00259 
00260     if ((nei == -1) || (own != -1 && own < nei))
00261     {
00262         // Ordering ok.
00263         if (debug & 2)
00264         {
00265             Pout<< "Adding face " << newFace
00266                 << " with new owner:" << own
00267                 << " with new neighbour:" << nei
00268                 << " patchID:" << patchID
00269                 << " anchor:" << masterPointI
00270                 << " zoneID:" << zoneID
00271                 << " zoneFlip:" << zoneFlip
00272                 << endl;
00273         }
00274 
00275         meshMod.setAction
00276         (
00277             polyAddFace
00278             (
00279                 newFace,                    // face
00280                 own,                        // owner
00281                 nei,                        // neighbour
00282                 masterPointI,               // master point
00283                 -1,                         // master edge
00284                 -1,                         // master face for addition
00285                 false,                      // flux flip
00286                 patchID,                    // patch for face
00287                 zoneID,                     // zone for face
00288                 zoneFlip                    // face zone flip
00289             )
00290         );
00291     }
00292     else
00293     {
00294         // Reverse owner/neighbour
00295         if (debug & 2)
00296         {
00297             Pout<< "Adding (reversed) face " << newFace.reverseFace()
00298                 << " with new owner:" << nei
00299                 << " with new neighbour:" << own
00300                 << " patchID:" << patchID
00301                 << " anchor:" << masterPointI
00302                 << " zoneID:" << zoneID
00303                 << " zoneFlip:" << zoneFlip
00304                 << endl;
00305         }
00306 
00307         meshMod.setAction
00308         (
00309             polyAddFace
00310             (
00311                 newFace.reverseFace(),      // face
00312                 nei,                        // owner
00313                 own,                        // neighbour
00314                 masterPointI,               // master point
00315                 -1,                         // master edge
00316                 -1,                         // master face for addition
00317                 false,                      // flux flip
00318                 patchID,                    // patch for face
00319                 zoneID,                     // zone for face
00320                 zoneFlip                    // face zone flip
00321             )
00322         );
00323     }
00324 }
00325 
00326 
00327 // Modifies existing faceI for either new owner/neighbour or new face points.
00328 void Foam::meshCutAndRemove::modFace
00329 (
00330     polyTopoChange& meshMod,
00331     const label faceI,
00332     const face& newFace,
00333     const label own,
00334     const label nei,
00335     const label patchID
00336 )
00337 {
00338     label zoneID;
00339     bool zoneFlip;
00340 
00341     getZoneInfo(faceI, zoneID, zoneFlip);
00342 
00343     if
00344     (
00345         (own != mesh().faceOwner()[faceI])
00346      || (
00347             mesh().isInternalFace(faceI)
00348          && (nei != mesh().faceNeighbour()[faceI])
00349         )
00350      || (newFace != mesh().faces()[faceI])
00351     )
00352     {
00353         if (debug & 2)
00354         {
00355             Pout<< "Modifying face " << faceI
00356                 << " old vertices:" << mesh().faces()[faceI]
00357                 << " new vertices:" << newFace
00358                 << " new owner:" << own
00359                 << " new neighbour:" << nei
00360                 << " new patch:" << patchID
00361                 << " new zoneID:" << zoneID
00362                 << " new zoneFlip:" << zoneFlip
00363                 << endl;
00364         }
00365 
00366         if ((nei == -1) || (own != -1 && own < nei))
00367         {
00368             meshMod.setAction
00369             (
00370                 polyModifyFace
00371                 (
00372                     newFace,            // modified face
00373                     faceI,              // label of face being modified
00374                     own,                // owner
00375                     nei,                // neighbour
00376                     false,              // face flip
00377                     patchID,            // patch for face
00378                     false,              // remove from zone
00379                     zoneID,             // zone for face
00380                     zoneFlip            // face flip in zone
00381                 )
00382             );
00383         }
00384         else
00385         {
00386             meshMod.setAction
00387             (
00388                 polyModifyFace
00389                 (
00390                     newFace.reverseFace(),  // modified face
00391                     faceI,                  // label of face being modified
00392                     nei,                    // owner
00393                     own,                    // neighbour
00394                     false,                  // face flip
00395                     patchID,                // patch for face
00396                     false,                  // remove from zone
00397                     zoneID,                 // zone for face
00398                     zoneFlip                // face flip in zone
00399                 )
00400             );
00401         }
00402     }
00403 }
00404 
00405 
00406 // Copies face starting from startFp up to and including endFp.
00407 void Foam::meshCutAndRemove::copyFace
00408 (
00409     const face& f,
00410     const label startFp,
00411     const label endFp,
00412     face& newFace
00413 ) const
00414 {
00415     label fp = startFp;
00416 
00417     label newFp = 0;
00418 
00419     while (fp != endFp)
00420     {
00421         newFace[newFp++] = f[fp];
00422 
00423         fp = (fp + 1) % f.size();
00424     }
00425     newFace[newFp] = f[fp];
00426 }
00427 
00428 
00429 // Actually split face in two along splitEdge v0, v1 (the two vertices in new
00430 // vertex numbering). Generates faces in same ordering
00431 // as original face. Replaces cutEdges by the points introduced on them
00432 // (addedPoints_).
00433 void Foam::meshCutAndRemove::splitFace
00434 (
00435     const face& f,
00436     const label v0,
00437     const label v1,
00438 
00439     face& f0,
00440     face& f1
00441 ) const
00442 {
00443     // Check if we find any new vertex which is part of the splitEdge.
00444     label startFp = findIndex(f, v0);
00445 
00446     if (startFp == -1)
00447     {
00448         FatalErrorIn
00449         (
00450             "meshCutAndRemove::splitFace"
00451             ", const face&, const label, const label, face&, face&)"
00452         )   << "Cannot find vertex (new numbering) " << v0
00453             << " on face " << f
00454             << abort(FatalError);
00455     }
00456 
00457     label endFp = findIndex(f, v1);
00458 
00459     if (endFp == -1)
00460     {
00461         FatalErrorIn
00462         (
00463             "meshCutAndRemove::splitFace("
00464             ", const face&, const label, const label, face&, face&)"
00465         )   << "Cannot find vertex (new numbering) " << v1
00466             << " on face " << f
00467             << abort(FatalError);
00468     }
00469 
00470 
00471     f0.setSize((endFp + 1 + f.size() - startFp) % f.size());
00472     f1.setSize(f.size() - f0.size() + 2);
00473 
00474     copyFace(f, startFp, endFp, f0);
00475     copyFace(f, endFp, startFp, f1);
00476 }
00477 
00478 
00479 // Adds additional vertices (from edge cutting) to face. Used for faces which
00480 // are not split but still might use edge that has been cut.
00481 Foam::face Foam::meshCutAndRemove::addEdgeCutsToFace(const label faceI) const
00482 {
00483     const face& f = mesh().faces()[faceI];
00484 
00485     face newFace(2 * f.size());
00486 
00487     label newFp = 0;
00488 
00489     forAll(f, fp)
00490     {
00491         // Duplicate face vertex.
00492         newFace[newFp++] = f[fp];
00493 
00494         // Check if edge has been cut.
00495         label fp1 = f.fcIndex(fp);
00496 
00497         HashTable<label, edge, Hash<edge> >::const_iterator fnd =
00498             addedPoints_.find(edge(f[fp], f[fp1]));
00499 
00500         if (fnd != addedPoints_.end())
00501         {
00502             // edge has been cut. Introduce new vertex.
00503             newFace[newFp++] = fnd();
00504         }
00505     }
00506 
00507     newFace.setSize(newFp);
00508 
00509     return newFace;
00510 }
00511 
00512 
00513 // Walk loop (loop of cuts) across circumference of cellI. Returns face in
00514 // new vertices.
00515 // Note: tricky bit is that it can use existing edges which have been split.
00516 Foam::face Foam::meshCutAndRemove::loopToFace
00517 (
00518     const label cellI,
00519     const labelList& loop
00520 ) const
00521 {
00522     face newFace(2*loop.size());
00523 
00524     label newFaceI = 0;
00525 
00526     forAll(loop, fp)
00527     {
00528         label cut = loop[fp];
00529 
00530         if (isEdge(cut))
00531         {
00532             label edgeI = getEdge(cut);
00533 
00534             const edge& e = mesh().edges()[edgeI];
00535 
00536             label vertI = addedPoints_[e];
00537 
00538             newFace[newFaceI++] = vertI;
00539         }
00540         else
00541         {
00542             // cut is vertex.
00543             label vertI = getVertex(cut);
00544 
00545             newFace[newFaceI++] = vertI;
00546 
00547             label nextCut = loop[loop.fcIndex(fp)];
00548 
00549             if (!isEdge(nextCut))
00550             {
00551                 // From vertex to vertex -> cross cut only if no existing edge.
00552                 label nextVertI = getVertex(nextCut);
00553 
00554                 label edgeI = meshTools::findEdge(mesh(), vertI, nextVertI);
00555 
00556                 if (edgeI != -1)
00557                 {
00558                     // Existing edge. Insert split-edge point if any.
00559                     HashTable<label, edge, Hash<edge> >::const_iterator fnd =
00560                         addedPoints_.find(mesh().edges()[edgeI]);
00561 
00562                     if (fnd != addedPoints_.end())
00563                     {
00564                         newFace[newFaceI++] = fnd();
00565                     }
00566                 }
00567             }
00568         }
00569     }
00570     newFace.setSize(newFaceI);
00571 
00572     return newFace;
00573 }
00574 
00575 
00576 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
00577 
00578 // Construct from components
00579 Foam::meshCutAndRemove::meshCutAndRemove(const polyMesh& mesh)
00580 :
00581     edgeVertex(mesh),
00582     addedFaces_(),
00583     addedPoints_()
00584 {}
00585 
00586 
00587 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
00588 
00589 void Foam::meshCutAndRemove::setRefinement
00590 (
00591     const label exposedPatchI,
00592     const cellCuts& cuts,
00593     const labelList& cutPatch,
00594     polyTopoChange& meshMod
00595 )
00596 {
00597     // Clear and size maps here since mesh size will change.
00598     addedFaces_.clear();
00599     addedFaces_.resize(cuts.nLoops());
00600 
00601     addedPoints_.clear();
00602     addedPoints_.resize(cuts.nLoops());
00603 
00604     if (cuts.nLoops() == 0)
00605     {
00606         return;
00607     }
00608 
00609     const labelListList& anchorPts = cuts.cellAnchorPoints();
00610     const labelListList& cellLoops = cuts.cellLoops();
00611     const polyBoundaryMesh& patches = mesh().boundaryMesh();
00612 
00613     if (exposedPatchI < 0 || exposedPatchI >= patches.size())
00614     {
00615         FatalErrorIn
00616         (
00617             "meshCutAndRemove::setRefinement("
00618             ", const label, const cellCuts&, const labelList&"
00619             ", polyTopoChange&)"
00620         )   << "Illegal exposed patch " << exposedPatchI
00621             << abort(FatalError);
00622     }
00623 
00624 
00625     //
00626     // Add new points along cut edges.
00627     //
00628 
00629     forAll(cuts.edgeIsCut(), edgeI)
00630     {
00631         if (cuts.edgeIsCut()[edgeI])
00632         {
00633             const edge& e = mesh().edges()[edgeI];
00634 
00635             // Check if there is any cell using this edge.
00636             if (debug && findCutCell(cuts, mesh().edgeCells()[edgeI]) == -1)
00637             {
00638                 FatalErrorIn
00639                 (
00640                     "meshCutAndRemove::setRefinement("
00641                     ", const label, const cellCuts&, const labelList&"
00642                     ", polyTopoChange&)"
00643                 )   << "Problem: cut edge but none of the cells using it is\n"
00644                     << "edge:" << edgeI << " verts:" << e
00645                     << abort(FatalError);
00646             }
00647 
00648             // One of the edge end points should be master point of nbCellI.
00649             label masterPointI = e.start();
00650 
00651             const point& v0 = mesh().points()[e.start()];
00652             const point& v1 = mesh().points()[e.end()];
00653 
00654             scalar weight = cuts.edgeWeight()[edgeI];
00655 
00656             point newPt = weight*v1 + (1.0-weight)*v0;
00657 
00658             label addedPointI =
00659                 meshMod.setAction
00660                 (
00661                     polyAddPoint
00662                     (
00663                         newPt,              // point
00664                         masterPointI,       // master point
00665                         -1,                 // zone for point
00666                         true                // supports a cell
00667                     )
00668                 );
00669 
00670             // Store on (hash of) edge.
00671             addedPoints_.insert(e, addedPointI);
00672 
00673             if (debug & 2)
00674             {
00675                 Pout<< "Added point " << addedPointI
00676                     << " to vertex "
00677                     << masterPointI << " of edge " << edgeI
00678                     << " vertices " << e << endl;
00679             }
00680         }
00681     }
00682 
00683 
00684     //
00685     // Remove all points that will not be used anymore
00686     //
00687     {
00688         boolList usedPoint(mesh().nPoints(), false);
00689 
00690         forAll(cellLoops, cellI)
00691         {
00692             const labelList& loop = cellLoops[cellI];
00693 
00694             if (loop.size())
00695             {
00696                 // Cell is cut. Uses only anchor points and loop itself.
00697                 forAll(loop, fp)
00698                 {
00699                     label cut = loop[fp];
00700 
00701                     if (!isEdge(cut))
00702                     {
00703                         usedPoint[getVertex(cut)] = true;
00704                     }
00705                 }
00706 
00707                 const labelList& anchors = anchorPts[cellI];
00708 
00709                 forAll(anchors, i)
00710                 {
00711                     usedPoint[anchors[i]] = true;
00712                 }
00713             }
00714             else
00715             {
00716                 // Cell is not cut so use all its points
00717                 const labelList& cPoints = mesh().cellPoints()[cellI];
00718 
00719                 forAll(cPoints, i)
00720                 {
00721                     usedPoint[cPoints[i]] = true;
00722                 }
00723             }
00724         }
00725 
00726 
00727         // Check
00728         const Map<edge>& faceSplitCut = cuts.faceSplitCut();
00729 
00730         forAllConstIter(Map<edge>, faceSplitCut, iter)
00731         {
00732             const edge& fCut = iter();
00733 
00734             forAll(fCut, i)
00735             {
00736                 label cut = fCut[i];
00737 
00738                 if (!isEdge(cut))
00739                 {
00740                     label pointI = getVertex(cut);
00741 
00742                     if (!usedPoint[pointI])
00743                     {
00744                         FatalErrorIn
00745                         (
00746                             "meshCutAndRemove::setRefinement("
00747                             ", const label, const cellCuts&, const labelList&"
00748                             ", polyTopoChange&)"
00749                         )   << "Problem: faceSplitCut not used by any loop"
00750                             << " or cell anchor point"
00751                             << "face:" << iter.key() << " point:" << pointI
00752                             << " coord:" << mesh().points()[pointI]
00753                             << abort(FatalError);
00754                     }
00755                 }
00756             }
00757         }
00758 
00759         forAll(cuts.pointIsCut(), pointI)
00760         {
00761             if (cuts.pointIsCut()[pointI])
00762             {
00763                 if (!usedPoint[pointI])
00764                 {
00765                     FatalErrorIn
00766                     (
00767                         "meshCutAndRemove::setRefinement("
00768                         ", const label, const cellCuts&, const labelList&"
00769                         ", polyTopoChange&)"
00770                     )   << "Problem: point is marked as cut but"
00771                         << " not used by any loop"
00772                         << " or cell anchor point"
00773                         << "point:" << pointI
00774                         << " coord:" << mesh().points()[pointI]
00775                         << abort(FatalError);
00776                 }
00777             }
00778         }
00779 
00780 
00781         // Remove unused points.
00782         forAll(usedPoint, pointI)
00783         {
00784             if (!usedPoint[pointI])
00785             {
00786                 meshMod.setAction(polyRemovePoint(pointI));
00787 
00788                 if (debug & 2)
00789                 {
00790                     Pout<< "Removing unused point " << pointI << endl;
00791                 }
00792             }
00793         }
00794     }
00795 
00796 
00797     //
00798     // For all cut cells add an internal or external face
00799     //
00800 
00801     forAll(cellLoops, cellI)
00802     {
00803         const labelList& loop = cellLoops[cellI];
00804 
00805         if (loop.size())
00806         {
00807             if (cutPatch[cellI] < 0 || cutPatch[cellI] >= patches.size())
00808             {
00809                 FatalErrorIn
00810                 (
00811                     "meshCutAndRemove::setRefinement("
00812                     ", const label, const cellCuts&, const labelList&"
00813                     ", polyTopoChange&)"
00814                 )   << "Illegal patch " << cutPatch[cellI]
00815                     << " provided for cut cell " << cellI
00816                     << abort(FatalError);
00817             }
00818 
00819             //
00820             // Convert loop (=list of cuts) into proper face.
00821             // cellCuts sets orientation is towards anchor side so reverse.
00822             //
00823             face newFace(loopToFace(cellI, loop));
00824 
00825             reverse(newFace);
00826 
00827             // Pick any anchor point on cell
00828             label masterPointI = findPatchFacePoint(newFace, exposedPatchI);
00829 
00830             label addedFaceI =
00831                 meshMod.setAction
00832                 (
00833                     polyAddFace
00834                     (
00835                         newFace,                // face
00836                         cellI,                  // owner
00837                         -1,                     // neighbour
00838                         masterPointI,           // master point
00839                         -1,                     // master edge
00840                         -1,                     // master face for addition
00841                         false,                  // flux flip
00842                         cutPatch[cellI],        // patch for face
00843                         -1,                     // zone for face
00844                         false                   // face zone flip
00845                     )
00846                 );
00847 
00848             addedFaces_.insert(cellI, addedFaceI);
00849 
00850             if (debug & 2)
00851             {
00852                 Pout<< "Added splitting face " << newFace << " index:"
00853                     << addedFaceI << " from masterPoint:" << masterPointI
00854                     << " to owner " << cellI << " with anchors:"
00855                     << anchorPts[cellI]
00856                     << " from Loop:";
00857 
00858                 // Gets edgeweights of loop
00859                 scalarField weights(loop.size());
00860                 forAll(loop, i)
00861                 {
00862                     label cut = loop[i];
00863 
00864                     weights[i] =
00865                     (
00866                         isEdge(cut)
00867                       ? cuts.edgeWeight()[getEdge(cut)]
00868                       : -GREAT
00869                     );
00870                 }
00871                 writeCuts(Pout, loop, weights);
00872                 Pout<< endl;
00873             }
00874         }
00875     }
00876 
00877 
00878     //
00879     // Modify faces to use only anchorpoints and loop points
00880     // (so throw away part without anchorpoints)
00881     //
00882 
00883 
00884     // Maintain whether face has been updated (for -split edges
00885     // -new owner/neighbour)
00886     boolList faceUptodate(mesh().nFaces(), false);
00887 
00888     const Map<edge>& faceSplitCuts = cuts.faceSplitCut();
00889 
00890     for
00891     (
00892         Map<edge>::const_iterator iter = faceSplitCuts.begin();
00893         iter != faceSplitCuts.end();
00894         ++iter
00895     )
00896     {
00897         label faceI = iter.key();
00898 
00899         // Renumber face to include split edges.
00900         face newFace(addEdgeCutsToFace(faceI));
00901 
00902         // Edge splitting the face. Convert edge to new vertex numbering.
00903         const edge& splitEdge = iter();
00904 
00905         label cut0 = splitEdge[0];
00906 
00907         label v0;
00908         if (isEdge(cut0))
00909         {
00910             label edgeI = getEdge(cut0);
00911             v0 = addedPoints_[mesh().edges()[edgeI]];
00912         }
00913         else
00914         {
00915             v0 = getVertex(cut0);
00916         }
00917 
00918         label cut1 = splitEdge[1];
00919         label v1;
00920         if (isEdge(cut1))
00921         {
00922             label edgeI = getEdge(cut1);
00923             v1 = addedPoints_[mesh().edges()[edgeI]];
00924         }
00925         else
00926         {
00927             v1 = getVertex(cut1);
00928         }
00929 
00930         // Split face along the elements of the splitEdge.
00931         face f0, f1;
00932         splitFace(newFace, v0, v1, f0, f1);
00933 
00934         label own = mesh().faceOwner()[faceI];
00935 
00936         label nei = -1;
00937 
00938         if (mesh().isInternalFace(faceI))
00939         {
00940             nei = mesh().faceNeighbour()[faceI];
00941         }
00942 
00943         if (debug & 2)
00944         {
00945             Pout<< "Split face " << mesh().faces()[faceI]
00946                 << " own:" << own << " nei:" << nei
00947                 << " into f0:" << f0
00948                 << " and f1:" << f1 << endl;
00949         }
00950 
00951 
00952         // Check which cell using face uses anchorPoints (so is kept)
00953         // and which one doesn't (gets removed)
00954 
00955         // Bit tricky. We have to know whether this faceSplit splits owner/
00956         // neighbour or both. Even if cell is cut we have to make sure this is
00957         // the one that cuts it (this face cut might not be the one splitting
00958         // the cell)
00959         // The face f gets split into two parts, f0 and f1.
00960         // Each of these can have a different owner and or neighbour.
00961 
00962         const face& f = mesh().faces()[faceI];
00963 
00964         label f0Own = -1;
00965         label f1Own = -1;
00966 
00967         if (cellLoops[own].empty())
00968         {
00969             // Owner side is not split so keep both halves.
00970             f0Own = own;
00971             f1Own = own;
00972         }
00973         else if (isIn(splitEdge, cellLoops[own]))
00974         {
00975             // Owner is cut by this splitCut. See which of f0, f1 gets
00976             // preserved and becomes owner, and which gets removed.
00977             if (firstCommon(f0, anchorPts[own]) != -1)
00978             {
00979                 // f0 preserved so f1 gets deleted
00980                 f0Own = own;
00981                 f1Own = -1;
00982             }
00983             else
00984             {
00985                 f0Own = -1;
00986                 f1Own = own;
00987             }
00988         }
00989         else
00990         {
00991             // Owner not cut by this splitCut but by another.
00992             // Check on original face whether
00993             // use anchorPts.
00994             if (firstCommon(f, anchorPts[own]) != -1)
00995             {
00996                 // both f0 and f1 owner side preserved
00997                 f0Own = own;
00998                 f1Own = own;
00999             }
01000             else
01001             {
01002                 // both f0 and f1 owner side removed
01003                 f0Own = -1;
01004                 f1Own = -1;
01005             }
01006         }
01007 
01008 
01009         label f0Nei = -1;
01010         label f1Nei = -1;
01011 
01012         if (nei != -1)
01013         {
01014             if (cellLoops[nei].empty())
01015             {
01016                 f0Nei = nei;
01017                 f1Nei = nei;
01018             }
01019             else if (isIn(splitEdge, cellLoops[nei]))
01020             {
01021                 // Neighbour is cut by this splitCut. So anchor part of it
01022                 // gets kept, non-anchor bit gets removed. See which of f0, f1
01023                 // connects to which part.
01024 
01025                 if (firstCommon(f0, anchorPts[nei]) != -1)
01026                 {
01027                     f0Nei = nei;
01028                     f1Nei = -1;
01029                 }
01030                 else
01031                 {
01032                     f0Nei = -1;
01033                     f1Nei = nei;
01034                 }
01035             }
01036             else
01037             {
01038                 // neighbour not cut by this splitCut. Check on original face
01039                 // whether use anchorPts.
01040 
01041                 if (firstCommon(f, anchorPts[nei]) != -1)
01042                 {
01043                     f0Nei = nei;
01044                     f1Nei = nei;
01045                 }
01046                 else
01047                 {
01048                     // both f0 and f1 on neighbour side removed
01049                     f0Nei = -1;
01050                     f1Nei = -1;
01051                 }
01052             }
01053         }
01054 
01055 
01056         if (debug & 2)
01057         {
01058             Pout<< "f0 own:" << f0Own << " nei:" << f0Nei
01059                 << "  f1 own:" << f1Own << " nei:" << f1Nei
01060                 << endl;
01061         }
01062 
01063 
01064         // If faces were internal but now become external set a patch.
01065         // If they were external already keep the patch.
01066         label patchID = patches.whichPatch(faceI);
01067 
01068         if (patchID == -1)
01069         {
01070             patchID = exposedPatchI;
01071         }
01072 
01073 
01074         // Do as much as possible by modifying faceI. Delay any remove
01075         // face. Keep track of whether faceI has been used.
01076 
01077         bool modifiedFaceI = false;
01078 
01079         if (f0Own == -1)
01080         {
01081             if (f0Nei != -1)
01082             {
01083                 // f0 becomes external face (note:modFace will reverse face)
01084                 modFace(meshMod, faceI, f0, f0Own, f0Nei, patchID);
01085                 modifiedFaceI = true;
01086             }
01087         }
01088         else
01089         {
01090             if (f0Nei == -1)
01091             {
01092                 // f0 becomes external face
01093                 modFace(meshMod, faceI, f0, f0Own, f0Nei, patchID);
01094                 modifiedFaceI = true;
01095             }
01096             else
01097             {
01098                 // f0 stays internal face.
01099                 modFace(meshMod, faceI, f0, f0Own, f0Nei, -1);
01100                 modifiedFaceI = true;
01101             }
01102         }
01103 
01104 
01105         // f1 is added face (if at all)
01106 
01107         if (f1Own == -1)
01108         {
01109             if (f1Nei == -1)
01110             {
01111                 // f1 not needed.
01112             }
01113             else
01114             {
01115                 // f1 becomes external face (note:modFace will reverse face)
01116                 if (!modifiedFaceI)
01117                 {
01118                     modFace(meshMod, faceI, f1, f1Own, f1Nei, patchID);
01119                     modifiedFaceI = true;
01120                 }
01121                 else
01122                 {
01123                     label masterPointI = findPatchFacePoint(f1, patchID);
01124 
01125                     addFace
01126                     (
01127                         meshMod,
01128                         faceI,          // face for zone info
01129                         masterPointI,   // inflation point
01130                         f1,             // vertices of face
01131                         f1Own,
01132                         f1Nei,
01133                         patchID         // patch for new face
01134                     );
01135                 }
01136             }
01137         }
01138         else
01139         {
01140             if (f1Nei == -1)
01141             {
01142                 // f1 becomes external face
01143                 if (!modifiedFaceI)
01144                 {
01145                     modFace(meshMod, faceI, f1, f1Own, f1Nei, patchID);
01146                     modifiedFaceI = true;
01147                 }
01148                 else
01149                 {
01150                     label masterPointI = findPatchFacePoint(f1, patchID);
01151 
01152                     addFace
01153                     (
01154                         meshMod,
01155                         faceI,
01156                         masterPointI,
01157                         f1,
01158                         f1Own,
01159                         f1Nei,
01160                         patchID
01161                     );
01162                 }
01163             }
01164             else
01165             {
01166                 // f1 is internal face.
01167                 if (!modifiedFaceI)
01168                 {
01169                     modFace(meshMod, faceI, f1, f1Own, f1Nei, -1);
01170                     modifiedFaceI = true;
01171                 }
01172                 else
01173                 {
01174                     label masterPointI = findPatchFacePoint(f1, -1);
01175 
01176                     addFace(meshMod, faceI, masterPointI, f1, f1Own, f1Nei, -1);
01177                 }
01178             }
01179         }
01180 
01181         if (f0Own == -1 && f0Nei == -1 && !modifiedFaceI)
01182         {
01183             meshMod.setAction(polyRemoveFace(faceI));
01184 
01185             if (debug & 2)
01186             {
01187                 Pout<< "Removed face " << faceI << endl;
01188             }
01189         }
01190 
01191         faceUptodate[faceI] = true;
01192     }
01193 
01194 
01195     //
01196     // Faces that have not been split but just appended to. Are guaranteed
01197     // to be reachable from an edgeCut.
01198     //
01199 
01200     const boolList& edgeIsCut = cuts.edgeIsCut();
01201 
01202     forAll(edgeIsCut, edgeI)
01203     {
01204         if (edgeIsCut[edgeI])
01205         {
01206             const labelList& eFaces = mesh().edgeFaces()[edgeI];
01207 
01208             forAll(eFaces, i)
01209             {
01210                 label faceI = eFaces[i];
01211 
01212                 if (!faceUptodate[faceI])
01213                 {
01214                     // So the face has not been split itself (i.e. its owner
01215                     // or neighbour have not been split) so it only
01216                     // borders by edge a cell which has been split.
01217 
01218                     // Get (new or original) owner and neighbour of faceI
01219                     label own, nei, patchID;
01220                     faceCells(cuts, exposedPatchI, faceI, own, nei, patchID);
01221 
01222 
01223                     if (own == -1 && nei == -1)
01224                     {
01225                         meshMod.setAction(polyRemoveFace(faceI));
01226 
01227                         if (debug & 2)
01228                         {
01229                             Pout<< "Removed face " << faceI << endl;
01230                         }
01231                     }
01232                     else
01233                     {
01234                         // Renumber face to include split edges.
01235                         face newFace(addEdgeCutsToFace(faceI));
01236 
01237                         if (debug & 2)
01238                         {
01239                             Pout<< "Added edge cuts to face " << faceI
01240                                 << " f:" << mesh().faces()[faceI]
01241                                 << " newFace:" << newFace << endl;
01242                         }
01243 
01244                         modFace
01245                         (
01246                             meshMod,
01247                             faceI,
01248                             newFace,
01249                             own,
01250                             nei,
01251                             patchID
01252                         );
01253                     }
01254 
01255                     faceUptodate[faceI] = true;
01256                 }
01257             }
01258         }
01259     }
01260 
01261 
01262     //
01263     // Remove any faces on the non-anchor side of a split cell.
01264     // Note: could loop through all cut cells only and check their faces but
01265     //       looping over all faces is cleaner and probably faster for dense
01266     //       cut patterns.
01267 
01268     const faceList& faces = mesh().faces();
01269 
01270     forAll(faces, faceI)
01271     {
01272         if (!faceUptodate[faceI])
01273         {
01274             // Get (new or original) owner and neighbour of faceI
01275             label own, nei, patchID;
01276             faceCells(cuts, exposedPatchI, faceI, own, nei, patchID);
01277 
01278             if (own == -1 && nei == -1)
01279             {
01280                 meshMod.setAction(polyRemoveFace(faceI));
01281 
01282                 if (debug & 2)
01283                 {
01284                     Pout<< "Removed face " << faceI << endl;
01285                 }
01286             }
01287             else
01288             {
01289                 modFace(meshMod, faceI, faces[faceI], own, nei, patchID);
01290             }
01291 
01292             faceUptodate[faceI] = true;
01293         }
01294     }
01295 
01296     if (debug)
01297     {
01298         Pout<< "meshCutAndRemove:" << nl
01299             << "    cells split:" << cuts.nLoops() << nl
01300             << "    faces added:" << addedFaces_.size() << nl
01301             << "    points added on edges:" << addedPoints_.size() << nl
01302             << endl;
01303     }
01304 }
01305 
01306 
01307 void Foam::meshCutAndRemove::updateMesh(const mapPolyMesh& map)
01308 {
01309     // Update stored labels for mesh change.
01310     {
01311         Map<label> newAddedFaces(addedFaces_.size());
01312 
01313         for
01314         (
01315             Map<label>::const_iterator iter = addedFaces_.begin();
01316             iter != addedFaces_.end();
01317             ++iter
01318         )
01319         {
01320             label cellI = iter.key();
01321 
01322             label newCellI = map.reverseCellMap()[cellI];
01323 
01324             label addedFaceI = iter();
01325 
01326             label newAddedFaceI = map.reverseFaceMap()[addedFaceI];
01327 
01328             if ((newCellI >= 0) && (newAddedFaceI >= 0))
01329             {
01330                 if
01331                 (
01332                     (debug & 2)
01333                  && (newCellI != cellI || newAddedFaceI != addedFaceI)
01334                 )
01335                 {
01336                     Pout<< "meshCutAndRemove::updateMesh :"
01337                         << " updating addedFace for cell " << cellI
01338                         << " from " << addedFaceI
01339                         << " to " << newAddedFaceI
01340                         << endl;
01341                 }
01342                 newAddedFaces.insert(newCellI, newAddedFaceI);
01343             }
01344         }
01345 
01346         // Copy
01347         addedFaces_.transfer(newAddedFaces);
01348     }
01349 
01350     {
01351         HashTable<label, edge, Hash<edge> > newAddedPoints(addedPoints_.size());
01352 
01353         for
01354         (
01355             HashTable<label, edge, Hash<edge> >::const_iterator iter =
01356                 addedPoints_.begin();
01357             iter != addedPoints_.end();
01358             ++iter
01359         )
01360         {
01361             const edge& e = iter.key();
01362 
01363             label newStart = map.reversePointMap()[e.start()];
01364 
01365             label newEnd = map.reversePointMap()[e.end()];
01366 
01367             label addedPointI = iter();
01368 
01369             label newAddedPointI = map.reversePointMap()[addedPointI];
01370 
01371             if ((newStart >= 0) && (newEnd >= 0) && (newAddedPointI >= 0))
01372             {
01373                 edge newE = edge(newStart, newEnd);
01374 
01375                 if
01376                 (
01377                     (debug & 2)
01378                  && (e != newE || newAddedPointI != addedPointI)
01379                 )
01380                 {
01381                     Pout<< "meshCutAndRemove::updateMesh :"
01382                         << " updating addedPoints for edge " << e
01383                         << " from " << addedPointI
01384                         << " to " << newAddedPointI
01385                         << endl;
01386                 }
01387 
01388                 newAddedPoints.insert(newE, newAddedPointI);
01389             }
01390         }
01391 
01392         // Copy
01393         addedPoints_.transfer(newAddedPoints);
01394     }
01395 }
01396 
01397 
01398 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines