basicChemistryModel.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) 2011-2021 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::basicChemistryModel
26 
27 Description
28  Base class for chemistry models
29 
30 SourceFiles
31  basicChemistryModelI.H
32  basicChemistryModel.C
33 
34 \*---------------------------------------------------------------------------*/
35 
36 #ifndef basicChemistryModel_H
37 #define basicChemistryModel_H
38 
39 #include "fluidReactionThermo.H"
40 
41 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
42 
43 namespace Foam
44 {
45 
46 // Forward declaration of classes
47 class fvMesh;
48 
49 /*---------------------------------------------------------------------------*\
50  class basicChemistryModel Declaration
51 \*---------------------------------------------------------------------------*/
52 
54 :
55  public IOdictionary
56 {
57 public:
58 
59  // Public Enumerations
60 
61  //- Enumeration for the type of Jacobian to be calculated by the
62  // chemistry model
63  enum class jacobianType
64  {
65  fast,
66  exact
67  };
68 
69  //- Jacobian type names
71 
72 
73 protected:
74 
75  // Protected data
76 
77  //- Reference to the mesh database
78  const fvMesh& mesh_;
79 
80  //- Reference to the thermo
82 
83  //- Chemistry activation switch
85 
86  //- Initial chemical time step
87  const scalar deltaTChemIni_;
88 
89  //- Maximum chemical time step
90  const scalar deltaTChemMax_;
91 
92  //- Latest estimation of integration step
94 
95 
96  // Protected Member Functions
97 
98  //- Correct function - updates due to mesh changes
99  void correct();
100 
101 
102 public:
103 
104  //- Runtime type information
105  TypeName("basicChemistryModel");
106 
107 
108  //- Declare run-time constructor selection tables
110  (
111  autoPtr,
113  thermo,
114  (const fluidReactionThermo& thermo),
115  (thermo)
116  );
117 
118 
119  // Constructors
120 
121  //- Construct from thermo
123 
124  //- Disallow default bitwise copy construction
125  basicChemistryModel(const basicChemistryModel&) = delete;
126 
127 
128  // Selectors
129 
130  //- Select based on fluid reaction thermo
132  (
133  const fluidReactionThermo& thermo
134  );
135 
136 
137  //- Destructor
138  virtual ~basicChemistryModel();
139 
140 
141  // Member Functions
142 
143  //- Return const access to the mesh database
144  inline const fvMesh& mesh() const;
145 
146  //- Return const access to the thermo
147  inline const fluidReactionThermo& thermo() const;
148 
149  //- Chemistry activation switch
150  inline Switch chemistry() const;
151 
152  //- The number of species
153  virtual label nSpecie() const = 0;
154 
155  //- The number of reactions
156  virtual label nReaction() const = 0;
157 
158  //- Return the latest estimation of integration step
159  inline const volScalarField::Internal& deltaTChem() const;
160 
161 
162  // Functions to be derived in derived classes
163 
164  // Fields
165 
166  //- Return const access to chemical source terms [kg/m^3/s]
167  virtual const volScalarField::Internal& RR
168  (
169  const label i
170  ) const = 0;
171 
172  //- Return access to chemical source terms [kg/m^3/s]
174  (
175  const label i
176  ) = 0;
177 
178  //- Return reaction rate of the speciei in reactioni
180  (
181  const label reactioni,
182  const label speciei
183  ) const = 0;
184 
185 
186  // Chemistry solution
187 
188  //- Calculates the reaction rates
189  virtual void calculate() = 0;
190 
191  //- Solve the reaction system for the given time step
192  // and return the characteristic time
193  virtual scalar solve(const scalar deltaT) = 0;
194 
195  //- Solve the reaction system for the given time step
196  // and return the characteristic time
197  virtual scalar solve(const scalarField& deltaT) = 0;
198 
199  //- Return the chemical time scale
200  virtual tmp<volScalarField> tc() const = 0;
201 
202  //- Return the heat release rate [kg/m/s^3]
203  virtual tmp<volScalarField> Qdot() const = 0;
204 
205 
206  // Member Operators
207 
208  //- Disallow default bitwise assignment
209  void operator=(const basicChemistryModel&) = delete;
210 };
211 
212 
213 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
214 
215 } // End namespace Foam
216 
217 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
218 
219 #include "basicChemistryModelI.H"
220 
221 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
222 
223 #endif
224 
225 // ************************************************************************* //
const fvMesh & mesh() const
Return const access to the mesh database.
virtual label nSpecie() const =0
The number of species.
const scalar deltaTChemMax_
Maximum chemical time step.
volScalarField::Internal deltaTChem_
Latest estimation of integration step.
declareRunTimeSelectionTable(autoPtr, basicChemistryModel, thermo,(const fluidReactionThermo &thermo),(thermo))
Declare run-time constructor selection tables.
jacobianType
Enumeration for the type of Jacobian to be calculated by the.
const fvMesh & mesh_
Reference to the mesh database.
A simple wrapper around bool so that it can be read as a word: true/false, on/off, yes/no, y/n, t/f, or none/any.
Definition: Switch.H:60
Base-class for multi-component fluid thermodynamic properties.
const fluidReactionThermo & thermo_
Reference to the thermo.
Switch chemistry() const
Chemistry activation switch.
void correct()
Correct function - updates due to mesh changes.
void operator=(const basicChemistryModel &)=delete
Disallow default bitwise assignment.
virtual tmp< volScalarField > tc() const =0
Return the chemical time scale.
static const NamedEnum< jacobianType, 2 > jacobianTypeNames_
Jacobian type names.
basicChemistryModel(const fluidReactionThermo &thermo)
Construct from thermo.
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:53
virtual ~basicChemistryModel()
Destructor.
TypeName("basicChemistryModel")
Runtime type information.
static autoPtr< basicChemistryModel > New(const fluidReactionThermo &thermo)
Select based on fluid reaction thermo.
Switch chemistry_
Chemistry activation switch.
virtual void calculate()=0
Calculates the reaction rates.
const fluidReactionThermo & thermo() const
Return const access to the thermo.
const scalar deltaTChemIni_
Initial chemical time step.
virtual label nReaction() const =0
The number of reactions.
virtual const volScalarField::Internal & RR(const label i) const =0
Return const access to chemical source terms [kg/m^3/s].
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:95
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
virtual tmp< volScalarField::Internal > calculateRR(const label reactioni, const label speciei) const =0
Return reaction rate of the speciei in reactioni.
const volScalarField::Internal & deltaTChem() const
Return the latest estimation of integration step.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
A class for managing temporary objects.
Definition: PtrList.H:53
virtual tmp< volScalarField > Qdot() const =0
Return the heat release rate [kg/m/s^3].
virtual scalar solve(const scalar deltaT)=0
Solve the reaction system for the given time step.
Namespace for OpenFOAM.
Base class for chemistry models.