chemistryReductionMethod.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2021-2022 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  chemistryReductionMethods.C
33  chemistryReductionMethodNew.C
34  chemistryReductionMethodI.H
35 
36 \*---------------------------------------------------------------------------*/
37 
38 #ifndef chemistryReductionMethod_H
39 #define chemistryReductionMethod_H
40 
41 #include "IOdictionary.H"
42 #include "DynamicField.H"
43 #include "Switch.H"
44 #include "cpuTime.H"
45 #include "OFstream.H"
46 
47 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
48 
49 namespace Foam
50 {
51 
52 template<class ThermoType>
53 class chemistryModel;
54 
55 /*---------------------------------------------------------------------------*\
56  Class chemistryReductionMethod Declaration
57 \*---------------------------------------------------------------------------*/
58 
59 template<class ThermoType>
61 {
62 protected:
63 
64  // Protected Data
65 
66  //- Dictionary that store the algorithm data
67  const dictionary coeffsDict_;
68 
69  //- Reference to the chemistry model
71 
72  //- Total number of species
73  const label nSpecie_;
74 
75  //- Number of active species
77 
78  //- List of disabled reactions (disabled = true)
80 
81  //- List of active species (active = true)
83 
84 
85  //- Protected Member Functions
86 
87  //- Initialise reduction of the mechanism
88  void initReduceMechanism();
89 
90  //- End reduction of the mechanism
92 
93 
94 private:
95 
96  // Private Data
97 
98  //- Switch to select performance logging
99  Switch log_;
100 
101  //- Tolerance for the mechanism reduction algorithm
102  scalar tolerance_;
103 
104  //- CPU time for logging
105  cpuTime cpuTime_;
106 
107  //- ...
108  int64_t sumnActiveSpecies_;
109 
110  //- ...
111  int64_t sumn_;
112 
113  //- ...
114  scalar reduceMechCpuTime_;
115 
116  // Log file for the average time spent reducing the chemistry
117  autoPtr<OFstream> cpuReduceFile_;
118 
119  // Write average number of species
120  autoPtr<OFstream> nActiveSpeciesFile_;
121 
122 
123 public:
124 
125  //- Runtime type information
126  TypeName("chemistryReductionMethod");
127 
128 
129  // Declare runtime constructor selection table
131  (
132  autoPtr,
134  dictionary,
135  (
136  const IOdictionary& dict,
138  ),
139  (dict, chemistry)
140  );
141 
142 
143  // Constructors
144 
145  //- Construct from components
147  (
149  );
150 
151  //- Construct from components
153  (
154  const IOdictionary& dict,
156  );
157 
158 
159  // Selector
160 
162  (
163  const IOdictionary& dict,
165  );
166 
167 
168  //- Destructor
169  virtual ~chemistryReductionMethod();
170 
171 
172  // Member Functions
173 
174  //- Return mechanism reduction active
175  virtual bool active() const
176  {
177  return true;
178  }
179 
180  //- Return the number of species
181  inline label nSpecie();
182 
183  //- Return the number of active species
184  inline label nActiveSpecies() const;
185 
186  //- Return the tolerance
187  inline scalar tolerance() const;
188 
189  //- Return the active species
190  inline const List<bool>& activeSpecies() const;
191 
192  //- Return whether or not a reaction is disabled
193  inline bool reactionDisabled(const label i) const;
194 
195  //- Reduce the mechanism
196  virtual void reduceMechanism
197  (
198  const scalar p,
199  const scalar T,
200  const scalarField& c,
201  List<label>& ctos,
202  DynamicList<label>& stoc,
203  const label li
204  ) = 0;
205 
206  //- ...
207  virtual void update();
208 };
209 
210 
211 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
212 
213 } // End namespace Foam
214 
215 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
216 
218 
219 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
220 
221 #ifdef NoRepository
222  #include "chemistryReductionMethod.C"
224 #endif
225 
226 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
227 
228 #endif
229 
230 // ************************************************************************* //
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:57
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.H:61
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
Extends base chemistry model by adding a thermo package, and ODE functions. Introduces chemistry equa...
An abstract class for methods of chemical mechanism reduction.
chemistryModel< ThermoType > & chemistry_
Reference to the chemistry model.
static autoPtr< chemistryReductionMethod< ThermoType > > New(const IOdictionary &dict, chemistryModel< ThermoType > &chemistry)
const List< bool > & activeSpecies() const
Return the active species.
TypeName("chemistryReductionMethod")
Runtime type information.
virtual ~chemistryReductionMethod()
Destructor.
void initReduceMechanism()
Protected Member Functions.
const label nSpecie_
Total number of species.
List< bool > activeSpecies_
List of active species (active = true)
const dictionary coeffsDict_
Dictionary that store the algorithm data.
Field< bool > reactionsDisabled_
List of disabled reactions (disabled = true)
bool reactionDisabled(const label i) const
Return whether or not a reaction is disabled.
scalar tolerance() const
Return the tolerance.
label nSpecie()
Return the number of species.
declareRunTimeSelectionTable(autoPtr, chemistryReductionMethod, dictionary,(const IOdictionary &dict, chemistryModel< ThermoType > &chemistry),(dict, chemistry))
void endReduceMechanism(List< label > &ctos, DynamicList< label > &stoc)
End reduction of the mechanism.
label nActiveSpecies() const
Return the number of active species.
label nActiveSpecies_
Number of active species.
chemistryReductionMethod(chemistryModel< ThermoType > &chemistry)
Construct from components.
virtual bool active() const
Return mechanism reduction active.
virtual void reduceMechanism(const scalar p, const scalar T, const scalarField &c, List< label > &ctos, DynamicList< label > &stoc, const label li)=0
Reduce the mechanism.
Starts timing CPU usage and return elapsed time from start.
Definition: cpuTime.H:55
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
const dimensionedScalar c
Speed of light in a vacuum.
Namespace for OpenFOAM.
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
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
dictionary dict
volScalarField & p
basicChemistryModel & chemistry