/*---------------------------------------------------------------------------*\ ========= | \\ / F ield | OpenFOAM: The Open Source CFD Toolbox \\ / O peration | \\ / A nd | www.openfoam.com \\/ M anipulation | ------------------------------------------------------------------------------- Copyright (C) 2016-2017 Wikki Ltd Copyright (C) 2021-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 "faFieldDecomposer.H" #include "GeometricField.H" #include "IOobjectList.H" #include "processorFaPatchField.H" #include "processorFaePatchField.H" // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * // template Foam::tmp> Foam::faFieldDecomposer::decomposeField ( const GeometricField& field ) const { // Create and map the internal field values Field internalField(field.internalField(), faceAddressing_); // Create and map the patch field values PtrList> patchFields(boundaryAddressing_.size()); forAll(boundaryAddressing_, patchi) { const label oldPatchi = boundaryAddressing_[patchi]; if (oldPatchi >= 0) { patchFields.set ( patchi, faPatchField::New ( field.boundaryField()[oldPatchi], procMesh_.boundary()[patchi], faPatchField::Internal::null(), patchFieldDecomposerPtrs_[patchi] ) ); } else { patchFields.set ( patchi, new processorFaPatchField ( procMesh_.boundary()[patchi], faPatchField::Internal::null(), Field ( field.internalField(), processorAreaPatchFieldDecomposerPtrs_[patchi] ) ) ); } } // Create the field for the processor return tmp>::New ( IOobject ( field.name(), procMesh_.thisDb().time().timeName(), procMesh_.thisDb(), IOobject::NO_READ, IOobject::NO_WRITE ), procMesh_, field.dimensions(), internalField, patchFields ); } template Foam::tmp> Foam::faFieldDecomposer::decomposeField ( const GeometricField& field ) const { labelList mapAddr ( labelList::subList ( edgeAddressing_, procMesh_.nInternalEdges() ) ); forAll(mapAddr, i) { mapAddr[i] -= 1; } // Create and map the internal field values Field internalField ( field.internalField(), mapAddr ); // Problem with addressing when a processor patch picks up both internal // edges and edges from cyclic boundaries. This is a bit of a hack, but // I cannot find a better solution without making the internal storage // mechanism for edgeFields correspond to the one of edges in polyMesh // (i.e. using slices) Field allEdgeField(field.mesh().nEdges()); forAll(field.internalField(), i) { allEdgeField[i] = field.internalField()[i]; } forAll(field.boundaryField(), patchi) { const Field& p = field.boundaryField()[patchi]; const label patchStart = field.mesh().boundary()[patchi].start(); forAll(p, i) { allEdgeField[patchStart + i] = p[i]; } } // Create and map the patch field values PtrList> patchFields(boundaryAddressing_.size()); forAll(boundaryAddressing_, patchi) { const label oldPatchi = boundaryAddressing_[patchi]; if (oldPatchi >= 0) { patchFields.set ( patchi, faePatchField::New ( field.boundaryField()[oldPatchi], procMesh_.boundary()[patchi], faePatchField::Internal::null(), patchFieldDecomposerPtrs_[patchi] ) ); } else { patchFields.set ( patchi, new processorFaePatchField ( procMesh_.boundary()[patchi], faePatchField::Internal::null(), Field ( allEdgeField, processorEdgePatchFieldDecomposerPtrs_[patchi] ) ) ); } } // Create the field for the processor return tmp>::New ( IOobject ( field.name(), procMesh_.thisDb().time().timeName(), procMesh_.thisDb(), IOobject::NO_READ, IOobject::NO_WRITE ), procMesh_, field.dimensions(), internalField, patchFields ); } template void Foam::faFieldDecomposer::decomposeFields ( const PtrList& fields ) const { forAll(fields, fieldi) { decomposeField(fields[fieldi])().write(); } } // ************************************************************************* //