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

boundaryMesh.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 "boundaryMesh.H"
00027 #include <OpenFOAM/Time.H>
00028 #include <OpenFOAM/polyMesh.H>
00029 #include <dynamicMesh/repatchPolyTopoChanger.H>
00030 #include <OpenFOAM/faceList.H>
00031 #include <meshTools/octree.H>
00032 #include "octreeDataFaceList.H"
00033 #include <triSurface/triSurface.H>
00034 #include <OpenFOAM/SortableList.H>
00035 #include <OpenFOAM/OFstream.H>
00036 
00037 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
00038 
00039 defineTypeNameAndDebug(Foam::boundaryMesh, 0);
00040 
00041 // Normal along which to divide faces into categories (used in getNearest)
00042 const Foam::vector Foam::boundaryMesh::splitNormal_(3, 2, 1);
00043 
00044 // Distance to face tolerance for getNearest
00045 const Foam::scalar Foam::boundaryMesh::distanceTol_ = 1E-2;
00046 
00047 
00048 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
00049 
00050 // Returns number of feature edges connected to pointI
00051 Foam::label Foam::boundaryMesh::nFeatureEdges(label pointI) const
00052 {
00053     label nFeats = 0;
00054 
00055     const labelList& pEdges = mesh().pointEdges()[pointI];
00056 
00057     forAll(pEdges, pEdgeI)
00058     {
00059         label edgeI = pEdges[pEdgeI];
00060 
00061         if (edgeToFeature_[edgeI] != -1)
00062         {
00063             nFeats++;
00064         }
00065     }
00066     return nFeats;
00067 }
00068 
00069 
00070 // Returns next feature edge connected to pointI
00071 Foam::label Foam::boundaryMesh::nextFeatureEdge
00072 (
00073     const label edgeI,
00074     const label vertI
00075 ) const
00076 {
00077     const labelList& pEdges = mesh().pointEdges()[vertI];
00078 
00079     forAll(pEdges, pEdgeI)
00080     {
00081         label nbrEdgeI = pEdges[pEdgeI];
00082 
00083         if (nbrEdgeI != edgeI)
00084         {
00085             label featI = edgeToFeature_[nbrEdgeI];
00086 
00087             if (featI != -1)
00088             {
00089                 return nbrEdgeI;
00090             }
00091         }
00092     }
00093 
00094     return -1;
00095 }
00096 
00097 
00098 // Finds connected feature edges, starting from startPointI and returns
00099 // feature labels (not edge labels). Marks feature edges handled in
00100 // featVisited.
00101 Foam::labelList Foam::boundaryMesh::collectSegment
00102 (
00103     const boolList& isFeaturePoint,
00104     const label startEdgeI,
00105     boolList& featVisited
00106 ) const
00107 {
00108     // Find starting feature point on edge.
00109 
00110     label edgeI = startEdgeI;
00111 
00112     const edge& e = mesh().edges()[edgeI];
00113 
00114     label vertI = e.start();
00115 
00116     while (!isFeaturePoint[vertI])
00117     {
00118         // Step to next feature edge
00119 
00120         edgeI = nextFeatureEdge(edgeI, vertI);
00121 
00122         if ((edgeI == -1) || (edgeI == startEdgeI))
00123         {
00124             break;
00125         }
00126 
00127         // Step to next vertex on edge
00128 
00129         const edge& e = mesh().edges()[edgeI];
00130 
00131         vertI = e.otherVertex(vertI);
00132     }
00133 
00134     //
00135     // Now we have:
00136     //    edgeI : first edge on this segment
00137     //    vertI : one of the endpoints of this segment
00138     //
00139     // Start walking other way and storing edges as we go along.
00140     //
00141 
00142     // Untrimmed storage for current segment
00143     labelList featLabels(featureEdges_.size());
00144 
00145     label featLabelI = 0;
00146 
00147     label initEdgeI = edgeI;
00148 
00149     do
00150     {
00151         // Mark edge as visited
00152         label featI = edgeToFeature_[edgeI];
00153 
00154         if (featI == -1)
00155         {
00156             FatalErrorIn("boundaryMesh::collectSegment")
00157                 << "Problem" << abort(FatalError);
00158         }
00159         featLabels[featLabelI++] = featI;
00160 
00161         featVisited[featI] = true;
00162 
00163         // Step to next vertex on edge
00164 
00165         const edge& e = mesh().edges()[edgeI];
00166 
00167         vertI = e.otherVertex(vertI);
00168 
00169         // Step to next feature edge
00170 
00171         edgeI = nextFeatureEdge(edgeI, vertI);
00172 
00173         if ((edgeI == -1) || (edgeI == initEdgeI))
00174         {
00175             break;
00176         }
00177     }
00178     while (!isFeaturePoint[vertI]);
00179 
00180 
00181     // Trim to size
00182     featLabels.setSize(featLabelI);
00183 
00184     return featLabels;
00185 }
00186 
00187 
00188 void Foam::boundaryMesh::markEdges
00189 (
00190     const label maxDistance,
00191     const label edgeI,
00192     const label distance,
00193     labelList& minDistance,
00194     DynamicList<label>& visited
00195 ) const
00196 {
00197     if (distance < maxDistance)
00198     {
00199         // Don't do anything if reached beyond maxDistance.
00200 
00201         if (minDistance[edgeI] == -1)
00202         {
00203             // First visit of edge. Store edge label.
00204             visited.append(edgeI);
00205         }
00206         else if (minDistance[edgeI] <= distance)
00207         {
00208             // Already done this edge
00209             return;
00210         }
00211 
00212         minDistance[edgeI] = distance;
00213 
00214         const edge& e = mesh().edges()[edgeI];
00215 
00216         // Do edges connected to e.start
00217         const labelList& startEdges = mesh().pointEdges()[e.start()];
00218 
00219         forAll(startEdges, pEdgeI)
00220         {
00221             markEdges
00222             (
00223                 maxDistance,
00224                 startEdges[pEdgeI],
00225                 distance+1,
00226                 minDistance,
00227                 visited
00228             );
00229         }
00230 
00231         // Do edges connected to e.end
00232         const labelList& endEdges = mesh().pointEdges()[e.end()];
00233 
00234         forAll(endEdges, pEdgeI)
00235         {
00236             markEdges
00237             (
00238                 maxDistance,
00239                 endEdges[pEdgeI],
00240                 distance+1,
00241                 minDistance,
00242                 visited
00243             );
00244         }
00245     }
00246 }
00247 
00248 
00249 Foam::label Foam::boundaryMesh::findPatchID
00250 (
00251     const polyPatchList& patches,
00252     const word& patchName
00253 ) const
00254 {
00255     forAll(patches, patchI)
00256     {
00257         if (patches[patchI].name() == patchName)
00258         {
00259             return patchI;
00260         }
00261     }
00262 
00263     return -1;
00264 }
00265 
00266 
00267 Foam::wordList Foam::boundaryMesh::patchNames() const
00268 {
00269     wordList names(patches_.size());
00270 
00271     forAll(patches_, patchI)
00272     {
00273         names[patchI] = patches_[patchI].name();
00274     }
00275     return names;
00276 }
00277 
00278 
00279 Foam::label Foam::boundaryMesh::whichPatch
00280 (
00281     const polyPatchList& patches,
00282     const label faceI
00283 ) const
00284 {
00285     forAll(patches, patchI)
00286     {
00287         const polyPatch& pp = patches[patchI];
00288 
00289         if ((faceI >= pp.start()) && (faceI < (pp.start() + pp.size())))
00290         {
00291             return patchI;
00292         }
00293     }
00294     return -1;
00295 }
00296 
00297 
00298 // Gets labels of changed faces and propagates them to the edges. Returns
00299 // labels of edges changed.
00300 Foam::labelList Foam::boundaryMesh::faceToEdge
00301 (
00302     const boolList& regionEdge,
00303     const label region,
00304     const labelList& changedFaces,
00305     labelList& edgeRegion
00306 ) const
00307 {
00308     labelList changedEdges(mesh().nEdges(), -1);
00309     label changedI = 0;
00310 
00311     forAll(changedFaces, i)
00312     {
00313         label faceI = changedFaces[i];
00314 
00315         const labelList& fEdges = mesh().faceEdges()[faceI];
00316 
00317         forAll(fEdges, fEdgeI)
00318         {
00319             label edgeI = fEdges[fEdgeI];
00320 
00321             if (!regionEdge[edgeI] && (edgeRegion[edgeI] == -1))
00322             {
00323                 edgeRegion[edgeI] = region;
00324 
00325                 changedEdges[changedI++] = edgeI;
00326             }
00327         }
00328     }
00329 
00330     changedEdges.setSize(changedI);
00331 
00332     return changedEdges;
00333 }
00334 
00335 
00336 // Reverse of faceToEdge: gets edges and returns faces
00337 Foam::labelList Foam::boundaryMesh::edgeToFace
00338 (
00339     const label region,
00340     const labelList& changedEdges,
00341     labelList& faceRegion
00342 ) const
00343 {
00344     labelList changedFaces(mesh().size(), -1);
00345     label changedI = 0;
00346 
00347     forAll(changedEdges, i)
00348     {
00349         label edgeI = changedEdges[i];
00350 
00351         const labelList& eFaces = mesh().edgeFaces()[edgeI];
00352 
00353         forAll(eFaces, eFaceI)
00354         {
00355             label faceI = eFaces[eFaceI];
00356 
00357             if (faceRegion[faceI] == -1)
00358             {
00359                 faceRegion[faceI] = region;
00360 
00361                 changedFaces[changedI++] = faceI;
00362             }
00363         }
00364     }
00365 
00366     changedFaces.setSize(changedI);
00367 
00368     return changedFaces;
00369 }
00370 
00371 
00372 // Finds area, starting at faceI, delimited by borderEdge
00373 void Foam::boundaryMesh::markZone
00374 (
00375     const boolList& borderEdge,
00376     label faceI,
00377     label currentZone,
00378     labelList& faceZone
00379 ) const
00380 {
00381     faceZone[faceI] = currentZone;
00382 
00383     // List of faces whose faceZone has been set.
00384     labelList changedFaces(1, faceI);
00385     // List of edges whose faceZone has been set.
00386     labelList changedEdges;
00387 
00388     // Zones on all edges.
00389     labelList edgeZone(mesh().nEdges(), -1);
00390 
00391     while(1)
00392     {
00393         changedEdges = faceToEdge
00394         (
00395             borderEdge,
00396             currentZone,
00397             changedFaces,
00398             edgeZone
00399         );
00400 
00401         if (debug)
00402         {
00403             Pout<< "From changedFaces:" << changedFaces.size()
00404                 << " to changedEdges:" << changedEdges.size()
00405                 << endl;
00406         }
00407 
00408         if (changedEdges.empty())
00409         {
00410             break;
00411         }
00412 
00413         changedFaces = edgeToFace(currentZone, changedEdges, faceZone);
00414 
00415         if (debug)
00416         {
00417             Pout<< "From changedEdges:" << changedEdges.size()
00418                 << " to changedFaces:" << changedFaces.size()
00419                 << endl;
00420         }
00421 
00422         if (changedFaces.empty())
00423         {
00424             break;
00425         }
00426     }
00427 }
00428 
00429 
00430 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
00431 
00432 // Null constructor
00433 Foam::boundaryMesh::boundaryMesh()
00434 :
00435     meshPtr_(NULL),
00436     patches_(),
00437     meshFace_(),
00438     featurePoints_(),
00439     featureEdges_(),
00440     featureToEdge_(),
00441     edgeToFeature_(),
00442     featureSegments_(),
00443     extraEdges_()
00444 {}
00445 
00446 
00447 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
00448 
00449 Foam::boundaryMesh::~boundaryMesh()
00450 {
00451     clearOut();
00452 }
00453 
00454 
00455 void Foam::boundaryMesh::clearOut()
00456 {
00457     if (meshPtr_)
00458     {
00459         delete meshPtr_;
00460 
00461         meshPtr_ = NULL;
00462     }
00463 }
00464 
00465 
00466 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
00467 
00468 void Foam::boundaryMesh::read(const polyMesh& mesh)
00469 {
00470     patches_.clear();
00471 
00472     patches_.setSize(mesh.boundaryMesh().size());
00473 
00474     // Number of boundary faces
00475     label nBFaces = mesh.nFaces() - mesh.nInternalFaces();
00476 
00477     faceList bFaces(nBFaces);
00478 
00479     meshFace_.setSize(nBFaces);
00480 
00481     label bFaceI = 0;
00482 
00483     // Collect all boundary faces.
00484     forAll(mesh.boundaryMesh(), patchI)
00485     {
00486         const polyPatch& pp = mesh.boundaryMesh()[patchI];
00487 
00488         patches_.set
00489         (
00490             patchI,
00491             new boundaryPatch
00492             (
00493                 pp.name(),
00494                 patchI,
00495                 pp.size(),
00496                 bFaceI,
00497                 pp.type()
00498             )
00499         );
00500 
00501         // Collect all faces in global numbering.
00502         forAll(pp, patchFaceI)
00503         {
00504             meshFace_[bFaceI] = pp.start() + patchFaceI;
00505 
00506             bFaces[bFaceI] = pp[patchFaceI];
00507 
00508             bFaceI++;
00509         }
00510     }
00511 
00512 
00513     if (debug)
00514     {
00515         Pout<< "read : patches now:" << endl;
00516 
00517         forAll(patches_, patchI)
00518         {
00519             const boundaryPatch& bp = patches_[patchI];
00520 
00521             Pout<< "    name  : " << bp.name() << endl
00522                 << "    size  : " << bp.size() << endl
00523                 << "    start : " << bp.start() << endl
00524                 << "    type  : " << bp.physicalType() << endl
00525                 << endl;
00526         }
00527     }
00528 
00529     //
00530     // Construct single patch for all of boundary
00531     //
00532 
00533     // Temporary primitivePatch to calculate compact points & faces.
00534     PrimitivePatch<face, List, const pointField&> globalPatch
00535     (
00536         bFaces,
00537         mesh.points()
00538     );
00539 
00540     // Store in local(compact) addressing
00541     clearOut();
00542 
00543     meshPtr_ = new bMesh(globalPatch.localFaces(), globalPatch.localPoints());
00544 
00545 
00546     if (debug & 2)
00547     {
00548         const bMesh& msh = *meshPtr_;
00549 
00550         Pout<< "** Start of Faces **" << endl;
00551 
00552         forAll(msh, faceI)
00553         {
00554             const face& f = msh[faceI];
00555 
00556             point ctr(vector::zero);
00557 
00558             forAll(f, fp)
00559             {
00560                 ctr += msh.points()[f[fp]];
00561             }
00562             ctr /= f.size();
00563 
00564             Pout<< "    " << faceI
00565                 << " ctr:" << ctr
00566                 << " verts:" << f
00567                 << endl;
00568         }
00569 
00570         Pout<< "** End of Faces **" << endl;
00571 
00572         Pout<< "** Start of Points **" << endl;
00573 
00574         forAll(msh.points(), pointI)
00575         {
00576             Pout<< "    " << pointI
00577                 << " coord:" << msh.points()[pointI]
00578                 << endl;
00579         }
00580 
00581         Pout<< "** End of Points **" << endl;
00582     }
00583 
00584     // Clear edge storage
00585     featurePoints_.setSize(0);
00586     featureEdges_.setSize(0);
00587 
00588     featureToEdge_.setSize(0);
00589     edgeToFeature_.setSize(meshPtr_->nEdges());
00590     edgeToFeature_ = -1;
00591 
00592     featureSegments_.setSize(0);
00593 
00594     extraEdges_.setSize(0);
00595 }
00596 
00597 
00598 void Foam::boundaryMesh::readTriSurface(const fileName& fName)
00599 {
00600     triSurface surf(fName);
00601 
00602     if (surf.empty())
00603     {
00604         return;
00605     }
00606 
00607     // Sort according to region
00608     SortableList<label> regions(surf.size());
00609 
00610     forAll(surf, triI)
00611     {
00612         regions[triI] = surf[triI].region();
00613     }
00614     regions.sort();
00615 
00616     // Determine region mapping.
00617     Map<label> regionToBoundaryPatch;
00618 
00619     label oldRegion = -1111;
00620     label boundPatch = 0;
00621 
00622     forAll(regions, i)
00623     {
00624         if (regions[i] != oldRegion)
00625         {
00626             regionToBoundaryPatch.insert(regions[i], boundPatch);
00627 
00628             oldRegion = regions[i];
00629             boundPatch++;
00630         }
00631     }
00632 
00633     const geometricSurfacePatchList& surfPatches = surf.patches();
00634 
00635     patches_.clear();
00636 
00637     if (surfPatches.size() == regionToBoundaryPatch.size())
00638     {
00639         // There are as many surface patches as region numbers in triangles
00640         // so use the surface patches
00641 
00642         patches_.setSize(surfPatches.size());
00643 
00644         // Take over patches, setting size to 0 for now.
00645         forAll(surfPatches, patchI)
00646         {
00647             const geometricSurfacePatch& surfPatch = surfPatches[patchI];
00648 
00649             patches_.set
00650             (
00651                 patchI,
00652                 new boundaryPatch
00653                 (
00654                     surfPatch.name(),
00655                     patchI,
00656                     0,
00657                     0,
00658                     surfPatch.geometricType()
00659                 )
00660             );
00661         }
00662     }
00663     else
00664     {
00665         // There are not enough surface patches. Make up my own.
00666 
00667         patches_.setSize(regionToBoundaryPatch.size());
00668 
00669         forAll(patches_, patchI)
00670         {
00671             patches_.set
00672             (
00673                 patchI,
00674                 new boundaryPatch
00675                 (
00676                     "patch" + name(patchI),
00677                     patchI,
00678                     0,
00679                     0,
00680                     "empty"
00681                 )
00682             );
00683         }
00684     }
00685 
00686     //
00687     // Copy according into bFaces according to regions
00688     //
00689 
00690     const labelList& indices = regions.indices();
00691 
00692     faceList bFaces(surf.size());
00693 
00694     meshFace_.setSize(surf.size());
00695 
00696     label bFaceI = 0;
00697 
00698     // Current region number
00699     label surfRegion = regions[0];
00700     label foamRegion = regionToBoundaryPatch[surfRegion];
00701 
00702     Pout<< "Surface region " << surfRegion << " becomes boundary patch "
00703         << foamRegion << " with name " << patches_[foamRegion].name() << endl;
00704 
00705 
00706     // Index in bFaces of start of current patch
00707     label startFaceI = 0;
00708 
00709     forAll(indices, indexI)
00710     {
00711         label triI = indices[indexI];
00712 
00713         const labelledTri& tri = surf.localFaces()[triI];
00714 
00715         if (tri.region() != surfRegion)
00716         {
00717             // Change of region. We now know the size of the previous one.
00718             boundaryPatch& bp = patches_[foamRegion];
00719 
00720             bp.size() = bFaceI - startFaceI;
00721             bp.start() = startFaceI;
00722 
00723             surfRegion = tri.region();
00724             foamRegion = regionToBoundaryPatch[surfRegion];
00725 
00726             Pout<< "Surface region " << surfRegion << " becomes boundary patch "
00727                 << foamRegion << " with name " << patches_[foamRegion].name()
00728                 << endl;
00729 
00730             startFaceI = bFaceI;
00731         }
00732 
00733         meshFace_[bFaceI] = triI;
00734 
00735         bFaces[bFaceI++] = face(tri);
00736     }
00737 
00738     // Final region
00739     boundaryPatch& bp = patches_[foamRegion];
00740 
00741     bp.size() = bFaceI - startFaceI;
00742     bp.start() = startFaceI;
00743 
00744     //
00745     // Construct single primitivePatch for all of boundary
00746     //
00747 
00748     clearOut();
00749 
00750     // Store compact.
00751     meshPtr_ = new bMesh(bFaces, surf.localPoints());
00752 
00753     // Clear edge storage
00754     featurePoints_.setSize(0);
00755     featureEdges_.setSize(0);
00756 
00757     featureToEdge_.setSize(0);
00758     edgeToFeature_.setSize(meshPtr_->nEdges());
00759     edgeToFeature_ = -1;
00760 
00761     featureSegments_.setSize(0);
00762 
00763     extraEdges_.setSize(0);
00764 }
00765 
00766 
00767 void Foam::boundaryMesh::writeTriSurface(const fileName& fName) const
00768 {
00769     geometricSurfacePatchList surfPatches(patches_.size());
00770 
00771     forAll(patches_, patchI)
00772     {
00773         const boundaryPatch& bp = patches_[patchI];
00774 
00775         surfPatches[patchI] =
00776             geometricSurfacePatch
00777             (
00778                 bp.physicalType(),
00779                 bp.name(),
00780                 patchI
00781             );
00782     }
00783 
00784     //
00785     // Simple triangulation.
00786     //
00787 
00788     // Get number of triangles per face
00789     labelList nTris(mesh().size());
00790 
00791     label totalNTris = getNTris(0, mesh().size(), nTris);
00792 
00793     // Determine per face the starting triangle.
00794     labelList startTri(mesh().size());
00795 
00796     label triI = 0;
00797 
00798     forAll(mesh(), faceI)
00799     {
00800         startTri[faceI] = triI;
00801 
00802         triI += nTris[faceI];
00803     }
00804 
00805     // Triangulate
00806     labelList triVerts(3*totalNTris);
00807 
00808     triangulate(0, mesh().size(), totalNTris, triVerts);
00809 
00810     // Convert to labelledTri
00811 
00812     List<labelledTri> tris(totalNTris);
00813 
00814     triI = 0;
00815 
00816     forAll(patches_, patchI)
00817     {
00818         const boundaryPatch& bp = patches_[patchI];
00819 
00820         forAll(bp, patchFaceI)
00821         {
00822             label faceI = bp.start() + patchFaceI;
00823 
00824             label triVertI = 3*startTri[faceI];
00825 
00826             for (label faceTriI = 0; faceTriI < nTris[faceI]; faceTriI++)
00827             {
00828                 label v0 = triVerts[triVertI++];
00829                 label v1 = triVerts[triVertI++];
00830                 label v2 = triVerts[triVertI++];
00831 
00832                 tris[triI++] = labelledTri(v0, v1, v2, patchI);
00833             }
00834         }
00835     }
00836 
00837     triSurface surf(tris, surfPatches, mesh().points());
00838 
00839     OFstream surfStream(fName);
00840 
00841     surf.write(surfStream);
00842 }
00843 
00844 
00845 // Get index in this (boundaryMesh) of face nearest to each boundary face in
00846 // pMesh.
00847 // Origininally all triangles/faces of boundaryMesh would be bunged into
00848 // one big octree. Problem was that faces on top of each other, differing
00849 // only in sign of normal, could not be found separately. It would always
00850 // find only one. We could detect that it was probably finding the wrong one
00851 // (based on normal) but could not 'tell' the octree to retrieve the other
00852 // one (since they occupy exactly the same space)
00853 // So now faces get put into different octrees depending on normal.
00854 // !It still will not be possible to differentiate between two faces on top
00855 // of each other having the same normal
00856 Foam::labelList Foam::boundaryMesh::getNearest
00857 (
00858     const primitiveMesh& pMesh,
00859     const vector& searchSpan
00860 ) const
00861 {
00862 
00863     // Divide faces into two bins acc. to normal
00864     // - left of splitNormal
00865     // - right ,,
00866     DynamicList<label> leftFaces(mesh().size()/2);
00867     DynamicList<label> rightFaces(mesh().size()/2);
00868 
00869     forAll(mesh(), bFaceI)
00870     {
00871         scalar sign = mesh().faceNormals()[bFaceI] & splitNormal_;
00872 
00873         if (sign > -1E-5)
00874         {
00875             rightFaces.append(bFaceI);
00876         }
00877         if (sign < 1E-5)
00878         {
00879             leftFaces.append(bFaceI);
00880         }
00881     }
00882 
00883     leftFaces.shrink();
00884     rightFaces.shrink();
00885 
00886     if (debug)
00887     {
00888         Pout<< "getNearest :"
00889             << " rightBin:" << rightFaces.size()
00890             << " leftBin:" << leftFaces.size()
00891             << endl;
00892     }
00893 
00894 
00895     // Overall bb
00896     treeBoundBox overallBb(mesh().localPoints());
00897 
00898     // Extend domain slightly (also makes it 3D if was 2D)
00899     // Note asymmetry to avoid having faces align with octree cubes.
00900     scalar tol = 1E-6 * overallBb.avgDim();
00901 
00902     point& bbMin = overallBb.min();
00903     bbMin.x() -= tol;
00904     bbMin.y() -= tol;
00905     bbMin.z() -= tol;
00906 
00907     point& bbMax = overallBb.max();
00908     bbMax.x() += 2*tol;
00909     bbMax.y() += 2*tol;
00910     bbMax.z() += 2*tol;
00911 
00912     // Create the octrees
00913     octree<octreeDataFaceList> leftTree
00914     (
00915         overallBb,
00916         octreeDataFaceList
00917         (
00918             mesh(),
00919             leftFaces
00920         ),
00921         1,
00922         20,
00923         10
00924     );
00925     octree<octreeDataFaceList> rightTree
00926     (
00927         overallBb,
00928         octreeDataFaceList
00929         (
00930             mesh(),
00931             rightFaces
00932         ),
00933         1,
00934         20,
00935         10
00936     );
00937 
00938     if (debug)
00939     {
00940         Pout<< "getNearest : built trees" << endl;
00941     }
00942 
00943 
00944     const vectorField& ns = mesh().faceNormals();
00945 
00946 
00947     //
00948     // Search nearest triangle centre for every polyMesh boundary face
00949     //
00950 
00951     labelList nearestBFaceI(pMesh.nFaces() - pMesh.nInternalFaces());
00952 
00953     treeBoundBox tightest;
00954 
00955     const scalar searchDim = mag(searchSpan);
00956 
00957     forAll(nearestBFaceI, patchFaceI)
00958     {
00959         label meshFaceI = pMesh.nInternalFaces() + patchFaceI;
00960 
00961         const point& ctr = pMesh.faceCentres()[meshFaceI];
00962 
00963         if (debug && (patchFaceI % 1000) == 0)
00964         {
00965             Pout<< "getNearest : patchFace:" << patchFaceI
00966                 << " meshFaceI:" << meshFaceI << " ctr:" << ctr << endl;
00967         }
00968 
00969 
00970         // Get normal from area vector
00971         vector n = pMesh.faceAreas()[meshFaceI];
00972         scalar area = mag(n);
00973         n /= area;
00974 
00975         scalar typDim = -GREAT;
00976         const face& f = pMesh.faces()[meshFaceI];
00977 
00978         forAll(f, fp)
00979         {
00980             typDim = max(typDim, mag(pMesh.points()[f[fp]] - ctr));
00981         }
00982 
00983         // Search right tree
00984         tightest.min() = ctr - searchSpan;
00985         tightest.max() = ctr + searchSpan;
00986         scalar rightDist = searchDim;
00987         label rightI = rightTree.findNearest(ctr, tightest, rightDist);
00988 
00989 
00990         // Search left tree. Note: could start from rightDist bounding box
00991         // instead of starting from top.
00992         tightest.min() = ctr - searchSpan;
00993         tightest.max() = ctr + searchSpan;
00994         scalar leftDist = searchDim;
00995         label leftI = leftTree.findNearest(ctr, tightest, leftDist);
00996 
00997 
00998         if (rightI == -1)
00999         {
01000             // No face found in right tree.
01001 
01002             if (leftI == -1)
01003             {
01004                 // No face found in left tree.
01005                 nearestBFaceI[patchFaceI] = -1;
01006             }
01007             else
01008             {
01009                 // Found in left but not in right. Choose left regardless if
01010                 // correct sign. Note: do we want this?
01011                 nearestBFaceI[patchFaceI] = leftFaces[leftI];
01012             }
01013         }
01014         else
01015         {
01016             if (leftI == -1)
01017             {
01018                 // Found in right but not in left. Choose right regardless if
01019                 // correct sign. Note: do we want this?
01020                 nearestBFaceI[patchFaceI] = rightFaces[rightI];
01021             }
01022             else
01023             {
01024                 // Found in both trees. Compare normals.
01025 
01026                 scalar rightSign = n & ns[rightFaces[rightI]];
01027                 scalar leftSign = n & ns[leftFaces[leftI]];
01028 
01029                 if
01030                 (
01031                     (rightSign > 0 && leftSign > 0)
01032                  || (rightSign < 0 && leftSign < 0)
01033                 )
01034                 {
01035                     // Both same sign. Choose nearest.
01036                     if (rightDist < leftDist)
01037                     {
01038                         nearestBFaceI[patchFaceI] = rightFaces[rightI];
01039                     }
01040                     else
01041                     {
01042                         nearestBFaceI[patchFaceI] = leftFaces[leftI];
01043                     }
01044                 }
01045                 else
01046                 {
01047                     // Differing sign.
01048                     // - if both near enough choose one with correct sign
01049                     // - otherwise choose nearest.
01050 
01051                     // Get local dimension as max of distance between ctr and
01052                     // any face vertex.
01053 
01054                     typDim *= distanceTol_;
01055 
01056                     if (rightDist < typDim && leftDist < typDim)
01057                     {
01058                         // Different sign and nearby. Choosing matching normal
01059                         if (rightSign > 0)
01060                         {
01061                             nearestBFaceI[patchFaceI] = rightFaces[rightI];
01062                         }
01063                         else
01064                         {
01065                             nearestBFaceI[patchFaceI] = leftFaces[leftI];
01066                         }
01067                     }
01068                     else
01069                     {
01070                         // Different sign but faraway. Choosing nearest.
01071                         if (rightDist < leftDist)
01072                         {
01073                             nearestBFaceI[patchFaceI] = rightFaces[rightI];
01074                         }
01075                         else
01076                         {
01077                             nearestBFaceI[patchFaceI] = leftFaces[leftI];
01078                         }
01079                     }
01080                 }
01081             }
01082         }
01083     }
01084 
01085     return nearestBFaceI;
01086 }
01087 
01088 
01089 void Foam::boundaryMesh::patchify
01090 (
01091     const labelList& nearest,
01092     const polyBoundaryMesh& oldPatches,
01093     polyMesh& newMesh
01094 ) const
01095 {
01096 
01097     // 2 cases to be handled:
01098     // A- patches in boundaryMesh patches_
01099     // B- patches not in boundaryMesh patches_ but in polyMesh
01100 
01101     // Create maps from patch name to new patch index.
01102     HashTable<label> nameToIndex(2*patches_.size());
01103 
01104     Map<word> indexToName(2*patches_.size());
01105 
01106 
01107     label nNewPatches = patches_.size();
01108 
01109     forAll(oldPatches, oldPatchI)
01110     {
01111         const polyPatch& patch = oldPatches[oldPatchI];
01112 
01113         label newPatchI = findPatchID(patch.name());
01114 
01115         if (newPatchI != -1)
01116         {
01117             nameToIndex.insert(patch.name(), newPatchI);
01118 
01119             indexToName.insert(newPatchI, patch.name());
01120         }
01121     }
01122 
01123     // Include all boundaryPatches not yet in nameToIndex (i.e. not in old
01124     // patches)
01125     forAll(patches_, bPatchI)
01126     {
01127         const boundaryPatch& bp = patches_[bPatchI];
01128 
01129         if (!nameToIndex.found(bp.name()))
01130         {
01131             nameToIndex.insert(bp.name(), bPatchI);
01132 
01133             indexToName.insert(bPatchI, bp.name());
01134         }
01135     }
01136 
01137     // Pass1:
01138     // Copy names&type of patches (with zero size) from old mesh as far as
01139     // possible. First patch created gets all boundary faces; all others get
01140     // zero faces (repatched later on). Exception is coupled patches which
01141     // keep their size.
01142 
01143     List<polyPatch*> newPatchPtrList(nNewPatches);
01144 
01145     label meshFaceI = newMesh.nInternalFaces();
01146 
01147     // First patch gets all non-coupled faces
01148     label facesToBeDone = newMesh.nFaces() - newMesh.nInternalFaces();
01149 
01150     forAll(patches_, bPatchI)
01151     {
01152         const boundaryPatch& bp = patches_[bPatchI];
01153 
01154         label newPatchI = nameToIndex[bp.name()];
01155 
01156         // Find corresponding patch in polyMesh
01157         label oldPatchI = findPatchID(oldPatches, bp.name());
01158 
01159         if (oldPatchI == -1)
01160         {
01161             // Newly created patch. Gets all or zero faces.
01162             if (debug)
01163             {
01164                 Pout<< "patchify : Creating new polyPatch:" << bp.name()
01165                     << " type:" << bp.physicalType() << endl;
01166             }
01167 
01168             newPatchPtrList[newPatchI] = polyPatch::New
01169             (
01170                 bp.physicalType(),
01171                 bp.name(),
01172                 facesToBeDone,
01173                 meshFaceI,
01174                 newPatchI,
01175                 newMesh.boundaryMesh()
01176             ).ptr();
01177 
01178             meshFaceI += facesToBeDone;
01179 
01180             // first patch gets all boundary faces; all others get 0.
01181             facesToBeDone = 0;
01182         }
01183         else
01184         {
01185             // Existing patch. Gets all or zero faces.
01186             const polyPatch& oldPatch = oldPatches[oldPatchI];
01187 
01188             if (debug)
01189             {
01190                 Pout<< "patchify : Cloning existing polyPatch:"
01191                     << oldPatch.name() << endl;
01192             }
01193 
01194             newPatchPtrList[newPatchI] = oldPatch.clone
01195             (
01196                 newMesh.boundaryMesh(),
01197                 newPatchI,
01198                 facesToBeDone,
01199                 meshFaceI
01200             ).ptr();
01201 
01202             meshFaceI += facesToBeDone;
01203 
01204             // first patch gets all boundary faces; all others get 0.
01205             facesToBeDone = 0;
01206         }
01207     }
01208 
01209 
01210     if (debug)
01211     {
01212         Pout<< "Patchify : new polyPatch list:" << endl;
01213 
01214         forAll(newPatchPtrList, patchI)
01215         {
01216             const polyPatch& newPatch = *newPatchPtrList[patchI];
01217 
01218             if (debug)
01219             {
01220                 Pout<< "polyPatch:" << newPatch.name() << endl
01221                     << "    type :" << newPatch.typeName << endl
01222                     << "    size :" << newPatch.size() << endl
01223                     << "    start:" << newPatch.start() << endl
01224                     << "    index:" << patchI << endl;
01225             }
01226         }
01227     }
01228 
01229     // Actually add new list of patches
01230     repatchPolyTopoChanger polyMeshRepatcher(newMesh);
01231     polyMeshRepatcher.changePatches(newPatchPtrList);
01232 
01233 
01234     // Pass2:
01235     // Change patch type for face
01236 
01237     if (newPatchPtrList.size())
01238     {
01239         List<DynamicList<label> > patchFaces(nNewPatches);
01240 
01241         // Give reasonable estimate for size of patches
01242         label nAvgFaces =
01243             (newMesh.nFaces() - newMesh.nInternalFaces())
01244           / nNewPatches;
01245 
01246         forAll(patchFaces, newPatchI)
01247         {
01248             patchFaces[newPatchI].setCapacity(nAvgFaces);
01249         }
01250 
01251         //
01252         // Sort faces acc. to new patch index. Can loop over all old patches
01253         // since will contain all faces.
01254         //
01255 
01256         forAll(oldPatches, oldPatchI)
01257         {
01258             const polyPatch& patch = oldPatches[oldPatchI];
01259 
01260             forAll(patch, patchFaceI)
01261             {
01262                 // Put face into region given by nearest boundary face
01263 
01264                 label meshFaceI = patch.start() + patchFaceI;
01265 
01266                 label bFaceI = meshFaceI - newMesh.nInternalFaces();
01267 
01268                 patchFaces[whichPatch(nearest[bFaceI])].append(meshFaceI);
01269             }
01270         }
01271 
01272         forAll(patchFaces, newPatchI)
01273         {
01274             patchFaces[newPatchI].shrink();
01275         }
01276 
01277 
01278         // Change patch > 0. (since above we put all faces into the zeroth
01279         // patch)
01280 
01281         for (label newPatchI = 1; newPatchI < patchFaces.size(); newPatchI++)
01282         {
01283             const labelList& pFaces = patchFaces[newPatchI];
01284 
01285             forAll(pFaces, pFaceI)
01286             {
01287                 polyMeshRepatcher.changePatchID(pFaces[pFaceI], newPatchI);
01288             }
01289         }
01290 
01291         polyMeshRepatcher.repatch();
01292     }
01293 }
01294 
01295 
01296 void Foam::boundaryMesh::setFeatureEdges(const scalar minCos)
01297 {
01298     edgeToFeature_.setSize(mesh().nEdges());
01299 
01300     edgeToFeature_ = -1;
01301 
01302     // 1. Mark feature edges
01303 
01304     // Storage for edge labels that are features. Trim later.
01305     featureToEdge_.setSize(mesh().nEdges());
01306 
01307     label featureI = 0;
01308 
01309     if (minCos >= 0.9999)
01310     {
01311         // Select everything
01312         forAll(mesh().edges(), edgeI)
01313         {
01314             edgeToFeature_[edgeI] = featureI;
01315             featureToEdge_[featureI++] = edgeI;
01316         }
01317     }
01318     else
01319     {
01320         forAll(mesh().edges(), edgeI)
01321         {
01322             const labelList& eFaces = mesh().edgeFaces()[edgeI];
01323 
01324             if (eFaces.size() == 2)
01325             {
01326                 label face0I = eFaces[0];
01327 
01328                 label face1I = eFaces[1];
01329 
01332                 //if (whichPatch(face0I) != whichPatch(face1I))
01333                 //{
01334                 //    edgeToFeature_[edgeI] = featureI;
01335                 //    featureToEdge_[featureI++] = edgeI;
01336                 //}
01337                 //else
01338                 {
01339                     const vector& n0 = mesh().faceNormals()[face0I];
01340 
01341                     const vector& n1 = mesh().faceNormals()[face1I];
01342 
01343                     float cosAng = n0 & n1;
01344 
01345                     if (cosAng < minCos)
01346                     {
01347                         edgeToFeature_[edgeI] = featureI;
01348                         featureToEdge_[featureI++] = edgeI;
01349                     }
01350                 }
01351             }
01352             else
01353             {
01354                 //Should not occur: 0 or more than two faces
01355                 edgeToFeature_[edgeI] = featureI;
01356                 featureToEdge_[featureI++] = edgeI;
01357             }
01358         }
01359     }
01360 
01361     // Trim featureToEdge_ to actual number of edges.
01362     featureToEdge_.setSize(featureI);
01363 
01364     //
01365     // Compact edges i.e. relabel vertices.
01366     //
01367 
01368     featureEdges_.setSize(featureI);
01369     featurePoints_.setSize(mesh().nPoints());
01370 
01371     labelList featToMeshPoint(mesh().nPoints(), -1);
01372 
01373     label featPtI = 0;
01374 
01375     forAll(featureToEdge_, fEdgeI)
01376     {
01377         label edgeI = featureToEdge_[fEdgeI];
01378 
01379         const edge& e = mesh().edges()[edgeI];
01380 
01381         label start = featToMeshPoint[e.start()];
01382 
01383         if (start == -1)
01384         {
01385             featToMeshPoint[e.start()] = featPtI;
01386 
01387             featurePoints_[featPtI] = mesh().points()[e.start()];
01388 
01389             start = featPtI;
01390 
01391             featPtI++;
01392         }
01393 
01394         label end = featToMeshPoint[e.end()];
01395 
01396         if (end == -1)
01397         {
01398             featToMeshPoint[e.end()] = featPtI;
01399 
01400             featurePoints_[featPtI] = mesh().points()[e.end()];
01401 
01402             end = featPtI;
01403 
01404             featPtI++;
01405         }
01406 
01407         // Store with renumbered vertices.
01408         featureEdges_[fEdgeI] = edge(start, end);
01409     }
01410 
01411     // Compact points
01412     featurePoints_.setSize(featPtI);
01413 
01414 
01415     //
01416     // 2. Mark endpoints of feature segments. These are points with
01417     // != 2 feature edges connected.
01418     // Note: can add geometric constraint here as well that if 2 feature
01419     // edges the angle between them should be less than xxx.
01420     //
01421 
01422     boolList isFeaturePoint(mesh().nPoints(), false);
01423 
01424     forAll(featureToEdge_, featI)
01425     {
01426         label edgeI = featureToEdge_[featI];
01427 
01428         const edge& e = mesh().edges()[edgeI];
01429 
01430         if (nFeatureEdges(e.start()) != 2)
01431         {
01432             isFeaturePoint[e.start()] = true;
01433         }
01434 
01435         if (nFeatureEdges(e.end()) != 2)
01436         {
01437             isFeaturePoint[e.end()] = true;
01438         }
01439     }
01440 
01441 
01442     //
01443     // 3: Split feature edges into segments:
01444     // find point with not 2 feature edges -> start of feature segment
01445     //
01446 
01447     DynamicList<labelList> segments;
01448 
01449 
01450     boolList featVisited(featureToEdge_.size(), false);
01451 
01452     do
01453     {
01454         label startFeatI = -1;
01455 
01456         forAll(featVisited, featI)
01457         {
01458             if (!featVisited[featI])
01459             {
01460                 startFeatI = featI;
01461 
01462                 break;
01463             }
01464         }
01465 
01466         if (startFeatI == -1)
01467         {
01468             // No feature lines left.
01469             break;
01470         }
01471 
01472         segments.append
01473         (
01474             collectSegment
01475             (
01476                 isFeaturePoint,
01477                 featureToEdge_[startFeatI],
01478                 featVisited
01479             )
01480         );
01481     }
01482     while (true);
01483 
01484 
01485     //
01486     // Store in *this
01487     //
01488     featureSegments_.setSize(segments.size());
01489 
01490     forAll(featureSegments_, segmentI)
01491     {
01492         featureSegments_[segmentI] = segments[segmentI];
01493     }
01494 }
01495 
01496 
01497 void Foam::boundaryMesh::setExtraEdges(const label edgeI)
01498 {
01499     labelList minDistance(mesh().nEdges(), -1);
01500 
01501     // All edge labels encountered
01502     DynamicList<label> visitedEdges;
01503 
01504     // Floodfill from edgeI starting from distance 0. Stop at distance.
01505     markEdges(8, edgeI, 0, minDistance, visitedEdges);
01506 
01507     // Set edge labels to display
01508     extraEdges_.transfer(visitedEdges);
01509 }
01510 
01511 
01512 Foam::label Foam::boundaryMesh::whichPatch(const label faceI) const
01513 {
01514     forAll(patches_, patchI)
01515     {
01516         const boundaryPatch& bp = patches_[patchI];
01517 
01518         if ((faceI >= bp.start()) && (faceI < (bp.start() + bp.size())))
01519         {
01520             return patchI;
01521         }
01522     }
01523 
01524     FatalErrorIn("boundaryMesh::whichPatch(const label)")
01525         << "Cannot find face " << faceI << " in list of boundaryPatches "
01526         << patches_
01527         << abort(FatalError);
01528 
01529     return -1;
01530 }
01531 
01532 
01533 Foam::label Foam::boundaryMesh::findPatchID(const word& patchName) const
01534 {
01535     forAll(patches_, patchI)
01536     {
01537         if (patches_[patchI].name() == patchName)
01538         {
01539             return patchI;
01540         }
01541     }
01542 
01543     return -1;
01544 }
01545 
01546 
01547 void Foam::boundaryMesh::addPatch(const word& patchName)
01548 {
01549     patches_.setSize(patches_.size() + 1);
01550 
01551     // Add empty patch at end of patch list.
01552 
01553     label patchI = patches_.size()-1;
01554 
01555     boundaryPatch* bpPtr = new boundaryPatch
01556     (
01557         patchName,
01558         patchI,
01559         0,
01560         mesh().size(),
01561         "empty"
01562     );
01563 
01564     patches_.set(patchI, bpPtr);
01565 
01566     if (debug)
01567     {
01568         Pout<< "addPatch : patches now:" << endl;
01569 
01570         forAll(patches_, patchI)
01571         {
01572             const boundaryPatch& bp = patches_[patchI];
01573 
01574             Pout<< "    name  : " << bp.name() << endl
01575                 << "    size  : " << bp.size() << endl
01576                 << "    start : " << bp.start() << endl
01577                 << "    type  : " << bp.physicalType() << endl
01578                 << endl;
01579         }
01580     }
01581 }
01582 
01583 
01584 void Foam::boundaryMesh::deletePatch(const word& patchName)
01585 {
01586     label delPatchI = findPatchID(patchName);
01587 
01588     if (delPatchI == -1)
01589     {
01590         FatalErrorIn("boundaryMesh::deletePatch(const word&")
01591             << "Can't find patch named " << patchName
01592             << abort(FatalError);
01593     }
01594 
01595     if (patches_[delPatchI].size())
01596     {
01597         FatalErrorIn("boundaryMesh::deletePatch(const word&")
01598             << "Trying to delete non-empty patch " << patchName
01599             << endl << "Current size:" << patches_[delPatchI].size()
01600             << abort(FatalError);
01601     }
01602 
01603     PtrList<boundaryPatch> newPatches(patches_.size() - 1);
01604 
01605     for (label patchI = 0; patchI < delPatchI; patchI++)
01606     {
01607         newPatches.set(patchI, patches_[patchI].clone());
01608     }
01609 
01610     // Move patches down, starting from delPatchI.
01611 
01612     for (label patchI = delPatchI + 1; patchI < patches_.size(); patchI++)
01613     {
01614         newPatches.set(patchI - 1, patches_[patchI].clone());
01615     }
01616 
01617     patches_.clear();
01618 
01619     patches_ = newPatches;
01620 
01621     if (debug)
01622     {
01623         Pout<< "deletePatch : patches now:" << endl;
01624 
01625         forAll(patches_, patchI)
01626         {
01627             const boundaryPatch& bp = patches_[patchI];
01628 
01629             Pout<< "    name  : " << bp.name() << endl
01630                 << "    size  : " << bp.size() << endl
01631                 << "    start : " << bp.start() << endl
01632                 << "    type  : " << bp.physicalType() << endl
01633                 << endl;
01634         }
01635     }
01636 }
01637 
01638 
01639 void Foam::boundaryMesh::changePatchType
01640 (
01641     const word& patchName,
01642     const word& patchType
01643 )
01644 {
01645     label changeI = findPatchID(patchName);
01646 
01647     if (changeI == -1)
01648     {
01649         FatalErrorIn("boundaryMesh::changePatchType(const word&, const word&)")
01650             << "Can't find patch named " << patchName
01651             << abort(FatalError);
01652     }
01653 
01654 
01655     // Cause we can't reassign to individual PtrList elems ;-(
01656     // work on copy
01657 
01658 
01659     PtrList<boundaryPatch> newPatches(patches_.size());
01660 
01661     forAll(patches_, patchI)
01662     {
01663         if (patchI == changeI)
01664         {
01665             // Create copy but for type
01666             const boundaryPatch& oldBp = patches_[patchI];
01667 
01668             boundaryPatch* bpPtr = new boundaryPatch
01669             (
01670                 oldBp.name(),
01671                 oldBp.index(),
01672                 oldBp.size(),
01673                 oldBp.start(),
01674                 patchType
01675             );
01676 
01677             newPatches.set(patchI, bpPtr);
01678         }
01679         else
01680         {
01681             // Create copy
01682             newPatches.set(patchI, patches_[patchI].clone());
01683         }
01684     }
01685 
01686     patches_ = newPatches;
01687 }
01688 
01689 
01690 void Foam::boundaryMesh::changeFaces
01691 (
01692     const labelList& patchIDs,
01693     labelList& oldToNew
01694 )
01695 {
01696     if (patchIDs.size() != mesh().size())
01697     {
01698         FatalErrorIn("boundaryMesh::changeFaces(const labelList& patchIDs)")
01699             << "List of patchIDs not equal to number of faces." << endl
01700             << "PatchIDs size:" << patchIDs.size()
01701             << " nFaces::" << mesh().size()
01702             << abort(FatalError);
01703     }
01704 
01705     // Count number of faces for each patch
01706 
01707     labelList nFaces(patches_.size(), 0);
01708 
01709     forAll(patchIDs, faceI)
01710     {
01711         label patchID = patchIDs[faceI];
01712 
01713         if (patchID < 0 || patchID >= patches_.size())
01714         {
01715             FatalErrorIn("boundaryMesh::changeFaces(const labelList&)")
01716                 << "PatchID " << patchID << " out of range"
01717                 << abort(FatalError);
01718         }
01719         nFaces[patchID]++;
01720     }
01721 
01722 
01723     // Determine position in faces_ for each patch
01724 
01725     labelList startFace(patches_.size());
01726 
01727     startFace[0] = 0;
01728 
01729     for (label patchI = 1; patchI < patches_.size(); patchI++)
01730     {
01731         startFace[patchI] = startFace[patchI-1] + nFaces[patchI-1];
01732     }
01733 
01734     // Update patch info
01735     PtrList<boundaryPatch> newPatches(patches_.size());
01736 
01737     forAll(patches_, patchI)
01738     {
01739         const boundaryPatch& bp = patches_[patchI];
01740 
01741         newPatches.set
01742         (
01743             patchI,
01744             new boundaryPatch
01745             (
01746                 bp.name(),
01747                 patchI,
01748                 nFaces[patchI],
01749                 startFace[patchI],
01750                 bp.physicalType()
01751             )
01752         );
01753     }
01754     patches_ = newPatches;
01755 
01756     if (debug)
01757     {
01758         Pout<< "changeFaces : patches now:" << endl;
01759 
01760         forAll(patches_, patchI)
01761         {
01762             const boundaryPatch& bp = patches_[patchI];
01763 
01764             Pout<< "    name  : " << bp.name() << endl
01765                 << "    size  : " << bp.size() << endl
01766                 << "    start : " << bp.start() << endl
01767                 << "    type  : " << bp.physicalType() << endl
01768                 << endl;
01769         }
01770     }
01771 
01772 
01773     // Construct face mapping array
01774     oldToNew.setSize(patchIDs.size());
01775 
01776     forAll(patchIDs, faceI)
01777     {
01778         int patchID = patchIDs[faceI];
01779 
01780         oldToNew[faceI] = startFace[patchID]++;
01781     }
01782 
01783     // Copy faces into correct position and maintain label of original face
01784     faceList newFaces(mesh().size());
01785 
01786     labelList newMeshFace(mesh().size());
01787 
01788     forAll(oldToNew, faceI)
01789     {
01790         newFaces[oldToNew[faceI]] = mesh()[faceI];
01791         newMeshFace[oldToNew[faceI]] = meshFace_[faceI];
01792     }
01793 
01794     // Reconstruct 'mesh' from new faces and (copy of) existing points.
01795     bMesh* newMeshPtr_ = new bMesh(newFaces, mesh().points());
01796 
01797     // Reset meshFace_ to new ordering.
01798     meshFace_.transfer(newMeshFace);
01799 
01800 
01801     // Remove old PrimitivePatch on meshPtr_.
01802     clearOut();
01803 
01804     // And insert new 'mesh'.
01805     meshPtr_ = newMeshPtr_;
01806 }
01807 
01808 
01809 Foam::label Foam::boundaryMesh::getNTris(const label faceI) const
01810 {
01811     const face& f = mesh()[faceI];
01812 
01813     return f.nTriangles(mesh().points());
01814 }
01815 
01816 
01817 Foam::label Foam::boundaryMesh::getNTris
01818 (
01819     const label startFaceI,
01820     const label nFaces,
01821     labelList& nTris
01822 ) const
01823 {
01824     label totalNTris = 0;
01825 
01826     nTris.setSize(nFaces);
01827 
01828     for (label i = 0; i < nFaces; i++)
01829     {
01830         label faceNTris = getNTris(startFaceI + i);
01831 
01832         nTris[i] = faceNTris;
01833 
01834         totalNTris += faceNTris;
01835     }
01836     return totalNTris;
01837 }
01838 
01839 
01840 // Simple triangulation of face subset. Stores vertices in tris[] as three
01841 // consecutive vertices per triangle.
01842 void Foam::boundaryMesh::triangulate
01843 (
01844     const label startFaceI,
01845     const label nFaces,
01846     const label totalNTris,
01847     labelList& triVerts
01848 ) const
01849 {
01850     // Triangulate faces.
01851     triVerts.setSize(3*totalNTris);
01852 
01853     label vertI = 0;
01854 
01855     for (label i = 0; i < nFaces; i++)
01856     {
01857         label faceI = startFaceI + i;
01858 
01859         const face& f = mesh()[faceI];
01860 
01861         // Have face triangulate itself (results in faceList)
01862         faceList triFaces(f.nTriangles(mesh().points()));
01863 
01864         label nTri = 0;
01865 
01866         f.triangles(mesh().points(), nTri, triFaces);
01867 
01868         // Copy into triVerts
01869 
01870         forAll(triFaces, triFaceI)
01871         {
01872             const face& triF = triFaces[triFaceI];
01873 
01874             triVerts[vertI++] = triF[0];
01875             triVerts[vertI++] = triF[1];
01876             triVerts[vertI++] = triF[2];
01877         }
01878     }
01879 }
01880 
01881 
01882 // Number of local points in subset.
01883 Foam::label Foam::boundaryMesh::getNPoints
01884 (
01885     const label startFaceI,
01886     const label nFaces
01887 ) const
01888 {
01889     SubList<face> patchFaces(mesh(), nFaces, startFaceI);
01890 
01891     primitivePatch patch(patchFaces, mesh().points());
01892 
01893     return patch.nPoints();
01894 }
01895 
01896 
01897 // Triangulation of face subset in local coords.
01898 void Foam::boundaryMesh::triangulateLocal
01899 (
01900     const label startFaceI,
01901     const label nFaces,
01902     const label totalNTris,
01903     labelList& triVerts,
01904     labelList& localToGlobal
01905 ) const
01906 {
01907     SubList<face> patchFaces(mesh(), nFaces, startFaceI);
01908 
01909     primitivePatch patch(patchFaces, mesh().points());
01910 
01911     // Map from local to mesh().points() addressing
01912     localToGlobal = patch.meshPoints();
01913 
01914     // Triangulate patch faces.
01915     triVerts.setSize(3*totalNTris);
01916 
01917     label vertI = 0;
01918 
01919     for (label i = 0; i < nFaces; i++)
01920     {
01921         // Local face
01922         const face& f = patch.localFaces()[i];
01923 
01924         // Have face triangulate itself (results in faceList)
01925         faceList triFaces(f.nTriangles(patch.localPoints()));
01926 
01927         label nTri = 0;
01928 
01929         f.triangles(patch.localPoints(), nTri, triFaces);
01930 
01931         // Copy into triVerts
01932 
01933         forAll(triFaces, triFaceI)
01934         {
01935             const face& triF = triFaces[triFaceI];
01936 
01937             triVerts[vertI++] = triF[0];
01938             triVerts[vertI++] = triF[1];
01939             triVerts[vertI++] = triF[2];
01940         }
01941     }
01942 }
01943 
01944 
01945 void Foam::boundaryMesh::markFaces
01946 (
01947     const labelList& protectedEdges,
01948     const label seedFaceI,
01949     boolList& visited
01950 ) const
01951 {
01952     boolList protectedEdge(mesh().nEdges(), false);
01953 
01954     forAll(protectedEdges, i)
01955     {
01956         protectedEdge[protectedEdges[i]] = true;
01957     }
01958 
01959 
01960     // Initialize zone for all faces to -1
01961     labelList currentZone(mesh().size(), -1);
01962 
01963     // Mark with 0 all faces reachable from seedFaceI
01964     markZone(protectedEdge, seedFaceI, 0, currentZone);
01965 
01966     // Set in visited all reached ones.
01967     visited.setSize(mesh().size());
01968 
01969     forAll(currentZone, faceI)
01970     {
01971         if (currentZone[faceI] == 0)
01972         {
01973             visited[faceI] = true;
01974         }
01975         else
01976         {
01977             visited[faceI] = false;
01978         }
01979     }
01980 }
01981 
01982 
01983 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines