pyrolysisChemistryModel.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-2017 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 Class
27  Foam::pyrolysisChemistryModel
28 
29 Description
30  Pyrolysis chemistry model. It includes gas phase in the solid
31  reaction.
32 
33 SourceFiles
34  pyrolysisChemistryModelI.H
35  pyrolysisChemistryModel.C
36 
37 \*---------------------------------------------------------------------------*/
38 
39 #ifndef pyrolysisChemistryModel_H
40 #define pyrolysisChemistryModel_H
41 
42 #include "volFieldsFwd.H"
43 #include "DimensionedField.H"
44 #include "solidChemistryModel.H"
45 
46 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
47 
48 namespace Foam
49 {
50 
51 // Forward Declarations
52 class fvMesh;
53 
54 /*---------------------------------------------------------------------------*\
55  Class pyrolysisChemistryModel Declaration
56 \*---------------------------------------------------------------------------*/
57 
58 template<class CompType, class SolidThermo, class GasThermo>
60 :
61  public solidChemistryModel<CompType, SolidThermo>
62 {
63 protected:
64 
65  // Protected Data
66 
67  //- List of gas species present in reaction system
69 
70  //- Thermodynamic data of gases
72 
73  //- Number of gas species
74  label nGases_;
75 
76  //- Number of components being solved by ODE
77  label nSpecie_;
78 
79  //- List of reaction rate per gas [kg/m3/s]
81 
82  //- List of accumulative solid concentrations
84 
85 
86  // Protected Member Functions
87 
88  //- Write access to source terms for gases
90 
91  //- No copy assignment
92  void operator=(const pyrolysisChemistryModel&) = delete;
93 
94 
95 private:
96 
97  //- Cell counter
98  label cellCounter_;
99 
100 
101 public:
102 
103  //- Runtime type information
104  TypeName("pyrolysis");
105 
106 
107  // Constructors
108 
109  //- Construct from thermo
110  pyrolysisChemistryModel(typename CompType::reactionThermo& thermo);
111 
112 
113  //- Destructor
114  virtual ~pyrolysisChemistryModel();
115 
116 
117  // Member Functions
118 
119  //- Thermodynamic data of gases
120  inline const PtrList<GasThermo>& gasThermo() const;
121 
122  //- Gases table
123  inline const speciesTable& gasTable() const;
124 
125  //- The number of solids
126  inline label nSpecie() const;
127 
128  //- The number of solids
129  inline label nGases() const;
130 
131 
132  //- dc/dt = omega, rate of change in concentration, for each species
133  virtual scalarField omega
134  (
135  const scalarField& c,
136  const scalar T,
137  const scalar p,
138  const bool updateC0 = false
139  ) const;
140 
141  //- Return the reaction rate for reaction r
142  // NOTE: Currently does not calculate reference specie and
143  // characteristic times (pf, cf,l Ref, etc.)
144  virtual scalar omega
145  (
146  const Reaction<SolidThermo>& r,
147  const scalarField& c,
148  const scalar T,
149  const scalar p,
150  scalar& pf,
151  scalar& cf,
152  label& lRef,
153  scalar& pr,
154  scalar& cr,
155  label& rRef
156  ) const;
157 
158 
159  //- Return the reaction rate for iReaction
160  // NOTE: Currently does not calculate reference specie and
161  // characteristic times (pf, cf,l Ref, etc.)
162  virtual scalar omegaI
163  (
164  label iReaction,
165  const scalarField& c,
166  const scalar T,
167  const scalar p,
168  scalar& pf,
169  scalar& cf,
170  label& lRef,
171  scalar& pr,
172  scalar& cr,
173  label& rRef
174  ) const;
175 
176 
177  //- Calculates the reaction rates
178  virtual void calculate();
179 
180 
181  // Chemistry model functions
182 
183  //- Return const access to the chemical source terms for gases
184  inline const volScalarField::Internal& RRg
185  (
186  const label i
187  ) const;
188 
189  //- Return total gas source term
190  inline tmp<volScalarField::Internal> RRg() const;
191 
192  //- Return sensible enthalpy for gas i [J/Kg]
193  virtual tmp<volScalarField> gasHs
194  (
195  const volScalarField& p,
196  const volScalarField& T,
197  const label i
198  ) const;
199 
200  //- Solve the reaction system for the given time step
201  // and return the characteristic time
202  virtual scalar solve(const scalar deltaT);
203 
204 
205  // ODE functions (overriding abstract functions in ODE.H)
206 
207  //- Number of ODE's to solve
208  virtual label nEqns() const;
209 
210  virtual void derivatives
211  (
212  const scalar t,
213  const scalarField& c,
214  scalarField& dcdt
215  ) const;
216 
217  virtual void jacobian
218  (
219  const scalar t,
220  const scalarField& c,
221  scalarField& dcdt,
222  scalarSquareMatrix& dfdc
223  ) const;
224 
225  virtual void solve
226  (
227  scalarField &c,
228  scalar& T,
229  scalar& p,
230  scalar& deltaT,
231  scalar& subDeltaT
232  ) const;
233 };
234 
235 
236 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
237 
238 } // End namespace Foam
239 
240 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
241 
242  #include "pyrolysisChemistryModelI.H"
243 
244 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
245 
246 #ifdef NoRepository
247  #include "pyrolysisChemistryModel.C"
248 #endif
249 
250 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
251 
252 #endif
253 
254 // ************************************************************************* //
virtual ~pyrolysisChemistryModel()
Destructor.
PtrList< GasThermo > gasThermo_
Thermodynamic data of gases.
const PtrList< GasThermo > & gasThermo() const
Thermodynamic data of gases.
Forwards and collection of common volume field types.
TypeName("pyrolysis")
Runtime type information.
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.
label nSpecie_
Number of components being solved by ODE.
virtual void derivatives(const scalar t, const scalarField &c, scalarField &dcdt) const
Calculate the derivatives in dydx.
PtrList< volScalarField::Internal > & RRg()
Write access to source terms for gases.
speciesTable pyrolisisGases_
List of gas species present in reaction system.
Pyrolysis chemistry model. It includes gas phase in the solid reaction.
Extends base solid chemistry model by adding a thermo package, and ODE functions. ...
Simple extension of ReactionThermo to handle reaction kinetics in addition to the equilibrium thermod...
Definition: Reaction.H:55
label nSpecie() const
The number of solids.
virtual scalarField omega(const scalarField &c, const scalar T, const scalar p, const bool updateC0=false) const
dc/dt = omega, rate of change in concentration, for each species
virtual label nEqns() const
Number of ODE&#39;s to solve.
virtual tmp< volScalarField > gasHs(const volScalarField &p, const volScalarField &T, const label i) const
Return sensible enthalpy for gas i [J/Kg].
label nGases() const
The number of solids.
label nGases_
Number of gas species.
virtual scalar solve(const scalar deltaT)
Solve the reaction system for the given time step.
PtrList< volScalarField::Internal > RRg_
List of reaction rate per gas [kg/m3/s].
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
void operator=(const pyrolysisChemistryModel &)=delete
No copy assignment.
PtrList< volScalarField > Ys0_
List of accumulative solid concentrations.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
pyrolysisChemistryModel(typename CompType::reactionThermo &thermo)
Construct from thermo.
A wordList with hashed named lookup, which can be faster in some situations than using the normal lis...
const dimensionedScalar c
Speed of light in a vacuum.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: areaFieldsFwd.H:42
volScalarField & p
A class for managing temporary objects.
Definition: HashPtrTable.H:50
virtual void calculate()
Calculates the reaction rates.
virtual void jacobian(const scalar t, const scalarField &c, scalarField &dcdt, scalarSquareMatrix &dfdc) const
Calculate the Jacobian of the system.
const speciesTable & gasTable() const
Gases table.
Namespace for OpenFOAM.