adjustTimeStepToCombustion.C
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) 2023 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 
27 #include "combustionModel.H"
28 #include "solver.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace functionObjects
36 {
38 
40  (
44  );
45 }
46 }
47 
48 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
49 
51 Foam::functionObjects::adjustTimeStepToCombustion::propsDictIo
52 (
53  const IOobject::readOption& r
54 ) const
55 {
56  return
57  typeIOobject<timeIOdictionary>
58  (
59  name() + "Properties",
60  obr_.time().name(),
61  "uniform",
62  obr_,
63  r,
65  false
66  );
67 }
68 
69 
70 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
71 
73 (
74  const word& name,
75  const Time& runTime,
76  const dictionary& dict
77 )
78 :
79  regionFunctionObject(name, runTime, dict),
80  phaseName_(word::null),
81  maxCo_(NaN),
82  extrapolate_(false),
83  haveCombustionDeltaT0_(false),
84  combustionDeltaT0_(NaN)
85 {
86  read(dict);
87 }
88 
89 
90 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
91 
93 {}
94 
95 
96 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
97 
99 (
100  const dictionary& dict
101 )
102 {
103  phaseName_ = dict.lookupOrDefault<word>("phase", word::null);
104  maxCo_ = dict.lookupOrDefault<scalar>("maxCo", 1);
105  extrapolate_ = dict.lookupOrDefault<bool>("extrapolate", false);
106 
108  (
109  this->propsDictIo(IOobject::MUST_READ_IF_MODIFIED)
110  );
111 
112  if (propsDictIo.headerOk())
113  {
114  const timeIOdictionary propsDict(propsDictIo);
115 
116  haveCombustionDeltaT0_ = true;
117  combustionDeltaT0_ = propsDict.lookup<scalar>("combustionDeltaT");
118  }
119  else
120  {
121  haveCombustionDeltaT0_ = false;
122  }
123 
124  return true;
125 }
126 
127 
129 {
130  return true;
131 }
132 
133 
135 {
136  if (extrapolate_ && obr_.time().writeTime())
137  {
139 
140  propsDict.add("combustionDeltaT", combustionDeltaT0_);
141 
142  propsDict.regIOobject::write();
143  }
144 
145  return true;
146 }
147 
148 
149 Foam::scalar
151 {
152  if (!time_.controlDict().lookupOrDefault("adjustTimeStep", false))
153  {
154  return vGreat;
155  }
156 
157  const combustionModel& combustion =
158  obr_.lookupObject<combustionModel>
159  (
161  (
163  phaseName_
164  )
165  );
166 
167  const fluidMulticomponentThermo& thermo = combustion.thermo();
168 
169  // Build a mass turnover rate
170  volScalarField::Internal rhoDotByRho
171  (
173  (
174  "rhoDotByRho",
175  combustion.mesh(),
177  )
178  );
179  forAll(thermo.Y(), i)
180  {
181  if (thermo.solveSpecie(i))
182  {
183  rhoDotByRho += mag(combustion.R(i))/2/thermo.rho()();
184  }
185  }
186 
187  // Convert to a time-scale
188  const scalar combustionDeltaT1 = maxCo_/max(gMax(rhoDotByRho), vSmall);
189 
190  // We want to clip the time-step to the time-scale, but also additionally
191  // reduce the time-step significantly if that time-scale is reducing
192  // rapidly. This helps us catch the onset of reactions. The following does
193  // this by passing the change ratio into a steeply non-linear function.
194  scalar deltaT;
195  if (extrapolate_ && haveCombustionDeltaT0_)
196  {
197  const scalar x = min(max(combustionDeltaT1/combustionDeltaT0_, 0), 1);
198  const scalar f = (1 - rootSmall)*(1 - sqrt(1 - sqr(x))) + rootSmall;
199  deltaT = f*combustionDeltaT1;
200  }
201  else
202  {
203  deltaT = combustionDeltaT1;
204  }
205 
206  // Store the latest combustion time-step
207  haveCombustionDeltaT0_ = true;
208  combustionDeltaT0_ = combustionDeltaT1;
209 
210  return deltaT;
211 }
212 
213 
214 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
static tmp< DimensionedField< Type, GeoMesh > > New(const word &name, const Mesh &mesh, const dimensionSet &, const Field< Type > &)
Return a temporary field constructed from name, mesh,.
@ MUST_READ_IF_MODIFIED
Definition: IOobject.H:119
static word groupName(Name name, const word &group)
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:76
Base class for combustion models.
const fluidMulticomponentThermo & thermo() const
Return const access to the thermo.
virtual tmp< volScalarField::Internal > R(const label speciei) const =0
Consumption rate, i.e. source term for specie equation.
static const word combustionPropertiesName
Default combustionProperties dictionary name.
const fvMesh & mesh() const
Return const access to the mesh database.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
const word & name() const
Return const reference to name.
Base-class for multi-component fluid thermodynamic properties.
Abstract base-class for Time/database functionObjects.
const word & name() const
Return the name of this functionObject.
Returns the minimum bulk reaction time scale.
virtual scalar maxDeltaT() const
Return the minimum combustion time-scale.
adjustTimeStepToCombustion(const word &name, const Time &runTime, const dictionary &dict)
Construct from components.
virtual bool read(const dictionary &)
Read and reset the timeStep Function1.
Specialisation of Foam::functionObject for a region and providing a reference to the region Foam::obj...
const objectRegistry & obr_
Reference to the region objectRegistry.
const Time & time() const
Return time.
timeIOdictionary derived from IOdictionary with globalFile set false to enable writing to processor t...
Templated form of IOobject providing type information for file reading and header type checking.
Definition: IOobject.H:531
A class for handling words, derived from string.
Definition: word.H:62
static const word null
An empty word.
Definition: word.H:77
defineTypeNameAndDebug(adjustTimeStepToCombustion, 0)
addToRunTimeSelectionTable(functionObject, adjustTimeStepToCombustion, dictionary)
Namespace for OpenFOAM.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
const dimensionSet dimless
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
const dimensionSet dimTime
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< scalar > mag(const dimensioned< Type > &)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
Type gMax(const FieldField< Field, Type > &f)
labelList f(nPoints)
dictionary dict
fluidMulticomponentThermo & thermo
Definition: createFields.H:31
IOdictionary propsDict(systemDict("particleTracksDict", args, runTime))