chemistryModel.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) 2016-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 Class
25  Foam::chemistryModel
26 
27 Description
28  Extends base chemistry model by adding a thermo package, and ODE functions.
29  Introduces chemistry equation system and evaluation of chemical source terms
30  with optional support for TDAC mechanism reduction and tabulation.
31 
32  References:
33  \verbatim
34  Contino, F., Jeanmart, H., Lucchini, T., & D’Errico, G. (2011).
35  Coupling of in situ adaptive tabulation and dynamic adaptive chemistry:
36  An effective method for solving combustion in engine simulations.
37  Proceedings of the Combustion Institute, 33(2), 3057-3064.
38 
39  Contino, F., Lucchini, T., D'Errico, G., Duynslaegher, C.,
40  Dias, V., & Jeanmart, H. (2012).
41  Simulations of advanced combustion modes using detailed chemistry
42  combined with tabulation and mechanism reduction techniques.
43  SAE International Journal of Engines,
44  5(2012-01-0145), 185-196.
45 
46  Contino, F., Foucher, F., Dagaut, P., Lucchini, T., D’Errico, G., &
47  Mounaïm-Rousselle, C. (2013).
48  Experimental and numerical analysis of nitric oxide effect on the
49  ignition of iso-octane in a single cylinder HCCI engine.
50  Combustion and Flame, 160(8), 1476-1483.
51 
52  Contino, F., Masurier, J. B., Foucher, F., Lucchini, T., D’Errico, G., &
53  Dagaut, P. (2014).
54  CFD simulations using the TDAC method to model iso-octane combustion
55  for a large range of ozone seeding and temperature conditions
56  in a single cylinder HCCI engine.
57  Fuel, 137, 179-184.
58  \endverbatim
59 
60 SourceFiles
61  chemistryModelI.H
62  chemistryModel.C
63 
64 \*---------------------------------------------------------------------------*/
65 
66 #ifndef chemistryModel_H
67 #define chemistryModel_H
68 
69 #include "odeChemistryModel.H"
70 #include "ReactionList.H"
71 #include "ODESystem.H"
72 #include "volFields.H"
73 #include "multiComponentMixture.H"
76 #include "DynamicField.H"
77 
78 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
79 
80 namespace Foam
81 {
82 
83 /*---------------------------------------------------------------------------*\
84  Class chemistryModel Declaration
85 \*---------------------------------------------------------------------------*/
86 
87 template<class ThermoType>
88 class chemistryModel
89 :
90  public odeChemistryModel
91 {
92  // Private classes
93 
94  //- Class to define scope of reaction evaluation. Runs pre-evaluate
95  // hook on all reactions on construction and post-evaluate on
96  // destruction.
97  class reactionEvaluationScope
98  {
100 
101  public:
102 
103  reactionEvaluationScope
104  (
106  )
107  :
108  chemistry_(chemistry)
109  {
110  forAll(chemistry_.reactions_, i)
111  {
112  chemistry_.reactions_[i].preEvaluate();
113  }
114  }
115 
116  ~reactionEvaluationScope()
117  {
118  forAll(chemistry_.reactions_, i)
119  {
120  chemistry_.reactions_[i].postEvaluate();
121  }
122  }
123  };
124 
125 
126  // Private data
127 
128  //- Switch to select performance logging
129  Switch log_;
130 
131  //- Switch to enable loadBalancing performance logging
132  Switch loadBalancing_;
133 
134  //- Type of the Jacobian to be calculated
135  const jacobianType jacobianType_;
136 
137  //- Reference to the multi component mixture
138  const multiComponentMixture<ThermoType>& mixture_;
139 
140  //- Thermodynamic data of the species
141  const PtrList<ThermoType>& specieThermos_;
142 
143  //- Reactions
144  const ReactionList<ThermoType> reactions_;
145 
146  //- List of reaction rate per specie [kg/m^3/s]
148 
149  //- Temporary mass fraction field
150  mutable scalarField Y_;
151 
152  //- Temporary simplified mechanism mass fraction field
154 
155  //- Temporary concentration field
156  mutable scalarField c_;
157 
158  //- Temporary simplified mechanism concentration field
160 
161  //- Specie-temperature-pressure workspace fields
162  mutable FixedList<scalarField, 5> YTpWork_;
163 
164  //- Specie-temperature-pressure workspace matrices
165  mutable FixedList<scalarSquareMatrix, 2> YTpYTpWork_;
166 
167  //- Mechanism reduction method
169 
170  //- Mechanism reduction method reference
172 
173  //- Tabulation method
174  autoPtr<chemistryTabulationMethod> tabulationPtr_;
175 
176  //- Tabulation method reference
177  chemistryTabulationMethod& tabulation_;
178 
179  //- Log file for average time spent solving the chemistry
180  autoPtr<OFstream> cpuSolveFile_;
181 
182 
183  // Private Member Functions
184 
185  //- Write access to chemical source terms
186  // (e.g. for multi-chemistry model)
188 
189  //- Solve the reaction system for the given time step
190  // of given type and return the characteristic time
191  // Variable number of species added
192  template<class DeltaTType>
193  scalar solve(const DeltaTType& deltaT);
194 
195 
196 public:
197 
198  //- Runtime type information
199  TypeName("chemistryModel");
200 
201 
202  // Constructors
203 
204  //- Construct from thermo
206 
207  //- Disallow default bitwise copy construction
208  chemistryModel(const chemistryModel&) = delete;
209 
210 
211  //- Destructor
212  virtual ~chemistryModel();
213 
214 
215  // Member Functions
216 
217  //- Return reference to the mixture
218  inline const multiComponentMixture<ThermoType>& mixture() const;
219 
220  //- The reactions
221  inline const PtrList<Reaction<ThermoType>>& reactions() const;
222 
223  //- Thermodynamic data of the species
224  inline const PtrList<ThermoType>& specieThermos() const;
225 
226  //- The number of reactions
227  virtual inline label nReaction() const;
228 
229  //- Calculates the reaction rates
230  virtual void calculate();
231 
232 
233  // Chemistry model functions (overriding abstract functions in
234  // basicChemistryModel.H)
235 
236  //- Return const access to the chemical source terms for specie, i
237  inline const volScalarField::Internal& RR
238  (
239  const label i
240  ) const;
241 
242  //- Return non const access to chemical source terms [kg/m^3/s]
243  virtual volScalarField::Internal& RR
244  (
245  const label i
246  );
247 
248  //- Return reaction rate of the speciei in reactionI
250  (
251  const label reactionI,
252  const label speciei
253  ) const;
254 
255  //- Solve the reaction system for the given time step
256  // and return the characteristic time
257  virtual scalar solve(const scalar deltaT);
258 
259  //- Solve the reaction system for the given time step
260  // and return the characteristic time
261  virtual scalar solve(const scalarField& deltaT);
262 
263  //- Return the chemical time scale
264  virtual tmp<volScalarField> tc() const;
265 
266  //- Return the heat release rate [kg/m/s^3]
267  virtual tmp<volScalarField> Qdot() const;
268 
269 
270  // ODE functions (overriding abstract functions in ODE.H)
271 
272  virtual void derivatives
273  (
274  const scalar t,
275  const scalarField& YTp,
276  const label li,
277  scalarField& dYTpdt
278  ) const;
279 
280  virtual void jacobian
281  (
282  const scalar t,
283  const scalarField& YTp,
284  const label li,
285  scalarField& dYTpdt,
287  ) const;
288 
289  virtual void solve
290  (
291  scalar& p,
292  scalar& T,
293  scalarField& Y,
294  const label li,
295  scalar& deltaT,
296  scalar& subDeltaT
297  ) const = 0;
298 
299 
300  // Mechanism reduction access functions
301 
302  //- Return true if specie i is active
303  inline bool active(const label i) const;
304 
305  //- Set specie i active
306  inline void setActive(const label i);
307 
308 
309  // Member Operators
310 
311  //- Disallow default bitwise assignment
312  void operator=(const chemistryModel&) = delete;
313 };
314 
315 
316 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
317 
318 } // End namespace Foam
319 
320 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
321 
322 #include "chemistryModelI.H"
323 
324 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
325 
326 #ifdef NoRepository
327  #include "chemistryModel.C"
328 #endif
329 
330 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
331 
332 #endif
333 
334 // ************************************************************************* //
virtual tmp< volScalarField::Internal > calculateRR(const label reactionI, const label speciei) const
Return reaction rate of the speciei in reactionI.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
virtual void derivatives(const scalar t, const scalarField &YTp, const label li, scalarField &dYTpdt) const
Calculate the derivatives in dydx.
virtual tmp< volScalarField > tc() const
Return the chemical time scale.
jacobianType
Enumeration for the type of Jacobian to be calculated by the.
bool active(const label i) const
Return true if specie i is active.
A 1D vector of objects of type <T> with a fixed size <Size>.
Definition: FixedList.H:54
const PtrList< volScalarField > & Y() const
Return a reference to the list of mass fraction fields.
virtual ~chemistryModel()
Destructor.
A simple wrapper around bool so that it can be read as a word: true/false, on/off, yes/no, y/n, t/f, or none/any.
Definition: Switch.H:60
Base-class for multi-component fluid thermodynamic properties.
Switch chemistry() const
Chemistry activation switch.
const PtrList< ThermoType > & specieThermos() const
Thermodynamic data of the species.
Switch chemistry_
Chemistry activation switch.
void operator=(const chemistryModel &)=delete
Disallow default bitwise assignment.
virtual void jacobian(const scalar t, const scalarField &YTp, const label li, scalarField &dYTpdt, scalarSquareMatrix &J) const
Calculate the Jacobian of the system.
Foam::multiComponentMixture.
chemistryModel(const fluidReactionThermo &thermo)
Construct from thermo.
TypeName("chemistryModel")
Runtime type information.
virtual label nReaction() const
The number of reactions.
const multiComponentMixture< ThermoType > & mixture() const
Return reference to the mixture.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
List of templated reactions.
Definition: ReactionList.H:51
virtual void calculate()
Calculates the reaction rates.
const fluidReactionThermo & thermo() const
Return const access to the thermo.
Extends base chemistry model by adding a thermo package, and ODE functions. Introduces chemistry equa...
rhoEqn solve()
An abstract class for chemistry tabulation.
void setActive(const label i)
Set specie i active.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
An abstract class for methods of chemical mechanism reduction.
virtual tmp< volScalarField > Qdot() const
Return the heat release rate [kg/m/s^3].
volScalarField & p
A class for managing temporary objects.
Definition: PtrList.H:53
Extends base chemistry model adding an ODESystem and the reduction maps needed for tabulation...
const PtrList< Reaction< ThermoType > > & reactions() const
The reactions.
Namespace for OpenFOAM.