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-2018 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  //- Disallow copy constructor
135  TDACChemistryModel(const TDACChemistryModel&);
136 
137  //- Disallow default bitwise assignment
138  void operator=(const TDACChemistryModel&);
139 
140  //- Solve the reaction system for the given time step
141  // of given type and return the characteristic time
142  // Variable number of species added
143  template<class DeltaTType>
144  scalar solve(const DeltaTType& deltaT);
145 
146 
147 public:
148 
149  //- Runtime type information
150  TypeName("TDAC");
151 
152 
153  // Constructors
154 
155  //- Construct from thermo
156  TDACChemistryModel(ReactionThermo& thermo);
157 
158 
159  //- Destructor
160  virtual ~TDACChemistryModel();
161 
162 
163  // Member Functions
164 
165  //- Return true if the time-step is variable and/or non-uniform
166  inline bool variableTimeStep() const;
167 
168  //- Return the number of chemistry evaluations (used by ISAT)
169  inline label timeSteps() const;
170 
171  //- Create and return a TDAC log file of the given name
172  inline autoPtr<OFstream> logFile(const word& name) const;
173 
174  inline PtrList<volScalarField>& Y();
175 
176  //- dc/dt = omega, rate of change in concentration, for each species
177  virtual void omega
178  (
179  const scalarField& c,
180  const scalar T,
181  const scalar p,
182  scalarField& dcdt
183  ) const;
184 
185  //- Return the reaction rate for reaction r and the reference
186  // species and charateristic times
187  virtual scalar omega
188  (
189  const Reaction<ThermoType>& r,
190  const scalarField& c,
191  const scalar T,
192  const scalar p,
193  scalar& pf,
194  scalar& cf,
195  label& lRef,
196  scalar& pr,
197  scalar& cr,
198  label& rRef
199  ) const;
200 
201 
202  // Chemistry model functions (overriding functions in
203  // StandardChemistryModel to use the private solve function)
204 
205  //- Solve the reaction system for the given time step
206  // and return the characteristic time
207  virtual scalar solve(const scalar deltaT);
208 
209  //- Solve the reaction system for the given time step
210  // and return the characteristic time
211  virtual scalar solve(const scalarField& deltaT);
212 
213 
214  // ODE functions (overriding functions in StandardChemistryModel to take
215  // into account the variable number of species)
216 
217  virtual void derivatives
218  (
219  const scalar t,
220  const scalarField& c,
221  scalarField& dcdt
222  ) const;
223 
224  virtual void jacobian
225  (
226  const scalar t,
227  const scalarField& c,
228  scalarField& dcdt,
230  ) const;
231 
232  virtual void solve
233  (
234  scalarField& c,
235  scalar& T,
236  scalar& p,
237  scalar& deltaT,
238  scalar& subDeltaT
239  ) const = 0;
240 
241 
242  // Mechanism reduction access functions
243 
244  inline void setNsDAC(const label newNsDAC);
245 
246  inline void setNSpecie(const label newNs);
247 
248  inline scalarField& completeC();
249 
250  inline scalarField& simplifiedC();
251 
252  inline Field<bool>& reactionsDisabled();
253 
254  inline bool active(const label i) const;
255 
256  inline void setActive(const label i);
257 
258  inline DynamicList<label>& simplifiedToCompleteIndex();
259 
260  inline Field<label>& completeToSimplifiedIndex();
261 
262  inline const Field<label>& completeToSimplifiedIndex() const;
263 
264  inline List<List<specieElement>>& specieComp();
265 
266  inline
267  autoPtr<chemistryReductionMethod<ReactionThermo, ThermoType>>&
268  mechRed();
271  {
272  return tabulationResults_;
273  }
274 
275  void setTabulationResultsAdd(const label celli);
276 
277  void setTabulationResultsGrow(const label celli);
278 
279  void setTabulationResultsRetrieve(const label celli);
280 
281  inline void resetTabulationResults();
282 };
283 
284 
285 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
286 
287 } // End namespace Foam
288 
289 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
290 
291 #include "TDACChemistryModelI.H"
292 
293 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
294 
295 #ifdef NoRepository
296  #include "TDACChemistryModel.C"
297 #endif
298 
299 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
300 
301 #endif
302 
303 // ************************************************************************* //
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()
ReactionThermo & thermo()
Return access to the thermo package.
void setTabulationResultsRetrieve(const label celli)
autoPtr< chemistryReductionMethod< ReactionThermo, ThermoType > > & mechRed()
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
TypeName("TDAC")
Runtime type information.
DynamicList< label > & simplifiedToCompleteIndex()
Field< bool > & reactionsDisabled()
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
void setNsDAC(const label newNsDAC)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
virtual void jacobian(const scalar t, const scalarField &c, scalarField &dcdt, scalarSquareMatrix &J) const
Calculate the Jacobian of the system.
label timeSteps() const
Return the number of chemistry evaluations (used by ISAT)
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)
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)
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.
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.