Loading...
Searching...
No Matches
setRDeltaT.H
Go to the documentation of this file.
1{
2 volScalarField& rDeltaT = trDeltaT.ref();
3
4 const dictionary& pimpleDict = pimple.dict();
5
6 scalar maxCo
7 (
8 pimpleDict.getOrDefault<scalar>("maxCo", 0.8)
9 );
10
12 (
13 pimpleDict.getOrDefault<scalar>("rDeltaTSmoothingCoeff", 0.02)
14 );
15
17 (
18 pimpleDict.getOrDefault<scalar>("rDeltaTDampingCoeff", 1.0)
19 );
20
21 scalar maxDeltaT
22 (
23 pimpleDict.getOrDefault<scalar>("maxDeltaT", GREAT)
24 );
25
26 volScalarField rDeltaT0("rDeltaT0", rDeltaT);
27
28 // Set the reciprocal time-step from the local Courant number
29 rDeltaT.ref() = max
30 (
31 1/dimensionedScalar("maxDeltaT", dimTime, maxDeltaT),
32 fvc::surfaceSum(mag(phi))()()
33 /((2*maxCo)*mesh.V())
34 );
35
36 // Update the boundary values of the reciprocal time-step
37 rDeltaT.correctBoundaryConditions();
38
39 {
40 auto limits = gMinMax(rDeltaT.primitiveField());
41 limits.reset(1/(limits.max()+VSMALL), 1/(limits.min()+VSMALL));
42
43 Info<< "Flow time scale min/max = "
44 << limits.min() << ", " << limits.max() << endl;
45 }
46
47 if (rDeltaTSmoothingCoeff < 1.0)
48 {
49 fvc::smooth(rDeltaT, rDeltaTSmoothingCoeff);
50 }
51
52 {
53 auto limits = gMinMax(rDeltaT.primitiveField());
54 limits.reset(1/(limits.max()+VSMALL), 1/(limits.min()+VSMALL));
55
56 Info<< "Smoothed flow time scale min/max = "
57 << limits.min() << ", " << limits.max() << endl;
58 }
59
60 // Limit rate of change of time scale
61 // - reduce as much as required
62 // - only increase at a fraction of old time scale
63 if
64 (
66 && runTime.timeIndex() > runTime.startTimeIndex() + 1
67 )
68 {
69 rDeltaT =
71 *max(rDeltaT/rDeltaT0, scalar(1) - rDeltaTDampingCoeff);
72
73 auto limits = gMinMax(rDeltaT.primitiveField());
74 limits.reset(1/(limits.max()+VSMALL), 1/(limits.min()+VSMALL));
75
76 Info<< "Damped flow time scale min/max = "
77 << limits.min() << ", " << limits.max() << endl;
78 }
79}
const dictionary & pimpleDict
pimpleControl & pimple
scalar rDeltaTDampingCoeff(pimpleDict.getOrDefault< scalar >("rDeltaTDampingCoeff", 1.0))
auto limits
Definition setRDeltaT.H:186
scalar rDeltaTSmoothingCoeff(pimpleDict.getOrDefault< scalar >("rDeltaTSmoothingCoeff", 0.1))
volScalarField rDeltaT0("rDeltaT0", rDeltaT)
dynamicFvMesh & mesh
engineTime & runTime
tmp< volScalarField > trDeltaT
scalar maxCo
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
GeometricField< scalar, fvPatchField, volMesh > volScalarField
messageStream Info
Information stream (stdout output on master, null elsewhere).
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
MinMax< Type > gMinMax(const FieldField< Field, Type > &f)