solidChemistryModel.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::solidChemistryModel
26 
27 Description
28  Extends base solid chemistry model by adding a thermo package, and ODE
29  functions.
30  Introduces chemistry equation system and evaluation of chemical source
31  terms.
32 
33 SourceFiles
34  solidChemistryModelI.H
35  solidChemistryModel.C
36 
37 \*---------------------------------------------------------------------------*/
38 
39 #ifndef solidChemistryModel_H
40 #define solidChemistryModel_H
41 
42 #include "Reaction.H"
43 #include "ODESystem.H"
44 #include "volFields.H"
45 #include "simpleMatrix.H"
46 
47 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
48 
49 namespace Foam
50 {
51 
52 // Forward declaration of classes
53 class fvMesh;
54 
55 /*---------------------------------------------------------------------------*\
56  Class solidChemistryModel Declaration
57 \*---------------------------------------------------------------------------*/
58 
59 template<class CompType, class SolidThermo>
61 :
62  public CompType,
63  public ODESystem
64 {
65 protected:
66 
67  //- Reference to solid mass fractions
69 
70  //- Reactions
72 
73  //- Thermodynamic data of solids
75 
76  //- Number of solid components
78 
79  //- Number of solid reactions
81 
82  //- List of reaction rate per solid [kg/m^3/s]
84 
85  //- List of active reacting cells
87 
88 
89  // Protected Member Functions
90 
91  //- Write access to source terms for solids
93 
94  //- Set reacting status of cell, celli
95  void setCellReacting(const label celli, const bool active);
96 
97 
98 public:
99 
100  //- Runtime type information
101  TypeName("solidChemistryModel");
102 
103 
104  // Constructors
105 
106  //- Construct from thermo
107  solidChemistryModel(typename CompType::reactionThermo& thermo);
108 
109  //- Disallow default bitwise copy construction
111 
112 
113  //- Destructor
114  virtual ~solidChemistryModel();
115 
116 
117  // Member Functions
118 
119  //- The reactions
120  inline const PtrList<Reaction<SolidThermo>>& reactions() const;
121 
122  //- The number of reactions
123  inline label nReaction() const;
124 
125 
126  //- dc/dt = omega, rate of change in concentration, for each species
127  virtual scalarField omega
128  (
129  const scalarField& c,
130  const scalar T,
131  const scalar p,
132  const bool updateC0 = false
133  ) const = 0;
134 
135  //- Return the reaction rate for reaction r and the reference
136  // species and charateristic times
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 = 0;
150 
151 
152  //- Return the reaction rate for iReaction and the reference
153  // species and charateristic times
154  virtual scalar omegaI
155  (
156  label iReaction,
157  const scalarField& c,
158  const scalar T,
159  const scalar p,
160  scalar& pf,
161  scalar& cf,
162  label& lRef,
163  scalar& pr,
164  scalar& cr,
165  label& rRef
166  ) const = 0;
167 
168  //- Calculates the reaction rates
169  virtual void calculate() = 0;
170 
171 
172  // Solid Chemistry model functions
173 
174  //- Return const access to the chemical source terms for solids
175  inline const volScalarField::Internal& RRs
176  (
177  const label i
178  ) const;
179 
180  //- Return total solid source term
181  inline tmp<volScalarField::Internal> RRs() const;
182 
183  //- Solve the reaction system for the given time step
184  // and return the characteristic time
185  virtual scalar solve(const scalar deltaT) = 0;
186 
187  //- Solve the reaction system for the given time step
188  // and return the characteristic time
189  virtual scalar solve(const scalarField& deltaT);
190 
191  //- Return the chemical time scale
192  virtual tmp<volScalarField> tc() const;
193 
194  //- Return the heat release rate [kg/m/s^3]
195  virtual tmp<volScalarField> Qdot() const;
196 
197 
198  // ODE functions (overriding abstract functions in ODE.H)
199 
200  //- Number of ODE's to solve
201  virtual label nEqns() const = 0;
202 
203  virtual void derivatives
204  (
205  const scalar t,
206  const scalarField& c,
207  scalarField& dcdt
208  ) const = 0;
209 
210  virtual void jacobian
211  (
212  const scalar t,
213  const scalarField& c,
214  scalarField& dcdt,
215  scalarSquareMatrix& dfdc
216  ) const = 0;
217 
218  virtual void solve
219  (
220  scalarField &c,
221  scalar& T,
222  scalar& p,
223  scalar& deltaT,
224  scalar& subDeltaT
225  ) const = 0;
226 
227 
228  // Member Operators
229 
230  //- Disallow default bitwise assignment
231  void operator=(const solidChemistryModel&) = delete;
232 };
233 
234 
235 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
236 
237 } // End namespace Foam
238 
239 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
240 
241  #include "solidChemistryModelI.H"
242 
243 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
244 
245 #ifdef NoRepository
246  #include "solidChemistryModel.C"
247 #endif
248 
249 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
250 
251 #endif
252 
253 // ************************************************************************* //
List< bool > reactingCells_
List of active reacting cells.
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 =0
Return the reaction rate for iReaction and the reference.
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
Abstract base class for the systems of ordinary differential equations.
Definition: ODESystem.H:46
void setCellReacting(const label celli, const bool active)
Set reacting status of cell, celli.
PtrList< volScalarField::Internal > & RRs()
Write access to source terms for solids.
const PtrList< Reaction< SolidThermo > > & reactions() const
The reactions.
label nReaction_
Number of solid reactions.
rhoReactionThermo & thermo
Definition: createFields.H:28
Extends base solid chemistry model by adding a thermo package, and ODE functions. Introduces chemistr...
solidChemistryModel(typename CompType::reactionThermo &thermo)
Construct from thermo.
Simple extension of ReactionThermo to handle reaction kinetics in addition to the equilibrium thermod...
Definition: Reaction.H:56
virtual ~solidChemistryModel()
Destructor.
PtrList< volScalarField::Internal > RRs_
List of reaction rate per solid [kg/m^3/s].
virtual tmp< volScalarField > tc() const
Return the chemical time scale.
TypeName("solidChemistryModel")
Runtime type information.
const PtrList< SolidThermo > & solidThermo_
Thermodynamic data of solids.
virtual scalar solve(const scalar deltaT)=0
Solve the reaction system for the given time step.
label nSolids_
Number of solid components.
virtual scalarField omega(const scalarField &c, const scalar T, const scalar p, const bool updateC0=false) const =0
dc/dt = omega, rate of change in concentration, for each species
virtual void jacobian(const scalar t, const scalarField &c, scalarField &dcdt, scalarSquareMatrix &dfdc) const =0
Calculate the Jacobian of the system.
const PtrList< Reaction< SolidThermo > > & reactions_
Reactions.
virtual tmp< volScalarField > Qdot() const
Return the heat release rate [kg/m/s^3].
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
PtrList< volScalarField > & Ys_
Reference to solid mass fractions.
void operator=(const solidChemistryModel &)=delete
Disallow default bitwise assignment.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:70
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...
label nReaction() const
The number of reactions.
virtual void calculate()=0
Calculates the reaction rates.
virtual void derivatives(const scalar t, const scalarField &c, scalarField &dcdt) const =0
Calculate the derivatives in dydx.
volScalarField & p
A class for managing temporary objects.
Definition: PtrList.H:53
virtual label nEqns() const =0
Number of ODE&#39;s to solve.
Namespace for OpenFOAM.