/*---------------------------------------------------------------------------*\ ========= | \\ / F ield | OpenFOAM: The Open Source CFD Toolbox \\ / O peration | \\ / A nd | www.openfoam.com \\/ M anipulation | ------------------------------------------------------------------------------- Copyright (C) 2012-2016 OpenFOAM Foundation Copyright (C) 2015-2022 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 "fvMeshTools.H" #include "pointSet.H" #include "faceSet.H" #include "cellSet.H" #include "IOobjectList.H" #include "basicFvGeometryScheme.H" #include "processorPolyPatch.H" #include "processorCyclicPolyPatch.H" #include "polyBoundaryMeshEntries.H" // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // // Adds patch if not yet there. Returns patchID. Foam::label Foam::fvMeshTools::addPatch ( fvMesh& mesh, const polyPatch& patch, const dictionary& patchFieldDict, const word& defaultPatchFieldType, const bool validBoundary ) { auto& polyPatches = const_cast(mesh.boundaryMesh()); label patchi = polyPatches.findPatchID(patch.name()); if (patchi != -1) { // Already there return patchi; } // Append at end unless there are processor patches label insertPatchi = polyPatches.size(); label startFacei = mesh.nFaces(); if (!isA(patch)) { forAll(polyPatches, patchi) { const polyPatch& pp = polyPatches[patchi]; if (isA(pp)) { insertPatchi = patchi; startFacei = pp.start(); break; } } } // Below is all quite a hack. Feel free to change once there is a better // mechanism to insert and reorder patches. // Clear local fields and e.g. polyMesh parallelInfo. mesh.clearOut(); label sz = polyPatches.size(); auto& fvPatches = const_cast(mesh.boundary()); // Add polyPatch at the end polyPatches.resize(sz+1); polyPatches.set ( sz, patch.clone ( polyPatches, insertPatchi, //index 0, //size startFacei //start ) ); fvPatches.resize(sz+1); fvPatches.set ( sz, fvPatch::New ( polyPatches[sz], // point to newly added polyPatch mesh.boundary() ) ); do { #undef doLocalCode #define doLocalCode(FieldType) \ { \ addPatchFields \ ( \ mesh, patchFieldDict, defaultPatchFieldType, Zero \ ); \ } // Volume fields doLocalCode(volScalarField); doLocalCode(volVectorField); doLocalCode(volSphericalTensorField); doLocalCode(volSymmTensorField); doLocalCode(volTensorField); // Surface fields doLocalCode(surfaceScalarField); doLocalCode(surfaceVectorField); doLocalCode(surfaceSphericalTensorField); doLocalCode(surfaceSymmTensorField); doLocalCode(surfaceTensorField); #undef doLocalCode } while (false); // Create reordering list // - patches before insert position stay as is // - patches after insert position move one up // - appended patch gets moved to insert position labelList oldToNew(identity(sz+1)); for (label i = insertPatchi; i < sz; ++i) { oldToNew[i] = i+1; } oldToNew[sz] = insertPatchi; // Shuffle into place polyPatches.reorder(oldToNew, validBoundary); fvPatches.reorder(oldToNew); do { #undef doLocalCode #define doLocalCode(FieldType) \ { \ reorderPatchFields(mesh, oldToNew); \ } // Volume fields doLocalCode(volScalarField); doLocalCode(volVectorField); doLocalCode(volSphericalTensorField); doLocalCode(volSymmTensorField); doLocalCode(volTensorField); // Surface fields doLocalCode(surfaceScalarField); doLocalCode(surfaceVectorField); doLocalCode(surfaceSphericalTensorField); doLocalCode(surfaceSymmTensorField); doLocalCode(surfaceTensorField); #undef doLocalCode } while (false); return insertPatchi; } void Foam::fvMeshTools::setPatchFields ( fvMesh& mesh, const label patchi, const dictionary& patchFieldDict ) { do { #undef doLocalCode #define doLocalCode(FieldType) \ { \ setPatchFields(mesh, patchi, patchFieldDict); \ } // Volume fields doLocalCode(volScalarField); doLocalCode(volVectorField); doLocalCode(volSphericalTensorField); doLocalCode(volSymmTensorField); doLocalCode(volTensorField); // Surface fields doLocalCode(surfaceScalarField); doLocalCode(surfaceVectorField); doLocalCode(surfaceSphericalTensorField); doLocalCode(surfaceSymmTensorField); doLocalCode(surfaceTensorField); #undef doLocalCode } while (false); } void Foam::fvMeshTools::zeroPatchFields(fvMesh& mesh, const label patchi) { do { #undef doLocalCode #define doLocalCode(FieldType) \ { \ setPatchFields(mesh, patchi, Zero); \ } // Volume fields doLocalCode(volScalarField); doLocalCode(volVectorField); doLocalCode(volSphericalTensorField); doLocalCode(volSymmTensorField); doLocalCode(volTensorField); // Surface fields doLocalCode(surfaceScalarField); doLocalCode(surfaceVectorField); doLocalCode(surfaceSphericalTensorField); doLocalCode(surfaceSymmTensorField); doLocalCode(surfaceTensorField); #undef doLocalCode } while (false); } // Deletes last patch void Foam::fvMeshTools::trimPatches(fvMesh& mesh, const label nPatches) { // Clear local fields and e.g. polyMesh globalMeshData. mesh.clearOut(); auto& polyPatches = const_cast(mesh.boundaryMesh()); auto& fvPatches = const_cast(mesh.boundary()); if (polyPatches.empty()) { FatalErrorInFunction << "No patches in mesh" << abort(FatalError); } label nFaces = 0; for (label patchi = nPatches; patchi < polyPatches.size(); patchi++) { nFaces += polyPatches[patchi].size(); } reduce(nFaces, sumOp