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 | Copyright (C) 2016-2017 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 chemistryModel 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 "chemistryModel.H"
70 #include "OFstream.H"
71 
72 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
73 
74 namespace Foam
75 {
76 
77 /*---------------------------------------------------------------------------*\
78  Class TDACChemistryModel Declaration
79 \*---------------------------------------------------------------------------*/
80 
81 template<class CompType, class ThermoType>
82 class TDACChemistryModel
83 :
84  public chemistryModel<CompType, 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<CompType, ThermoType>> mechRed_;
101 
102  // Tabulation
103  autoPtr<chemistryTabulationMethod<CompType, ThermoType>> tabulation_;
104 
105  // Log file for the average time spent reducing the chemistry
106  autoPtr<OFstream> cpuReduceFile_;
107 
108  // Write average number of species
109  autoPtr<OFstream> nActiveSpeciesFile_;
110 
111  //- Log file for the average time spent adding tabulated data
112  autoPtr<OFstream> cpuAddFile_;
113 
114  //- Log file for the average time spent growing tabulated data
115  autoPtr<OFstream> cpuGrowFile_;
116 
117  //- Log file for the average time spent retrieving tabulated data
118  autoPtr<OFstream> cpuRetrieveFile_;
119 
120  //- Log file for average time spent solving the chemistry
121  autoPtr<OFstream> cpuSolveFile_;
122 
123  // Field containing information about tabulation:
124  // 0 -> add (direct integration)
125  // 1 -> grow
126  // 2 -> retrieve
127  volScalarField tabulationResults_;
128 
129 
130  // Private Member Functions
131 
132  //- Disallow copy constructor
133  TDACChemistryModel(const TDACChemistryModel&);
134 
135  //- Disallow default bitwise assignment
136  void operator=(const TDACChemistryModel&);
137 
138  //- Solve the reaction system for the given time step
139  // of given type and return the characteristic time
140  // Variable number of species added
141  template<class DeltaTType>
142  scalar solve(const DeltaTType& deltaT);
143 
144 
145 public:
146 
147  //- Runtime type information
148  TypeName("TDACChemistryModel");
149 
150 
151  // Constructors
152 
153  //- Construct from mesh
154  TDACChemistryModel
155  (
156  const fvMesh& mesh,
157  const word& phaseName
158  );
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 charateristic 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  // chemistryModel 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 chemistryModel to take into
217  // 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 autoPtr<chemistryReductionMethod<CompType, ThermoType>>&
277  mechRed();
280  {
281  return tabulationResults_;
282  }
283 
284  void setTabulationResultsAdd(const label celli);
285 
286  void setTabulationResultsGrow(const label celli);
287 
288  void setTabulationResultsRetrieve(const label celli);
289 
290  inline void resetTabulationResults();
291 };
292 
293 
294 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
295 
296 } // End namespace Foam
297 
298 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
299 
300 #include "TDACChemistryModelI.H"
301 
302 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
303 
304 #ifdef NoRepository
305  #include "TDACChemistryModel.C"
306 #endif
307 
308 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
309 
310 #endif
311 
312 // ************************************************************************* //
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)
PtrList< volScalarField > & Y()
void setTabulationResultsRetrieve(const label celli)
Field< label > & completeToSimplifiedIndex()
virtual ~TDACChemistryModel()
Destructor.
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
DynamicList< label > & simplifiedToCompleteIndex()
Field< bool > & reactionsDisabled()
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
void setNsDAC(const label newNsDAC)
dynamicFvMesh & mesh
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
label timeSteps() const
Return the number of chemistry evaluations (used by ISAT)
void setTabulationResultsAdd(const label celli)
void jacobian(const scalar t, const scalarField &c, scalarSquareMatrix &dfdc) const
Pure jacobian function for tabulation.
autoPtr< chemistryReductionMethod< CompType, ThermoType > > & mechRed()
List< List< specieElement > > & specieComp()
void setNSpecie(const label newNs)
virtual void derivatives(const scalar t, const scalarField &c, scalarField &dcdt) const
Calculate the derivatives in dydx.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
bool variableTimeStep() const
Return true if the time-step is variable and/or non-uniform.
bool active(const label i) const
const dimensionedScalar c
Speed of light in a vacuum.
TypeName("TDACChemistryModel")
Runtime type information.
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)
tmp< volScalarField > tabulationResults() const
SquareMatrix< scalar > scalarSquareMatrix
Namespace for OpenFOAM.