solidChemistryModel.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  Copyright (C) 2016 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 Class
28  Foam::solidChemistryModel
29 
30 Description
31  Extends base solid chemistry model by adding a thermo package, and ODE
32  functions.
33 
34  Introduces chemistry equation system and evaluation of chemical source
35  terms.
36 
37 SourceFiles
38  solidChemistryModelI.H
39  solidChemistryModel.C
40 
41 \*---------------------------------------------------------------------------*/
42 
43 #ifndef solidChemistryModel_H
44 #define solidChemistryModel_H
45 
46 #include "Reaction.H"
47 #include "ODESystem.H"
48 #include "volFields.H"
49 #include "simpleMatrix.H"
50 
51 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
52 
53 namespace Foam
54 {
55 
56 // Forward declaration of classes
57 class fvMesh;
58 
59 /*---------------------------------------------------------------------------*\
60  Class solidChemistryModel Declaration
61 \*---------------------------------------------------------------------------*/
62 
63 template<class CompType, class SolidThermo>
65 :
66  public CompType,
67  public ODESystem
68 {
69  // Private Member Functions
70 
71  //- No copy construct
73 
74  //- No copy assignment
75  void operator=(const solidChemistryModel&) = delete;
76 
77 
78 protected:
79 
80  //- Reference to solid mass fractions
82 
83  //- Reactions
85 
86  //- Thermodynamic data of solids
88 
89  //- Number of solid components
90  label nSolids_;
91 
92  //- Number of solid reactions
93  label nReaction_;
94 
95  //- List of reaction rate per solid [kg/m3/s]
97 
98  //- List of active reacting cells
100 
101 
102  // Protected Member Functions
103 
104  //- Write access to source terms for solids
106 
107  //- Set reacting status of cell, celli
108  void setCellReacting(const label celli, const bool active);
109 
110 
111 public:
112 
113  //- Runtime type information
114  TypeName("solidChemistryModel");
115 
116 
117  // Constructors
118 
119  //- Construct from thermo
120  solidChemistryModel(typename CompType::reactionThermo& thermo);
121 
122 
123  //- Destructor
124  virtual ~solidChemistryModel();
125 
126 
127  // Member Functions
128 
129  //- The reactions
130  inline const PtrList<Reaction<SolidThermo>>& reactions() const;
131 
132  //- The number of reactions
133  inline label nReaction() const;
134 
135 
136  //- dc/dt = omega, rate of change in concentration, for each species
137  virtual scalarField omega
138  (
139  const scalarField& c,
140  const scalar T,
141  const scalar p,
142  const bool updateC0 = false
143  ) const = 0;
144 
145  //- Return the reaction rate for reaction r and the reference
146  // species and characteristic times
147  virtual scalar omega
148  (
149  const Reaction<SolidThermo>& r,
150  const scalarField& c,
151  const scalar T,
152  const scalar p,
153  scalar& pf,
154  scalar& cf,
155  label& lRef,
156  scalar& pr,
157  scalar& cr,
158  label& rRef
159  ) const = 0;
160 
161 
162  //- Return the reaction rate for iReaction and the reference
163  // species and characteristic times
164  virtual scalar omegaI
165  (
166  label iReaction,
167  const scalarField& c,
168  const scalar T,
169  const scalar p,
170  scalar& pf,
171  scalar& cf,
172  label& lRef,
173  scalar& pr,
174  scalar& cr,
175  label& rRef
176  ) const = 0;
177 
178  //- Calculates the reaction rates
179  virtual void calculate() = 0;
180 
181 
182  // Solid Chemistry model functions
183 
184  //- Return const access to the chemical source terms for solids
185  inline const volScalarField::Internal& RRs
186  (
187  const label i
188  ) const;
189 
190  //- Return total solid source term
191  inline tmp<volScalarField::Internal> RRs() const;
192 
193  //- Return net solid sensible enthalpy
195 
196  //- Solve the reaction system for the given time step
197  // and return the characteristic time
198  virtual scalar solve(const scalar deltaT) = 0;
199 
200  //- Solve the reaction system for the given time step
201  // and return the characteristic time
202  virtual scalar solve(const scalarField& deltaT);
203 
204  //- Return the chemical time scale
205  virtual tmp<volScalarField> tc() const;
206 
207  //- Return the heat release rate [kg/m/s3]
208  virtual tmp<volScalarField> Qdot() const;
209 
210 
211  // ODE functions (overriding abstract functions in ODE.H)
212 
213  //- Number of ODE's to solve
214  virtual label nEqns() const = 0;
215 
216  virtual void derivatives
217  (
218  const scalar t,
219  const scalarField& c,
220  scalarField& dcdt
221  ) const = 0;
222 
223  virtual void jacobian
224  (
225  const scalar t,
226  const scalarField& c,
227  scalarField& dcdt,
228  scalarSquareMatrix& dfdc
229  ) const = 0;
230 
231  virtual void solve
232  (
233  scalarField &c,
234  scalar& T,
235  scalar& p,
236  scalar& deltaT,
237  scalar& subDeltaT
238  ) const = 0;
239 };
240 
241 
242 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
243 
244 } // End namespace Foam
245 
246 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
247 
248 #include "solidChemistryModelI.H"
249 
250 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
251 
252 #ifdef NoRepository
253  #include "solidChemistryModel.C"
254 #endif
255 
256 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
257 
258 #endif
259 
260 // ************************************************************************* //
volFields.H
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::solidChemistryModel::reactingCells_
List< bool > reactingCells_
List of active reacting cells.
Definition: solidChemistryModel.H:98
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::solidChemistryModel::jacobian
virtual void jacobian(const scalar t, const scalarField &c, scalarField &dcdt, scalarSquareMatrix &dfdc) const =0
Calculate the Jacobian of the system.
Foam::solidChemistryModel::Ys_
PtrList< volScalarField > & Ys_
Reference to solid mass fractions.
Definition: solidChemistryModel.H:80
solidChemistryModelI.H
Foam::solidChemistryModel::nReaction
label nReaction() const
The number of reactions.
Definition: solidChemistryModelI.H:52
Foam::solidChemistryModel::~solidChemistryModel
virtual ~solidChemistryModel()
Destructor.
Definition: solidChemistryModel.C:89
thermo
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
ODESystem.H
Foam::solidChemistryModel::nSolids_
label nSolids_
Number of solid components.
Definition: solidChemistryModel.H:89
Foam::solidChemistryModel::TypeName
TypeName("solidChemistryModel")
Runtime type information.
Foam::solidChemistryModel::RRsHs
tmp< DimensionedField< scalar, volMesh > > RRsHs() const
Return net solid sensible enthalpy.
Definition: solidChemistryModelI.H:104
Foam::solidChemistryModel::reactions_
const PtrList< Reaction< SolidThermo > > & reactions_
Reactions.
Definition: solidChemistryModel.H:83
Foam::solidChemistryModel::calculate
virtual void calculate()=0
Calculates the reaction rates.
Foam::Field< scalar >
Foam::solidChemistryModel
Extends base solid chemistry model by adding a thermo package, and ODE functions.
Definition: solidChemistryModel.H:63
Foam::solidChemistryModel::setCellReacting
void setCellReacting(const label celli, const bool active)
Set reacting status of cell, celli.
Definition: solidChemistryModel.C:157
Foam::T
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Definition: FieldFieldFunctions.C:58
Foam::PtrList
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: List.H:59
Foam::solidChemistryModel::solidThermo_
const PtrList< SolidThermo > & solidThermo_
Thermodynamic data of solids.
Definition: solidChemistryModel.H:86
solidChemistryModel.C
Foam::solidChemistryModel::RRs
PtrList< volScalarField::Internal > & RRs()
Write access to source terms for solids.
Definition: solidChemistryModelI.H:35
Reaction.H
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::solidChemistryModel::omega
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
Foam::solidChemistryModel::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 =0
Return the reaction rate for iReaction and the reference.
Foam::solidChemistryModel::derivatives
virtual void derivatives(const scalar t, const scalarField &c, scalarField &dcdt) const =0
Calculate the derivatives in dydx.
Foam::SquareMatrix< scalar >
Foam::solidChemistryModel::Qdot
virtual tmp< volScalarField > Qdot() const
Return the heat release rate [kg/m/s3].
Definition: solidChemistryModel.C:117
Foam::solidChemistryModel::tc
virtual tmp< volScalarField > tc() const
Return the chemical time scale.
Definition: solidChemistryModel.C:108
Foam::ODESystem
Abstract base class for the systems of ordinary differential equations.
Definition: ODESystem.H:49
Foam::List< bool >
Foam::solidChemistryModel::RRs_
PtrList< volScalarField::Internal > RRs_
List of reaction rate per solid [kg/m3/s].
Definition: solidChemistryModel.H:95
Foam::solidChemistryModel::nEqns
virtual label nEqns() const =0
Number of ODE's to solve.
Foam::constant::universal::c
const dimensionedScalar c
Speed of light in a vacuum.
Foam::solidChemistryModel::nReaction_
label nReaction_
Number of solid reactions.
Definition: solidChemistryModel.H:92
simpleMatrix.H
Foam::Reaction
Simple extension of ReactionThermo to handle reaction kinetics in addition to the equilibrium thermod...
Definition: Reaction.H:59
Foam::solidChemistryModel::solve
virtual scalar solve(const scalar deltaT)=0
Solve the reaction system for the given time step.
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54
Foam::solidChemistryModel::reactions
const PtrList< Reaction< SolidThermo > > & reactions() const
The reactions.
Definition: solidChemistryModelI.H:43