pimpleFoam: Added LTS capability for demonstration and testing

For most steady cases simpleFoam is likely to converge faster than pimpleFoam
with LTS but this capability may be useful for testing meshes, BCs etc. for more
complex solver for which SIMPLE is not stable and LTS is provided instead.
This commit is contained in:
Henry Weller
2019-03-28 11:54:55 +00:00
parent 87d3d2d3d0
commit e1e3e2a333
15 changed files with 748 additions and 4 deletions

View File

@ -1,3 +1,5 @@
#include "createRDeltaT.H"
Info<< "Reading field p\n" << endl;
volScalarField p
(

View File

@ -39,6 +39,8 @@ Description
#include "pimpleControl.H"
#include "CorrectPhi.H"
#include "fvOptions.H"
#include "localEulerDdtScheme.H"
#include "fvcSmooth.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
@ -53,11 +55,15 @@ int main(int argc, char *argv[])
#include "createDyMControls.H"
#include "createFields.H"
#include "createUfIfPresent.H"
#include "CourantNo.H"
#include "setInitialDeltaT.H"
turbulence->validate();
if (!LTS)
{
#include "CourantNo.H"
#include "setInitialDeltaT.H"
}
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
Info<< "\nStarting time loop\n" << endl;
@ -65,8 +71,16 @@ int main(int argc, char *argv[])
while (runTime.run())
{
#include "readDyMControls.H"
#include "CourantNo.H"
#include "setDeltaT.H"
if (LTS)
{
#include "setRDeltaT.H"
}
else
{
#include "CourantNo.H"
#include "setDeltaT.H"
}
runTime++;

View File

@ -0,0 +1,69 @@
{
volScalarField& rDeltaT = trDeltaT.ref();
const dictionary& pimpleDict = pimple.dict();
scalar maxCo
(
pimpleDict.lookupOrDefault<scalar>("maxCo", 0.8)
);
scalar rDeltaTSmoothingCoeff
(
pimpleDict.lookupOrDefault<scalar>("rDeltaTSmoothingCoeff", 0.02)
);
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(dimTime, maxDeltaT),
fvc::surfaceSum(mag(phi))()()
/((2*maxCo)*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);
}
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 =
rDeltaT0
*max(rDeltaT/rDeltaT0, scalar(1) - rDeltaTDampingCoeff);
Info<< "Damped flow time scale min/max = "
<< gMin(1/rDeltaT.primitiveField())
<< ", " << gMax(1/rDeltaT.primitiveField()) << endl;
}
}