Loading...
Searching...
No Matches
setRDeltaT.H
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2013-2016 OpenFOAM Foundation
9 Copyright (C) 2020,2025 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
29{
30 volScalarField& rDeltaT = trDeltaT.ref();
32 const dictionary& pimpleDict = pimple.dict();
33
34 // Maximum flow Courant number
35 scalar maxCo(pimpleDict.get<scalar>("maxCo"));
36
37 // Maximum time scale
38 scalar maxDeltaT(pimpleDict.getOrDefault<scalar>("maxDeltaT", GREAT));
39
40 // Smoothing parameter (0-1) when smoothing iterations > 0
42 (
43 pimpleDict.getOrDefault<scalar>("rDeltaTSmoothingCoeff", 0.1)
44 );
45
46 // Damping coefficient (1-0)
48 (
49 pimpleDict.getOrDefault<scalar>("rDeltaTDampingCoeff", 1.0)
50 );
52 // Maximum change in cell temperature per iteration
53 // (relative to previous value)
54 scalar alphaTemp(pimpleDict.getOrDefault<scalar>("alphaTemp", 0.05));
56 // Maximum change in cell concentration per iteration
57 // (relative to reference value)
58 scalar alphaY(pimpleDict.getOrDefault<scalar>("alphaY", 1.0));
59
60
61 // The old reciprocal time scale field, with any damping factor
62 tmp<volScalarField> rDeltaT0_damped;
63
64 // Calculate damped value before applying any other changes
65 if
66 (
68 && runTime.timeIndex() > runTime.startTimeIndex() + 1
69 )
70 {
71 rDeltaT0_damped = (scalar(1) - rDeltaTDampingCoeff)*(rDeltaT);
72 }
73
74
75 Info<< "Time scales min/max:" << endl;
76
77 // Flow time scale
78 {
79 rDeltaT.ref() =
80 (
81 fvc::surfaceSum(mag(phi))()()
82 /((2*maxCo)*mesh.V()*rho())
83 );
84
85 // Limit the largest time scale (=> smallest reciprocal time)
86 rDeltaT.clamp_min(1/maxDeltaT);
87
88 auto limits = gMinMax(rDeltaT.primitiveField());
89 limits.reset(1/(limits.max()+VSMALL), 1/(limits.min()+VSMALL));
90
91 Info<< " Flow = "
92 << limits.min() << ", " << limits.max() << endl;
93 }
94
95 // Heat release rate time scale
96 if (alphaTemp < 1)
97 {
98 volScalarField::Internal rDeltaTT
99 (
100 mag(Qdot)/(alphaTemp*rho*thermo.Cp()*T)
101 );
102
103 rDeltaT.primitiveFieldRef().clamp_min(rDeltaTT);
104
105 auto limits = gMinMax(rDeltaTT.field());
106 limits.reset(1/(limits.max()+VSMALL), 1/(limits.min()+VSMALL));
107
108 Info<< " Temperature = "
109 << limits.min() << ", " << limits.max() << endl;
110 }
111
112 // Reaction rate time scale
113 if (alphaY < 1)
114 {
115 dictionary Yref(pimpleDict.subDict("Yref"));
116
117 volScalarField::Internal rDeltaTY
118 (
119 IOobject
120 (
121 "rDeltaTY",
122 runTime.timeName(),
123 mesh
124 ),
125 mesh,
126 dimensionedScalar(rDeltaT.dimensions(), Zero)
127 );
128
129 bool foundY = false;
130
131 forAll(Y, i)
132 {
133 if (i != inertIndex && composition.active(i))
134 {
135 volScalarField& Yi = Y[i];
136
137 if (Yref.found(Yi.name()))
138 {
139 foundY = true;
140 const scalar Yrefi = Yref.get<scalar>(Yi.name());
141
142 rDeltaTY.field() = max
143 (
144 mag
145 (
146 reaction->R(Yi)().source()
147 /((Yrefi*alphaY)*(rho*mesh.V()))
148 ),
149 rDeltaTY
150 );
151 }
152 }
153 }
154
155 if (foundY)
156 {
157 rDeltaT.primitiveFieldRef().clamp_min(rDeltaTY);
159 auto limits = gMinMax(rDeltaTY.field());
160 limits.reset(1/(limits.max()+VSMALL), 1/(limits.min()+VSMALL));
161
162 Info<< " Composition = "
163 << limits.min() << ", " << limits.max() << endl;
164 }
165 else
166 {
167 IOWarningIn(args.executable().c_str(), Yref)
168 << "Cannot find any active species in Yref " << Yref
169 << endl;
171 }
172
173 // Update tho boundary values of the reciprocal time-step
174 rDeltaT.correctBoundaryConditions();
175
176 // Spatially smooth the time scale field
177 if (rDeltaTSmoothingCoeff < 1)
179 fvc::smooth(rDeltaT, rDeltaTSmoothingCoeff);
180 }
181
182 // Limit rate of change of time scale (=> smallest reciprocal time)
183 // - reduce as much as required
184 // - only increase at a fraction of old time scale
185 if (rDeltaT0_damped)
187 rDeltaT.clamp_min(rDeltaT0_damped());
188 }
189
190 // Update tho boundary values of the reciprocal time-step
191 rDeltaT.correctBoundaryConditions();
192
193 auto limits = gMinMax(rDeltaT.field());
194 limits.reset(1/(limits.max()+VSMALL), 1/(limits.min()+VSMALL));
195
196 Info<< " Overall = "
197 << limits.min() << ", " << limits.max() << endl;
198}
199
200
201// ************************************************************************* //
const dictionary & pimpleDict
pimpleControl & pimple
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
basicSpecieMixture & composition
PtrList< volScalarField > & Y
const volScalarField & T
scalar rDeltaTDampingCoeff(pimpleDict.getOrDefault< scalar >("rDeltaTDampingCoeff", 1.0))
scalar alphaY(pimpleDict.getOrDefault< scalar >("alphaY", 1.0))
auto limits
Definition setRDeltaT.H:186
scalar rDeltaTSmoothingCoeff(pimpleDict.getOrDefault< scalar >("rDeltaTSmoothingCoeff", 0.1))
tmp< volScalarField > rDeltaT0_damped
Definition setRDeltaT.H:55
scalar alphaTemp(pimpleDict.getOrDefault< scalar >("alphaTemp", 0.05))
dynamicFvMesh & mesh
engineTime & runTime
tmp< volScalarField > trDeltaT
scalar maxCo
label inertIndex
CombustionModel< rhoReactionThermo > & reaction
#define IOWarningIn(functionName, ios)
Report an IO warning using Foam::Warning.
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
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
MinMax< Type > gMinMax(const FieldField< Field, Type > &f)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Foam::argList args(argc, argv)
scalar Qdot
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299