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) 2011-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", 0.2)
50 );
52 // Maximum change in cell temperature per iteration
53 // (relative to previous value)
54 scalar alphaTemp(pimpleDict.getOrDefault("alphaTemp", 0.05));
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 Info<< "Time scales min/max:" << endl;
72
73 // Flow time scale
74 {
75 rDeltaT.ref() =
76 (
77 fvc::surfaceSum(mag(phi))()()
78 /((2*maxCo)*mesh.V()*rho())
79 );
80
81 // Limit the largest time scale (=> smallest reciprocal time)
82 rDeltaT.clamp_min(1/maxDeltaT);
83
84 auto limits = gMinMax(rDeltaT.primitiveField());
85 limits.reset(1/(limits.max()+VSMALL), 1/(limits.min()+VSMALL));
86
87 Info<< " Flow = "
88 << limits.min() << ", " << limits.max() << endl;
89 }
90
91 // Reaction source time scale
92 {
93 volScalarField::Internal rDeltaTT
94 (
95 mag
96 (
97 parcels.hsTrans()/(mesh.V()*runTime.deltaT())
98 + Qdot
99 )
100 /(
102 *rho()
103 *thermo.Cp()()()
104 *T()
105 )
106 );
107
108 rDeltaT.primitiveFieldRef().clamp_min(rDeltaTT);
109
110 auto limits = gMinMax(rDeltaTT.field());
111 limits.reset(1/(limits.max()+VSMALL), 1/(limits.min()+VSMALL));
112
113 Info<< " Temperature = "
114 << limits.min() << ", " << limits.max() << endl;
115 }
116
117 // Update the boundary values of the reciprocal time-step
118 rDeltaT.correctBoundaryConditions();
119
120 // Spatially smooth the time scale field
121 if (rDeltaTSmoothingCoeff < 1.0)
122 {
123 fvc::smooth(rDeltaT, rDeltaTSmoothingCoeff);
124 }
125
126 // Limit rate of change of time scale (=> smallest reciprocal time)
127 // - reduce as much as required
128 // - only increase at a fraction of old time scale
129 if (rDeltaT0_damped)
130 {
131 rDeltaT.clamp_min(rDeltaT0_damped());
132 }
133
134 // Update the boundary values of the reciprocal time-step
135 rDeltaT.correctBoundaryConditions();
136
137 auto limits = gMinMax(rDeltaT.primitiveField());
138 limits.reset(1/(limits.max()+VSMALL), 1/(limits.min()+VSMALL));
139
140 Info<< " Overall = "
141 << limits.min() << ", " << limits.max() << endl;
142}
143
144
145// ************************************************************************* //
const dictionary & pimpleDict
pimpleControl & pimple
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
const volScalarField & T
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
scalar alphaTemp(pimpleDict.getOrDefault< scalar >("alphaTemp", 0.05))
dynamicFvMesh & mesh
engineTime & runTime
tmp< volScalarField > trDeltaT
scalar maxCo
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)
scalar Qdot