chemistryReductionMethod.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) 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::chemistryReductionMethod
26 
27 Description
28  An abstract class for methods of chemical mechanism reduction
29 
30 SourceFiles
31  chemistryReductionMethod.C
32 
33 \*---------------------------------------------------------------------------*/
34 
35 #ifndef chemistryReductionMethod_H
36 #define chemistryReductionMethod_H
37 
38 #include "IOdictionary.H"
39 #include "Switch.H"
40 #include "scalarField.H"
41 
42 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
43 
44 namespace Foam
45 {
46 
47 template<class CompType, class ThermoType>
48 class TDACChemistryModel;
49 
50 /*---------------------------------------------------------------------------*\
51  Class chemistrySolver Declaration
52 \*---------------------------------------------------------------------------*/
53 
54 template<class CompType, class ThermoType>
56 {
57 
58 protected:
59 
60  const IOdictionary& dict_;
61 
62  //- Dictionary that store the algorithm data
63  const dictionary coeffsDict_;
64 
65  //- Is mechanism reduction active?
67 
68  //- Switch to select performance logging
69  Switch log_;
70 
72 
73  //- List of active species (active = true)
75 
76  //- Number of active species
77  label NsSimp_;
78 
79  //- Number of species
80  const label nSpecie_;
81 
82  //- Tolerance for the mechanism reduction algorithm
83  scalar tolerance_;
84 
85 
86 public:
87 
88  //- Runtime type information
89  TypeName("chemistryReductionMethod");
90 
91 
92  // Declare runtime constructor selection table
94  (
95  autoPtr,
97  dictionary,
98  (
99  const IOdictionary& dict,
101  ),
102  (dict, chemistry)
103  );
104 
105 
106  // Constructors
107 
108  //- Construct from components
110  (
111  const IOdictionary& dict,
113  );
114 
115 
116  // Selector
117 
119  (
120  const IOdictionary& dict,
122  );
123 
124 
125  //- Destructor
126  virtual ~chemistryReductionMethod();
127 
128 
129  // Member Functions
130 
131  //- Is mechanism reduction active?
132  inline bool active() const;
133 
134  //- Is performance data logging enabled?
135  inline bool log() const;
136 
137  //- Return the active species
138  inline const List<bool>& activeSpecies() const;
139 
140  //- Return the number of active species
141  inline label NsSimp();
142 
143  //- Return the initial number of species
144  inline label nSpecie();
145 
146  //- Return the tolerance
147  inline scalar tolerance() const;
148 
149  //- Reduce the mechanism
150  virtual void reduceMechanism
151  (
152  const scalarField &c,
153  const scalar T,
154  const scalar p
155  ) = 0;
156 };
157 
158 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
159 
160 } // End namespace Foam
161 
162 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
163 
165 
166 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
167 
168 #ifdef NoRepository
169  #include "chemistryReductionMethod.C"
171 #endif
172 
173 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
174 
175 #endif
176 
177 // ************************************************************************* //
Extends chemistryModel by adding the TDAC method.
dictionary dict
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
TDACChemistryModel< CompType, ThermoType > & chemistry_
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
const List< bool > & activeSpecies() const
Return the active species.
static autoPtr< chemistryReductionMethod< CompType, ThermoType > > New(const IOdictionary &dict, TDACChemistryModel< CompType, ThermoType > &chemistry)
scalar tolerance() const
Return the tolerance.
bool active() const
Is mechanism reduction active?
chemistryReductionMethod(const IOdictionary &dict, TDACChemistryModel< CompType, ThermoType > &chemistry)
Construct from components.
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
virtual ~chemistryReductionMethod()
Destructor.
List< bool > activeSpecies_
List of active species (active = true)
scalar tolerance_
Tolerance for the mechanism reduction algorithm.
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:52
const label nSpecie_
Number of species.
bool log() const
Is performance data logging enabled?
Switch active_
Is mechanism reduction active?
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Switch log_
Switch to select performance logging.
const dictionary coeffsDict_
Dictionary that store the algorithm data.
label NsSimp()
Return the number of active species.
virtual void reduceMechanism(const scalarField &c, const scalar T, const scalar p)=0
Reduce the mechanism.
declareRunTimeSelectionTable(autoPtr, chemistryReductionMethod, dictionary,(const IOdictionary &dict, TDACChemistryModel< CompType, ThermoType > &chemistry),(dict, chemistry))
const dimensionedScalar c
Speed of light in a vacuum.
label nSpecie()
Return the initial number of species.
psiChemistryModel & chemistry
Definition: createFields.H:29
TypeName("chemistryReductionMethod")
Runtime type information.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
An abstract class for methods of chemical mechanism reduction.
volScalarField & p
label NsSimp_
Number of active species.
Namespace for OpenFOAM.