/*---------------------------------------------------------------------------*\ ========= | \\ / F ield | foam-extend: Open Source CFD \\ / O peration | Version: 4.1 \\ / A nd | Web: http://www.foam-extend.org \\/ M anipulation | For copyright notice see file Copyright ------------------------------------------------------------------------------- License This file is part of foam-extend. foam-extend is free software: you can redistribute it and/or modify it under the terms of the GNU General Public License as published by the Free Software Foundation, either version 3 of the License, or (at your option) any later version. foam-extend is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details. You should have received a copy of the GNU General Public License along with foam-extend. If not, see . \*---------------------------------------------------------------------------*/ #include "boundaryMesh.H" #include "objectRegistry.H" #include "foamTime.H" #include "polyMesh.H" #include "repatchPolyTopoChanger.H" #include "faceList.H" #include "indexedOctree.H" #include "treeDataPrimitivePatch.H" #include "triSurface.H" #include "SortableList.H" #include "OFstream.H" #include "uindirectPrimitivePatch.H" // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * // defineTypeNameAndDebug(Foam::boundaryMesh, 0); // Normal along which to divide faces into categories (used in getNearest) const Foam::vector Foam::boundaryMesh::splitNormal_(3, 2, 1); // Distance to face tolerance for getNearest const Foam::scalar Foam::boundaryMesh::distanceTol_ = 1E-2; // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * // // Returns number of feature edges connected to pointI Foam::label Foam::boundaryMesh::nFeatureEdges(label pointI) const { label nFeats = 0; const labelList& pEdges = mesh().pointEdges()[pointI]; forAll(pEdges, pEdgeI) { label edgeI = pEdges[pEdgeI]; if (edgeToFeature_[edgeI] != -1) { nFeats++; } } return nFeats; } // Returns next feature edge connected to pointI Foam::label Foam::boundaryMesh::nextFeatureEdge ( const label edgeI, const label vertI ) const { const labelList& pEdges = mesh().pointEdges()[vertI]; forAll(pEdges, pEdgeI) { label nbrEdgeI = pEdges[pEdgeI]; if (nbrEdgeI != edgeI) { label featI = edgeToFeature_[nbrEdgeI]; if (featI != -1) { return nbrEdgeI; } } } return -1; } // Finds connected feature edges, starting from startPointI and returns // feature labels (not edge labels). Marks feature edges handled in // featVisited. Foam::labelList Foam::boundaryMesh::collectSegment ( const boolList& isFeaturePoint, const label startEdgeI, boolList& featVisited ) const { // Find starting feature point on edge. label edgeI = startEdgeI; const edge& e = mesh().edges()[edgeI]; label vertI = e.start(); while (!isFeaturePoint[vertI]) { // Step to next feature edge edgeI = nextFeatureEdge(edgeI, vertI); if ((edgeI == -1) || (edgeI == startEdgeI)) { break; } // Step to next vertex on edge const edge& e = mesh().edges()[edgeI]; vertI = e.otherVertex(vertI); } // // Now we have: // edgeI : first edge on this segment // vertI : one of the endpoints of this segment // // Start walking other way and storing edges as we go along. // // Untrimmed storage for current segment labelList featLabels(featureEdges_.size()); label featLabelI = 0; label initEdgeI = edgeI; do { // Mark edge as visited label featI = edgeToFeature_[edgeI]; if (featI == -1) { FatalErrorIn("boundaryMesh::collectSegment") << "Problem" << abort(FatalError); } featLabels[featLabelI++] = featI; featVisited[featI] = true; // Step to next vertex on edge const edge& e = mesh().edges()[edgeI]; vertI = e.otherVertex(vertI); // Step to next feature edge edgeI = nextFeatureEdge(edgeI, vertI); if ((edgeI == -1) || (edgeI == initEdgeI)) { break; } } while (!isFeaturePoint[vertI]); // Trim to size featLabels.setSize(featLabelI); return featLabels; } void Foam::boundaryMesh::markEdges ( const label maxDistance, const label edgeI, const label distance, labelList& minDistance, dynamicLabelList& visited ) const { if (distance < maxDistance) { // Don't do anything if reached beyond maxDistance. if (minDistance[edgeI] == -1) { // First visit of edge. Store edge label. visited.append(edgeI); } else if (minDistance[edgeI] <= distance) { // Already done this edge return; } minDistance[edgeI] = distance; const edge& e = mesh().edges()[edgeI]; // Do edges connected to e.start const labelList& startEdges = mesh().pointEdges()[e.start()]; forAll(startEdges, pEdgeI) { markEdges ( maxDistance, startEdges[pEdgeI], distance+1, minDistance, visited ); } // Do edges connected to e.end const labelList& endEdges = mesh().pointEdges()[e.end()]; forAll(endEdges, pEdgeI) { markEdges ( maxDistance, endEdges[pEdgeI], distance+1, minDistance, visited ); } } } Foam::label Foam::boundaryMesh::findPatchID ( const polyPatchList& patches, const word& patchName ) const { forAll(patches, patchI) { if (patches[patchI].name() == patchName) { return patchI; } } return -1; } Foam::wordList Foam::boundaryMesh::patchNames() const { wordList names(patches_.size()); forAll(patches_, patchI) { names[patchI] = patches_[patchI].name(); } return names; } Foam::label Foam::boundaryMesh::whichPatch ( const polyPatchList& patches, const label faceI ) const { forAll(patches, patchI) { const polyPatch& pp = patches[patchI]; if ((faceI >= pp.start()) && (faceI < (pp.start() + pp.size()))) { return patchI; } } return -1; } // Gets labels of changed faces and propagates them to the edges. Returns // labels of edges changed. Foam::labelList Foam::boundaryMesh::faceToEdge ( const boolList& regionEdge, const label region, const labelList& changedFaces, labelList& edgeRegion ) const { labelList changedEdges(mesh().nEdges(), -1); label changedI = 0; forAll(changedFaces, i) { label faceI = changedFaces[i]; const labelList& fEdges = mesh().faceEdges()[faceI]; forAll(fEdges, fEdgeI) { label edgeI = fEdges[fEdgeI]; if (!regionEdge[edgeI] && (edgeRegion[edgeI] == -1)) { edgeRegion[edgeI] = region; changedEdges[changedI++] = edgeI; } } } changedEdges.setSize(changedI); return changedEdges; } // Reverse of faceToEdge: gets edges and returns faces Foam::labelList Foam::boundaryMesh::edgeToFace ( const label region, const labelList& changedEdges, labelList& faceRegion ) const { labelList changedFaces(mesh().size(), -1); label changedI = 0; forAll(changedEdges, i) { label edgeI = changedEdges[i]; const labelList& eFaces = mesh().edgeFaces()[edgeI]; forAll(eFaces, eFaceI) { label faceI = eFaces[eFaceI]; if (faceRegion[faceI] == -1) { faceRegion[faceI] = region; changedFaces[changedI++] = faceI; } } } changedFaces.setSize(changedI); return changedFaces; } // Finds area, starting at faceI, delimited by borderEdge void Foam::boundaryMesh::markZone ( const boolList& borderEdge, label faceI, label currentZone, labelList& faceZone ) const { faceZone[faceI] = currentZone; // List of faces whose faceZone has been set. labelList changedFaces(1, faceI); // List of edges whose faceZone has been set. labelList changedEdges; // Zones on all edges. labelList edgeZone(mesh().nEdges(), -1); while(1) { changedEdges = faceToEdge ( borderEdge, currentZone, changedFaces, edgeZone ); if (debug) { Pout<< "From changedFaces:" << changedFaces.size() << " to changedEdges:" << changedEdges.size() << endl; } if (changedEdges.empty()) { break; } changedFaces = edgeToFace(currentZone, changedEdges, faceZone); if (debug) { Pout<< "From changedEdges:" << changedEdges.size() << " to changedFaces:" << changedFaces.size() << endl; } if (changedFaces.empty()) { break; } } } // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * // // Null constructor Foam::boundaryMesh::boundaryMesh() : meshPtr_(nullptr), patches_(), meshFace_(), featurePoints_(), featureEdges_(), featureToEdge_(), edgeToFeature_(), featureSegments_(), extraEdges_() {} // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * // Foam::boundaryMesh::~boundaryMesh() { clearOut(); } void Foam::boundaryMesh::clearOut() { if (meshPtr_) { delete meshPtr_; meshPtr_ = nullptr; } } // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * // void Foam::boundaryMesh::read(const polyMesh& mesh) { patches_.clear(); patches_.setSize(mesh.boundaryMesh().size()); // Number of boundary faces label nBFaces = mesh.nFaces() - mesh.nInternalFaces(); faceList bFaces(nBFaces); meshFace_.setSize(nBFaces); label bFaceI = 0; // Collect all boundary faces. forAll(mesh.boundaryMesh(), patchI) { const polyPatch& pp = mesh.boundaryMesh()[patchI]; patches_.set ( patchI, new boundaryPatch ( pp.name(), patchI, pp.size(), bFaceI, pp.type() ) ); // Collect all faces in global numbering. forAll(pp, patchFaceI) { meshFace_[bFaceI] = pp.start() + patchFaceI; bFaces[bFaceI] = pp[patchFaceI]; bFaceI++; } } if (debug) { Pout<< "read : patches now:" << endl; forAll(patches_, patchI) { const boundaryPatch& bp = patches_[patchI]; Pout<< " name : " << bp.name() << endl << " size : " << bp.size() << endl << " start : " << bp.start() << endl << " type : " << bp.physicalType() << endl << endl; } } // // Construct single patch for all of boundary // // Temporary primitivePatch to calculate compact points & faces. PrimitivePatch globalPatch ( bFaces, mesh.points() ); // Store in local(compact) addressing clearOut(); meshPtr_ = new bMesh(globalPatch.localFaces(), globalPatch.localPoints()); if (debug & 2) { const bMesh& msh = *meshPtr_; Pout<< "** Start of Faces **" << endl; forAll(msh, faceI) { const face& f = msh[faceI]; point ctr(vector::zero); forAll(f, fp) { ctr += msh.points()[f[fp]]; } ctr /= f.size(); Pout<< " " << faceI << " ctr:" << ctr << " verts:" << f << endl; } Pout<< "** End of Faces **" << endl; Pout<< "** Start of Points **" << endl; forAll(msh.points(), pointI) { Pout<< " " << pointI << " coord:" << msh.points()[pointI] << endl; } Pout<< "** End of Points **" << endl; } // Clear edge storage featurePoints_.setSize(0); featureEdges_.setSize(0); featureToEdge_.setSize(0); edgeToFeature_.setSize(meshPtr_->nEdges()); edgeToFeature_ = -1; featureSegments_.setSize(0); extraEdges_.setSize(0); } void Foam::boundaryMesh::readTriSurface(const fileName& fName) { triSurface surf(fName); if (surf.empty()) { return; } // Sort according to region SortableList