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