ISAT.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::chemistryTabulationMethods::ISAT
28
29Description
30 Implementation of the ISAT (In-situ adaptive tabulation), for chemistry
31 calculation.
32
33 Reference:
34 \verbatim
35 Pope, S. B. (1997).
36 Computationally efficient implementation of combustion chemistry using
37 in situ adaptive tabulation.
38 Combustion Theory and Modelling, 1, 41-63.
39 \endverbatim
40
41\*---------------------------------------------------------------------------*/
42
43#ifndef ISAT_H
44#define ISAT_H
45
46#include "binaryTree.H"
47
48// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
49
50namespace Foam
52namespace chemistryTabulationMethods
53{
54
55/*---------------------------------------------------------------------------*\
56 Class ISAT Declaration
57\*---------------------------------------------------------------------------*/
58
59template<class CompType, class ThermoType>
60class ISAT
61:
62 public chemistryTabulationMethod<CompType, ThermoType>
63{
64 // Private data
65
66 //- List of the stored 'points' organized in a binary tree
68
69 //- List of scale factors for species, temperature and pressure
70 scalarField scaleFactor_;
71
72 const Time& runTime_;
73
74 //- Lifetime (number of time steps) of a stored point
75 label chPMaxLifeTime_;
76
77 //- Maximum number of growths before removing from the tree
78 label maxGrowth_;
79
80 //- Check the binary tree for leafs to remove every interval
81 label checkEntireTreeInterval_;
82
83 //- Factor that multiply the ideal depth of a binary tree to decide
84 //- whether to try to balance of not
85 scalar maxDepthFactor_;
86
87 //- Minimal size before trying to balance the tree
88 label minBalanceThreshold_;
89
90 //- After a failed primary retrieve, look in the MRU list
91 Switch MRURetrieve_;
92
93 //- Most Recently Used (MRU) list of chemPoint
95
96 //- Maximum size of the MRU list
97 label maxMRUSize_;
98
99 //- Store a pointer to the last chemPointISAT found
101
102 //- Switch to allow growth (on by default)
103 Switch growPoints_;
104
105 // Statistics on ISAT usage
106 label nRetrieved_;
107 label nGrowth_;
108 label nAdd_;
109
110 autoPtr<OFstream> nRetrievedFile_;
111 autoPtr<OFstream> nGrowthFile_;
112 autoPtr<OFstream> nAddFile_;
113 autoPtr<OFstream> sizeFile_;
114
115 bool cleaningRequired_;
116
117 //- Number of equations in addition to the species eqs.
118 label nAdditionalEqns_;
119
120
121 // Private Member Functions
122
123 //- No copy construct
124 ISAT(const ISAT&) = delete;
125
126 //- Add a chemPoint to the MRU list
127 void addToMRU(chemPointISAT<CompType, ThermoType>* phi0);
128
129 //- Compute and return the mapping of the composition phiq
130 // Input : phi0 the nearest chemPoint used in the linear interpolation
131 // phiq the composition of the query point for which we want to
132 // compute the mapping
133 // Rphiq the mapping of the new composition point (given as empty)
134 // Output: void (the mapping is stored in the Rphiq array)
135 // Rphiq = Rphi0 + A * (phiq-phi0)
136 void calcNewC
137 (
139 const scalarField& phiq,
140 scalarField& Rphiq
141 );
142
143 //- Check if the composition of the query point phiq lies in the
144 //- ellipsoid of accuracy approximating the region of accuracy of the
145 //- stored chemPoint phi0
146 // Input : phi0 the nearest chemPoint used in the linear interpolation
147 // phiq the composition of the query point for which we want to
148 // compute the mapping
149 // Output: true if phiq is in the EOA, false if not
150 bool grow
151 (
153 const scalarField& phiq,
154 const scalarField& Rphiq
155 );
156
157 //- Clean and balance the tree
158 bool cleanAndBalance();
159
160 //- Functions to construct the gradients matrix
161 // When mechanism reduction is active, the A matrix is given by
162 // Aaa Aad
163 // A = ( Ada Add ), where the sub gradient matrices are:
164 // (Aaa) active species according to active species, (Aad) active
165 // species according to disabled species, (Ada) disabled species
166 // according to active species, and (Add) disabled species according to
167 // disabled species.
168 // The current implementation computes Aaa with the Jacobian of the
169 // reduced set of species. Aad = 0, Ada = 0, and Add = I.
170 // To be implemented: add options to compute the A matrix for different
171 // strategies
172 void computeA
173 (
175 const scalarField& Rphiq,
176 const scalar rho,
177 const scalar dt
178 );
179
180
181public:
182
183 //- Runtime type information
184 TypeName("ISAT");
185
186
187 // Constructors
188
189 //- Construct from dictionary
190 ISAT
191 (
192 const dictionary& chemistryProperties,
194 );
195
196
197 // Destructor
198 virtual ~ISAT();
199
200
201 // Member Functions
204 {
205 return chemisTree_;
206 }
208 inline const scalarField& scaleFactor() const
209 {
210 return scaleFactor_;
211 }
212
213 //- Return the size of the binary tree
214 virtual inline label size()
215 {
216 return chemisTree_.size();
217 }
218
219 virtual void writePerformance();
220
221 //- Find the closest stored leaf of phiQ and store the result in
222 //- RphiQ or return false.
223 virtual bool retrieve
224 (
225 const Foam::scalarField& phiq,
226 scalarField& Rphiq
227 );
228
229 //- Add information to the tabulation.
230 // This function can grow an existing point or add a new leaf to the
231 // binary tree Input : phiq the new composition to store Rphiq the
232 // mapping of the new composition point
233 virtual label add
234 (
235 const scalarField& phiq,
236 const scalarField& Rphiq,
237 const scalar rho,
238 const scalar deltaT
239 );
241 virtual bool update()
242 {
243 return cleanAndBalance();
244 }
245};
246
247
248// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
249
250} // End namespace chemistryTabulationMethods
251} // End namespace Foam
252
253// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
254
255#ifdef NoRepository
256 #include "ISAT.C"
257#endif
258
259// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
260
261#endif
262
263// ************************************************************************* //
static const Foam::dimensionedScalar A("", Foam::dimPressure, 611.21)
Template class for non-intrusive linked lists.
Definition: LList.H:79
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.H:78
Extends StandardChemistryModel by adding the TDAC method.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:80
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
Data storage of the chemistryOnLineLibrary according to a binary tree structure.
Definition: binaryTree.H:60
Leaf of the binary tree. The chemPoint stores the composition 'phi', the mapping of this composition ...
An abstract class for chemistry tabulation.
Implementation of the ISAT (In-situ adaptive tabulation), for chemistry calculation.
Definition: ISAT.H:62
virtual label add(const scalarField &phiq, const scalarField &Rphiq, const scalar rho, const scalar deltaT)
Add information to the tabulation.
Definition: ISAT.C:510
virtual bool retrieve(const Foam::scalarField &phiq, scalarField &Rphiq)
Definition: ISAT.C:435
virtual label size()
Return the size of the binary tree.
Definition: ISAT.H:213
TypeName("ISAT")
Runtime type information.
binaryTree< CompType, ThermoType > & chemisTree()
Definition: ISAT.H:202
const scalarField & scaleFactor() const
Definition: ISAT.H:207
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
BasicChemistryModel< psiReactionThermo > & chemistry
Namespace for OpenFOAM.
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition: typeInfo.H:73