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-2022 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  Info<< "Time scales min/max:" << endl;
32 
33  // Cache old reciprocal time scale field
34  const volScalarField rDeltaT0("rDeltaT0", rDeltaT);
35 
36  // Flow time scale
37  {
38  // Maximum flow Courant number
39  const scalar maxCo(pimpleDict.lookup<scalar>("maxCo"));
40 
41  rDeltaT.ref() =
42  (
43  fvc::surfaceSum(mag(phi))()()
44  /((2*maxCo)*mesh.V()*rho())
45  );
46 
47  // Limit the largest time step
48  if (pimpleDict.found("maxDeltaT"))
49  {
50  rDeltaT.max(1/pimpleDict.lookup<scalar>("maxDeltaT"));
51  }
52 
53  // Limit the smallest time step
54  if (pimpleDict.found("minDeltaT"))
55  {
56  rDeltaT.min(1/pimpleDict.lookup<scalar>("minDeltaT"));
57  }
58 
59  Info<< " Flow = "
60  << 1/gMax(rDeltaT.primitiveField()) << ", "
61  << 1/gMin(rDeltaT.primitiveField()) << endl;
62  }
63 
64  // Maximum change in cell temperature per iteration
65  // (relative to previous value)
66  const scalar alphaTemp(pimpleDict.lookupOrDefault("alphaTemp", 0.05));
67 
68  // Heat release rate time scale
69  if (alphaTemp < 1)
70  {
71  volScalarField::Internal rDeltaTT
72  (
73  mag(reaction->Qdot())/(alphaTemp*rho*thermo.Cp()*T)
74  );
75 
76  Info<< " Temperature = "
77  << 1/(gMax(rDeltaTT.field()) + vSmall) << ", "
78  << 1/(gMin(rDeltaTT.field()) + vSmall) << endl;
79 
80  rDeltaT.ref() = max(rDeltaT(), rDeltaTT);
81  }
82 
83  // Reaction rate time scale
84  if (pimpleDict.found("alphaY"))
85  {
86  // Maximum change in cell concentration per iteration
87  // (relative to reference value)
88  const scalar alphaY(pimpleDict.lookup<scalar>("alphaY"));
89 
90  const dictionary Yref(pimpleDict.subDict("Yref"));
91 
92  volScalarField::Internal rDeltaTY
93  (
94  IOobject
95  (
96  "rDeltaTY",
97  runTime.timeName(),
98  mesh
99  ),
100  mesh,
101  dimensionedScalar(rDeltaT.dimensions(), 0)
102  );
103 
104  bool foundY = false;
105 
106  forAll(Y, i)
107  {
108  if (composition.solve(i))
109  {
110  volScalarField& Yi = Y[i];
111 
112  if (Yref.found(Yi.name()))
113  {
114  foundY = true;
115  scalar Yrefi = Yref.lookup<scalar>(Yi.name());
116 
117  rDeltaTY.field() = max
118  (
119  mag
120  (
121  reaction->R(Yi)().source()
122  /((Yrefi*alphaY)*(rho*mesh.V()))
123  ),
124  rDeltaTY
125  );
126  }
127  }
128  }
129 
130  if (foundY)
131  {
132  Info<< " Composition = "
133  << 1/(gMax(rDeltaTY.field()) + vSmall) << ", "
134  << 1/(gMin(rDeltaTY.field()) + vSmall) << endl;
135 
136  rDeltaT.ref() = max(rDeltaT(), rDeltaTY);
137  }
138  else
139  {
140  IOWarningIn(args.executable().c_str(), Yref)
141  << "Cannot find any active species in Yref " << Yref
142  << endl;
143  }
144  }
145 
146  // Update the boundary values of the reciprocal time-step
147  rDeltaT.correctBoundaryConditions();
148 
149  // Smoothing parameter (0-1) when smoothing iterations > 0
150  const scalar rDeltaTSmoothingCoeff
151  (
152  pimpleDict.lookupOrDefault<scalar>("rDeltaTSmoothingCoeff", 0.1)
153  );
154 
155  // Spatially smooth the time scale field
156  if (rDeltaTSmoothingCoeff < 1)
157  {
158  fvc::smooth(rDeltaT, rDeltaTSmoothingCoeff);
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  (
166  pimpleDict.found("rDeltaTDampingCoeff")
167  && runTime.timeIndex() > runTime.startTimeIndex() + 1
168  )
169  {
170  // Damping coefficient (1-0)
171  const scalar rDeltaTDampingCoeff
172  (
173  pimpleDict.lookup<scalar>("rDeltaTDampingCoeff")
174  );
175 
176  rDeltaT = max
177  (
178  rDeltaT,
179  (scalar(1) - rDeltaTDampingCoeff)*rDeltaT0
180  );
181  }
182 
183  // Update the boundary values of the reciprocal time-step
184  rDeltaT.correctBoundaryConditions();
185 
186  Info<< " Overall = "
187  << 1/gMax(rDeltaT.primitiveField())
188  << ", " << 1/gMin(rDeltaT.primitiveField()) << endl;
189 }
190 
191 
192 // ************************************************************************* //
const word & executable() const
Name of executable without the path.
Definition: argListI.H:36
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
fluidReactionThermo & thermo
Definition: createFields.H:28
pimpleNoLoopControl & pimple
Type gMin(const FieldField< Field, Type > &f)
basicSpecieMixture & composition
combustionModel & reaction
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
fvMesh & mesh
void smooth(volScalarField &field, const scalar coeff)
Definition: fvcSmooth.C:35
const volScalarField rDeltaT0("rDeltaT0", rDeltaT)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:58
forAll(phases, phasei)
Definition: setRDeltaT.H:28
phi
Definition: correctPhi.H:3
#define IOWarningIn(functionName, ios)
Report an IO warning using Foam::Warning.
const scalar maxCo(pimpleDict.lookupOrDefault< scalar >("maxCo", 0.9))
Type gMax(const FieldField< Field, Type > &f)
const volScalarField & T
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
PtrList< volScalarField > & Y
const scalar rDeltaTSmoothingCoeff(pimpleDict.lookupOrDefault< scalar >("rDeltaTSmoothingCoeff", 0.1))
messageStream Info
dimensioned< scalar > mag(const dimensioned< Type > &)
const scalar rDeltaTDampingCoeff(pimpleDict.lookupOrDefault< scalar >("rDeltaTDampingCoeff", 1.0))
tmp< GeometricField< Type, fvPatchField, volMesh > > surfaceSum(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
tmp< volScalarField > trDeltaT
Definition: createRDeltaT.H:3
Foam::argList args(argc, argv)