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 -------------------------------------------------------------------------------
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::TDACChemistryModel
28 
29 Description
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 
60 SourceFiles
61  TDACChemistryModelI.H
62  TDACChemistryModel.C
63 
64 \*---------------------------------------------------------------------------*/
65 
66 #ifndef TDACChemistryModel_H
67 #define TDACChemistryModel_H
68 
69 #include "StandardChemistryModel.H"
72 #include "OFstream.H"
73 
74 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
75 
76 namespace Foam
77 {
78 
79 /*---------------------------------------------------------------------------*\
80  Class TDACChemistryModel Declaration
81 \*---------------------------------------------------------------------------*/
82 
83 template<class ReactionThermo, class ThermoType>
84 class TDACChemistryModel
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_;
102  autoPtr<chemistryReductionMethod<ReactionThermo, ThermoType>>
103  mechRed_;
104 
105  // Tabulation
106  autoPtr<chemistryTabulationMethod<ReactionThermo, ThermoType>>
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 
149 public:
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,
231  scalarSquareMatrix& dfdc
232  ) const;
233 
234  virtual void jacobian
235  (
236  const scalar t,
237  const scalarField& c,
238  scalarField& dcdt,
239  scalarSquareMatrix& dfdc
240  ) const;
241 
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 
262  inline Field<bool>& reactionsDisabled();
263 
264  inline bool active(const label i) const;
265 
266  inline void setActive(const label i);
267 
268  inline DynamicList<label>& simplifiedToCompleteIndex();
269 
270  inline Field<label>& completeToSimplifiedIndex();
271 
272  inline const Field<label>& completeToSimplifiedIndex() const;
273 
274  inline List<List<specieElement>>& specieComp();
275 
276  inline
277  autoPtr<chemistryReductionMethod<ReactionThermo, ThermoType>>&
278  mechRed();
279 
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 // ************************************************************************* //
Foam::TDACChemistryModel::completeC
scalarField & completeC()
Definition: TDACChemistryModelI.H:152
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
p
volScalarField & p
Definition: createFieldRefs.H:8
chemistryTabulationMethod.H
Foam::TDACChemistryModel::~TDACChemistryModel
virtual ~TDACChemistryModel()
Destructor.
Definition: TDACChemistryModel.C:148
Foam::TDACChemistryModel::timeSteps
label timeSteps() const
Return the number of chemistry evaluations (used by ISAT)
Definition: TDACChemistryModelI.H:40
Foam::TDACChemistryModel::tabulationResults
tmp< volScalarField > tabulationResults() const
Definition: TDACChemistryModel.H:279
Foam::TDACChemistryModel::derivatives
virtual void derivatives(const scalar t, const scalarField &c, scalarField &dcdt) const
Definition: TDACChemistryModel.C:319
TDACChemistryModelI.H
Foam::TDACChemistryModel::Y
PtrList< volScalarField > & Y()
Definition: TDACChemistryModelI.H:64
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::TDACChemistryModel::setNsDAC
void setNsDAC(const label newNsDAC)
Definition: TDACChemistryModelI.H:101
StandardChemistryModel.H
thermo
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
TDACChemistryModel.C
Foam::TDACChemistryModel::active
bool active(const label i) const
Definition: TDACChemistryModelI.H:91
Foam::TDACChemistryModel::reactionsDisabled
Field< bool > & reactionsDisabled()
Definition: TDACChemistryModelI.H:144
Foam::TDACChemistryModel::setTabulationResultsGrow
void setTabulationResultsGrow(const label celli)
Definition: TDACChemistryModel.C:932
Foam::TDACChemistryModel::completeToSimplifiedIndex
Field< label > & completeToSimplifiedIndex()
Definition: TDACChemistryModelI.H:127
OFstream.H
Foam::TDACChemistryModel::resetTabulationResults
void resetTabulationResults()
Definition: TDACChemistryModelI.H:176
Foam::TDACChemistryModel::simplifiedToCompleteIndex
DynamicList< label > & simplifiedToCompleteIndex()
Definition: TDACChemistryModelI.H:118
Foam::T
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Definition: FieldFieldFunctions.C:58
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
Foam::TDACChemistryModel::setActive
void setActive(const label i)
Definition: TDACChemistryModelI.H:81
Foam::scalarSquareMatrix
SquareMatrix< scalar > scalarSquareMatrix
Definition: scalarMatrices.H:57
Foam::TDACChemistryModel::setNSpecie
void setNSpecie(const label newNs)
Definition: TDACChemistryModelI.H:109
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::TDACChemistryModel::omega
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
Definition: TDACChemistryModel.C:156
Foam::TDACChemistryModel::logFile
autoPtr< OFstream > logFile(const word &name) const
Create and return a TDAC log file of the given name.
Definition: TDACChemistryModelI.H:49
Foam::TDACChemistryModel::variableTimeStep
bool variableTimeStep() const
Return true if the time-step is variable and/or non-uniform.
Definition: TDACChemistryModelI.H:32
Foam::TDACChemistryModel::simplifiedC
scalarField & simplifiedC()
Definition: TDACChemistryModelI.H:160
Foam::TDACChemistryModel::setTabulationResultsAdd
void setTabulationResultsAdd(const label celli)
Definition: TDACChemistryModel.C:922
Foam::TDACChemistryModel::jacobian
void jacobian(const scalar t, const scalarField &c, scalarSquareMatrix &dfdc) const
Pure jacobian function for tabulation.
Definition: TDACChemistryModel.C:403
Foam::TDACChemistryModel::mechRed
autoPtr< chemistryReductionMethod< ReactionThermo, ThermoType > > & mechRed()
Definition: TDACChemistryModelI.H:73
chemistryReductionMethod.H
Foam::constant::universal::c
const dimensionedScalar c
Speed of light in a vacuum.
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
Foam::TDACChemistryModel::specieComp
List< List< specieElement > > & specieComp()
Definition: TDACChemistryModelI.H:168
Foam::TDACChemistryModel::TypeName
TypeName("TDAC")
Runtime type information.
Foam::TDACChemistryModel::setTabulationResultsRetrieve
void setTabulationResultsRetrieve(const label celli)
Definition: TDACChemistryModel.C:941