chemistryModel.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 | Copyright (C) 2011-2015 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
30  terms.
31 
32 SourceFiles
33  chemistryModelI.H
34  chemistryModel.C
35 
36 \*---------------------------------------------------------------------------*/
37 
38 #ifndef chemistryModel_H
39 #define chemistryModel_H
40 
41 #include "Reaction.H"
42 #include "ODESystem.H"
43 #include "volFieldsFwd.H"
44 #include "simpleMatrix.H"
45 #include "DimensionedField.H"
46 
47 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
48 
49 namespace Foam
50 {
51 
52 // Forward declaration of classes
53 class fvMesh;
54 
55 /*---------------------------------------------------------------------------*\
56  Class chemistryModel Declaration
57 \*---------------------------------------------------------------------------*/
58 
59 template<class CompType, class ThermoType>
60 class chemistryModel
61 :
62  public CompType,
63  public ODESystem
64 {
65  // Private Member Functions
66 
67  //- Disallow copy constructor
69 
70  //- Disallow default bitwise assignment
71  void operator=(const chemistryModel&);
72 
73  //- Solve the reaction system for the given time step
74  // of given type and return the characteristic time
75  template<class DeltaTType>
76  scalar solve(const DeltaTType& deltaT);
77 
78 
79 protected:
80 
81  typedef ThermoType thermoType;
82 
83  // Private data
84 
85  //- Reference to the field of specie mass fractions
87 
88  //- Reactions
90 
91  //- Thermodynamic data of the species
93 
94  //- Number of species
96 
97  //- Number of reactions
99 
100  //- Temperature below which the reaction rates are assumed 0
101  scalar Treact_;
102 
103  //- List of reaction rate per specie [kg/m3/s]
105 
106 
107  // Protected Member Functions
108 
109  //- Write access to chemical source terms
110  // (e.g. for multi-chemistry model)
112 
113 
114 public:
115 
116  //- Runtime type information
117  TypeName("chemistryModel");
118 
119 
120  // Constructors
121 
122  //- Construct from mesh
123  chemistryModel(const fvMesh& mesh, const word& phaseName);
124 
125 
126  //- Destructor
127  virtual ~chemistryModel();
128 
129 
130  // Member Functions
131 
132  //- The reactions
133  inline const PtrList<Reaction<ThermoType> >& reactions() const;
134 
135  //- Thermodynamic data of the species
136  inline const PtrList<ThermoType>& specieThermo() const;
137 
138  //- The number of species
139  inline label nSpecie() const;
140 
141  //- The number of reactions
142  inline label nReaction() const;
143 
144  //- Temperature below which the reaction rates are assumed 0
145  inline scalar Treact() const;
146 
147  //- Temperature below which the reaction rates are assumed 0
148  inline scalar& Treact();
149 
150  //- dc/dt = omega, rate of change in concentration, for each species
151  virtual tmp<scalarField> omega
152  (
153  const scalarField& c,
154  const scalar T,
155  const scalar p
156  ) const;
157 
158  //- Return the reaction rate for reaction r and the reference
159  // species and charateristic times
160  virtual scalar omega
161  (
162  const Reaction<ThermoType>& r,
163  const scalarField& c,
164  const scalar T,
165  const scalar p,
166  scalar& pf,
167  scalar& cf,
168  label& lRef,
169  scalar& pr,
170  scalar& cr,
171  label& rRef
172  ) const;
173 
174 
175  //- Return the reaction rate for iReaction and the reference
176  // species and charateristic times
177  virtual scalar omegaI
178  (
179  label iReaction,
180  const scalarField& c,
181  const scalar T,
182  const scalar p,
183  scalar& pf,
184  scalar& cf,
185  label& lRef,
186  scalar& pr,
187  scalar& cr,
188  label& rRef
189  ) const;
190 
191  //- Calculates the reaction rates
192  virtual void calculate();
193 
194 
195  // Chemistry model functions (overriding abstract functions in
196  // basicChemistryModel.H)
197 
198  //- Return const access to the chemical source terms for specie, i
200  (
201  const label i
202  ) const;
203 
204  //- Return non const access to chemical source terms [kg/m3/s]
206  (
207  const label i
208  );
209 
210  //- Return reaction rate of the speciei in reactionI
212  (
213  const label reactionI,
214  const label speciei
215  ) const;
216 
217  //- Solve the reaction system for the given time step
218  // and return the characteristic time
219  virtual scalar solve(const scalar deltaT);
220 
221  //- Solve the reaction system for the given time step
222  // and return the characteristic time
223  virtual scalar solve(const scalarField& deltaT);
224 
225  //- Return the chemical time scale
226  virtual tmp<volScalarField> tc() const;
227 
228  //- Return source for enthalpy equation [kg/m/s3]
229  virtual tmp<volScalarField> Sh() const;
230 
231  //- Return the heat release, i.e. enthalpy/sec [kg/m2/s3]
232  virtual tmp<volScalarField> dQ() const;
233 
234 
235  // ODE functions (overriding abstract functions in ODE.H)
236 
237  //- Number of ODE's to solve
238  virtual label nEqns() const;
239 
240  virtual void derivatives
241  (
242  const scalar t,
243  const scalarField& c,
244  scalarField& dcdt
245  ) const;
246 
247  virtual void jacobian
248  (
249  const scalar t,
250  const scalarField& c,
251  scalarField& dcdt,
252  scalarSquareMatrix& dfdc
253  ) const;
254 
255  virtual void solve
256  (
257  scalarField &c,
258  scalar& T,
259  scalar& p,
260  scalar& deltaT,
261  scalar& subDeltaT
262  ) const;
263 };
264 
265 
266 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
267 
268 } // End namespace Foam
269 
270 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
271 
272 #include "chemistryModelI.H"
273 
274 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
275 
276 #ifdef NoRepository
277 # include "chemistryModel.C"
278 #endif
279 
280 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
281 
282 #endif
283 
284 // ************************************************************************* //
Foam::chemistryModel::Treact
scalar Treact() const
Temperature below which the reaction rates are assumed 0.
Definition: chemistryModelI.H:73
Foam::chemistryModel::nSpecie_
label nSpecie_
Number of species.
Definition: chemistryModel.H:94
volFieldsFwd.H
Foam::chemistryModel::calculate
virtual void calculate()
Calculates the reaction rates.
Definition: chemistryModel.C:714
Foam::chemistryModel::dQ
virtual tmp< volScalarField > dQ() const
Return the heat release, i.e. enthalpy/sec [kg/m2/s3].
Definition: chemistryModel.C:592
Foam::chemistryModel::chemistryModel
chemistryModel(const chemistryModel &)
Disallow copy constructor.
p
p
Definition: pEqn.H:62
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
chemistryModel.C
Foam::chemistryModel::derivatives
virtual void derivatives(const scalar t, const scalarField &c, scalarField &dcdt) const
Calculate the derivatives in dydx.
Definition: chemistryModel.C:278
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
DimensionedField.H
Foam::chemistryModel::reactions_
const PtrList< Reaction< ThermoType > > & reactions_
Reactions.
Definition: chemistryModel.H:88
Foam::chemistryModel::nReaction
label nReaction() const
The number of reactions.
Definition: chemistryModelI.H:65
Foam::chemistryModel::reactions
const PtrList< Reaction< ThermoType > > & reactions() const
The reactions.
Definition: chemistryModelI.H:41
Foam::chemistryModel::calculateRR
virtual tmp< DimensionedField< scalar, volMesh > > calculateRR(const label reactionI, const label speciei) const
Return reaction rate of the speciei in reactionI.
Definition: chemistryModel.C:634
ODESystem.H
Foam::chemistryModel::TypeName
TypeName("chemistryModel")
Runtime type information.
Foam::chemistryModel::RR_
PtrList< DimensionedField< scalar, volMesh > > RR_
List of reaction rate per specie [kg/m3/s].
Definition: chemistryModel.H:103
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::chemistryModel::Y_
PtrList< volScalarField > & Y_
Reference to the field of specie mass fractions.
Definition: chemistryModel.H:85
Foam::chemistryModel::tc
virtual tmp< volScalarField > tc() const
Return the chemical time scale.
Definition: chemistryModel.C:464
Foam::T
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Definition: FieldFieldFunctions.C:55
Foam::chemistryModel::specieThermo
const PtrList< ThermoType > & specieThermo() const
Thermodynamic data of the species.
Definition: chemistryModelI.H:49
Foam::PtrList
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:61
chemistryModelI.H
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Reaction.H
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::chemistryModel::operator=
void operator=(const chemistryModel &)
Disallow default bitwise assignment.
Foam::chemistryModel::nSpecie
label nSpecie() const
The number of species.
Definition: chemistryModelI.H:57
Foam::chemistryModel::omega
virtual tmp< scalarField > omega(const scalarField &c, const scalar T, const scalar p) const
dc/dt = omega, rate of change in concentration, for each species
Definition: chemistryModel.C:96
Foam::chemistryModel::omegaI
virtual scalar omegaI(label iReaction, const scalarField &c, const scalar T, const scalar p, scalar &pf, scalar &cf, label &lRef, scalar &pr, scalar &cr, label &rRef) const
Return the reaction rate for iReaction and the reference.
Definition: chemistryModel.C:138
Foam::SquareMatrix< scalar >
Foam::ODESystem
Abstract base class for the systems of ordinary differential equations.
Definition: ODESystem.H:46
Foam::chemistryModel::solve
scalar solve(const DeltaTType &deltaT)
Solve the reaction system for the given time step.
Foam::chemistryModel::nReaction_
label nReaction_
Number of reactions.
Definition: chemistryModel.H:97
Foam::chemistryModel::jacobian
virtual void jacobian(const scalar t, const scalarField &c, scalarField &dcdt, scalarSquareMatrix &dfdc) const
Calculate the Jacobian of the system.
Definition: chemistryModel.C:323
Foam::chemistryModel::RR
PtrList< DimensionedField< scalar, volMesh > > & RR()
Write access to chemical source terms.
Definition: chemistryModelI.H:33
Foam::chemistryModel::thermoType
ThermoType thermoType
Definition: chemistryModel.H:80
Foam::chemistryModel::Treact_
scalar Treact_
Temperature below which the reaction rates are assumed 0.
Definition: chemistryModel.H:100
Foam::constant::universal::c
const dimensionedScalar c
Speed of light in a vacuum.
simpleMatrix.H
Foam::chemistryModel
Extends base chemistry model by adding a thermo package, and ODE functions. Introduces chemistry equa...
Definition: chemistryModel.H:59
Foam::Reaction
Simple extension of ReactionThermo to handle reaction kinetics in addition to the equilibrium thermod...
Definition: Reaction.H:53
Foam::chemistryModel::nEqns
virtual label nEqns() const
Number of ODE's to solve.
Definition: chemistryModel.C:624
Foam::chemistryModel::Sh
virtual tmp< volScalarField > Sh() const
Return source for enthalpy equation [kg/m/s3].
Definition: chemistryModel.C:550
Foam::chemistryModel::specieThermo_
const PtrList< ThermoType > & specieThermo_
Thermodynamic data of the species.
Definition: chemistryModel.H:91
Foam::chemistryModel::~chemistryModel
virtual ~chemistryModel()
Destructor.
Definition: chemistryModel.C:87
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:51