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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2016 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::chemistryReductionMethod
28 
29 Description
30  An abstract class for methods of chemical mechanism reduction
31 
32 SourceFiles
33  chemistryReductionMethod.C
34 
35 \*---------------------------------------------------------------------------*/
36 
37 #ifndef chemistryReductionMethod_H
38 #define chemistryReductionMethod_H
39 
40 #include "IOdictionary.H"
41 #include "Switch.H"
42 #include "scalarField.H"
43 
44 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
45 
46 namespace Foam
47 {
48 
49 template<class CompType, class ThermoType>
50 class TDACChemistryModel;
51 
52 /*---------------------------------------------------------------------------*\
53  Class chemistrySolver Declaration
54 \*---------------------------------------------------------------------------*/
55 
56 template<class CompType, class ThermoType>
58 {
59 
60 protected:
61 
63 
64  //- Dictionary that store the algorithm data
65  const dictionary coeffsDict_;
66 
67  //- Is mechanism reduction active?
69 
70  //- Switch to select performance logging
71  Switch log_;
72 
74 
75  //- List of active species (active = true)
77 
78  //- Number of active species
79  label NsSimp_;
80 
81  //- Number of species
82  const label nSpecie_;
83 
84  //- Tolerance for the mechanism reduction algorithm
85  scalar tolerance_;
86 
87 
88 public:
89 
90  //- Runtime type information
91  TypeName("chemistryReductionMethod");
92 
93 
94  // Declare runtime constructor selection table
96  (
97  autoPtr,
99  dictionary,
100  (
101  const IOdictionary& dict,
103  ),
104  (dict, chemistry)
105  );
106 
107 
108  // Constructors
109 
110  //- Construct from components
112  (
113  const IOdictionary& dict,
115  );
116 
117 
118  // Selector
119 
121  (
122  const IOdictionary& dict,
124  );
125 
126 
127  //- Destructor
128  virtual ~chemistryReductionMethod();
129 
130 
131  // Member Functions
132 
133  //- Is mechanism reduction active?
134  inline bool active() const;
135 
136  //- Is performance data logging enabled?
137  inline bool log() const;
138 
139  //- Return the active species
140  inline const List<bool>& activeSpecies() const;
141 
142  //- Return the number of active species
143  inline label NsSimp();
144 
145  //- Return the initial number of species
146  inline label nSpecie();
147 
148  //- Return the tolerance
149  inline scalar tolerance() const;
150 
151  //- Reduce the mechanism
152  virtual void reduceMechanism
153  (
154  const scalarField &c,
155  const scalar T,
156  const scalar p
157  ) = 0;
158 };
159 
160 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
161 
162 } // End namespace Foam
163 
164 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
165 
167 
168 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
169 
170 #ifdef NoRepository
171  #include "chemistryReductionMethod.C"
173 #endif
174 
175 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
176 
177 #endif
178 
179 // ************************************************************************* //
Extends StandardChemistryModel by adding the TDAC method.
dictionary dict
TDACChemistryModel< CompType, ThermoType > & chemistry_
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
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, any/none. Also accepts 0/1 as a string and shortcuts t/f, y/n.
Definition: Switch.H:77
virtual ~chemistryReductionMethod()
Destructor.
BasicChemistryModel< psiReactionThermo > & chemistry
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:50
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.
TypeName("chemistryReductionMethod")
Runtime type information.
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
An abstract class for methods of chemical mechanism reduction.
volScalarField & p
label NsSimp_
Number of active species.
Namespace for OpenFOAM.