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