setRDeltaT.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2013-2020 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 {
27  volScalarField& rDeltaT = trDeltaT.ref();
28 
29  const dictionary& pimpleDict = pimple.dict();
30 
31  // Maximum flow Courant number
32  scalar maxCo(pimpleDict.lookup<scalar>("maxCo"));
33 
34  // Maximum time scale
35  scalar maxDeltaT(pimpleDict.lookupOrDefault<scalar>("maxDeltaT", great));
36 
37  // Smoothing parameter (0-1) when smoothing iterations > 0
39  (
40  pimpleDict.lookupOrDefault<scalar>("rDeltaTSmoothingCoeff", 0.1)
41  );
42 
43  // Damping coefficient (1-0)
44  scalar rDeltaTDampingCoeff
45  (
46  pimpleDict.lookupOrDefault<scalar>("rDeltaTDampingCoeff", 1.0)
47  );
48 
49  // Maximum change in cell temperature per iteration
50  // (relative to previous value)
51  scalar alphaTemp(pimpleDict.lookupOrDefault("alphaTemp", 0.05));
52 
53  // Maximum change in cell concentration per iteration
54  // (relative to reference value)
55  scalar alphaY(pimpleDict.lookupOrDefault("alphaY", 1.0));
56 
57  Info<< "Time scales min/max:" << endl;
58 
59  // Cache old reciprocal time scale field
60  volScalarField rDeltaT0("rDeltaT0", rDeltaT);
61 
62  // Flow time scale
63  {
64  rDeltaT.ref() =
65  (
66  fvc::surfaceSum(mag(phi))()()
67  /((2*maxCo)*mesh.V()*rho())
68  );
69 
70  // Limit the largest time scale
71  rDeltaT.max(1/maxDeltaT);
72 
73  Info<< " Flow = "
74  << 1/gMax(rDeltaT.primitiveField()) << ", "
75  << 1/gMin(rDeltaT.primitiveField()) << endl;
76  }
77 
78  // Heat release rate time scale
79  if (alphaTemp < 1)
80  {
81  volScalarField::Internal rDeltaTT
82  (
83  mag(reaction->Qdot())/(alphaTemp*rho*thermo.Cp()*T)
84  );
85 
86  Info<< " Temperature = "
87  << 1/(gMax(rDeltaTT.field()) + vSmall) << ", "
88  << 1/(gMin(rDeltaTT.field()) + vSmall) << endl;
89 
90  rDeltaT.ref() = max(rDeltaT(), rDeltaTT);
91  }
92 
93  // Reaction rate time scale
94  if (alphaY < 1)
95  {
96  dictionary Yref(pimpleDict.subDict("Yref"));
97 
98  volScalarField::Internal rDeltaTY
99  (
100  IOobject
101  (
102  "rDeltaTY",
103  runTime.timeName(),
104  mesh
105  ),
106  mesh,
107  dimensionedScalar(rDeltaT.dimensions(), 0)
108  );
109 
110  bool foundY = false;
111 
112  forAll(Y, i)
113  {
114  if (i != inertIndex && composition.active(i))
115  {
116  volScalarField& Yi = Y[i];
117 
118  if (Yref.found(Yi.name()))
119  {
120  foundY = true;
121  scalar Yrefi = Yref.lookup<scalar>(Yi.name());
122 
123  rDeltaTY.field() = max
124  (
125  mag
126  (
127  reaction->R(Yi)().source()
128  /((Yrefi*alphaY)*(rho*mesh.V()))
129  ),
130  rDeltaTY
131  );
132  }
133  }
134  }
135 
136  if (foundY)
137  {
138  Info<< " Composition = "
139  << 1/(gMax(rDeltaTY.field()) + vSmall) << ", "
140  << 1/(gMin(rDeltaTY.field()) + vSmall) << endl;
141 
142  rDeltaT.ref() = max(rDeltaT(), rDeltaTY);
143  }
144  else
145  {
146  IOWarningIn(args.executable().c_str(), Yref)
147  << "Cannot find any active species in Yref " << Yref
148  << endl;
149  }
150  }
151 
152  // Update the boundary values of the reciprocal time-step
153  rDeltaT.correctBoundaryConditions();
154 
155  // Spatially smooth the time scale field
157  {
159  }
160 
161  // Limit rate of change of time scale
162  // - reduce as much as required
163  // - only increase at a fraction of old time scale
164  if
165  (
167  && runTime.timeIndex() > runTime.startTimeIndex() + 1
168  )
169  {
170  rDeltaT = max
171  (
172  rDeltaT,
173  (scalar(1) - rDeltaTDampingCoeff)*rDeltaT0
174  );
175  }
176 
177  // Update the boundary values of the reciprocal time-step
178  rDeltaT.correctBoundaryConditions();
179 
180  Info<< " Overall = "
181  << 1/gMax(rDeltaT.primitiveField())
182  << ", " << 1/gMin(rDeltaT.primitiveField()) << endl;
183 }
184 
185 
186 // ************************************************************************* //
const word & executable() const
Name of executable without the path.
Definition: argListI.H:36
pimpleNoLoopControl & pimple
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
Type gMin(const FieldField< Field, Type > &f)
basicSpecieMixture & composition
scalar alphaY(pimpleDict.lookupOrDefault("alphaY", 1.0))
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
scalar rDeltaTDampingCoeff(pimpleDict.lookupOrDefault< scalar >("rDeltaTDampingCoeff", 1.0))
rhoReactionThermo & thermo
Definition: createFields.H:28
void smooth(volScalarField &field, const scalar coeff)
Definition: fvcSmooth.C:35
phi
Definition: pEqn.H:104
label inertIndex
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
dynamicFvMesh & mesh
scalar maxCo(pimpleDict.lookupOrDefault< scalar >("maxCo", 0.9))
scalar rDeltaTSmoothingCoeff(pimpleDict.lookupOrDefault< scalar >("rDeltaTSmoothingCoeff", 0.1))
forAll(phases, phasei)
Definition: setRDeltaT.H:23
CombustionModel< rhoReactionThermo > & reaction
scalar alphaTemp(pimpleDict.lookupOrDefault("alphaTemp", 0.05))
#define IOWarningIn(functionName, ios)
Report an IO warning using Foam::Warning.
Type gMax(const FieldField< Field, Type > &f)
const volScalarField & T
volScalarField rDeltaT0("rDeltaT0", rDeltaT)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
PtrList< volScalarField > & Y
scalar maxDeltaT(pimpleDict.lookupOrDefault< scalar >("maxDeltaT", great))
messageStream Info
dimensioned< scalar > mag(const dimensioned< Type > &)
tmp< GeometricField< Type, fvPatchField, volMesh > > surfaceSum(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
tmp< volScalarField > trDeltaT
Definition: createRDeltaT.H:3
Foam::argList args(argc, argv)