Files
openfoam/applications/solvers/multiphase/VoF/setRDeltaT.H
Mark Olesen 3e43edf056 ENH: unify use of dictionary method names
- previously introduced `getOrDefault` as a dictionary _get_ method,
  now complete the transition and use it everywhere instead of
  `lookupOrDefault`. This avoids mixed usage of the two methods that
  are identical in behaviour, makes for shorter names, and promotes
  the distinction between "lookup" access (ie, return a token stream,
  locate and return an entry) and "get" access (ie, the above with
  conversion to concrete types such as scalar, label etc).
2020-06-02 17:26:03 +02:00

137 lines
3.2 KiB
C

{
volScalarField& rDeltaT = trDeltaT.ref();
const dictionary& pimpleDict = pimple.dict();
scalar maxCo
(
pimpleDict.getOrDefault<scalar>("maxCo", 0.9)
);
scalar maxAlphaCo
(
pimpleDict.getOrDefault<scalar>("maxAlphaCo", 0.2)
);
scalar rDeltaTSmoothingCoeff
(
pimpleDict.getOrDefault<scalar>("rDeltaTSmoothingCoeff", 0.1)
);
label nAlphaSpreadIter
(
pimpleDict.getOrDefault<label>("nAlphaSpreadIter", 1)
);
scalar alphaSpreadDiff
(
pimpleDict.getOrDefault<scalar>("alphaSpreadDiff", 0.2)
);
scalar alphaSpreadMax
(
pimpleDict.getOrDefault<scalar>("alphaSpreadMax", 0.99)
);
scalar alphaSpreadMin
(
pimpleDict.getOrDefault<scalar>("alphaSpreadMin", 0.01)
);
label nAlphaSweepIter
(
pimpleDict.getOrDefault<label>("nAlphaSweepIter", 5)
);
scalar rDeltaTDampingCoeff
(
pimpleDict.getOrDefault<scalar>("rDeltaTDampingCoeff", 1.0)
);
scalar maxDeltaT
(
pimpleDict.getOrDefault<scalar>("maxDeltaT", GREAT)
);
volScalarField rDeltaT0("rDeltaT0", rDeltaT);
// Set the reciprocal time-step from the local Courant number
rDeltaT.ref() = max
(
1/dimensionedScalar("maxDeltaT", dimTime, maxDeltaT),
fvc::surfaceSum(mag(rhoPhi))()()
/((2*maxCo)*mesh.V()*rho())
);
if (maxAlphaCo < maxCo)
{
// Further limit the reciprocal time-step
// in the vicinity of the interface
volScalarField alpha1Bar(fvc::average(alpha1));
rDeltaT.ref() = max
(
rDeltaT(),
pos0(alpha1Bar() - alphaSpreadMin)
*pos0(alphaSpreadMax - alpha1Bar())
*fvc::surfaceSum(mag(phi))()()
/((2*maxAlphaCo)*mesh.V())
);
}
// Update tho boundary values of the reciprocal time-step
rDeltaT.correctBoundaryConditions();
Info<< "Flow time scale min/max = "
<< gMin(1/rDeltaT.primitiveField())
<< ", " << gMax(1/rDeltaT.primitiveField()) << endl;
if (rDeltaTSmoothingCoeff < 1.0)
{
fvc::smooth(rDeltaT, rDeltaTSmoothingCoeff);
}
if (nAlphaSpreadIter > 0)
{
fvc::spread
(
rDeltaT,
alpha1,
nAlphaSpreadIter,
alphaSpreadDiff,
alphaSpreadMax,
alphaSpreadMin
);
}
if (nAlphaSweepIter > 0)
{
fvc::sweep(rDeltaT, alpha1, nAlphaSweepIter, alphaSpreadDiff);
}
Info<< "Smoothed flow time scale min/max = "
<< gMin(1/rDeltaT.primitiveField())
<< ", " << gMax(1/rDeltaT.primitiveField()) << endl;
// Limit rate of change of time scale
// - reduce as much as required
// - only increase at a fraction of old time scale
if
(
rDeltaTDampingCoeff < 1.0
&& runTime.timeIndex() > runTime.startTimeIndex() + 1
)
{
rDeltaT = max
(
rDeltaT,
(scalar(1) - rDeltaTDampingCoeff)*rDeltaT0
);
Info<< "Damped flow time scale min/max = "
<< gMin(1/rDeltaT.primitiveField())
<< ", " << gMax(1/rDeltaT.primitiveField()) << endl;
}
}