mirror of
https://develop.openfoam.com/Development/openfoam.git
synced 2025-11-28 03:28:01 +00:00
BUG: setTurbulenceFields: update processor boundaries (fixes #2527)
This commit is contained in:
@ -109,6 +109,7 @@ Note
|
||||
#include "singlePhaseTransportModel.H"
|
||||
#include "turbulentTransportModel.H"
|
||||
#include "turbulentFluidThermoModel.H"
|
||||
#include "processorFvPatchField.H"
|
||||
#include "wallFvPatch.H"
|
||||
#include "fixedValueFvPatchFields.H"
|
||||
|
||||
@ -122,6 +123,40 @@ void InfoField(const word& fldName)
|
||||
}
|
||||
|
||||
|
||||
template<class Type>
|
||||
void correctProcessorPatches
|
||||
(
|
||||
GeometricField<Type, fvPatchField, volMesh>& vf
|
||||
)
|
||||
{
|
||||
if (!Pstream::parRun())
|
||||
{
|
||||
return;
|
||||
}
|
||||
|
||||
// Not possible to use correctBoundaryConditions on fields as they may
|
||||
// use local info as opposed to the constraint values employed here,
|
||||
// but still need to update processor patches
|
||||
auto& bf = vf.boundaryFieldRef();
|
||||
|
||||
forAll(bf, patchi)
|
||||
{
|
||||
if (isA<processorFvPatchField<Type>>(bf[patchi]))
|
||||
{
|
||||
bf[patchi].initEvaluate();
|
||||
}
|
||||
}
|
||||
|
||||
forAll(bf, patchi)
|
||||
{
|
||||
if (isA<processorFvPatchField<Type>>(bf[patchi]))
|
||||
{
|
||||
bf[patchi].evaluate();
|
||||
}
|
||||
}
|
||||
}
|
||||
|
||||
|
||||
IOobject createIOobject
|
||||
(
|
||||
const fvMesh& mesh,
|
||||
@ -447,22 +482,26 @@ int main(int argc, char *argv[])
|
||||
// (M:Eq. 9)
|
||||
const dimensionedScalar maxU(dimVelocity, SMALL);
|
||||
U *= min(scalar(1), fRei*uTau/max(mag(U), maxU));
|
||||
correctProcessorPatches<vector>(U);
|
||||
}
|
||||
|
||||
if (tepsilon.valid())
|
||||
{
|
||||
tepsilon.ref() = epsilon;
|
||||
correctProcessorPatches<scalar>(tepsilon.ref());
|
||||
}
|
||||
|
||||
if (tk.valid())
|
||||
{
|
||||
tk.ref() = k;
|
||||
correctProcessorPatches<scalar>(tk.ref());
|
||||
}
|
||||
|
||||
if (tomega.valid())
|
||||
{
|
||||
const dimensionedScalar k0(sqr(dimLength/dimTime), SMALL);
|
||||
tomega.ref() = Cmu*epsilon/(k + k0);
|
||||
correctProcessorPatches<scalar>(tomega.ref());
|
||||
}
|
||||
|
||||
if (tR.valid())
|
||||
@ -475,6 +514,7 @@ int main(int argc, char *argv[])
|
||||
{
|
||||
R[celli] = Rdiag[celli];
|
||||
}
|
||||
correctProcessorPatches<symmTensor>(R);
|
||||
}
|
||||
|
||||
|
||||
|
||||
Reference in New Issue
Block a user