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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2018 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 Class
27  Foam::basicChemistryModel
28 
29 Description
30  Base class for chemistry models
31 
32 SourceFiles
33  basicChemistryModelI.H
34  basicChemistryModel.C
35 
36 \*---------------------------------------------------------------------------*/
37 
38 #ifndef basicChemistryModel_H
39 #define basicChemistryModel_H
40 
41 #include "IOdictionary.H"
42 #include "Switch.H"
43 #include "scalarField.H"
44 #include "volFields.H"
45 #include "basicThermo.H"
46 
47 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
48 
49 namespace Foam
50 {
51 
52 // Forward declaration of classes
53 class fvMesh;
54 
55 /*---------------------------------------------------------------------------*\
56  class basicChemistryModel Declaration
57 \*---------------------------------------------------------------------------*/
58 
60 :
61  public IOdictionary
62 {
63  // Private Member Functions
64 
65  //- Construct as copy (not implemented)
67 
68  //- No copy assignment
69  void operator=(const basicChemistryModel&) = delete;
70 
71 
72 protected:
73 
74  // Protected data
75 
76  //- Reference to the mesh database
77  const fvMesh& mesh_;
78 
79  //- Chemistry activation switch
81 
82  //- Initial chemical time step
83  const scalar deltaTChemIni_;
84 
85  //- Maximum chemical time step
86  const scalar deltaTChemMax_;
87 
88  //- Latest estimation of integration step
90 
91 
92  // Protected Member Functions
93 
94  //- Return non-const access to the latest estimation of integration
95  // step, e.g. for multi-chemistry model
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  // Constructors
109 
110  //- Construct from thermo
112 
113 
114  // Selectors
115 
116  //- Generic New for each of the related chemistry model
117  template<class ChemistryModel>
119  (
120  typename ChemistryModel::reactionThermo& thermo
121  );
122 
123 
124  //- Destructor
125  virtual ~basicChemistryModel();
126 
127 
128  // Member Functions
129 
130  //- Return const access to the mesh database
131  inline const fvMesh& mesh() const;
132 
133  //- Chemistry activation switch
134  inline Switch chemistry() const;
135 
136  //- The number of species
137  virtual label nSpecie() const = 0;
138 
139  //- The number of reactions
140  virtual label nReaction() const = 0;
141 
142  //- Return the latest estimation of integration step
143  inline const volScalarField::Internal& deltaTChem() const;
144 
145 
146  // Functions to be derived in derived classes
147 
148  // Fields
149 
150  //- Return const access to chemical source terms [kg/m3/s]
151  virtual const volScalarField::Internal& RR
152  (
153  const label i
154  ) const = 0;
155 
156  //- Return access to chemical source terms [kg/m3/s]
158  (
159  const label i
160  ) = 0;
161 
162  //- Return reaction rate of the speciei in reactioni
164  (
165  const label reactioni,
166  const label speciei
167  ) const = 0;
168 
169 
170  // Chemistry solution
171 
172  //- Calculates the reaction rates
173  virtual void calculate() = 0;
174 
175  //- Solve the reaction system for the given time step
176  // and return the characteristic time
177  virtual scalar solve(const scalar deltaT) = 0;
178 
179  //- Solve the reaction system for the given time step
180  // and return the characteristic time
181  virtual scalar solve(const scalarField& deltaT) = 0;
182 
183  //- Return the chemical time scale
184  virtual tmp<volScalarField> tc() const = 0;
185 
186  //- Return the heat release rate [kg/m/s3]
187  virtual tmp<volScalarField> Qdot() const = 0;
188 };
189 
190 
191 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
192 
193 } // End namespace Foam
194 
195 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
196 
197 #include "basicChemistryModelI.H"
198 
199 #ifdef NoRepository
201 #endif
202 
203 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
204 
205 #endif
206 
207 // ************************************************************************* //
volFields.H
Foam::IOdictionary
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:54
Foam::Switch
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.H:77
basicThermo.H
Foam::basicChemistryModel::calculate
virtual void calculate()=0
Calculates the reaction rates.
Foam::basicChemistryModel::deltaTChemIni_
const scalar deltaTChemIni_
Initial chemical time step.
Definition: basicChemistryModel.H:82
Foam::basicChemistryModel::New
static autoPtr< ChemistryModel > New(typename ChemistryModel::reactionThermo &thermo)
Generic New for each of the related chemistry model.
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::basicChemistryModel
Base class for chemistry models.
Definition: basicChemistryModel.H:58
scalarField.H
thermo
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
Foam::basicThermo
Abstract base-class for fluid and solid thermodynamic properties.
Definition: basicThermo.H:63
Foam::basicChemistryModel::chemistry
Switch chemistry() const
Chemistry activation switch.
Definition: basicChemistryModelI.H:36
Foam::basicChemistryModel::nReaction
virtual label nReaction() const =0
The number of reactions.
basicChemistryModelI.H
Foam::basicChemistryModel::calculateRR
virtual tmp< volScalarField::Internal > calculateRR(const label reactioni, const label speciei) const =0
Return reaction rate of the speciei in reactioni.
Foam::Field< scalar >
Foam::basicChemistryModel::mesh
const fvMesh & mesh() const
Return const access to the mesh database.
Definition: basicChemistryModelI.H:30
Foam::basicChemistryModel::Qdot
virtual tmp< volScalarField > Qdot() const =0
Return the heat release rate [kg/m/s3].
Switch.H
Foam::basicChemistryModel::solve
virtual scalar solve(const scalar deltaT)=0
Solve the reaction system for the given time step.
Foam::basicChemistryModel::nSpecie
virtual label nSpecie() const =0
The number of species.
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:85
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::basicChemistryModel::correct
void correct()
Correct function - updates due to mesh changes.
Definition: basicChemistryModel.C:42
Foam::basicChemistryModel::mesh_
const fvMesh & mesh_
Reference to the mesh database.
Definition: basicChemistryModel.H:76
Foam::basicChemistryModel::deltaTChem_
volScalarField::Internal deltaTChem_
Latest estimation of integration step.
Definition: basicChemistryModel.H:88
IOdictionary.H
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:53
Foam::basicChemistryModel::TypeName
TypeName("basicChemistryModel")
Runtime type information.
Foam::basicChemistryModel::~basicChemistryModel
virtual ~basicChemistryModel()
Destructor.
Definition: basicChemistryModel.C:83
Foam::basicChemistryModel::deltaTChem
volScalarField::Internal & deltaTChem()
Return non-const access to the latest estimation of integration.
Definition: basicChemistryModelI.H:50
Foam::basicChemistryModel::tc
virtual tmp< volScalarField > tc() const =0
Return the chemical time scale.
basicChemistryModelTemplates.C
Foam::basicChemistryModel::RR
virtual const volScalarField::Internal & RR(const label i) const =0
Return const access to chemical source terms [kg/m3/s].
Foam::basicChemistryModel::deltaTChemMax_
const scalar deltaTChemMax_
Maximum chemical time step.
Definition: basicChemistryModel.H:85
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::basicChemistryModel::chemistry_
Switch chemistry_
Chemistry activation switch.
Definition: basicChemistryModel.H:79