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()*rho())
34 );
35
36 if (pimple.transonic())
37 {
38 surfaceScalarField phid
39 (
40 "phid",
41 fvc::interpolate(psi)*fvc::flux(U)
42 );
43
44 rDeltaT.ref() = max
45 (
46 rDeltaT(),
47 fvc::surfaceSum(mag(phid))()()
48 /((2*maxCo)*mesh.V()*psi())
49 );
50 }
51
52 // Update the boundary values of the reciprocal time-step
53 rDeltaT.correctBoundaryConditions();
54
55 {
56 auto limits = gMinMax(rDeltaT.primitiveField());
57 limits.reset(1/(limits.max()+VSMALL), 1/(limits.min()+VSMALL));
58
59 Info<< "Flow time scale min/max = "
60 << limits.min() << ", " << limits.max() << endl;
61 }
62
63 if (rDeltaTSmoothingCoeff < 1.0)
64 {
65 fvc::smooth(rDeltaT, rDeltaTSmoothingCoeff);
66 }
67
68 {
69 auto limits = gMinMax(rDeltaT.primitiveField());
70 limits.reset(1/(limits.max()+VSMALL), 1/(limits.min()+VSMALL));
71
72 Info<< "Smoothed flow time scale min/max = "
73 << limits.min() << ", " << limits.max() << endl;
74 }
75
76 // Limit rate of change of time scale
77 // - reduce as much as required
78 // - only increase at a fraction of old time scale
79 if
80 (
82 && runTime.timeIndex() > runTime.startTimeIndex() + 1
83 )
84 {
85 rDeltaT =
87 *max(rDeltaT/rDeltaT0, scalar(1) - rDeltaTDampingCoeff);
88
89 auto limits = gMinMax(rDeltaT.primitiveField());
90 limits.reset(1/(limits.max()+VSMALL), 1/(limits.min()+VSMALL));
91
92 Info<< "Damped flow time scale min/max = "
93 << limits.min() << ", " << limits.max() << endl;
94 }
95}
const dictionary & pimpleDict
pimpleControl & pimple
U
Definition pEqn.H:72
const volScalarField & psi
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)
surfaceScalarField phid("phid", fvc::interpolate(psi) *(fvc::flux(HbyA)+MRF.zeroFilter(rhorAUf *fvc::ddtCorr(rho, U, phi)/fvc::interpolate(rho))))
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)