chemPointISAT.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::chemPointISAT
28
29Description
30 Leaf of the binary tree.
31 The chemPoint stores the composition 'phi', the mapping of this
32 composition Rphi, the mapping gradient matrix A and the matrix describing
33 the Ellipsoid Of Accuracy (EOA).
34
35 1)When the chemPoint is created the region of accuracy is approximated by
36 an ellipsoid E centered in 'phi' (obtained with the constant):
37 E = {x| ||L^T.(x-phi)|| <= 1},
38 with x a point in the composition space and L^T the transpose of an upper
39 triangular matrix describing the EOA (see below: "Computation of L" ).
40
41 2)To RETRIEVE the mapping from the chemPoint phi, the query point phiq has to
42 be in the EOA of phi. It follows that, dphi=phiq-phi and to test if phiq
43 is in the ellipsoid there are two methods. First, compare r=||dphi|| with
44 rmin and rmax. If r < rmin, phiq is in the EOA. If r > rmax, phiq is out of
45 the EOA. This operations is O(completeSpaceSize) and is performed first.
46 If rmin < r < rmax, then the second method is used:
47 ||L^T.dphi|| <= 1 to be in the EOA.
48
49 If phiq is in the EOA, Rphiq is obtained by linear interpolation:
50 Rphiq= Rphi + A.dphi.
51
52 3)If phiq is not in the EOA, then the mapping is computed. But as the EOA
53 is a conservative approximation of the region of accuracy surrounding the
54 point phi, we could expand it by comparing the computed results with the
55 one obtained by linear interpolation. The error epsGrow is calculated:
56 epsGrow = ||B.(dR - dRl)||,
57 with dR = Rphiq - Rphi, dRl = A.dphi and B the diagonal scale factor
58 matrix.
59 If epsGrow <= tolerance, the EOA is too conservative and a GROW is perforned
60 otherwise, the newly computed mapping is associated to the initial
61 composition and added to the tree.
62
63 4)To GROW the EOA, we expand it to include the previous EOA and the query
64 point phiq. The rank-one matrix method is used. The EOA is transformed
65 to a hypersphere centered at the origin. Then it is expanded to include
66 the transformed point phiq' on its boundary. Then the inverse transformation
67 give the modified matrix L' (see below: "Grow the EOA").
68
69
70 Computation of L :
71 In [1], the EOA of the constant approximation is given by
72 E = {x| ||B.A/tolerance.(x-phi)|| <= 1},
73 with B a scale factor diagonal matrix, A the mapping gradient matrix and
74 tolerance the absolute tolerance. If we take the QR decomposition of
75 (B.A)/tolerance= Q.R, with Q an orthogonal matrix and R an upper triangular
76 matrix such that the EOA is described by
77 (phiq-phi0)^T.R^T.R.(phiq-phi0) <= 1
78 L^T = R, both Cholesky decomposition of A^T.B^T.B.A/tolerance^2
79 This representation of the ellipsoid is used in [2] and in order to avoid
80 large value of semi-axe length in certain direction, a Singular Value
81 Decomposition (SVD) is performed on the L matrix:
82 L = UDV^T,
83 with the orthogonal matrix U giving the directions of the principal axes
84 and 1/di the inverse of the element of the diagonal matrix D giving the
85 length of the principal semi-axes. To avoid very large value of those
86 length,
87 di' = max(di, 1/(alphaEOA*sqrt(tolerance))), with alphaEOA = 0.1 (see [2])
88 di' = max(di, 1/2), see [1]. The latter will be used in this implementation.
89 And L' = UD'V^T, with D' the diagonal matrix with the modified di'.
90
91 Grow the EOA :
92 More details about the minimum-volume ellipsoid covering an ellispoid E and
93 a point p are found in [3]. Here is the main steps to obtain the modified
94 matrix L' describing the new ellipsoid.
95 1) calculate the point p' in the transformed space :
96 p' = L^T.(p-phi)
97 2) compute the rank-one decomposition:
98 G = I + gamma.p'.p'^T,
99 with gamma = (1/|p'|-1)*1/|p'|^2
100 3) compute L':
101 L' = L.G.
102
103 References:
104 \verbatim
105 [1] Pope, S. B. (1997).
106 Computationally efficient implementation of combustion chemistry using
107 in situ adaptive tabulation.
108 Combustion Theory and Modelling, 1, 41-63.
109
110 [2] Lu, L., & Pope, S. B. (2009).
111 An improved algorithm for in situ adaptive tabulation.
112 Journal of Computational Physics, 228(2), 361-386.
113
114 [3] Pope, S. B. (2008).
115 Algorithms for ellipsoids.
116 Cornell University Report No. FDA, 08-01.
117 \endverbatim
118
119\*---------------------------------------------------------------------------*/
120
121#ifndef chemPointISAT_H
122#define chemPointISAT_H
123
124// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
125
126namespace Foam
127{
128
129// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
130
131template<class CompType, class ThermoType>
132class binaryNode;
133
134template<class CompType, class ThermoType>
135class TDACChemistryModel;
136
137
138/*---------------------------------------------------------------------------*\
139 Class chemPointISAT Declaration
140\*---------------------------------------------------------------------------*/
141
142template<class CompType, class ThermoType>
143class chemPointISAT
144{
145 // Private data
146
147 //- Reference to the chemistryModel object
149
150 //- Vector storing the composition, temperature and pressure
151 // and deltaT if a variable time step is set on
152 scalarField phi_;
153
154 //- Vector storing the mapping of the composition phi
155 scalarField Rphi_;
156
157 //- LT the transpose of the L matrix describing the Ellipsoid Of
158 // Accuracy use List of Lists to be able to change size if DAC is used
160
161 //- A the mapping gradient matrix
163
164 //- Vector storing the scale factor
165 scalarField scaleFactor_;
166
167 //- Reference to the node in the binary tree holding this chemPoint
169
170 //- The size of the composition space (size of the vector phi)
171 label completeSpaceSize_;
172
173 //- Number of times the element has been grown
174 label nGrowth_;
175
176 //- Tolerance for the Ellipsoid of accuracy
177 static scalar tolerance_;
178
179 //- Number of active species stored in the chemPoint
180 label nActiveSpecies_;
181
182 //- Vectors that store the index conversion between the simplified
183 // and the complete chemical mechanism
184 List<label> simplifiedToCompleteIndex_;
185
186 label timeTag_;
187 label lastTimeUsed_;
188
189 bool toRemove_;
190
191 label maxNumNewDim_;
192
193 Switch printProportion_;
194
195 //- Variable to store the number of retrieves the chemPoint
196 // will generate at each time step
197 label numRetrieve_;
198
199 //- Variable to store the number of time steps the chempoint is allowed
200 // to still live according to the maxChPLifeTime_ parameter
201 label nLifeTime_;
202
203 List<label> completeToSimplifiedIndex_;
204
205 //- Number of equations in addition to the species eqs.
206 label nAdditionalEqns_;
207
208 label idT_;
209 label idp_;
210 label iddeltaT_;
211
212 //- QR decomposition of a matrix
213 // Input : nCols cols number
214 // R the matrix to decompose
215 // QT an empty matrix that stores the transpose of the Q matrix
216 // R is returned in the given R matrix
217 // which is used to store the ellipsoid of accuracy
218 void qrDecompose
219 (
220 const label nCols,
222 );
223
224 //- QR update of the matrix A
225 void qrUpdate
226 (
228 const label n,
229 const scalarField& u,
230 const scalarField& v
231 );
232
233 void rotate
234 (
236 const label i,
237 const scalar a,
238 const scalar b,
239 label n
240 );
241
242
243public:
244
245 // Constructors
246
247 //- Construct from components
249 (
251 const scalarField& phi,
252 const scalarField& Rphi,
253 const scalarSquareMatrix& A,
255 const scalar& tolerance,
256 const label& completeSpaceSize,
257 const dictionary& coeffsDict,
259 );
260
261 //- Construct from another chemPoint and reference to a binary node
263 (
266 );
267
268 //- Construct from another chemPoint
270 (
272 );
273
274
275 // Member functions
276
277 //- Access to the TDACChemistryModel
279 {
280 return chemistry_;
281 }
283 inline label nGrowth()
284 {
285 return nGrowth_;
286 }
288 inline label& completeSpaceSize()
289 {
290 return completeSpaceSize_;
291 }
293 inline const scalarField& phi() const
294 {
295 return phi_;
296 }
298 inline const scalarField& Rphi() const
299 {
300 return Rphi_;
301 }
303 inline const scalarField& scaleFactor()
304 {
305 return scaleFactor_;
306 }
308 inline const scalar& tolerance()
309 {
310 return tolerance_;
311 }
313 static void changeTolerance(scalar newTol)
314 {
315 tolerance_ = newTol;
316 }
319 {
320 return node_;
321 }
323 inline const scalarSquareMatrix& A() const
324 {
325 return A_;
326 }
328 inline scalarSquareMatrix& A()
329 {
330 return A_;
331 }
333 inline const scalarSquareMatrix& LT() const
334 {
335 return LT_;
336 }
338 inline scalarSquareMatrix& LT()
339 {
340 return LT_;
341 }
343 inline label nActiveSpecies()
344 {
345 return nActiveSpecies_;
346 }
349 {
350 return completeToSimplifiedIndex_;
351 }
354 {
355 return simplifiedToCompleteIndex_;
356 }
357
358 //- Increases the number of retrieves the chempoint has generated
359 void increaseNumRetrieve();
360
361 //- Resets the number of retrieves at each time step
362 void resetNumRetrieve();
363
364 //- Increases the "counter" of the chP life
365 void increaseNLifeTime();
366
367 label simplifiedToCompleteIndex(const label i);
369 inline const label& timeTag()
370 {
371 return timeTag_;
372 }
374 inline label& lastTimeUsed()
375 {
376 return lastTimeUsed_;
377 }
379 inline bool& toRemove()
380 {
381 return toRemove_;
382 }
384 inline label& maxNumNewDim()
385 {
386 return maxNumNewDim_;
387 }
389 inline const label& numRetrieve()
390 {
391 return numRetrieve_;
392 }
394 inline const label& nLifeTime()
395 {
396 return nLifeTime_;
397 }
399 inline bool variableTimeStep() const
400 {
401 return chemistry_.variableTimeStep();
402 }
403
404 // ISAT functions
405
406 //- To RETRIEVE the mapping from the stored chemPoint phi, the query
407 // point phiq has to be in the EOA of phi.
408 // To test if phiq is in the ellipsoid:
409 // ||L^T.dphi|| <= 1
410 bool inEOA(const scalarField& phiq);
411
412 //- More details about the minimum-volume ellipsoid covering an
413 // ellipsoid E and a point p are found in [1].
414 // Here is the main steps to obtain the
415 // modified matrix L' describind the new ellipsoid.
416 // 1) calculate the point p' in the transformed space :
417 // p' = L^T.(p-phi)
418 // 2) compute the rank-one decomposition:
419 // G = I + gamma.p'.p'^T,
420 // with gamma = (1/|p'|-1)*1/|p'|^2
421 // 3) compute L':
422 // L'L'^T = (L.G)(L.G)^T,
423 // L'^T is then obtained by QR decomposition of (L.G)^T = G^T.L^T
424 // [1] Stephen B. Pope, "Algorithms for ellipsoids", FDA 08-01,
425 // Cornell University, 2008
426 bool grow(const scalarField& phiq);
427
428 //- If phiq is not in the EOA, then the mapping is computed.
429 // But as the EOA is a conservative approximation of the region of
430 // accuracy surrounding the point phi, we could expand it by
431 // comparing the computed results with the one obtained by linear
432 // interpolation. The error eps is calculated:
433 // eps = ||B.(dR - dRl)||,
434 // with dR = Rphiq - Rphi, dRl = A.dphi and B the diagonal scale
435 // factor matrix.
436 // If eps <= tolerance, the EOA is too conservative and a GROW is
437 // performed,
438 // otherwise, the newly computed mapping is associated to the
439 // initial composition and added to the tree.
440 bool checkSolution
441 (
442 const scalarField& phiq,
443 const scalarField& Rphiq
444 );
445};
446
447
448// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
449
450} // End namespace Foam
451
452// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
453
454#ifdef NoRepository
455 #include "chemPointISAT.C"
456#endif
457
458// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
459
460#endif
461
462// ************************************************************************* //
#define R(A, B, C, D, E, F, K, M)
label n
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: List.H:77
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.
bool variableTimeStep() const
Return true if the time-step is variable and/or non-uniform.
Node of the binary tree.
Definition: binaryNode.H:52
Leaf of the binary tree. The chemPoint stores the composition 'phi', the mapping of this composition ...
void resetNumRetrieve()
Resets the number of retrieves at each time step.
chemPointISAT(const chemPointISAT< CompType, ThermoType > &p, binaryNode< CompType, ThermoType > *node)
Construct from another chemPoint and reference to a binary node.
bool inEOA(const scalarField &phiq)
To RETRIEVE the mapping from the stored chemPoint phi, the query.
void increaseNumRetrieve()
Increases the number of retrieves the chempoint has generated.
const scalarField & scaleFactor()
static void changeTolerance(scalar newTol)
scalarSquareMatrix & LT()
label & completeSpaceSize()
const scalar & tolerance()
const label & numRetrieve()
bool grow(const scalarField &phiq)
More details about the minimum-volume ellipsoid covering an.
const label & timeTag()
void increaseNLifeTime()
Increases the "counter" of the chP life.
const scalarSquareMatrix & LT() const
bool checkSolution(const scalarField &phiq, const scalarField &Rphiq)
If phiq is not in the EOA, then the mapping is computed.
const scalarField & Rphi() const
bool variableTimeStep() const
const label & nLifeTime()
const scalarField & phi() const
binaryNode< CompType, ThermoType > *& node()
List< label > & completeToSimplifiedIndex()
scalarSquareMatrix & A()
List< label > & simplifiedToCompleteIndex()
TDACChemistryModel< CompType, ThermoType > & chemistry()
Access to the TDACChemistryModel.
const scalarSquareMatrix & A() const
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
volScalarField & p
Namespace for OpenFOAM.
volScalarField & b
Definition: createFields.H:27