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

meshDualiser.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 "meshDualiser.H"
00027 #include <meshTools/meshTools.H>
00028 #include <OpenFOAM/polyMesh.H>
00029 #include <dynamicMesh/polyTopoChange.H>
00030 #include <OpenFOAM/mapPolyMesh.H>
00031 #include <meshTools/edgeFaceCirculator.H>
00032 #include <OpenFOAM/mergePoints.H>
00033 #include <OpenFOAM/OFstream.H>
00034 
00035 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
00036 
00037 defineTypeNameAndDebug(Foam::meshDualiser, 0);
00038 
00039 
00040 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
00041 
00042 void Foam::meshDualiser::checkPolyTopoChange(const polyTopoChange& meshMod)
00043 {
00044     // Assume no removed points
00045     pointField points(meshMod.points().size());
00046     forAll(meshMod.points(), i)
00047     {
00048         points[i] = meshMod.points()[i];
00049     }
00050 
00051     labelList oldToNew;
00052     pointField newPoints;
00053     bool hasMerged = mergePoints
00054     (
00055         points,
00056         1E-6,
00057         false,
00058         oldToNew,
00059         newPoints
00060     );
00061 
00062     if (hasMerged)
00063     {
00064         labelListList newToOld(invertOneToMany(newPoints.size(), oldToNew));
00065 
00066         forAll(newToOld, newI)
00067         {
00068             if (newToOld[newI].size() != 1)
00069             {
00070                 FatalErrorIn
00071                 (
00072                     "meshDualiser::checkPolyTopoChange(const polyTopoChange&)"
00073                 )   << "duplicate verts:" << newToOld[newI]
00074                     << " coords:"
00075                     << UIndirectList<point>(points, newToOld[newI])()
00076                     << abort(FatalError);
00077             }
00078         }
00079     }
00080 }
00081 
00082 
00083 // Dump state so far.
00084 void Foam::meshDualiser::dumpPolyTopoChange
00085 (
00086     const polyTopoChange& meshMod,
00087     const fileName& prefix
00088 )
00089 {
00090     OFstream str1(prefix + "Faces.obj");
00091     OFstream str2(prefix + "Edges.obj");
00092 
00093     Info<< "Dumping current polyTopoChange. Faces to " << str1.name()
00094         << " , points and edges to " << str2.name() << endl;
00095 
00096     const DynamicList<point>& points = meshMod.points();
00097 
00098     forAll(points, pointI)
00099     {
00100         meshTools::writeOBJ(str1, points[pointI]);
00101         meshTools::writeOBJ(str2, points[pointI]);
00102     }
00103 
00104     const DynamicList<face>& faces = meshMod.faces();
00105 
00106     forAll(faces, faceI)
00107     {
00108         const face& f = faces[faceI];
00109 
00110         str1<< 'f';
00111         forAll(f, fp)
00112         {
00113             str1<< ' ' << f[fp]+1;
00114         }
00115         str1<< nl;
00116 
00117         str2<< 'l';
00118         forAll(f, fp)
00119         {
00120             str2<< ' ' << f[fp]+1;
00121         }
00122         str2<< ' ' << f[0]+1 << nl;
00123     }
00124 }
00125 
00126 
00127 //- Given cell and point on mesh finds the corresponding dualCell. Most
00128 //  points become only one cell but the feature points might be split.
00129 Foam::label Foam::meshDualiser::findDualCell
00130 (
00131     const label cellI,
00132     const label pointI
00133 ) const
00134 {
00135     const labelList& dualCells = pointToDualCells_[pointI];
00136 
00137     if (dualCells.size() == 1)
00138     {
00139         return dualCells[0];
00140     }
00141     else
00142     {
00143         label index = findIndex(mesh_.pointCells()[pointI], cellI);
00144 
00145         return dualCells[index];
00146     }
00147 }
00148 
00149 
00150 // Helper function to generate dualpoints on all boundary edges emanating
00151 // from (boundary & feature) point
00152 void Foam::meshDualiser::generateDualBoundaryEdges
00153 (
00154     const PackedBoolList& isBoundaryEdge,
00155     const label pointI,
00156     polyTopoChange& meshMod
00157 )
00158 {
00159     const labelList& pEdges = mesh_.pointEdges()[pointI];
00160 
00161     forAll(pEdges, pEdgeI)
00162     {
00163         label edgeI = pEdges[pEdgeI];
00164 
00165         if (edgeToDualPoint_[edgeI] == -1 && isBoundaryEdge.get(edgeI) == 1)
00166         {
00167             const edge& e = mesh_.edges()[edgeI];
00168 
00169             edgeToDualPoint_[edgeI] = meshMod.addPoint
00170             (
00171                 e.centre(mesh_.points()),
00172                 pointI, // masterPoint
00173                 -1,     // zoneID
00174                 true    // inCell
00175             );
00176         }
00177     }
00178 }
00179 
00180 
00181 // Return true if point on face has same dual cells on both owner and neighbour
00182 // sides.
00183 bool Foam::meshDualiser::sameDualCell
00184 (
00185     const label faceI,
00186     const label pointI
00187 ) const
00188 {
00189     if (!mesh_.isInternalFace(faceI))
00190     {
00191         FatalErrorIn("sameDualCell(const label, const label)")
00192             << "face:" << faceI << " is not internal face."
00193             << abort(FatalError);
00194     }
00195 
00196     label own = mesh_.faceOwner()[faceI];
00197     label nei = mesh_.faceNeighbour()[faceI];
00198 
00199     return findDualCell(own, pointI) == findDualCell(nei, pointI);
00200 }
00201 
00202 
00203 Foam::label Foam::meshDualiser::addInternalFace
00204 (
00205     const label masterPointI,
00206     const label masterEdgeI,
00207     const label masterFaceI,
00208 
00209     const bool edgeOrder,
00210     const label dualCell0,
00211     const label dualCell1,
00212     const DynamicList<label>& verts,
00213     polyTopoChange& meshMod
00214 ) const
00215 {
00216     face newFace(verts);
00217 
00218     if (edgeOrder != (dualCell0 < dualCell1))
00219     {
00220         reverse(newFace);
00221     }
00222 
00223     if (debug)
00224     {
00225         pointField facePoints(meshMod.points(), newFace);
00226 
00227         labelList oldToNew;
00228         pointField newPoints;
00229         bool hasMerged = mergePoints
00230         (
00231             facePoints,
00232             1E-6,
00233             false,
00234             oldToNew,
00235             newPoints
00236         );
00237 
00238         if (hasMerged)
00239         {
00240             FatalErrorIn("addInternalFace(..)")
00241                 << "verts:" << verts << " newFace:" << newFace
00242                 << " face points:" << facePoints
00243                 << abort(FatalError);
00244         }
00245     }
00246 
00247 
00248     label zoneID = -1;
00249     bool zoneFlip = false;
00250     if (masterFaceI != -1)
00251     {
00252         zoneID = mesh_.faceZones().whichZone(masterFaceI);
00253 
00254         if (zoneID != -1)
00255         {
00256             const faceZone& fZone = mesh_.faceZones()[zoneID];
00257 
00258             zoneFlip = fZone.flipMap()[fZone.whichFace(masterFaceI)];
00259         }
00260     }
00261 
00262     label dualFaceI;
00263 
00264     if (dualCell0 < dualCell1)
00265     {
00266         dualFaceI = meshMod.addFace
00267         (
00268             newFace,
00269             dualCell0,      // own
00270             dualCell1,      // nei
00271             masterPointI,   // masterPointID
00272             masterEdgeI,    // masterEdgeID
00273             masterFaceI,    // masterFaceID
00274             false,          // flipFaceFlux
00275             -1,             // patchID
00276             zoneID,         // zoneID
00277             zoneFlip        // zoneFlip
00278         );
00279 
00280         //pointField dualPoints(meshMod.points());
00281         //vector n(newFace.normal(dualPoints));
00282         //n /= mag(n);
00283         //Pout<< "Generated internal dualFace:" << dualFaceI
00284         //    << " verts:" << newFace
00285         //    << " points:" << UIndirectList<point>(meshMod.points(), newFace)()
00286         //    << " n:" << n
00287         //    << " between dualowner:" << dualCell0
00288         //    << " dualneigbour:" << dualCell1
00289         //    << endl;
00290     }
00291     else
00292     {
00293         dualFaceI = meshMod.addFace
00294         (
00295             newFace,
00296             dualCell1,      // own
00297             dualCell0,      // nei
00298             masterPointI,   // masterPointID
00299             masterEdgeI,    // masterEdgeID
00300             masterFaceI,    // masterFaceID
00301             false,          // flipFaceFlux
00302             -1,             // patchID
00303             zoneID,         // zoneID
00304             zoneFlip        // zoneFlip
00305         );
00306 
00307         //pointField dualPoints(meshMod.points());
00308         //vector n(newFace.normal(dualPoints));
00309         //n /= mag(n);
00310         //Pout<< "Generated internal dualFace:" << dualFaceI
00311         //    << " verts:" << newFace
00312         //    << " points:" << UIndirectList<point>(meshMod.points(), newFace)()
00313         //    << " n:" << n
00314         //    << " between dualowner:" << dualCell1
00315         //    << " dualneigbour:" << dualCell0
00316         //    << endl;
00317     }
00318     return dualFaceI;
00319 }
00320 
00321 
00322 Foam::label Foam::meshDualiser::addBoundaryFace
00323 (
00324     const label masterPointI,
00325     const label masterEdgeI,
00326     const label masterFaceI,
00327 
00328     const label dualCellI,
00329     const label patchI,
00330     const DynamicList<label>& verts,
00331     polyTopoChange& meshMod
00332 ) const
00333 {
00334     face newFace(verts);
00335 
00336     label zoneID = -1;
00337     bool zoneFlip = false;
00338     if (masterFaceI != -1)
00339     {
00340         zoneID = mesh_.faceZones().whichZone(masterFaceI);
00341 
00342         if (zoneID != -1)
00343         {
00344             const faceZone& fZone = mesh_.faceZones()[zoneID];
00345 
00346             zoneFlip = fZone.flipMap()[fZone.whichFace(masterFaceI)];
00347         }
00348     }
00349 
00350     label dualFaceI = meshMod.addFace
00351     (
00352         newFace,
00353         dualCellI,      // own
00354         -1,             // nei
00355         masterPointI,   // masterPointID
00356         masterEdgeI,    // masterEdgeID
00357         masterFaceI,    // masterFaceID
00358         false,          // flipFaceFlux
00359         patchI,         // patchID
00360         zoneID,         // zoneID
00361         zoneFlip        // zoneFlip
00362     );
00363 
00364     //pointField dualPoints(meshMod.points());
00365     //vector n(newFace.normal(dualPoints));
00366     //n /= mag(n);
00367     //Pout<< "Generated boundary dualFace:" << dualFaceI
00368     //    << " verts:" << newFace
00369     //    << " points:" << UIndirectList<point>(meshMod.points(), newFace)()
00370     //    << " n:" << n
00371     //    << " on dualowner:" << dualCellI
00372     //    << endl;
00373     return dualFaceI;
00374 }
00375 
00376 
00377 // Walks around edgeI.
00378 // splitFace=true : creates multiple faces
00379 // splitFace=false: creates single face if same dual cells on both sides,
00380 //                  multiple faces otherwise.
00381 void Foam::meshDualiser::createFacesAroundEdge
00382 (
00383     const bool splitFace,
00384     const PackedBoolList& isBoundaryEdge,
00385     const label edgeI,
00386     const label startFaceI,
00387     polyTopoChange& meshMod,
00388     boolList& doneEFaces
00389 ) const
00390 {
00391     const edge& e = mesh_.edges()[edgeI];
00392     const labelList& eFaces = mesh_.edgeFaces()[edgeI];
00393 
00394     label fp = edgeFaceCirculator::getMinIndex
00395     (
00396         mesh_.faces()[startFaceI],
00397         e[0],
00398         e[1]
00399     );
00400 
00401     edgeFaceCirculator ie
00402     (
00403         mesh_,
00404         startFaceI, // face
00405         true,       // ownerSide
00406         fp,         // fp
00407         isBoundaryEdge.get(edgeI) == 1  // isBoundaryEdge
00408     );
00409     ie.setCanonical();
00410 
00411     bool edgeOrder = ie.sameOrder(e[0], e[1]);
00412     label startFaceLabel = ie.faceLabel();
00413 
00414     //Pout<< "At edge:" << edgeI << " verts:" << e
00415     //    << " points:" << mesh_.points()[e[0]] << mesh_.points()[e[1]]
00416     //    << " started walking at face:" << ie.faceLabel()
00417     //    << " verts:" << mesh_.faces()[ie.faceLabel()]
00418     //    << " edgeOrder:" << edgeOrder
00419     //    << " in direction of cell:" << ie.cellLabel()
00420     //    << endl;
00421 
00422     // Walk and collect face.
00423     DynamicList<label> verts(100);
00424 
00425     if (edgeToDualPoint_[edgeI] != -1)
00426     {
00427         verts.append(edgeToDualPoint_[edgeI]);
00428     }
00429     if (faceToDualPoint_[ie.faceLabel()] != -1)
00430     {
00431         doneEFaces[findIndex(eFaces, ie.faceLabel())] = true;
00432         verts.append(faceToDualPoint_[ie.faceLabel()]);
00433     }
00434     if (cellToDualPoint_[ie.cellLabel()] != -1)
00435     {
00436         verts.append(cellToDualPoint_[ie.cellLabel()]);
00437     }
00438 
00439     label currentDualCell0 = findDualCell(ie.cellLabel(), e[0]);
00440     label currentDualCell1 = findDualCell(ie.cellLabel(), e[1]);
00441 
00442     ++ie;
00443 
00444     while (true)
00445     {
00446         label faceI = ie.faceLabel();
00447 
00448         // Mark face as visited.
00449         doneEFaces[findIndex(eFaces, faceI)] = true;
00450 
00451         if (faceToDualPoint_[faceI] != -1)
00452         {
00453             verts.append(faceToDualPoint_[faceI]);
00454         }
00455 
00456         label cellI = ie.cellLabel();
00457 
00458         if (cellI == -1)
00459         {
00460             // At ending boundary face. We've stored the face point above
00461             // so this is the whole face.
00462             break;
00463         }
00464 
00465 
00466         label dualCell0 = findDualCell(cellI, e[0]);
00467         label dualCell1 = findDualCell(cellI, e[1]);
00468 
00469         // Generate face. (always if splitFace=true; only if needed to
00470         // separate cells otherwise)
00471         if
00472         (
00473             splitFace
00474          || (
00475                 dualCell0 != currentDualCell0
00476              || dualCell1 != currentDualCell1
00477             )
00478         )
00479         {
00480             // Close current face.
00481             addInternalFace
00482             (
00483                 -1,         // masterPointI
00484                 edgeI,      // masterEdgeI
00485                 -1,         // masterFaceI
00486                 edgeOrder,
00487                 currentDualCell0,
00488                 currentDualCell1,
00489                 verts.shrink(),
00490                 meshMod
00491             );
00492 
00493             // Restart
00494             currentDualCell0 = dualCell0;
00495             currentDualCell1 = dualCell1;
00496 
00497             verts.clear();
00498             if (edgeToDualPoint_[edgeI] != -1)
00499             {
00500                 verts.append(edgeToDualPoint_[edgeI]);
00501             }
00502             if (faceToDualPoint_[faceI] != -1)
00503             {
00504                 verts.append(faceToDualPoint_[faceI]);
00505             }
00506         }
00507 
00508         if (cellToDualPoint_[cellI] != -1)
00509         {
00510             verts.append(cellToDualPoint_[cellI]);
00511         }
00512 
00513         ++ie;
00514 
00515         if (ie == ie.end())
00516         {
00517             // Back at start face (for internal edge only). See if this needs
00518             // adding.
00519             if (isBoundaryEdge.get(edgeI) == 0)
00520             {
00521                 label startDual = faceToDualPoint_[startFaceLabel];
00522 
00523                 if (startDual != -1 && findIndex(verts, startDual) == -1)
00524                 {
00525                     verts.append(startDual);
00526                 }
00527             }
00528             break;
00529         }
00530     }
00531 
00532     verts.shrink();
00533     addInternalFace
00534     (
00535         -1,         // masterPointI
00536         edgeI,      // masterEdgeI
00537         -1,         // masterFaceI
00538         edgeOrder,
00539         currentDualCell0,
00540         currentDualCell1,
00541         verts,
00542         meshMod
00543     );
00544 }
00545 
00546 
00547 // Walks around circumference of faceI. Creates single face. Gets given
00548 // starting (feature) edge to start from. Returns ending edge. (all edges
00549 // in form of index in faceEdges)
00550 void Foam::meshDualiser::createFaceFromInternalFace
00551 (
00552     const label faceI,
00553     label& fp,
00554     polyTopoChange& meshMod
00555 ) const
00556 {
00557     const face& f = mesh_.faces()[faceI];
00558     const labelList& fEdges = mesh_.faceEdges()[faceI];
00559     label own = mesh_.faceOwner()[faceI];
00560     label nei = mesh_.faceNeighbour()[faceI];
00561 
00562     //Pout<< "createFaceFromInternalFace : At face:" << faceI
00563     //    << " verts:" << f
00564     //    << " points:" << UIndirectList<point>(mesh_.points(), f)()
00565     //    << " started walking at edge:" << fEdges[fp]
00566     //    << " verts:" << mesh_.edges()[fEdges[fp]]
00567     //    << endl;
00568 
00569 
00570     // Walk and collect face.
00571     DynamicList<label> verts(100);
00572 
00573     verts.append(faceToDualPoint_[faceI]);
00574     verts.append(edgeToDualPoint_[fEdges[fp]]);
00575 
00576     // Step to vertex after edge mid
00577     fp = f.fcIndex(fp);
00578 
00579     // Get cells on either side of face at that point
00580     label currentDualCell0 = findDualCell(own, f[fp]);
00581     label currentDualCell1 = findDualCell(nei, f[fp]);
00582 
00583     forAll(f, i)
00584     {
00585         // Check vertex
00586         if (pointToDualPoint_[f[fp]] != -1)
00587         {
00588             verts.append(pointToDualPoint_[f[fp]]);
00589         }
00590 
00591         // Edge between fp and fp+1
00592         label edgeI = fEdges[fp];
00593 
00594         if (edgeToDualPoint_[edgeI] != -1)
00595         {
00596             verts.append(edgeToDualPoint_[edgeI]);
00597         }
00598 
00599         // Next vertex on edge
00600         label nextFp = f.fcIndex(fp);
00601 
00602         // Get dual cells on nextFp to check whether face needs closing.
00603         label dualCell0 = findDualCell(own, f[nextFp]);
00604         label dualCell1 = findDualCell(nei, f[nextFp]);
00605 
00606         if (dualCell0 != currentDualCell0 || dualCell1 != currentDualCell1)
00607         {
00608             // Check: make sure that there is a midpoint on the edge.
00609             if (edgeToDualPoint_[edgeI] == -1)
00610             {
00611                 FatalErrorIn("createFacesFromInternalFace(..)")
00612                     << "face:" << faceI << " verts:" << f
00613                     << " points:" << UIndirectList<point>(mesh_.points(), f)()
00614                     << " no feature edge between " << f[fp]
00615                     << " and " << f[nextFp] << " although have different"
00616                     << " dual cells." << endl
00617                     << "point " << f[fp] << " has dual cells "
00618                     << currentDualCell0 << " and " << currentDualCell1
00619                     << " ; point "<< f[nextFp] << " has dual cells "
00620                     << dualCell0 << " and " << dualCell1
00621                     << abort(FatalError);
00622             }
00623 
00624 
00625             // Close current face.
00626             verts.shrink();
00627             addInternalFace
00628             (
00629                 -1,         // masterPointI
00630                 -1,         // masterEdgeI
00631                 faceI,      // masterFaceI
00632                 true,       // edgeOrder,
00633                 currentDualCell0,
00634                 currentDualCell1,
00635                 verts,
00636                 meshMod
00637             );
00638             break;
00639         }
00640 
00641         fp = nextFp;
00642     }
00643 }
00644 
00645 
00646 // Given a point on a face converts the faces around the point.
00647 // (pointFaces()). Gets starting face and marks off visited faces in donePFaces.
00648 void Foam::meshDualiser::createFacesAroundBoundaryPoint
00649 (
00650     const label patchI,
00651     const label patchPointI,
00652     const label startFaceI,
00653     polyTopoChange& meshMod,
00654     boolList& donePFaces            // pFaces visited
00655 ) const
00656 {
00657     const polyBoundaryMesh& patches = mesh_.boundaryMesh();
00658     const polyPatch& pp = patches[patchI];
00659     const labelList& pFaces = pp.pointFaces()[patchPointI];
00660     const labelList& own = mesh_.faceOwner();
00661 
00662     label pointI = pp.meshPoints()[patchPointI];
00663 
00664     if (pointToDualPoint_[pointI] == -1)
00665     {
00666         // Not a feature point. Loop over all connected
00667         // pointFaces.
00668 
00669         // Starting face
00670         label faceI = startFaceI;
00671 
00672         DynamicList<label> verts(4);
00673 
00674         while (true)
00675         {
00676             label index = findIndex(pFaces, faceI-pp.start());
00677 
00678             // Has face been visited already?
00679             if (donePFaces[index])
00680             {
00681                 break;
00682             }
00683             donePFaces[index] = true;
00684 
00685             // Insert face centre
00686             verts.append(faceToDualPoint_[faceI]);
00687 
00688             label dualCellI = findDualCell(own[faceI], pointI);
00689 
00690             // Get the edge before the patchPointI
00691             const face& f = mesh_.faces()[faceI];
00692             label fp = findIndex(f, pointI);
00693             label prevFp = f.rcIndex(fp);
00694             label edgeI = mesh_.faceEdges()[faceI][prevFp];
00695 
00696             if (edgeToDualPoint_[edgeI] != -1)
00697             {
00698                 verts.append(edgeToDualPoint_[edgeI]);
00699             }
00700 
00701             // Get next boundary face (whilst staying on edge).
00702             edgeFaceCirculator circ
00703             (
00704                 mesh_,
00705                 faceI,
00706                 true,   // ownerSide
00707                 prevFp, // index of edge in face
00708                 true    // isBoundaryEdge
00709             );
00710 
00711             do
00712             {
00713                 ++circ;
00714             }
00715             while (mesh_.isInternalFace(circ.faceLabel()));
00716 
00717             // Step to next face
00718             faceI = circ.faceLabel();
00719 
00720             if (faceI < pp.start() || faceI >= pp.start()+pp.size())
00721             {
00722                 FatalErrorIn
00723                 (
00724                     "meshDualiser::createFacesAroundBoundaryPoint(..)"
00725                 )   << "Walked from face on patch:" << patchI
00726                     << " to face:" << faceI
00727                     << " fc:" << mesh_.faceCentres()[faceI]
00728                     << " on patch:" << patches.whichPatch(faceI)
00729                     << abort(FatalError);
00730             }
00731 
00732             // Check if different cell.
00733             if (dualCellI != findDualCell(own[faceI], pointI))
00734             {
00735                 FatalErrorIn
00736                 (
00737                     "meshDualiser::createFacesAroundBoundaryPoint(..)"
00738                 )   << "Different dual cells but no feature edge"
00739                     << " inbetween point:" << pointI
00740                     << " coord:" << mesh_.points()[pointI]
00741                     << abort(FatalError);
00742             }
00743         }
00744 
00745         verts.shrink();
00746 
00747         label dualCellI = findDualCell(own[faceI], pointI);
00748 
00749         //Bit dodgy: create dualface from the last face (instead of from
00750         // the central point). This will also use the original faceZone to
00751         // put the new face (which might span multiple original faces) in.
00752 
00753         addBoundaryFace
00754         (
00755             //pointI,     // masterPointI
00756             -1,         // masterPointI
00757             -1,         // masterEdgeI
00758             faceI,      // masterFaceI
00759             dualCellI,
00760             patchI,
00761             verts,
00762             meshMod
00763         );
00764     }
00765     else
00766     {
00767         label faceI = startFaceI;
00768 
00769         // Storage for face
00770         DynamicList<label> verts(mesh_.faces()[faceI].size());
00771 
00772         // Starting point.
00773         verts.append(pointToDualPoint_[pointI]);
00774 
00775         // Find edge between pointI and next point on face.
00776         const labelList& fEdges = mesh_.faceEdges()[faceI];
00777         label nextEdgeI = fEdges[findIndex(mesh_.faces()[faceI], pointI)];
00778         if (edgeToDualPoint_[nextEdgeI] != -1)
00779         {
00780             verts.append(edgeToDualPoint_[nextEdgeI]);
00781         }
00782 
00783         do
00784         {
00785             label index = findIndex(pFaces, faceI-pp.start());
00786 
00787             // Has face been visited already?
00788             if (donePFaces[index])
00789             {
00790                 break;
00791             }
00792             donePFaces[index] = true;
00793 
00794             // Face centre
00795             verts.append(faceToDualPoint_[faceI]);
00796 
00797             // Find edge before pointI on faceI
00798             const labelList& fEdges = mesh_.faceEdges()[faceI];
00799             const face& f = mesh_.faces()[faceI];
00800             label prevFp = f.rcIndex(findIndex(f, pointI));
00801             label edgeI = fEdges[prevFp];
00802 
00803             if (edgeToDualPoint_[edgeI] != -1)
00804             {
00805                 // Feature edge. Close any face so far. Note: uses face to
00806                 // create dualFace from. Could use pointI instead.
00807                 verts.append(edgeToDualPoint_[edgeI]);
00808                 addBoundaryFace
00809                 (
00810                     -1,     // masterPointI
00811                     -1,     // masterEdgeI
00812                     faceI,  // masterFaceI
00813                     findDualCell(own[faceI], pointI),
00814                     patchI,
00815                     verts.shrink(),
00816                     meshMod
00817                 );
00818                 verts.clear();
00819 
00820                 verts.append(pointToDualPoint_[pointI]);
00821                 verts.append(edgeToDualPoint_[edgeI]);
00822             }
00823 
00824             // Cross edgeI to next boundary face
00825             edgeFaceCirculator circ
00826             (
00827                 mesh_,
00828                 faceI,
00829                 true,   // ownerSide
00830                 prevFp, // index of edge in face
00831                 true    // isBoundaryEdge
00832             );
00833 
00834             do
00835             {
00836                 ++circ;
00837             }
00838             while (mesh_.isInternalFace(circ.faceLabel()));
00839 
00840             // Step to next face. Quit if not on same patch.
00841             faceI = circ.faceLabel();
00842         }
00843         while
00844         (
00845             faceI != startFaceI
00846          && faceI >= pp.start()
00847          && faceI < pp.start()+pp.size()
00848         );
00849 
00850         if (verts.size() > 2)
00851         {
00852             // Note: face created from face, not from pointI
00853             addBoundaryFace
00854             (
00855                 -1,             // masterPointI
00856                 -1,             // masterEdgeI
00857                 startFaceI,     // masterFaceI
00858                 findDualCell(own[faceI], pointI),
00859                 patchI,
00860                 verts.shrink(),
00861                 meshMod
00862             );
00863         }
00864     }
00865 }
00866 
00867 
00868 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
00869 
00870 Foam::meshDualiser::meshDualiser(const polyMesh& mesh)
00871 :
00872     mesh_(mesh),
00873     pointToDualCells_(mesh_.nPoints()),
00874     pointToDualPoint_(mesh_.nPoints(), -1),
00875     cellToDualPoint_(mesh_.nCells()),
00876     faceToDualPoint_(mesh_.nFaces(), -1),
00877     edgeToDualPoint_(mesh_.nEdges(), -1)
00878 {}
00879 
00880 
00881 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
00882 
00883 
00884 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
00885 
00886 void Foam::meshDualiser::setRefinement
00887 (
00888     const bool splitFace,
00889     const labelList& featureFaces,
00890     const labelList& featureEdges,
00891     const labelList& singleCellFeaturePoints,
00892     const labelList& multiCellFeaturePoints,
00893     polyTopoChange& meshMod
00894 )
00895 {
00896     const labelList& own = mesh_.faceOwner();
00897     const labelList& nei = mesh_.faceNeighbour();
00898     const vectorField& cellCentres = mesh_.cellCentres();
00899 
00900     // Mark boundary edges and points.
00901     // (Note: in 1.4.2 we can use the built-in mesh point ordering
00902     //  facility instead)
00903     PackedBoolList isBoundaryEdge(mesh_.nEdges());
00904     for (label faceI = mesh_.nInternalFaces(); faceI < mesh_.nFaces(); faceI++)
00905     {
00906         const labelList& fEdges = mesh_.faceEdges()[faceI];
00907 
00908         forAll(fEdges, i)
00909         {
00910             isBoundaryEdge.set(fEdges[i], 1);
00911         }
00912     }
00913 
00914 
00915     if (splitFace)
00916     {
00917         // This is a special mode where whenever we are walking around an edge
00918         // every area through a cell becomes a separate dualface. So two
00919         // dual cells will probably have more than one dualface between them!
00920         // This mode implies that
00921         // - all faces have to be feature faces since there has to be a
00922         //   dualpoint at the face centre.
00923         // - all edges have to be feature edges ,,
00924         boolList featureFaceSet(mesh_.nFaces(), false);
00925         forAll(featureFaces, i)
00926         {
00927             featureFaceSet[featureFaces[i]] = true;
00928         }
00929         label faceI = findIndex(featureFaceSet, false);
00930 
00931         if (faceI != -1)
00932         {
00933             FatalErrorIn
00934             (
00935                 "meshDualiser::setRefinement"
00936                 "(const labelList&, const labelList&, const labelList&"
00937                 ", const labelList, polyTopoChange&)"
00938             )   << "In split-face-mode (splitFace=true) but not all faces"
00939                 << " marked as feature faces." << endl
00940                 << "First conflicting face:" << faceI
00941                 << " centre:" << mesh_.faceCentres()[faceI]
00942                 << abort(FatalError);
00943         }
00944 
00945         boolList featureEdgeSet(mesh_.nEdges(), false);
00946         forAll(featureEdges, i)
00947         {
00948             featureEdgeSet[featureEdges[i]] = true;
00949         }
00950         label edgeI = findIndex(featureEdgeSet, false);
00951 
00952         if (edgeI != -1)
00953         {
00954             const edge& e = mesh_.edges()[edgeI];
00955             FatalErrorIn
00956             (
00957                 "meshDualiser::setRefinement"
00958                 "(const labelList&, const labelList&, const labelList&"
00959                 ", const labelList, polyTopoChange&)"
00960             )   << "In split-face-mode (splitFace=true) but not all edges"
00961                 << " marked as feature edges." << endl
00962                 << "First conflicting edge:" << edgeI
00963                 << " verts:" << e
00964                 << " coords:" << mesh_.points()[e[0]] << mesh_.points()[e[1]]
00965                 << abort(FatalError);
00966         }
00967     }
00968     else
00969     {
00970         // Check that all boundary faces are feature faces.
00971 
00972         boolList featureFaceSet(mesh_.nFaces(), false);
00973         forAll(featureFaces, i)
00974         {
00975             featureFaceSet[featureFaces[i]] = true;
00976         }
00977         for
00978         (
00979             label faceI = mesh_.nInternalFaces();
00980             faceI < mesh_.nFaces();
00981             faceI++
00982         )
00983         {
00984             if (!featureFaceSet[faceI])
00985             {
00986                 FatalErrorIn
00987                 (
00988                     "meshDualiser::setRefinement"
00989                     "(const labelList&, const labelList&, const labelList&"
00990                     ", const labelList, polyTopoChange&)"
00991                 )   << "Not all boundary faces marked as feature faces."
00992                     << endl
00993                     << "First conflicting face:" << faceI
00994                     << " centre:" << mesh_.faceCentres()[faceI]
00995                     << abort(FatalError);
00996             }
00997         }
00998     }
00999 
01000 
01001 
01002 
01003     // Start creating cells, points, and faces (in that order)
01004 
01005 
01006     // 1. Mark which cells to create
01007     // Mostly every point becomes one cell but sometimes (for feature points)
01008     // all cells surrounding a feature point become cells. Also a non-manifold
01009     // point can create two cells! So a dual cell is uniquely defined by a
01010     // mesh point + cell (as in pointCells index)
01011 
01012     // 2. Mark which face centres to create
01013 
01014     // 3. Internal faces can now consist of
01015     //      - only cell centres of walk around edge
01016     //      - cell centres + face centres of walk around edge
01017     //      - same but now other side is not a single cell
01018 
01019     // 4. Boundary faces (or internal faces between cell zones!) now consist of
01020     //      - walk around boundary point.
01021 
01022 
01023 
01024     autoPtr<OFstream> dualCcStr;
01025     if (debug)
01026     {
01027         dualCcStr.reset(new OFstream("dualCc.obj"));
01028         Pout<< "Dumping centres of dual cells to " << dualCcStr().name()
01029             << endl;
01030     }
01031 
01032 
01033     // Dual cells (from points)
01034     // ~~~~~~~~~~~~~~~~~~~~~~~~
01035 
01036     // pointToDualCells_[pointI]
01037     // - single entry : all cells surrounding point all become the same
01038     //                  cell.
01039     // - multiple entries: in order of pointCells.
01040 
01041 
01042     // feature points that become single cell
01043     forAll(singleCellFeaturePoints, i)
01044     {
01045         label pointI = singleCellFeaturePoints[i];
01046 
01047         pointToDualPoint_[pointI] = meshMod.addPoint
01048         (
01049             mesh_.points()[pointI],
01050             pointI,                                 // masterPoint
01051             mesh_.pointZones().whichZone(pointI),   // zoneID
01052             true                                    // inCell
01053         );
01054 
01055         // Generate single cell
01056         pointToDualCells_[pointI].setSize(1);
01057         pointToDualCells_[pointI][0] = meshMod.addCell
01058         (
01059             pointI, //masterPointID,
01060             -1,     //masterEdgeID,
01061             -1,     //masterFaceID,
01062             -1,     //masterCellID,
01063             -1      //zoneID
01064         );
01065         if (dualCcStr.valid())
01066         {
01067             meshTools::writeOBJ(dualCcStr(), mesh_.points()[pointI]);
01068         }
01069     }
01070 
01071     // feature points that become multiple cells
01072     forAll(multiCellFeaturePoints, i)
01073     {
01074         label pointI = multiCellFeaturePoints[i];
01075 
01076         if (pointToDualCells_[pointI].size() > 0)
01077         {
01078             FatalErrorIn
01079             (
01080                 "meshDualiser::setRefinement"
01081                 "(const labelList&, const labelList&, const labelList&"
01082                 ", const labelList, polyTopoChange&)"
01083             )   << "Point " << pointI << " at:" << mesh_.points()[pointI]
01084                 << " is both in singleCellFeaturePoints"
01085                 << " and multiCellFeaturePoints."
01086                 << abort(FatalError);
01087         }
01088 
01089         pointToDualPoint_[pointI] = meshMod.addPoint
01090         (
01091             mesh_.points()[pointI],
01092             pointI,                                 // masterPoint
01093             mesh_.pointZones().whichZone(pointI),   // zoneID
01094             true                                    // inCell
01095         );
01096 
01097         // Create dualcell for every cell connected to dual point
01098 
01099         const labelList& pCells = mesh_.pointCells()[pointI];
01100 
01101         pointToDualCells_[pointI].setSize(pCells.size());
01102 
01103         forAll(pCells, pCellI)
01104         {
01105             pointToDualCells_[pointI][pCellI] = meshMod.addCell
01106             (
01107                 pointI,                                     //masterPointID
01108                 -1,                                         //masterEdgeID
01109                 -1,                                         //masterFaceID
01110                 -1,                                         //masterCellID
01111                 mesh_.cellZones().whichZone(pCells[pCellI]) //zoneID
01112             );
01113             if (dualCcStr.valid())
01114             {
01115                 meshTools::writeOBJ
01116                 (
01117                     dualCcStr(),
01118                     0.5*(mesh_.points()[pointI]+cellCentres[pCells[pCellI]])
01119                 );
01120             }
01121         }
01122     }
01123     // Normal points
01124     forAll(mesh_.points(), pointI)
01125     {
01126         if (pointToDualCells_[pointI].empty())
01127         {
01128             pointToDualCells_[pointI].setSize(1);
01129             pointToDualCells_[pointI][0] = meshMod.addCell
01130             (
01131                 pointI, //masterPointID,
01132                 -1,     //masterEdgeID,
01133                 -1,     //masterFaceID,
01134                 -1,     //masterCellID,
01135                 -1      //zoneID
01136             );
01137 
01138             if (dualCcStr.valid())
01139             {
01140                 meshTools::writeOBJ(dualCcStr(), mesh_.points()[pointI]);
01141             }
01142         }
01143     }
01144 
01145 
01146     // Dual points (from cell centres, feature faces, feature edges)
01147     // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
01148 
01149     forAll(cellToDualPoint_, cellI)
01150     {
01151         cellToDualPoint_[cellI] = meshMod.addPoint
01152         (
01153             cellCentres[cellI],
01154             mesh_.faces()[mesh_.cells()[cellI][0]][0],     // masterPoint
01155             -1,     // zoneID
01156             true    // inCell
01157         );
01158     }
01159 
01160     // From face to dual point
01161 
01162     forAll(featureFaces, i)
01163     {
01164         label faceI = featureFaces[i];
01165 
01166         faceToDualPoint_[faceI] = meshMod.addPoint
01167         (
01168             mesh_.faceCentres()[faceI],
01169             mesh_.faces()[faceI][0],     // masterPoint
01170             -1,     // zoneID
01171             true    // inCell
01172         );
01173     }
01174     // Detect whether different dual cells on either side of a face. This
01175     // would neccesitate having a dual face built from the face and thus a
01176     // dual point at the face centre.
01177     for (label faceI = 0; faceI < mesh_.nInternalFaces(); faceI++)
01178     {
01179         if (faceToDualPoint_[faceI] == -1)
01180         {
01181             const face& f = mesh_.faces()[faceI];
01182 
01183             forAll(f, fp)
01184             {
01185                 label ownDualCell = findDualCell(own[faceI], f[fp]);
01186                 label neiDualCell = findDualCell(nei[faceI], f[fp]);
01187 
01188                 if (ownDualCell != neiDualCell)
01189                 {
01190                     faceToDualPoint_[faceI] = meshMod.addPoint
01191                     (
01192                         mesh_.faceCentres()[faceI],
01193                         f[fp],  // masterPoint
01194                         -1,     // zoneID
01195                         true    // inCell
01196                     );
01197 
01198                     break;
01199                 }
01200             }
01201         }
01202     }
01203 
01204     // From edge to dual point
01205 
01206     forAll(featureEdges, i)
01207     {
01208         label edgeI = featureEdges[i];
01209 
01210         const edge& e = mesh_.edges()[edgeI];
01211 
01212         edgeToDualPoint_[edgeI] = meshMod.addPoint
01213         (
01214             e.centre(mesh_.points()),
01215             e[0],   // masterPoint
01216             -1,     // zoneID
01217             true    // inCell
01218         );
01219     }
01220 
01221     // Detect whether different dual cells on either side of an edge. This
01222     // would neccesitate having a dual face built perpendicular to the edge
01223     // and thus a dual point at the mid of the edge.
01224     // Note: not really true - the face can be built without the edge centre!
01225     const labelListList& edgeCells = mesh_.edgeCells();
01226 
01227     forAll(edgeCells, edgeI)
01228     {
01229        if (edgeToDualPoint_[edgeI] == -1)
01230        {
01231             const edge& e = mesh_.edges()[edgeI];
01232 
01233             // We need a point on the edge if not all cells on both sides
01234             // are the same.
01235 
01236             const labelList& eCells = mesh_.edgeCells()[edgeI];
01237 
01238             label dualE0 = findDualCell(eCells[0], e[0]);
01239             label dualE1 = findDualCell(eCells[0], e[1]);
01240 
01241             for (label i = 1; i < eCells.size(); i++)
01242             {
01243                 label newDualE0 = findDualCell(eCells[i], e[0]);
01244 
01245                 if (dualE0 != newDualE0)
01246                 {
01247                     edgeToDualPoint_[edgeI] = meshMod.addPoint
01248                     (
01249                         e.centre(mesh_.points()),
01250                         e[0],                               // masterPoint
01251                         mesh_.pointZones().whichZone(e[0]), // zoneID
01252                         true                                // inCell
01253                     );
01254 
01255                     break;
01256                 }
01257 
01258                 label newDualE1 = findDualCell(eCells[i], e[1]);
01259 
01260                 if (dualE1 != newDualE1)
01261                 {
01262                     edgeToDualPoint_[edgeI] = meshMod.addPoint
01263                     (
01264                         e.centre(mesh_.points()),
01265                         e[1],   // masterPoint
01266                         mesh_.pointZones().whichZone(e[1]), // zoneID
01267                         true    // inCell
01268                     );
01269 
01270                     break;
01271                 }
01272             }
01273         }
01274     }
01275 
01276     // Make sure all boundary edges emanating from feature points are
01277     // feature edges as well.
01278     forAll(singleCellFeaturePoints, i)
01279     {
01280         generateDualBoundaryEdges
01281         (
01282             isBoundaryEdge,
01283             singleCellFeaturePoints[i],
01284             meshMod
01285         );
01286     }
01287     forAll(multiCellFeaturePoints, i)
01288     {
01289         generateDualBoundaryEdges
01290         (
01291             isBoundaryEdge,
01292             multiCellFeaturePoints[i],
01293             meshMod
01294         );
01295     }
01296 
01297 
01298     // Check for duplicate points
01299     if (debug)
01300     {
01301         dumpPolyTopoChange(meshMod, "generatedPoints_");
01302         checkPolyTopoChange(meshMod);
01303     }
01304 
01305 
01306     // Now we have all points and cells
01307     // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
01308     //  - pointToDualCells_ : per point a single dualCell or multiple dualCells
01309     //  - pointToDualPoint_ : per point -1 or the dual point at the coordinate
01310     //  - edgeToDualPoint_  : per edge -1 or the edge centre
01311     //  - faceToDualPoint_  : per face -1 or the face centre
01312     //  - cellToDualPoint_  : per cell the cell centre
01313     // Now we have to walk all edges and construct faces. Either single face
01314     // per edge or multiple (-if nonmanifold edge -if different dualcells)
01315 
01316     const edgeList& edges = mesh_.edges();
01317 
01318     forAll(edges, edgeI)
01319     {
01320         const labelList& eFaces = mesh_.edgeFaces()[edgeI];
01321 
01322         boolList doneEFaces(eFaces.size(), false);
01323 
01324         forAll(eFaces, i)
01325         {
01326             if (!doneEFaces[i])
01327             {
01328                 // We found a face that hasn't yet been visited. This might
01329                 // happen for non-manifold edges where a single edge can
01330                 // become multiple faces.
01331 
01332                 label startFaceI = eFaces[i];
01333 
01334                 //Pout<< "Walking edge:" << edgeI
01335                 //    << " points:" << mesh_.points()[e[0]]
01336                 //    << mesh_.points()[e[1]]
01337                 //    << " startFace:" << startFaceI
01338                 //    << " at:" << mesh_.faceCentres()[startFaceI]
01339                 //    << endl;
01340 
01341                 createFacesAroundEdge
01342                 (
01343                     splitFace,
01344                     isBoundaryEdge,
01345                     edgeI,
01346                     startFaceI,
01347                     meshMod,
01348                     doneEFaces
01349                 );
01350             }
01351         }
01352     }
01353 
01354     if (debug)
01355     {
01356         dumpPolyTopoChange(meshMod, "generatedFacesFromEdges_");
01357     }
01358 
01359     // Create faces from feature faces. These can be internal or external faces.
01360     // - feature face : centre needs to be included.
01361     //      - single cells on either side: triangulate
01362     //      - multiple cells: create single face between unique cell pair. Only
01363     //                        create face where cells differ on either side.
01364     // - non-feature face : inbetween cell zones.
01365     forAll(faceToDualPoint_, faceI)
01366     {
01367         if (faceToDualPoint_[faceI] != -1 && mesh_.isInternalFace(faceI))
01368         {
01369             const face& f = mesh_.faces()[faceI];
01370             const labelList& fEdges = mesh_.faceEdges()[faceI];
01371 
01372             // Starting edge
01373             label fp = 0;
01374 
01375             do
01376             {
01377                 // Find edge that is in dual mesh and where the
01378                 // next point (fp+1) has different dual cells on either side.
01379                 bool foundStart = false;
01380 
01381                 do
01382                 {
01383                     if
01384                     (
01385                         edgeToDualPoint_[fEdges[fp]] != -1
01386                     && !sameDualCell(faceI, f.nextLabel(fp))
01387                     )
01388                     {
01389                         foundStart = true;
01390                         break;
01391                     }
01392                     fp = f.fcIndex(fp);
01393                 }
01394                 while (fp != 0);
01395 
01396                 if (!foundStart)
01397                 {
01398                     break;
01399                 }
01400 
01401                 // Walk from edge fp and generate a face.
01402                 createFaceFromInternalFace
01403                 (
01404                     faceI,
01405                     fp,
01406                     meshMod
01407                 );
01408             }
01409             while(fp != 0);
01410         }
01411     }
01412 
01413     if (debug)
01414     {
01415         dumpPolyTopoChange(meshMod, "generatedFacesFromFeatFaces_");
01416     }
01417 
01418 
01419     // Create boundary faces. Every boundary point has one or more dualcells.
01420     // These need to be closed.
01421     const polyBoundaryMesh& patches = mesh_.boundaryMesh();
01422 
01423     forAll(patches, patchI)
01424     {
01425         const polyPatch& pp = patches[patchI];
01426 
01427         const labelListList& pointFaces = pp.pointFaces();
01428 
01429         forAll(pointFaces, patchPointI)
01430         {
01431             const labelList& pFaces = pointFaces[patchPointI];
01432 
01433             boolList donePFaces(pFaces.size(), false);
01434 
01435             forAll(pFaces, i)
01436             {
01437                 if (!donePFaces[i])
01438                 {
01439                     // Starting face
01440                     label startFaceI = pp.start()+pFaces[i];
01441 
01442                     //Pout<< "Walking around point:" << pointI
01443                     //    << " coord:" << mesh_.points()[pointI]
01444                     //    << " on patch:" << patchI
01445                     //    << " startFace:" << startFaceI
01446                     //    << " at:" << mesh_.faceCentres()[startFaceI]
01447                     //    << endl;
01448 
01449                     createFacesAroundBoundaryPoint
01450                     (
01451                         patchI,
01452                         patchPointI,
01453                         startFaceI,
01454                         meshMod,
01455                         donePFaces            // pFaces visited
01456                     );
01457                 }
01458             }
01459         }
01460     }
01461 
01462     if (debug)
01463     {
01464         dumpPolyTopoChange(meshMod, "generatedFacesFromBndFaces_");
01465     }
01466 }
01467 
01468 
01469 
01470 // * * * * * * * * * * * * * * * Member Operators  * * * * * * * * * * * * * //
01471 
01472 
01473 // * * * * * * * * * * * * * * * Friend Functions  * * * * * * * * * * * * * //
01474 
01475 
01476 // * * * * * * * * * * * * * * * Friend Operators  * * * * * * * * * * * * * //
01477 
01478 
01479 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines