StandardChemistryModel.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) 2011-2017 OpenFOAM Foundation
9-------------------------------------------------------------------------------
10License
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
26Class
27 Foam::StandardChemistryModel
28
29Description
30 Extends base chemistry model by adding a thermo package, and ODE functions.
31 Introduces chemistry equation system and evaluation of chemical source
32 terms.
33
34SourceFiles
35 StandardChemistryModelI.H
36 StandardChemistryModel.C
37
38\*---------------------------------------------------------------------------*/
39
40#ifndef StandardChemistryModel_H
41#define StandardChemistryModel_H
42
43#include "BasicChemistryModel.H"
44#include "Reaction.H"
45#include "ODESystem.H"
46#include "volFields.H"
47#include "simpleMatrix.H"
48
49// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
50
51namespace Foam
52{
53
54// Forward declaration of classes
55class fvMesh;
56
57/*---------------------------------------------------------------------------*\
58 Class StandardChemistryModel Declaration
59\*---------------------------------------------------------------------------*/
60
61template<class ReactionThermo, class ThermoType>
63:
64 public BasicChemistryModel<ReactionThermo>,
65 public ODESystem
66{
67 // Private Member Functions
68
69 //- Solve the reaction system for the given time step
70 // of given type and return the characteristic time
71 template<class DeltaTType>
72 scalar solve(const DeltaTType& deltaT);
73
74 //- No copy construct
76 (
78 ) = delete;
79
80 //- No copy assignment
81 void operator=
82 (
84 ) = delete;
85
86
87protected:
89 typedef ThermoType thermoType;
90
91
92 // Protected data
93
94 //- Reference to the field of specie mass fractions
96
97 //- Reactions
99
100 //- Thermodynamic data of the species
102
103 //- Number of species
104 label nSpecie_;
105
106 //- Number of reactions
107 label nReaction_;
108
109 //- Temperature below which the reaction rates are assumed 0
110 scalar Treact_;
111
112 //- List of reaction rate per specie [kg/m3/s]
114
115 //- Temporary concentration field
116 mutable scalarField c_;
117
118 //- Temporary rate-of-change of concentration field
119 mutable scalarField dcdt_;
120
121
122 // Protected Member Functions
123
124 //- Write access to chemical source terms
125 // (e.g. for multi-chemistry model)
127
128
129public:
130
131 //- Runtime type information
132 TypeName("standard");
133
134
135 // Constructors
136
137 //- Construct from thermo
138 StandardChemistryModel(ReactionThermo& thermo);
139
140
141 //- Destructor
142 virtual ~StandardChemistryModel();
143
144
145 // Member Functions
146
147 //- The reactions
148 inline const PtrList<Reaction<ThermoType>>& reactions() const;
149
150 //- Thermodynamic data of the species
151 inline const PtrList<ThermoType>& specieThermo() const;
152
153 //- The number of species
154 virtual inline label nSpecie() const;
155
156 //- The number of reactions
157 virtual inline label nReaction() const;
158
159 //- Temperature below which the reaction rates are assumed 0
160 inline scalar Treact() const;
161
162 //- Temperature below which the reaction rates are assumed 0
163 inline scalar& Treact();
164
165 //- dc/dt = omega, rate of change in concentration, for each species
166 virtual void omega
167 (
168 const scalarField& c,
169 const scalar T,
170 const scalar p,
171 scalarField& dcdt
172 ) const;
173
174 //- Return the reaction rate for reaction r and the reference
175 // species and characteristic times
176 virtual scalar omega
177 (
178 const Reaction<ThermoType>& r,
179 const scalarField& c,
180 const scalar T,
181 const scalar p,
182 scalar& pf,
183 scalar& cf,
184 label& lRef,
185 scalar& pr,
186 scalar& cr,
187 label& rRef
188 ) const;
189
190
191 //- Return the reaction rate for iReaction and the reference
192 // species and characteristic times
193 virtual scalar omegaI
194 (
195 label iReaction,
196 const scalarField& c,
197 const scalar T,
198 const scalar p,
199 scalar& pf,
200 scalar& cf,
201 label& lRef,
202 scalar& pr,
203 scalar& cr,
204 label& rRef
205 ) const;
206
207 //- Calculates the reaction rates
208 virtual void calculate();
209
210
211 // Chemistry model functions (overriding abstract functions in
212 // basicChemistryModel.H)
213
214 //- Return const access to the chemical source terms for specie, i
215 inline const volScalarField::Internal& RR
216 (
217 const label i
218 ) const;
219
220 //- Return non const access to chemical source terms [kg/m3/s]
222 (
223 const label i
224 );
225
226 //- Return reaction rate of the speciei in reactionI
228 (
229 const label reactionI,
230 const label speciei
231 ) const;
232
233 //- Solve the reaction system for the given time step
234 // and return the characteristic time
235 virtual scalar solve(const scalar deltaT);
236
237 //- Solve the reaction system for the given time step
238 // and return the characteristic time
239 virtual scalar solve(const scalarField& deltaT);
240
241 //- Return the chemical time scale
242 virtual tmp<volScalarField> tc() const;
243
244 //- Return the heat release rate [kg/m/s3]
245 virtual tmp<volScalarField> Qdot() const;
246
247
248 // ODE functions (overriding abstract functions in ODE.H)
249
250 //- Number of ODE's to solve
251 inline virtual label nEqns() const;
252
253 virtual void derivatives
254 (
255 const scalar t,
256 const scalarField& c,
257 scalarField& dcdt
258 ) const;
259
260 virtual void jacobian
261 (
262 const scalar t,
263 const scalarField& c,
264 scalarField& dcdt,
266 ) const;
268 virtual void solve
269 (
270 scalarField &c,
271 scalar& T,
272 scalar& p,
273 scalar& deltaT,
274 scalar& subDeltaT
275 ) const = 0;
276};
277
278
279// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
280
281} // End namespace Foam
282
283// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
284
286
287// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
288
289#ifdef NoRepository
290 #include "StandardChemistryModel.C"
291#endif
292
293// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
294
295#endif
296
297// ************************************************************************* //
Basic chemistry model templated on thermodynamics.
Abstract base class for the systems of ordinary differential equations.
Definition: ODESystem.H:50
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: PtrList.H:73
Simple extension of ReactionThermo to handle reaction kinetics in addition to the equilibrium thermod...
Definition: Reaction.H:73
Extends base chemistry model by adding a thermo package, and ODE functions. Introduces chemistry equa...
scalarField c_
Temporary concentration field.
virtual tmp< volScalarField > tc() const
Return the chemical time scale.
scalarField dcdt_
Temporary rate-of-change of concentration field.
const PtrList< ThermoType > & specieThermo() const
Thermodynamic data of the species.
virtual tmp< volScalarField::Internal > calculateRR(const label reactionI, const label speciei) const
Return reaction rate of the speciei in reactionI.
PtrList< volScalarField::Internal > & RR()
Write access to chemical source terms.
const PtrList< ThermoType > & specieThermo_
Thermodynamic data of the species.
label nSpecie_
Number of species.
PtrList< volScalarField > & Y_
Reference to the field of specie mass fractions.
virtual void omega(const scalarField &c, const scalar T, const scalar p, scalarField &dcdt) const
dc/dt = omega, rate of change in concentration, for each species
virtual void derivatives(const scalar t, const scalarField &c, scalarField &dcdt) const
Calculate the derivatives in dydx.
virtual ~StandardChemistryModel()
Destructor.
const PtrList< Reaction< ThermoType > > & reactions_
Reactions.
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.
virtual label nReaction() const
The number of reactions.
virtual void solve(scalarField &c, scalar &T, scalar &p, scalar &deltaT, scalar &subDeltaT) const =0
scalar Treact() const
Temperature below which the reaction rates are assumed 0.
scalar Treact_
Temperature below which the reaction rates are assumed 0.
TypeName("standard")
Runtime type information.
virtual label nEqns() const
Number of ODE's to solve.
PtrList< volScalarField::Internal > RR_
List of reaction rate per specie [kg/m3/s].
virtual label nSpecie() const
The number of species.
label nReaction_
Number of reactions.
const PtrList< Reaction< ThermoType > > & reactions() const
The reactions.
virtual tmp< volScalarField > Qdot() const
Return the heat release rate [kg/m/s3].
virtual void jacobian(const scalar t, const scalarField &c, scalarField &dcdt, scalarSquareMatrix &dfdc) const
Calculate the Jacobian of the system.
virtual void calculate()
Calculates the reaction rates.
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
A class for managing temporary objects.
Definition: tmp.H:65
volScalarField & p
const volScalarField & T
Namespace for OpenFOAM.
CEqn solve()
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition: typeInfo.H:73