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