Files
openfoam/applications/solvers/multiphase/VoF/setRDeltaT.H
Henry Weller b167c95f19 compressibleInterFoam: Completed LTS and semi-implicit MULES support
Now the interFoam and compressibleInterFoam families of solvers use the same
alphaEqn formulation and supporting all of the MULES options without
code-duplication.

The semi-implicit MULES support allows running with significantly larger
time-steps but this does reduce the interface sharpness.
2017-02-09 17:31:57 +00:00

137 lines
3.2 KiB
C

{
volScalarField& rDeltaT = trDeltaT.ref();
const dictionary& pimpleDict = pimple.dict();
scalar maxCo
(
pimpleDict.lookupOrDefault<scalar>("maxCo", 0.9)
);
scalar maxAlphaCo
(
pimpleDict.lookupOrDefault<scalar>("maxAlphaCo", 0.2)
);
scalar rDeltaTSmoothingCoeff
(
pimpleDict.lookupOrDefault<scalar>("rDeltaTSmoothingCoeff", 0.1)
);
label nAlphaSpreadIter
(
pimpleDict.lookupOrDefault<label>("nAlphaSpreadIter", 1)
);
scalar alphaSpreadDiff
(
pimpleDict.lookupOrDefault<scalar>("alphaSpreadDiff", 0.2)
);
scalar alphaSpreadMax
(
pimpleDict.lookupOrDefault<scalar>("alphaSpreadMax", 0.99)
);
scalar alphaSpreadMin
(
pimpleDict.lookupOrDefault<scalar>("alphaSpreadMin", 0.01)
);
label nAlphaSweepIter
(
pimpleDict.lookupOrDefault<label>("nAlphaSweepIter", 5)
);
scalar rDeltaTDampingCoeff
(
pimpleDict.lookupOrDefault<scalar>("rDeltaTDampingCoeff", 1.0)
);
scalar maxDeltaT
(
pimpleDict.lookupOrDefault<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(),
pos(alpha1Bar() - alphaSpreadMin)
*pos(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.0) - rDeltaTDampingCoeff)*rDeltaT0
);
Info<< "Damped flow time scale min/max = "
<< gMin(1/rDeltaT.primitiveField())
<< ", " << gMax(1/rDeltaT.primitiveField()) << endl;
}
}