TDACChemistryModel.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-2017 OpenFOAM Foundation
9-------------------------------------------------------------------------------
10License
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
26Class
27 Foam::TDACChemistryModel
28
29Description
30 Extends StandardChemistryModel by adding the TDAC method.
31
32 References:
33 \verbatim
34 Contino, F., Jeanmart, H., Lucchini, T., & D’Errico, G. (2011).
35 Coupling of in situ adaptive tabulation and dynamic adaptive chemistry:
36 An effective method for solving combustion in engine simulations.
37 Proceedings of the Combustion Institute, 33(2), 3057-3064.
38
39 Contino, F., Lucchini, T., D'Errico, G., Duynslaegher, C.,
40 Dias, V., & Jeanmart, H. (2012).
41 Simulations of advanced combustion modes using detailed chemistry
42 combined with tabulation and mechanism reduction techniques.
43 SAE International Journal of Engines,
44 5(2012-01-0145), 185-196.
45
46 Contino, F., Foucher, F., Dagaut, P., Lucchini, T., D’Errico, G., &
47 Mounaïm-Rousselle, C. (2013).
48 Experimental and numerical analysis of nitric oxide effect on the
49 ignition of iso-octane in a single cylinder HCCI engine.
50 Combustion and Flame, 160(8), 1476-1483.
51
52 Contino, F., Masurier, J. B., Foucher, F., Lucchini, T., D’Errico, G., &
53 Dagaut, P. (2014).
54 CFD simulations using the TDAC method to model iso-octane combustion
55 for a large range of ozone seeding and temperature conditions
56 in a single cylinder HCCI engine.
57 Fuel, 137, 179-184.
58 \endverbatim
59
60SourceFiles
61 TDACChemistryModelI.H
62 TDACChemistryModel.C
63
64\*---------------------------------------------------------------------------*/
65
66#ifndef TDACChemistryModel_H
67#define TDACChemistryModel_H
68
72#include "OFstream.H"
73
74// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
75
76namespace Foam
77{
78
79/*---------------------------------------------------------------------------*\
80 Class TDACChemistryModel Declaration
81\*---------------------------------------------------------------------------*/
82
83template<class ReactionThermo, class ThermoType>
85:
86 public StandardChemistryModel<ReactionThermo, ThermoType>
87{
88 // Private Member Data
89
90 bool variableTimeStep_;
91
92 label timeSteps_;
93
94 // Mechanism reduction
95 label NsDAC_;
96 scalarField completeC_;
97 scalarField simplifiedC_;
98 Field<bool> reactionsDisabled_;
99 List<List<specieElement>> specieComp_;
100 Field<label> completeToSimplifiedIndex_;
101 DynamicList<label> simplifiedToCompleteIndex_;
103 mechRed_;
104
105 // Tabulation
107 tabulation_;
108
109 // Log file for the average time spent reducing the chemistry
110 autoPtr<OFstream> cpuReduceFile_;
111
112 // Write average number of species
113 autoPtr<OFstream> nActiveSpeciesFile_;
114
115 //- Log file for the average time spent adding tabulated data
116 autoPtr<OFstream> cpuAddFile_;
117
118 //- Log file for the average time spent growing tabulated data
119 autoPtr<OFstream> cpuGrowFile_;
120
121 //- Log file for the average time spent retrieving tabulated data
122 autoPtr<OFstream> cpuRetrieveFile_;
123
124 //- Log file for average time spent solving the chemistry
125 autoPtr<OFstream> cpuSolveFile_;
126
127 // Field containing information about tabulation:
128 // 0 -> add (direct integration)
129 // 1 -> grow
130 // 2 -> retrieve
131 volScalarField tabulationResults_;
132
133
134 // Private Member Functions
135
136 //- No copy construct
137 TDACChemistryModel(const TDACChemistryModel&) = delete;
138
139 //- No copy assignment
140 void operator=(const TDACChemistryModel&) = delete;
141
142 //- Solve the reaction system for the given time step
143 // of given type and return the characteristic time
144 // Variable number of species added
145 template<class DeltaTType>
146 scalar solve(const DeltaTType& deltaT);
147
148
149public:
150
151 //- Runtime type information
152 TypeName("TDAC");
153
154
155 // Constructors
156
157 //- Construct from thermo
158 TDACChemistryModel(ReactionThermo& thermo);
159
160
161 //- Destructor
162 virtual ~TDACChemistryModel();
163
164
165 // Member Functions
166
167 //- Return true if the time-step is variable and/or non-uniform
168 inline bool variableTimeStep() const;
169
170 //- Return the number of chemistry evaluations (used by ISAT)
171 inline label timeSteps() const;
172
173 //- Create and return a TDAC log file of the given name
174 inline autoPtr<OFstream> logFile(const word& name) const;
175
176 inline PtrList<volScalarField>& Y();
177
178 //- dc/dt = omega, rate of change in concentration, for each species
179 virtual void omega
180 (
181 const scalarField& c,
182 const scalar T,
183 const scalar p,
184 scalarField& dcdt
185 ) const;
186
187 //- Return the reaction rate for reaction r and the reference
188 // species and characteristic times
189 virtual scalar omega
190 (
191 const Reaction<ThermoType>& r,
192 const scalarField& c,
193 const scalar T,
194 const scalar p,
195 scalar& pf,
196 scalar& cf,
197 label& lRef,
198 scalar& pr,
199 scalar& cr,
200 label& rRef
201 ) const;
202
203
204 // Chemistry model functions (overriding functions in
205 // StandardChemistryModel to use the private solve function)
206
207 //- Solve the reaction system for the given time step
208 // and return the characteristic time
209 virtual scalar solve(const scalar deltaT);
210
211 //- Solve the reaction system for the given time step
212 // and return the characteristic time
213 virtual scalar solve(const scalarField& deltaT);
214
215
216 // ODE functions (overriding functions in StandardChemistryModel to take
217 // into account the variable number of species)
218
219 virtual void derivatives
220 (
221 const scalar t,
222 const scalarField& c,
223 scalarField& dcdt
224 ) const;
225
226 //- Pure jacobian function for tabulation
227 void jacobian
228 (
229 const scalar t,
230 const scalarField& c,
232 ) const;
233
234 virtual void jacobian
235 (
236 const scalar t,
237 const scalarField& c,
238 scalarField& dcdt,
240 ) const;
242 virtual void solve
243 (
244 scalarField& c,
245 scalar& T,
246 scalar& p,
247 scalar& deltaT,
248 scalar& subDeltaT
249 ) const = 0;
250
251
252 // Mechanism reduction access functions
253
254 inline void setNsDAC(const label newNsDAC);
255
256 inline void setNSpecie(const label newNs);
257
258 inline scalarField& completeC();
259
260 inline scalarField& simplifiedC();
261
263
264 inline bool active(const label i) const;
265
266 inline void setActive(const label i);
267
269
271
272 inline const Field<label>& completeToSimplifiedIndex() const;
273
275
276 inline
278 mechRed();
281 {
282 return tabulationResults_;
283 }
284
285 void setTabulationResultsAdd(const label celli);
286
287 void setTabulationResultsGrow(const label celli);
288
289 void setTabulationResultsRetrieve(const label celli);
290
291 inline void resetTabulationResults();
292};
293
294
295// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
296
297} // End namespace Foam
298
299// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
300
301#include "TDACChemistryModelI.H"
302
303// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
304
305#ifdef NoRepository
306 #include "TDACChemistryModel.C"
307#endif
308
309// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
310
311#endif
312
313// ************************************************************************* //
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:72
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: List.H:77
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: PtrList.H:73
Simple extension of ReactionThermo to handle reaction kinetics in addition to the equilibrium thermod...
Definition: Reaction.H:73
Extends base chemistry model by adding a thermo package, and ODE functions. Introduces chemistry equa...
Extends StandardChemistryModel by adding the TDAC method.
bool active(const label i) const
virtual ~TDACChemistryModel()
Destructor.
void setNsDAC(const label newNsDAC)
tmp< volScalarField > tabulationResults() const
TypeName("TDAC")
Runtime type information.
autoPtr< chemistryReductionMethod< ReactionThermo, ThermoType > > & mechRed()
autoPtr< OFstream > logFile(const word &name) const
Create and return a TDAC log file of the given name.
void setNSpecie(const label newNs)
virtual void omega(const scalarField &c, const scalar T, const scalar p, scalarField &dcdt) const
dc/dt = omega, rate of change in concentration, for each species
virtual void derivatives(const scalar t, const scalarField &c, scalarField &dcdt) const
Calculate the derivatives in dydx.
DynamicList< label > & simplifiedToCompleteIndex()
List< List< specieElement > > & specieComp()
virtual void solve(scalarField &c, scalar &T, scalar &p, scalar &deltaT, scalar &subDeltaT) const =0
void setTabulationResultsAdd(const label celli)
PtrList< volScalarField > & Y()
label timeSteps() const
Return the number of chemistry evaluations (used by ISAT)
Field< label > & completeToSimplifiedIndex()
bool variableTimeStep() const
Return true if the time-step is variable and/or non-uniform.
Field< bool > & reactionsDisabled()
void setActive(const label i)
void setTabulationResultsRetrieve(const label celli)
void jacobian(const scalar t, const scalarField &c, scalarSquareMatrix &dfdc) const
Pure jacobian function for tabulation.
void setTabulationResultsGrow(const label celli)
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
const word & name() const
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
A class for managing temporary objects.
Definition: tmp.H:65
A class for handling words, derived from Foam::string.
Definition: word.H:68
volScalarField & p
const volScalarField & T
Namespace for OpenFOAM.
CEqn solve()
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition: typeInfo.H:73