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