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