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.9)
9 );
10
12 (
13 pimpleDict.getOrDefault<scalar>("maxAlphaCo", 0.2)
14 );
15
17 (
18 pimpleDict.getOrDefault<scalar>("rDeltaTSmoothingCoeff", 0.1)
19 );
20
22 (
23 pimpleDict.getOrDefault<label>("nAlphaSpreadIter", 1)
24 );
25
27 (
28 pimpleDict.getOrDefault<scalar>("alphaSpreadDiff", 0.2)
29 );
30
32 (
33 pimpleDict.getOrDefault<scalar>("alphaSpreadMax", 0.99)
34 );
35
37 (
38 pimpleDict.getOrDefault<scalar>("alphaSpreadMin", 0.01)
39 );
40
42 (
43 pimpleDict.getOrDefault<label>("nAlphaSweepIter", 5)
44 );
45
47 (
48 pimpleDict.getOrDefault<scalar>("rDeltaTDampingCoeff", 1.0)
49 );
50
51 scalar maxDeltaT
52 (
53 pimpleDict.getOrDefault<scalar>("maxDeltaT", GREAT)
54 );
55
56
57 // The old reciprocal time scale field, with any damping factor
58 tmp<volScalarField> rDeltaT0_damped;
59
60 // Calculate damped value before applying any other changes
61 if
62 (
64 && runTime.timeIndex() > runTime.startTimeIndex() + 1
65 )
66 {
67 rDeltaT0_damped = (scalar(1) - rDeltaTDampingCoeff)*(rDeltaT);
68 }
69
70
71 volScalarField rDeltaT0("rDeltaT0", rDeltaT);
72
73 // Set the reciprocal time-step from the local Courant number
74 rDeltaT.ref() = max
75 (
76 1/dimensionedScalar("maxDeltaT", dimTime, maxDeltaT),
77 fvc::surfaceSum(mag(rhoPhi))()()
78 /((2*maxCo)*mesh.V()*rho())
79 );
80
82 {
83 // Further limit the reciprocal time-step
84 // in the vicinity of the interface
85
86 volScalarField alpha1Bar(fvc::average(alpha1));
87
88 rDeltaT.ref() = max
89 (
90 rDeltaT(),
91 pos0(alpha1Bar() - alphaSpreadMin)
92 *pos0(alphaSpreadMax - alpha1Bar())
93 *fvc::surfaceSum(mag(phi))()()
94 /((2*maxAlphaCo)*mesh.V())
95 );
96 }
97
98 // Update tho boundary values of the reciprocal time-step
99 rDeltaT.correctBoundaryConditions();
100
101 {
102 auto limits = gMinMax(rDeltaT.primitiveField());
103 limits.reset(1/(limits.max()+VSMALL), 1/(limits.min()+VSMALL));
104
105 Info<< "Flow time scale min/max = "
106 << limits.min() << ", " << limits.max() << endl;
107 }
108 if (rDeltaTSmoothingCoeff < 1.0)
109 {
110 fvc::smooth(rDeltaT, rDeltaTSmoothingCoeff);
111 }
112
113 if (nAlphaSpreadIter > 0)
114 {
115 fvc::spread
116 (
117 rDeltaT,
118 alpha1,
123 );
124 }
125
126 if (nAlphaSweepIter > 0)
127 {
128 fvc::sweep(rDeltaT, alpha1, nAlphaSweepIter, alphaSpreadDiff);
129 }
130
131 {
132 auto limits = gMinMax(rDeltaT.primitiveField());
133 limits.reset(1/(limits.max()+VSMALL), 1/(limits.min()+VSMALL));
134
135 Info<< "Smoothed flow time scale min/max = "
136 << limits.min() << ", " << limits.max() << endl;
137 }
138
139 // Limit rate of change of time scale (=> smallest reciprocal time)
140 // - reduce as much as required
141 // - only increase at a fraction of old time scale
142 if (rDeltaT0_damped)
143 {
144 rDeltaT.clamp_min(rDeltaT0_damped());
145
146 auto limits = gMinMax(rDeltaT.primitiveField());
147 limits.reset(1/(limits.max()+VSMALL), 1/(limits.min()+VSMALL));
148
149 Info<< "Damped flow time scale min/max = "
150 << limits.min() << ", " << limits.max() << endl;
151 }
152}
const dictionary & pimpleDict
rhoPhi
Definition rhoEqn.H:10
pimpleControl & pimple
const volScalarField & alpha1
scalar rDeltaTDampingCoeff(pimpleDict.getOrDefault< scalar >("rDeltaTDampingCoeff", 1.0))
auto limits
Definition setRDeltaT.H:186
scalar rDeltaTSmoothingCoeff(pimpleDict.getOrDefault< scalar >("rDeltaTSmoothingCoeff", 0.1))
tmp< volScalarField > rDeltaT0_damped
Definition setRDeltaT.H:55
volScalarField rDeltaT0("rDeltaT0", rDeltaT)
dynamicFvMesh & mesh
engineTime & runTime
tmp< volScalarField > trDeltaT
scalar maxCo
scalar alphaSpreadMax(pimpleDict.getOrDefault< scalar >("alphaSpreadMax", 0.99))
scalar alphaSpreadDiff(pimpleDict.getOrDefault< scalar >("alphaSpreadDiff", 0.2))
scalar alphaSpreadMin(pimpleDict.getOrDefault< scalar >("alphaSpreadMin", 0.01))
label nAlphaSpreadIter(pimpleDict.getOrDefault< label >("nAlphaSpreadIter", 1))
label nAlphaSweepIter(pimpleDict.getOrDefault< label >("nAlphaSweepIter", 5))
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)
scalar maxAlphaCo