/*---------------------------------------------------------------------------*\ ========= | \\ / F ield | OpenFOAM: The Open Source CFD Toolbox \\ / O peration | \\ / A nd | Copyright (C) 2015-2019 OpenCFD Ltd. \\/ M anipulation | ------------------------------------------------------------------------------- | Copyright (C) 2011-2017 OpenFOAM Foundation ------------------------------------------------------------------------------- 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 "meshRefinement.H" #include "volMesh.H" #include "volFields.H" #include "surfaceMesh.H" #include "syncTools.H" #include "Time.H" #include "refinementSurfaces.H" #include "refinementFeatures.H" #include "decompositionMethod.H" #include "regionSplit.H" #include "fvMeshDistribute.H" #include "indirectPrimitivePatch.H" #include "polyTopoChange.H" #include "removeCells.H" #include "mapDistributePolyMesh.H" #include "localPointRegion.H" #include "pointMesh.H" #include "pointFields.H" #include "slipPointPatchFields.H" #include "fixedValuePointPatchFields.H" #include "calculatedPointPatchFields.H" #include "cyclicSlipPointPatchFields.H" #include "processorPointPatch.H" #include "globalIndex.H" #include "meshTools.H" #include "OFstream.H" #include "Random.H" #include "searchableSurfaces.H" #include "treeBoundBox.H" #include "zeroGradientFvPatchFields.H" #include "fvMeshTools.H" #include "motionSmoother.H" #include "faceSet.H" #include "topoDistanceData.H" #include "FaceCellWave.H" // Leak path #include "shortestPathSet.H" #include "meshSearch.H" // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * // namespace Foam { defineTypeNameAndDebug(meshRefinement, 0); } const Foam::Enum < Foam::meshRefinement::debugType > Foam::meshRefinement::debugTypeNames ({ { debugType::MESH, "mesh" }, { debugType::OBJINTERSECTIONS, "intersections" }, { debugType::FEATURESEEDS, "featureSeeds" }, { debugType::ATTRACTION, "attraction" }, { debugType::LAYERINFO, "layerInfo" }, }); //const Foam::Enum //< // Foam::meshRefinement::outputType //> //Foam::meshRefinement::outputTypeNames //({ // { outputType::OUTPUTLAYERINFO, "layerInfo" } //}); const Foam::Enum < Foam::meshRefinement::writeType > Foam::meshRefinement::writeTypeNames ({ { writeType::WRITEMESH, "mesh" }, { writeType::NOWRITEREFINEMENT, "noRefinement" }, { writeType::WRITELEVELS, "scalarLevels" }, { writeType::WRITELAYERSETS, "layerSets" }, { writeType::WRITELAYERFIELDS, "layerFields" }, }); Foam::meshRefinement::writeType Foam::meshRefinement::writeLevel_; //Foam::meshRefinement::outputType Foam::meshRefinement::outputLevel_; // Inside/outside test for polyMesh:.findCell() // 2.4.x : default = polyMesh::FACE_DIAG_TRIS // 1712 : default = polyMesh::CELL_TETS // // - CELL_TETS is better with concave cells, but much slower. // - use faster method (FACE_DIAG_TRIS) here static const Foam::polyMesh::cellDecomposition findCellMode(Foam::polyMesh::FACE_DIAG_TRIS); // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * // void Foam::meshRefinement::calcNeighbourData ( labelList& neiLevel, pointField& neiCc ) const { const labelList& cellLevel = meshCutter_.cellLevel(); const pointField& cellCentres = mesh_.cellCentres(); const label nBoundaryFaces = mesh_.nBoundaryFaces(); if (neiLevel.size() != nBoundaryFaces || neiCc.size() != nBoundaryFaces) { FatalErrorInFunction << nBoundaryFaces << " neiLevel:" << neiLevel.size() << abort(FatalError); } const polyBoundaryMesh& patches = mesh_.boundaryMesh(); labelHashSet addedPatchIDSet(meshedPatches()); forAll(patches, patchi) { const polyPatch& pp = patches[patchi]; const labelUList& faceCells = pp.faceCells(); const vectorField::subField faceCentres = pp.faceCentres(); const vectorField::subField faceAreas = pp.faceAreas(); label bFacei = pp.start()-mesh_.nInternalFaces(); if (pp.coupled()) { forAll(faceCells, i) { neiLevel[bFacei] = cellLevel[faceCells[i]]; neiCc[bFacei] = cellCentres[faceCells[i]]; bFacei++; } } else if (addedPatchIDSet.found(patchi)) { // Face was introduced from cell-cell intersection. Try to // reconstruct other side cell(centre). Three possibilities: // - cells same size. // - preserved cell smaller. Not handled. // - preserved cell larger. forAll(faceCells, i) { // Extrapolate the face centre. const vector fn = normalised(faceAreas[i]); label own = faceCells[i]; label ownLevel = cellLevel[own]; label faceLevel = meshCutter_.faceLevel(pp.start()+i); if (faceLevel < 0) { // Due to e.g. face merging no longer a consistent // refinementlevel of face. Assume same as cell faceLevel = ownLevel; } // Normal distance from face centre to cell centre scalar d = ((faceCentres[i] - cellCentres[own]) & fn); if (faceLevel > ownLevel) { // Other cell more refined. Adjust normal distance d *= 0.5; } neiLevel[bFacei] = faceLevel; // Calculate other cell centre by extrapolation neiCc[bFacei] = faceCentres[i] + d*fn; bFacei++; } } else { forAll(faceCells, i) { neiLevel[bFacei] = cellLevel[faceCells[i]]; neiCc[bFacei] = faceCentres[i]; bFacei++; } } } // Swap coupled boundaries. Apply separation to cc since is coordinate. syncTools::swapBoundaryFacePositions(mesh_, neiCc); syncTools::swapBoundaryFaceList(mesh_, neiLevel); } void Foam::meshRefinement::calcCellCellRays ( const pointField& neiCc, const labelList& neiLevel, const labelList& testFaces, pointField& start, pointField& end, labelList& minLevel ) const { const labelList& cellLevel = meshCutter_.cellLevel(); const pointField& cellCentres = mesh_.cellCentres(); start.setSize(testFaces.size()); end.setSize(testFaces.size()); minLevel.setSize(testFaces.size()); forAll(testFaces, i) { const label facei = testFaces[i]; const label own = mesh_.faceOwner()[facei]; if (mesh_.isInternalFace(facei)) { const label nei = mesh_.faceNeighbour()[facei]; start[i] = cellCentres[own]; end[i] = cellCentres[nei]; minLevel[i] = min(cellLevel[own], cellLevel[nei]); } else { const label bFacei = facei - mesh_.nInternalFaces(); start[i] = cellCentres[own]; end[i] = neiCc[bFacei]; minLevel[i] = min(cellLevel[own], neiLevel[bFacei]); } } // Extend segments a bit { const vectorField smallVec(ROOTSMALL*(end-start)); start -= smallVec; end += smallVec; } } void Foam::meshRefinement::updateIntersections(const labelList& changedFaces) { // Stats on edges to test. Count proc faces only once. bitSet isMasterFace(syncTools::getMasterFaces(mesh_)); { label nMasterFaces = isMasterFace.count(); reduce(nMasterFaces, sumOp