237 lines
6.7 KiB
C
237 lines
6.7 KiB
C
|
/*---------------------------------------------------------------------------*\
|
||
|
========= |
|
||
|
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
||
|
\\ / O peration |
|
||
|
\\ / A nd | Copyright held by original author
|
||
|
\\/ M anipulation |
|
||
|
-------------------------------------------------------------------------------
|
||
|
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 2 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, write to the Free Software Foundation,
|
||
|
Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
|
||
|
|
||
|
\*---------------------------------------------------------------------------*/
|
||
|
|
||
|
#include "fvFieldDecomposer.H"
|
||
|
#include "processorFvPatchField.H"
|
||
|
#include "processorFvsPatchField.H"
|
||
|
|
||
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||
|
|
||
|
namespace Foam
|
||
|
{
|
||
|
|
||
|
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
|
||
|
|
||
|
template<class Type>
|
||
|
tmp<GeometricField<Type, fvPatchField, volMesh> >
|
||
|
fvFieldDecomposer::decomposeField
|
||
|
(
|
||
|
const GeometricField<Type, fvPatchField, volMesh>& field
|
||
|
) const
|
||
|
{
|
||
|
// Create and map the internal field values
|
||
|
Field<Type> internalField(field.internalField(), cellAddressing_);
|
||
|
|
||
|
// Create and map the patch field values
|
||
|
PtrList<fvPatchField<Type> > patchFields(boundaryAddressing_.size());
|
||
|
|
||
|
forAll (boundaryAddressing_, patchi)
|
||
|
{
|
||
|
if (boundaryAddressing_[patchi] >= 0)
|
||
|
{
|
||
|
patchFields.set
|
||
|
(
|
||
|
patchi,
|
||
|
fvPatchField<Type>::New
|
||
|
(
|
||
|
field.boundaryField()[boundaryAddressing_[patchi]],
|
||
|
procMesh_.boundary()[patchi],
|
||
|
DimensionedField<Type, volMesh>::null(),
|
||
|
*patchFieldDecomposerPtrs_[patchi]
|
||
|
)
|
||
|
);
|
||
|
}
|
||
|
else
|
||
|
{
|
||
|
patchFields.set
|
||
|
(
|
||
|
patchi,
|
||
|
new processorFvPatchField<Type>
|
||
|
(
|
||
|
procMesh_.boundary()[patchi],
|
||
|
DimensionedField<Type, volMesh>::null(),
|
||
|
Field<Type>
|
||
|
(
|
||
|
field.internalField(),
|
||
|
*processorVolPatchFieldDecomposerPtrs_[patchi]
|
||
|
)
|
||
|
)
|
||
|
);
|
||
|
}
|
||
|
}
|
||
|
|
||
|
// Create the field for the processor
|
||
|
return tmp<GeometricField<Type, fvPatchField, volMesh> >
|
||
|
(
|
||
|
new GeometricField<Type, fvPatchField, volMesh>
|
||
|
(
|
||
|
IOobject
|
||
|
(
|
||
|
field.name(),
|
||
|
procMesh_.time().timeName(),
|
||
|
procMesh_,
|
||
|
IOobject::NO_READ,
|
||
|
IOobject::NO_WRITE
|
||
|
),
|
||
|
procMesh_,
|
||
|
field.dimensions(),
|
||
|
internalField,
|
||
|
patchFields
|
||
|
)
|
||
|
);
|
||
|
}
|
||
|
|
||
|
|
||
|
template<class Type>
|
||
|
tmp<GeometricField<Type, fvsPatchField, surfaceMesh> >
|
||
|
fvFieldDecomposer::decomposeField
|
||
|
(
|
||
|
const GeometricField<Type, fvsPatchField, surfaceMesh>& field
|
||
|
) const
|
||
|
{
|
||
|
labelList mapAddr
|
||
|
(
|
||
|
labelList::subList
|
||
|
(
|
||
|
faceAddressing_,
|
||
|
procMesh_.nInternalFaces()
|
||
|
)
|
||
|
);
|
||
|
forAll (mapAddr, i)
|
||
|
{
|
||
|
mapAddr[i] -= 1;
|
||
|
}
|
||
|
|
||
|
// Create and map the internal field values
|
||
|
Field<Type> internalField
|
||
|
(
|
||
|
field.internalField(),
|
||
|
mapAddr
|
||
|
);
|
||
|
|
||
|
// Problem with addressing when a processor patch picks up both internal
|
||
|
// faces and faces from cyclic boundaries. This is a bit of a hack, but
|
||
|
// I cannot find a better solution without making the internal storage
|
||
|
// mechanism for surfaceFields correspond to the one of faces in polyMesh
|
||
|
// (i.e. using slices)
|
||
|
Field<Type> allFaceField(field.mesh().nFaces());
|
||
|
|
||
|
forAll (field.internalField(), i)
|
||
|
{
|
||
|
allFaceField[i] = field.internalField()[i];
|
||
|
}
|
||
|
|
||
|
forAll (field.boundaryField(), patchi)
|
||
|
{
|
||
|
const Field<Type> & p = field.boundaryField()[patchi];
|
||
|
|
||
|
const label patchStart = field.mesh().boundaryMesh()[patchi].start();
|
||
|
|
||
|
forAll (p, i)
|
||
|
{
|
||
|
allFaceField[patchStart + i] = p[i];
|
||
|
}
|
||
|
}
|
||
|
|
||
|
// Create and map the patch field values
|
||
|
PtrList<fvsPatchField<Type> > patchFields(boundaryAddressing_.size());
|
||
|
|
||
|
forAll (boundaryAddressing_, patchi)
|
||
|
{
|
||
|
if (boundaryAddressing_[patchi] >= 0)
|
||
|
{
|
||
|
patchFields.set
|
||
|
(
|
||
|
patchi,
|
||
|
fvsPatchField<Type>::New
|
||
|
(
|
||
|
field.boundaryField()[boundaryAddressing_[patchi]],
|
||
|
procMesh_.boundary()[patchi],
|
||
|
DimensionedField<Type, surfaceMesh>::null(),
|
||
|
*patchFieldDecomposerPtrs_[patchi]
|
||
|
)
|
||
|
);
|
||
|
}
|
||
|
else
|
||
|
{
|
||
|
patchFields.set
|
||
|
(
|
||
|
patchi,
|
||
|
new processorFvsPatchField<Type>
|
||
|
(
|
||
|
procMesh_.boundary()[patchi],
|
||
|
DimensionedField<Type, surfaceMesh>::null(),
|
||
|
Field<Type>
|
||
|
(
|
||
|
allFaceField,
|
||
|
*processorSurfacePatchFieldDecomposerPtrs_[patchi]
|
||
|
)
|
||
|
)
|
||
|
);
|
||
|
}
|
||
|
}
|
||
|
|
||
|
// Create the field for the processor
|
||
|
return tmp<GeometricField<Type, fvsPatchField, surfaceMesh> >
|
||
|
(
|
||
|
new GeometricField<Type, fvsPatchField, surfaceMesh>
|
||
|
(
|
||
|
IOobject
|
||
|
(
|
||
|
field.name(),
|
||
|
procMesh_.time().timeName(),
|
||
|
procMesh_,
|
||
|
IOobject::NO_READ,
|
||
|
IOobject::NO_WRITE
|
||
|
),
|
||
|
procMesh_,
|
||
|
field.dimensions(),
|
||
|
internalField,
|
||
|
patchFields
|
||
|
)
|
||
|
);
|
||
|
}
|
||
|
|
||
|
|
||
|
template<class GeoField>
|
||
|
void fvFieldDecomposer::decomposeFields
|
||
|
(
|
||
|
const PtrList<GeoField>& fields
|
||
|
) const
|
||
|
{
|
||
|
forAll (fields, fieldI)
|
||
|
{
|
||
|
decomposeField(fields[fieldI])().write();
|
||
|
}
|
||
|
}
|
||
|
|
||
|
|
||
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||
|
|
||
|
} // End namespace Foam
|
||
|
|
||
|
// ************************************************************************* //
|