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-2018 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 "IOdictionary.H"
40 #include "Switch.H"
41 #include "scalarField.H"
42 #include "volFields.H"
43 #include "basicThermo.H"
44 
45 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
46 
47 namespace Foam
48 {
49 
50 // Forward declaration of classes
51 class fvMesh;
52 
53 /*---------------------------------------------------------------------------*\
54  class basicChemistryModel Declaration
55 \*---------------------------------------------------------------------------*/
56 
58 :
59  public IOdictionary
60 {
61  // Private Member Functions
62 
63  //- Construct as copy (not implemented)
65 
66  //- Disallow default bitwise assignment
67  void operator=(const basicChemistryModel&);
68 
69 
70 protected:
71 
72  // Protected data
73 
74  //- Reference to the mesh database
75  const fvMesh& mesh_;
76 
77  //- Chemistry activation switch
79 
80  //- Initial chemical time step
81  const scalar deltaTChemIni_;
82 
83  //- Maximum chemical time step
84  const scalar deltaTChemMax_;
85 
86  //- Latest estimation of integration step
88 
89 
90  // Protected Member Functions
91 
92  //- Return non-const access to the latest estimation of integration
93  // step, e.g. for multi-chemistry model
95 
96  //- Correct function - updates due to mesh changes
97  void correct();
98 
99 
100 public:
101 
102  //- Runtime type information
103  TypeName("basicChemistryModel");
104 
105 
106  // Constructors
107 
108  //- Construct from thermo
110 
111 
112  // Selectors
113 
114  //- Generic New for each of the related chemistry model
115  template<class ChemistryModel>
117  (
118  typename ChemistryModel::reactionThermo& thermo
119  );
120 
121 
122  //- Destructor
123  virtual ~basicChemistryModel();
124 
125 
126  // Member Functions
127 
128  //- Return const access to the mesh database
129  inline const fvMesh& mesh() const;
130 
131  //- Chemistry activation switch
132  inline Switch chemistry() const;
133 
134  //- The number of species
135  virtual label nSpecie() const = 0;
136 
137  //- The number of reactions
138  virtual label nReaction() const = 0;
139 
140  //- Return the latest estimation of integration step
141  inline const volScalarField::Internal& deltaTChem() const;
142 
143 
144  // Functions to be derived in derived classes
145 
146  // Fields
147 
148  //- Return const access to chemical source terms [kg/m3/s]
149  virtual const volScalarField::Internal& RR
150  (
151  const label i
152  ) const = 0;
153 
154  //- Return access to chemical source terms [kg/m3/s]
156  (
157  const label i
158  ) = 0;
159 
160  //- Return reaction rate of the speciei in reactioni
162  (
163  const label reactioni,
164  const label speciei
165  ) const = 0;
166 
167 
168  // Chemistry solution
169 
170  //- Calculates the reaction rates
171  virtual void calculate() = 0;
172 
173  //- Solve the reaction system for the given time step
174  // and return the characteristic time
175  virtual scalar solve(const scalar deltaT) = 0;
176 
177  //- Solve the reaction system for the given time step
178  // and return the characteristic time
179  virtual scalar solve(const scalarField& deltaT) = 0;
180 
181  //- Return the chemical time scale
182  virtual tmp<volScalarField> tc() const = 0;
183 
184  //- Return the heat release rate [kg/m/s3]
185  virtual tmp<volScalarField> Qdot() const = 0;
186 };
187 
188 
189 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
190 
191 } // End namespace Foam
192 
193 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
194 
195 #include "basicChemistryModelI.H"
196 
197 #ifdef NoRepository
199 #endif
200 
201 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
202 
203 #endif
204 
205 // ************************************************************************* //
const fvMesh & mesh() const
Return const access to the mesh database.
static autoPtr< ChemistryModel > New(typename ChemistryModel::reactionThermo &thermo)
Generic New for each of the related chemistry model.
virtual label nSpecie() const =0
The number of species.
const scalar deltaTChemMax_
Maximum chemical time step.
volScalarField::Internal deltaTChem_
Latest estimation of integration step.
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 fluid and solid thermodynamic properties.
Definition: basicThermo.H:52
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.
Definition: Switch.H:60
Switch chemistry() const
Chemistry activation switch.
void correct()
Correct function - updates due to mesh changes.
virtual tmp< volScalarField > tc() const =0
Return the chemical time scale.
rhoReactionThermo & thermo
Definition: createFields.H:28
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:52
virtual ~basicChemistryModel()
Destructor.
TypeName("basicChemistryModel")
Runtime type information.
Switch chemistry_
Chemistry activation switch.
virtual void calculate()=0
Calculates the reaction rates.
volScalarField::Internal & deltaTChem()
Return non-const access to the latest estimation of integration.
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/m3/s].
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
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.
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/s3].
virtual scalar solve(const scalar deltaT)=0
Solve the reaction system for the given time step.
Namespace for OpenFOAM.
Base class for chemistry models.