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