mirror of
https://develop.openfoam.com/Development/openfoam.git
synced 2025-11-28 03:28:01 +00:00
reactingFoam::setRDeltaT: Add support for limiting the local time-step by the reaction rates
e.g. in the reactingFoam/laminar/counterFlowFlame2DLTS tutorial:
PIMPLE
{
momentumPredictor no;
nOuterCorrectors 1;
nCorrectors 1;
nNonOrthogonalCorrectors 0;
maxDeltaT 1e-2;
maxCo 1;
alphaTemp 0.05;
alphaY 0.05;
Yref
{
O2 0.1;
".*" 1;
}
rDeltaTSmoothingCoeff 1;
rDeltaTDampingCoeff 1;
}
will limit the LTS time-step according to the rate of consumption of 'O2'
normalized by the reference mass-fraction of 0.1 and all other species
normalized by the reference mass-fraction of 1. Additionally the time-step
factor of 'alphaY' is applied to all species. Only the species specified in the
'Yref' sub-dictionary are included in the LTS limiter and if 'alphaY' is omitted
or set to 1 the reaction rates are not included in the LTS limiter.
This commit is contained in:
@ -43,13 +43,16 @@ License
|
||||
// Damping coefficient (1-0)
|
||||
scalar rDeltaTDampingCoeff
|
||||
(
|
||||
pimpleDict.lookupOrDefault<scalar>("rDeltaTDampingCoeff", 1)
|
||||
pimpleDict.lookupOrDefault<scalar>("rDeltaTDampingCoeff", 1.0)
|
||||
);
|
||||
|
||||
// Maximum change in cell temperature per iteration
|
||||
// (relative to previous value)
|
||||
scalar alphaTemp(pimpleDict.lookupOrDefault("alphaTemp", 0.05));
|
||||
|
||||
// Maximum change in cell concentration per iteration
|
||||
// (relative to reference value)
|
||||
scalar alphaY(pimpleDict.lookupOrDefault("alphaY", 1.0));
|
||||
|
||||
Info<< "Time scales min/max:" << endl;
|
||||
|
||||
@ -68,12 +71,12 @@ License
|
||||
rDeltaT.max(1/maxDeltaT);
|
||||
|
||||
Info<< " Flow = "
|
||||
<< gMin(1/rDeltaT.primitiveField()) << ", "
|
||||
<< gMax(1/rDeltaT.primitiveField()) << endl;
|
||||
<< 1/gMax(rDeltaT.primitiveField()) << ", "
|
||||
<< 1/gMin(rDeltaT.primitiveField()) << endl;
|
||||
}
|
||||
|
||||
// Reaction source time scale
|
||||
if (alphaTemp < 1.0)
|
||||
// Heat release rate time scale
|
||||
if (alphaTemp < 1)
|
||||
{
|
||||
volScalarField::Internal rDeltaTT
|
||||
(
|
||||
@ -81,21 +84,76 @@ License
|
||||
);
|
||||
|
||||
Info<< " Temperature = "
|
||||
<< gMin(1/(rDeltaTT.field() + VSMALL)) << ", "
|
||||
<< gMax(1/(rDeltaTT.field() + VSMALL)) << endl;
|
||||
<< 1/(gMax(rDeltaTT.field()) + VSMALL) << ", "
|
||||
<< 1/(gMin(rDeltaTT.field()) + VSMALL) << endl;
|
||||
|
||||
rDeltaT.ref() = max
|
||||
rDeltaT.ref() = max(rDeltaT(), rDeltaTT);
|
||||
}
|
||||
|
||||
// Reaction rate time scale
|
||||
if (alphaY < 1)
|
||||
{
|
||||
dictionary Yref(pimpleDict.subDict("Yref"));
|
||||
|
||||
volScalarField::Internal rDeltaTY
|
||||
(
|
||||
rDeltaT(),
|
||||
rDeltaTT
|
||||
IOobject
|
||||
(
|
||||
"rDeltaTY",
|
||||
runTime.timeName(),
|
||||
mesh
|
||||
),
|
||||
mesh,
|
||||
dimensionedScalar("rDeltaTY", rDeltaT.dimensions(), 0)
|
||||
);
|
||||
|
||||
bool foundY = false;
|
||||
|
||||
forAll(Y, i)
|
||||
{
|
||||
if (i != inertIndex && composition.active(i))
|
||||
{
|
||||
volScalarField& Yi = Y[i];
|
||||
|
||||
if (Yref.found(Yi.name()))
|
||||
{
|
||||
foundY = true;
|
||||
scalar Yrefi = readScalar(Yref.lookup(Yi.name()));
|
||||
|
||||
rDeltaTY.field() = max
|
||||
(
|
||||
mag
|
||||
(
|
||||
reaction->R(Yi)().source()
|
||||
/((Yrefi*alphaY)*(rho*mesh.V()))
|
||||
),
|
||||
rDeltaTY
|
||||
);
|
||||
}
|
||||
}
|
||||
}
|
||||
|
||||
if (foundY)
|
||||
{
|
||||
Info<< " Composition = "
|
||||
<< 1/(gMax(rDeltaTY.field()) + VSMALL) << ", "
|
||||
<< 1/(gMin(rDeltaTY.field()) + VSMALL) << endl;
|
||||
|
||||
rDeltaT.ref() = max(rDeltaT(), rDeltaTY);
|
||||
}
|
||||
else
|
||||
{
|
||||
IOWarningIn(args.executable().c_str(), Yref)
|
||||
<< "Cannot find any active species in Yref " << Yref
|
||||
<< endl;
|
||||
}
|
||||
}
|
||||
|
||||
// Update tho boundary values of the reciprocal time-step
|
||||
rDeltaT.correctBoundaryConditions();
|
||||
|
||||
// Spatially smooth the time scale field
|
||||
if (rDeltaTSmoothingCoeff < 1.0)
|
||||
if (rDeltaTSmoothingCoeff < 1)
|
||||
{
|
||||
fvc::smooth(rDeltaT, rDeltaTSmoothingCoeff);
|
||||
}
|
||||
@ -105,7 +163,7 @@ License
|
||||
// - only increase at a fraction of old time scale
|
||||
if
|
||||
(
|
||||
rDeltaTDampingCoeff < 1.0
|
||||
rDeltaTDampingCoeff < 1
|
||||
&& runTime.timeIndex() > runTime.startTimeIndex() + 1
|
||||
)
|
||||
{
|
||||
@ -120,8 +178,8 @@ License
|
||||
rDeltaT.correctBoundaryConditions();
|
||||
|
||||
Info<< " Overall = "
|
||||
<< gMin(1/rDeltaT.primitiveField())
|
||||
<< ", " << gMax(1/rDeltaT.primitiveField()) << endl;
|
||||
<< 1/gMax(rDeltaT.primitiveField())
|
||||
<< ", " << 1/gMin(rDeltaT.primitiveField()) << endl;
|
||||
}
|
||||
|
||||
|
||||
|
||||
Reference in New Issue
Block a user