- patchFields now get mapped (instead of created) - with -consistent it now maps all patches except for processor ones (they are the only ones that are processor-local) - all constraint patches get evaluated after mapping to bring them up to date. Patch contributed by Mattijs Janssens
194 lines
5.6 KiB
C++
194 lines
5.6 KiB
C++
/*---------------------------------------------------------------------------*\
|
|
========= |
|
|
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
|
\\ / O peration |
|
|
\\ / A nd | Copyright (C) 2011-2016 OpenFOAM Foundation
|
|
\\/ 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 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 <http://www.gnu.org/licenses/>.
|
|
|
|
\*---------------------------------------------------------------------------*/
|
|
|
|
#ifndef MapConsistentVolFields_H
|
|
#define MapConsistentVolFields_H
|
|
|
|
#include "GeometricField.H"
|
|
#include "meshToMesh.H"
|
|
#include "IOobjectList.H"
|
|
|
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
|
|
|
namespace Foam
|
|
{
|
|
|
|
template<class Type>
|
|
void evaluateConstraintTypes(GeometricField<Type, fvPatchField, volMesh>& fld)
|
|
{
|
|
typename GeometricField<Type, fvPatchField, volMesh>::
|
|
Boundary& fldBf = fld.boundaryFieldRef();
|
|
|
|
if
|
|
(
|
|
Pstream::defaultCommsType == Pstream::blocking
|
|
|| Pstream::defaultCommsType == Pstream::nonBlocking
|
|
)
|
|
{
|
|
label nReq = Pstream::nRequests();
|
|
|
|
forAll(fldBf, patchi)
|
|
{
|
|
fvPatchField<Type>& tgtField = fldBf[patchi];
|
|
|
|
if
|
|
(
|
|
tgtField.type() == tgtField.patch().patch().type()
|
|
&& polyPatch::constraintType(tgtField.patch().patch().type())
|
|
)
|
|
{
|
|
tgtField.initEvaluate(Pstream::defaultCommsType);
|
|
}
|
|
}
|
|
|
|
// Block for any outstanding requests
|
|
if
|
|
(
|
|
Pstream::parRun()
|
|
&& Pstream::defaultCommsType == Pstream::nonBlocking
|
|
)
|
|
{
|
|
Pstream::waitRequests(nReq);
|
|
}
|
|
|
|
forAll(fldBf, patchi)
|
|
{
|
|
fvPatchField<Type>& tgtField = fldBf[patchi];
|
|
|
|
if
|
|
(
|
|
tgtField.type() == tgtField.patch().patch().type()
|
|
&& polyPatch::constraintType(tgtField.patch().patch().type())
|
|
)
|
|
{
|
|
tgtField.evaluate(Pstream::defaultCommsType);
|
|
}
|
|
}
|
|
}
|
|
else if (Pstream::defaultCommsType == Pstream::scheduled)
|
|
{
|
|
const lduSchedule& patchSchedule =
|
|
fld.mesh().globalData().patchSchedule();
|
|
|
|
forAll(patchSchedule, patchEvali)
|
|
{
|
|
label patchi = patchSchedule[patchEvali].patch;
|
|
fvPatchField<Type>& tgtField = fldBf[patchi];
|
|
|
|
if
|
|
(
|
|
tgtField.type() == tgtField.patch().patch().type()
|
|
&& polyPatch::constraintType(tgtField.patch().patch().type())
|
|
)
|
|
{
|
|
if (patchSchedule[patchEvali].init)
|
|
{
|
|
tgtField.initEvaluate(Pstream::scheduled);
|
|
}
|
|
else
|
|
{
|
|
tgtField.evaluate(Pstream::scheduled);
|
|
}
|
|
}
|
|
}
|
|
}
|
|
}
|
|
|
|
|
|
template<class Type, class CombineOp>
|
|
void MapVolFields
|
|
(
|
|
const IOobjectList& objects,
|
|
const HashSet<word>& selectedFields,
|
|
const meshToMesh& interp,
|
|
const CombineOp& cop
|
|
)
|
|
{
|
|
typedef GeometricField<Type, fvPatchField, volMesh> fieldType;
|
|
|
|
const fvMesh& meshSource = static_cast<const fvMesh&>(interp.srcRegion());
|
|
const fvMesh& meshTarget = static_cast<const fvMesh&>(interp.tgtRegion());
|
|
|
|
IOobjectList fields = objects.lookupClass(fieldType::typeName);
|
|
|
|
forAllIter(IOobjectList, fields, fieldIter)
|
|
{
|
|
const word& fieldName = fieldIter()->name();
|
|
|
|
if (selectedFields.empty() || selectedFields.found(fieldName))
|
|
{
|
|
const fieldType fieldSource(*fieldIter(), meshSource);
|
|
|
|
IOobject targetIO
|
|
(
|
|
fieldName,
|
|
meshTarget.time().timeName(),
|
|
meshTarget,
|
|
IOobject::MUST_READ
|
|
);
|
|
|
|
if (targetIO.headerOk())
|
|
{
|
|
Info<< " interpolating onto existing field "
|
|
<< fieldName << endl;
|
|
fieldType fieldTarget(targetIO, meshTarget);
|
|
|
|
interp.mapSrcToTgt(fieldSource, cop, fieldTarget);
|
|
|
|
evaluateConstraintTypes(fieldTarget);
|
|
|
|
fieldTarget.write();
|
|
}
|
|
else
|
|
{
|
|
Info<< " creating new field "
|
|
<< fieldName << endl;
|
|
|
|
targetIO.readOpt() = IOobject::NO_READ;
|
|
|
|
tmp<fieldType>
|
|
tfieldTarget(interp.mapSrcToTgt(fieldSource, cop));
|
|
|
|
fieldType fieldTarget(targetIO, tfieldTarget);
|
|
|
|
evaluateConstraintTypes(fieldTarget);
|
|
|
|
fieldTarget.write();
|
|
}
|
|
}
|
|
}
|
|
}
|
|
|
|
|
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
|
|
|
} // End namespace Foam
|
|
|
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
|
|
|
#endif
|
|
|
|
// ************************************************************************* //
|