/*---------------------------------------------------------------------------*\ ========= | \\ / F ield | foam-extend: Open Source CFD \\ / O peration | Version: 3.2 \\ / 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 "Time.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, DynamicList