/*---------------------------------------------------------------------------*\ ========= | \\ / F ield | OpenFOAM: The Open Source CFD Toolbox \\ / O peration | \\ / A nd | www.openfoam.com \\/ M anipulation | ------------------------------------------------------------------------------- Copyright (C) 2012-2017 OpenFOAM Foundation Copyright (C) 2019 OpenCFD Ltd. ------------------------------------------------------------------------------- License This file is part of OpenFOAM. OpenFOAM 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. OpenFOAM 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 OpenFOAM. If not, see . \*---------------------------------------------------------------------------*/ #include "DelaunayMesh.H" #include "fvMesh.H" #include "pointConversion.H" #include "wallPolyPatch.H" #include "processorPolyPatch.H" #include "labelIOField.H" // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * // template void Foam::DelaunayMesh::sortFaces ( faceList& faces, labelList& owner, labelList& neighbour ) const { // Upper triangular order: // + owner is sorted in ascending cell order // + within each block of equal value for owner, neighbour is sorted in // ascending cell order. // + faces sorted to correspond // e.g. // owner | neighbour // 0 | 2 // 0 | 23 // 0 | 71 // 1 | 23 // 1 | 24 // 1 | 91 List ownerNeighbourPair(owner.size()); forAll(ownerNeighbourPair, oNI) { ownerNeighbourPair[oNI] = labelPair(owner[oNI], neighbour[oNI]); } Info<< nl << "Sorting faces, owner and neighbour into upper triangular order" << endl; labelList oldToNew(sortedOrder(ownerNeighbourPair)); oldToNew = invert(oldToNew.size(), oldToNew); inplaceReorder(oldToNew, faces); inplaceReorder(oldToNew, owner); inplaceReorder(oldToNew, neighbour); } template void Foam::DelaunayMesh::addPatches ( const label nInternalFaces, faceList& faces, labelList& owner, PtrList& patchDicts, const List>& patchFaces, const List>& patchOwners ) const { label nPatches = patchFaces.size(); patchDicts.setSize(nPatches); forAll(patchDicts, patchi) { patchDicts.set(patchi, new dictionary()); } label nBoundaryFaces = 0; forAll(patchFaces, p) { patchDicts[p].set("nFaces", patchFaces[p].size()); patchDicts[p].set("startFace", nInternalFaces + nBoundaryFaces); nBoundaryFaces += patchFaces[p].size(); } faces.setSize(nInternalFaces + nBoundaryFaces); owner.setSize(nInternalFaces + nBoundaryFaces); label facei = nInternalFaces; forAll(patchFaces, p) { forAll(patchFaces[p], f) { faces[facei] = patchFaces[p][f]; owner[facei] = patchOwners[p][f]; facei++; } } } // * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * * // template void Foam::DelaunayMesh::printInfo(Ostream& os) const { PrintTable triInfoTable("Mesh Statistics"); triInfoTable.add("Points", Triangulation::number_of_vertices()); triInfoTable.add("Edges", Triangulation::number_of_finite_edges()); triInfoTable.add("Faces", Triangulation::number_of_finite_facets()); triInfoTable.add("Cells", Triangulation::number_of_finite_cells()); scalar minSize = GREAT; scalar maxSize = 0; for ( Finite_vertices_iterator vit = Triangulation::finite_vertices_begin(); vit != Triangulation::finite_vertices_end(); ++vit ) { // Only internal or boundary vertices have a size if (vit->internalOrBoundaryPoint()) { minSize = min(vit->targetCellSize(), minSize); maxSize = max(vit->targetCellSize(), maxSize); } } Info<< incrIndent; triInfoTable.print(Info, true, true); Info<< "Size (Min/Max) = " << returnReduce(minSize, minOp()) << " " << returnReduce(maxSize, maxOp()) << endl; Info<< decrIndent; } template void Foam::DelaunayMesh::printVertexInfo(Ostream& os) const { label nInternal = 0; //label nInternalRef = 0; label nUnassigned = 0; //label nUnassignedRef = 0; label nInternalNearBoundary = 0; //label nInternalNearBoundaryRef = 0; label nInternalSurface = 0; //label nInternalSurfaceRef = 0; label nInternalFeatureEdge = 0; //label nInternalFeatureEdgeRef = 0; label nInternalFeaturePoint = 0; //label nInternalFeaturePointRef = 0; label nExternalSurface = 0; //label nExternalSurfaceRef = 0; label nExternalFeatureEdge = 0; //label nExternalFeatureEdgeRef = 0; label nExternalFeaturePoint = 0; //label nExternalFeaturePointRef = 0; label nFar = 0; label nReferred = 0; for ( Finite_vertices_iterator vit = Triangulation::finite_vertices_begin(); vit != Triangulation::finite_vertices_end(); ++vit ) { if (vit->type() == Vb::vtInternal) { if (vit->referred()) { nReferred++; //++nInternalRef; } nInternal++; } else if (vit->type() == Vb::vtUnassigned) { if (vit->referred()) { nReferred++; //++nUnassignedRef; } nUnassigned++; } else if (vit->type() == Vb::vtInternalNearBoundary) { if (vit->referred()) { nReferred++; //++nInternalNearBoundaryRef; } nInternalNearBoundary++; } else if (vit->type() == Vb::vtInternalSurface) { if (vit->referred()) { nReferred++; //++nInternalSurfaceRef; } nInternalSurface++; } else if (vit->type() == Vb::vtInternalFeatureEdge) { if (vit->referred()) { nReferred++; //++nInternalFeatureEdgeRef; } nInternalFeatureEdge++; } else if (vit->type() == Vb::vtInternalFeaturePoint) { if (vit->referred()) { nReferred++; //++nInternalFeaturePointRef; } nInternalFeaturePoint++; } else if (vit->type() == Vb::vtExternalSurface) { if (vit->referred()) { nReferred++; //++nExternalSurfaceRef; } nExternalSurface++; } else if (vit->type() == Vb::vtExternalFeatureEdge) { if (vit->referred()) { nReferred++; //++nExternalFeatureEdgeRef; } nExternalFeatureEdge++; } else if (vit->type() == Vb::vtExternalFeaturePoint) { if (vit->referred()) { nReferred++; //++nExternalFeaturePointRef; } nExternalFeaturePoint++; } else if (vit->type() == Vb::vtFar) { nFar++; } } label nTotalVertices = nUnassigned + nInternal + nInternalNearBoundary + nInternalSurface + nInternalFeatureEdge + nInternalFeaturePoint + nExternalSurface + nExternalFeatureEdge + nExternalFeaturePoint + nFar; if (nTotalVertices != label(Triangulation::number_of_vertices())) { WarningInFunction << nTotalVertices << " does not equal " << Triangulation::number_of_vertices() << endl; } PrintTable vertexTable("Vertex Type Information"); vertexTable.add("Total", nTotalVertices); vertexTable.add("Unassigned", nUnassigned); vertexTable.add("nInternal", nInternal); vertexTable.add("nInternalNearBoundary", nInternalNearBoundary); vertexTable.add("nInternalSurface", nInternalSurface); vertexTable.add("nInternalFeatureEdge", nInternalFeatureEdge); vertexTable.add("nInternalFeaturePoint", nInternalFeaturePoint); vertexTable.add("nExternalSurface", nExternalSurface); vertexTable.add("nExternalFeatureEdge", nExternalFeatureEdge); vertexTable.add("nExternalFeaturePoint", nExternalFeaturePoint); vertexTable.add("nFar", nFar); vertexTable.add("nReferred", nReferred); os << endl; vertexTable.print(os); } template Foam::autoPtr Foam::DelaunayMesh::createMesh ( const fileName& name, labelPairLookup& vertexMap, labelList& cellMap, const bool writeDelaunayData ) const { pointField points(Triangulation::number_of_vertices()); faceList faces(Triangulation::number_of_finite_facets()); labelList owner(Triangulation::number_of_finite_facets()); labelList neighbour(Triangulation::number_of_finite_facets()); wordList patchNames(1, "foamyHexMesh_defaultPatch"); wordList patchTypes(1, wallPolyPatch::typeName); PtrList patchDicts(1); patchDicts.set(0, new dictionary()); List> patchFaces(1, DynamicList()); List> patchOwners(1, DynamicList