ISAT.C
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-2021 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 \*---------------------------------------------------------------------------*/
25 
26 #include "ISAT.H"
27 #include "LUscalarMatrix.H"
28 
29 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
30 
31 template<class ThermoType>
33 (
34  const dictionary& chemistryProperties,
36 )
37 :
39  (
40  chemistryProperties,
41  chemistry
42  ),
43  chemisTree_(chemistry, this->coeffsDict_),
44  scaleFactor_(chemistry.nEqns() + 1, 1),
45  runTime_(chemistry.time()),
46  chPMaxLifeTime_
47  (
48  this->coeffsDict_.lookupOrDefault("chPMaxLifeTime", INT_MAX)
49  ),
50  maxGrowth_(this->coeffsDict_.lookupOrDefault("maxGrowth", INT_MAX)),
51  checkEntireTreeInterval_
52  (
53  this->coeffsDict_.lookupOrDefault("checkEntireTreeInterval", INT_MAX)
54  ),
55  maxDepthFactor_
56  (
57  this->coeffsDict_.lookupOrDefault
58  (
59  "maxDepthFactor",
60  (chemisTree_.maxNLeafs() - 1)
61  /(log(scalar(chemisTree_.maxNLeafs()))/log(2.0))
62  )
63  ),
64  minBalanceThreshold_
65  (
66  this->coeffsDict_.lookupOrDefault
67  (
68  "minBalanceThreshold", 0.1*chemisTree_.maxNLeafs()
69  )
70  ),
71  MRURetrieve_(this->coeffsDict_.lookupOrDefault("MRURetrieve", false)),
72  maxMRUSize_(this->coeffsDict_.lookupOrDefault("maxMRUSize", 0)),
73  lastSearch_(nullptr),
74  growPoints_(this->coeffsDict_.lookupOrDefault("growPoints", true)),
75  nRetrieved_(0),
76  nGrowth_(0),
77  nAdd_(0),
78  cleaningRequired_(false)
79 {
80  if (this->active_)
81  {
82  dictionary scaleDict(this->coeffsDict_.subDict("scaleFactor"));
83  label Ysize = this->chemistry_.Y().size();
84  scalar otherScaleFactor = scaleDict.lookup<scalar>("otherSpecies");
85  for (label i=0; i<Ysize; i++)
86  {
87  if (!scaleDict.found(this->chemistry_.Y()[i].member()))
88  {
89  scaleFactor_[i] = otherScaleFactor;
90  }
91  else
92  {
93  scaleFactor_[i] =
94  scaleDict.lookup<scalar>(this->chemistry_.Y()[i].member());
95  }
96  }
97  scaleFactor_[Ysize] = scaleDict.lookup<scalar>("Temperature");
98  scaleFactor_[Ysize + 1] = scaleDict.lookup<scalar>("Pressure");
99  scaleFactor_[Ysize + 2] = scaleDict.lookup<scalar>("deltaT");
100  }
101 
102  if (this->log())
103  {
104  nRetrievedFile_ = chemistry.logFile("found_isat.out");
105  nGrowthFile_ = chemistry.logFile("growth_isat.out");
106  nAddFile_ = chemistry.logFile("add_isat.out");
107  sizeFile_ = chemistry.logFile("size_isat.out");
108  }
109 }
110 
111 
112 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
113 
114 template<class ThermoType>
116 {}
117 
118 
119 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
120 
121 template<class ThermoType>
123 (
125 )
126 {
127  if (maxMRUSize_ > 0 && MRURetrieve_)
128  {
129  // First search if the chemPoint is already in the list
130  bool isInList = false;
131  typename SLList <chemPointISAT<ThermoType>*>::iterator iter =
132  MRUList_.begin();
133  for ( ; iter != MRUList_.end(); ++iter)
134  {
135  if (iter() == phi0)
136  {
137  isInList = true;
138  break;
139  }
140  }
141  // If it is in the list, then move it to front
142  if (isInList)
143  {
144  if (iter() != MRUList_.first())
145  {
146  // iter hold the position of the element to move
147  MRUList_.remove(iter);
148 
149  // Insert the element in front of the list
150  MRUList_.insert(phi0);
151  }
152  }
153  else // chemPoint not yet in the list, iter is last
154  {
155  if (MRUList_.size() == maxMRUSize_)
156  {
157  MRUList_.remove(iter);
158  MRUList_.insert(phi0);
159  }
160  else
161  {
162  MRUList_.insert(phi0);
163  }
164  }
165  }
166 }
167 
168 
169 template<class ThermoType>
171 (
173  const scalarField& phiq,
174  scalarField& Rphiq
175 )
176 {
177  label nEqns = this->chemistry_.nEqns(); // Species, T, p
178  bool mechRedActive = this->chemistry_.mechRed()->active();
179  Rphiq = phi0->Rphi();
180  scalarField dphi(phiq-phi0->phi());
181  const scalarSquareMatrix& gradientsMatrix = phi0->A();
182  List<label>& completeToSimplified(phi0->completeToSimplifiedIndex());
183 
184  // Rphiq[i]=Rphi0[i]+A(i, j)dphi[j]
185  // where Aij is dRi/dphi_j
186  for (label i=0; i<nEqns-3; i++)
187  {
188  if (mechRedActive)
189  {
190  label si = completeToSimplified[i];
191  // The species is active
192  if (si != -1)
193  {
194  for (label j=0; j<nEqns-2; j++)
195  {
196  label sj = completeToSimplified[j];
197  if (sj != -1)
198  {
199  Rphiq[i] += gradientsMatrix(si, sj)*dphi[j];
200  }
201  }
202  Rphiq[i] +=
203  gradientsMatrix(si, phi0->nActiveSpecies())*dphi[nEqns - 2];
204  Rphiq[i] +=
205  gradientsMatrix(si, phi0->nActiveSpecies() + 1)
206  *dphi[nEqns - 1];
207  Rphiq[i] +=
208  gradientsMatrix(si, phi0->nActiveSpecies() + 2)
209  *dphi[nEqns];
210 
211  // As we use an approximation of A, Rphiq should be checked for
212  // negative values
213  Rphiq[i] = max(0, Rphiq[i]);
214  }
215  // The species is not active A(i, j) = I(i, j)
216  else
217  {
218  Rphiq[i] += dphi[i];
219  Rphiq[i] = max(0, Rphiq[i]);
220  }
221  }
222  else // Mechanism reduction is not active
223  {
224  for (label j=0; j<nEqns; j++)
225  {
226  Rphiq[i] += gradientsMatrix(i, j)*dphi[j];
227  }
228  // As we use a first order gradient matrix, Rphiq should be checked
229  // for negative values
230  Rphiq[i] = max(0, Rphiq[i]);
231  }
232  }
233 }
234 
235 
236 template<class ThermoType>
238 (
240  const scalarField& phiq,
241  const scalarField& Rphiq
242 )
243 {
244  // If the pointer to the chemPoint is nullptr, the function stops
245  if (!phi0)
246  {
247  return false;
248  }
249 
250  // Raise a flag when the chemPoint used has been grown more than the
251  // allowed number of time
252  if (phi0->nGrowth() > maxGrowth_)
253  {
254  cleaningRequired_ = true;
255  phi0->toRemove() = true;
256  return false;
257  }
258 
259  // If the solution RphiQ is still within the tolerance we try to grow it
260  // in some cases this might result in a failure and the grow function of
261  // the chemPoint returns false
262  if (phi0->checkSolution(phiq, Rphiq))
263  {
264  return phi0->grow(phiq);
265  }
266  // The actual solution and the approximation given by ISAT are too different
267  else
268  {
269  return false;
270  }
271 }
272 
273 
274 template<class ThermoType>
275 bool
277 {
278  bool treeModified(false);
279 
280  // Check all chemPoints to see if we need to delete some of the chemPoints
281  // according to the elapsed time and number of growths
282  chemPointISAT<ThermoType>* x = chemisTree_.treeMin();
283  while(x != nullptr)
284  {
286  chemisTree_.treeSuccessor(x);
287 
288  scalar elapsedTimeSteps = this->chemistry_.timeSteps() - x->timeTag();
289 
290  if ((elapsedTimeSteps > chPMaxLifeTime_) || (x->nGrowth() > maxGrowth_))
291  {
292  chemisTree_.deleteLeaf(x);
293  treeModified = true;
294  }
295  x = xtmp;
296  }
297 
298  MRUList_.clear();
299 
300  // Check if the tree should be balanced according to criterion:
301  // -the depth of the tree bigger than a*log2(size), log2(size) being the
302  // ideal depth (e.g. 4 leafs can be stored in a tree of depth 2)
303  if
304  (
305  chemisTree_.size() > minBalanceThreshold_
306  && chemisTree_.depth() >
307  maxDepthFactor_*log(scalar(chemisTree_.size()))/log(2.0)
308  )
309  {
310  chemisTree_.balance();
311  treeModified = true;
312  }
313 
314  // Return a bool to specify if the tree structure has been modified and is
315  // now below the user specified limit (true if not full)
316  return (treeModified && !chemisTree_.isFull());
317 }
318 
319 
320 template<class ThermoType>
322 (
324  const scalarField& Rphiq,
325  const label li,
326  const scalar rhoi,
327  const scalar dt
328 )
329 {
330  bool mechRedActive = this->chemistry_.mechRed()->active();
331  label speciesNumber = this->chemistry_.nSpecie();
332  scalarField Rcq(this->chemistry_.nEqns() + 1);
333  for (label i=0; i<speciesNumber; i++)
334  {
335  label s2c = i;
336  if (mechRedActive)
337  {
338  s2c = this->chemistry_.simplifiedToCompleteIndex()[i];
339  }
340  Rcq[i] = rhoi*Rphiq[s2c]/this->chemistry_.specieThermos()[s2c].W();
341  }
342  Rcq[speciesNumber] = Rphiq[Rphiq.size() - 3];
343  Rcq[speciesNumber + 1] = Rphiq[Rphiq.size() - 2];
344  Rcq[speciesNumber + 2] = Rphiq[Rphiq.size() - 1];
345 
346  // Aaa is computed implicitly,
347  // A is given by A = C(psi0, t0+dt), where C is obtained through solving
348  // d/dt C(psi0, t) = J(psi(t))C(psi0, t)
349  // If we solve it implicitly:
350  // (C(psi0, t0+dt) - C(psi0, t0))/dt = J(psi(t0+dt))C(psi0, t0+dt)
351  // The Jacobian is thus computed according to the mapping
352  // C(psi0,t0+dt)*(I-dt*J(psi(t0+dt))) = C(psi0, t0)
353  // A = C(psi0,t0)/(I-dt*J(psi(t0+dt)))
354  // where C(psi0,t0) = I
355  scalarField dcdt(speciesNumber + 2, Zero);
356  this->chemistry_.jacobian(runTime_.value(), Rcq, li, dcdt, A);
357 
358  // The jacobian is computed according to the molar concentration
359  // the following conversion allows the code to use A with mass fraction
360  for (label i=0; i<speciesNumber; i++)
361  {
362  label si = i;
363 
364  if (mechRedActive)
365  {
366  si = this->chemistry_.simplifiedToCompleteIndex()[i];
367  }
368 
369  for (label j=0; j<speciesNumber; j++)
370  {
371  label sj = j;
372  if (mechRedActive)
373  {
374  sj = this->chemistry_.simplifiedToCompleteIndex()[j];
375  }
376  A(i, j) *=
377  -dt*this->chemistry_.specieThermos()[si].W()
378  /this->chemistry_.specieThermos()[sj].W();
379  }
380 
381  A(i, i) += 1;
382  // Columns for pressure and temperature
383  A(i, speciesNumber) *=
384  -dt*this->chemistry_.specieThermos()[si].W()/rhoi;
385  A(i, speciesNumber + 1) *=
386  -dt*this->chemistry_.specieThermos()[si].W()/rhoi;
387  }
388 
389  // For the temperature and pressure lines, ddc(dTdt)
390  // should be converted in ddY(dTdt)
391  for (label i=0; i<speciesNumber; i++)
392  {
393  label si = i;
394  if (mechRedActive)
395  {
396  si = this->chemistry_.simplifiedToCompleteIndex()[i];
397  }
398 
399  A(speciesNumber, i) *=
400  -dt*rhoi/this->chemistry_.specieThermos()[si].W();
401  A(speciesNumber + 1, i) *=
402  -dt*rhoi/this->chemistry_.specieThermos()[si].W();
403  }
404 
405  A(speciesNumber, speciesNumber) = -dt*A(speciesNumber, speciesNumber) + 1;
406 
407  A(speciesNumber + 1, speciesNumber + 1) =
408  -dt*A(speciesNumber + 1, speciesNumber + 1) + 1;
409 
410  A(speciesNumber + 2, speciesNumber + 2) = 1;
411 
412  // Inverse of (I-dt*J(psi(t0+dt)))
413  LUscalarMatrix LUA(A);
414  LUA.inv(A);
415 
416  // After inversion, lines of p and T are set to 0 except diagonal. This
417  // avoid skewness of the ellipsoid of accuracy and potential issues in the
418  // binary tree.
419  for (label i=0; i<speciesNumber; i++)
420  {
421  A(speciesNumber, i) = 0;
422  A(speciesNumber + 1, i) = 0;
423  }
424 }
425 
426 
427 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
428 
429 template<class ThermoType>
431 (
432  const Foam::scalarField& phiq,
433  scalarField& Rphiq
434 )
435 {
436  bool retrieved(false);
438 
439  // If the tree is not empty
440  if (chemisTree_.size())
441  {
442  chemisTree_.binaryTreeSearch(phiq, chemisTree_.root(), phi0);
443 
444  // lastSearch keeps track of the chemPoint we obtain by the regular
445  // binary tree search
446  lastSearch_ = phi0;
447  if (phi0->inEOA(phiq))
448  {
449  retrieved = true;
450  }
451  // After a successful secondarySearch, phi0 store a pointer to the
452  // found chemPoint
453  else if (chemisTree_.secondaryBTSearch(phiq, phi0))
454  {
455  retrieved = true;
456  }
457  else if (MRURetrieve_)
458  {
459  typename SLList
460  <
462  >::iterator iter = MRUList_.begin();
463 
464  for ( ; iter != MRUList_.end(); ++iter)
465  {
466  phi0 = iter();
467  if (phi0->inEOA(phiq))
468  {
469  retrieved = true;
470  break;
471  }
472  }
473  }
474  }
475  // The tree is empty, retrieved is still false
476  else
477  {
478  // There is no chempoints that we can try to grow
479  lastSearch_ = nullptr;
480  }
481 
482  if (retrieved)
483  {
484  phi0->increaseNumRetrieve();
485  scalar elapsedTimeSteps =
486  this->chemistry_.timeSteps() - phi0->timeTag();
487 
488  // Raise a flag when the chemPoint has been used more than the allowed
489  // number of time steps
490  if (elapsedTimeSteps > chPMaxLifeTime_ && !phi0->toRemove())
491  {
492  cleaningRequired_ = true;
493  phi0->toRemove() = true;
494  }
495  lastSearch_->lastTimeUsed() = this->chemistry_.timeSteps();
496  addToMRU(phi0);
497  calcNewC(phi0, phiq, Rphiq);
498  nRetrieved_++;
499  return true;
500  }
501  else
502  {
503  // This point is reached when every retrieve trials have failed
504  // or if the tree is empty
505  return false;
506  }
507 }
508 
509 
510 template<class ThermoType>
512 (
513  const scalarField& phiq,
514  const scalarField& Rphiq,
515  const label li,
516  const scalar rho,
517  const scalar deltaT
518 )
519 {
520  label growthOrAddFlag = 1;
521  // If lastSearch_ holds a valid pointer to a chemPoint AND the growPoints_
522  // option is on, the code first tries to grow the point hold by lastSearch_
523  if (lastSearch_ && growPoints_)
524  {
525  if (grow(lastSearch_, phiq, Rphiq))
526  {
527  nGrowth_++;
528  growthOrAddFlag = 0;
529  addToMRU(lastSearch_);
530  //the structure of the tree is not modified, return false
531  return growthOrAddFlag;
532  }
533  }
534 
535  // If the code reach this point, it is either because lastSearch_ is not
536  // valid, OR because growPoints_ is not on, OR because the grow operation
537  // has failed. In the three cases, a new point is added to the tree.
538  if (chemisTree().isFull())
539  {
540  // If cleanAndBalance operation do not result in a reduction of the tree
541  // size, the last possibility is to delete completely the tree.
542  // It can be partially rebuild with the MRU list if this is used.
543  if (!cleanAndBalance())
544  {
546  if (maxMRUSize_>0)
547  {
548  // Create a copy of each chemPointISAT of the MRUList_ before
549  // they are deleted
550  typename SLList
551  <
553  >::iterator iter = MRUList_.begin();
554  for ( ; iter != MRUList_.end(); ++iter)
555  {
556  tempList.append
557  (
558  new chemPointISAT<ThermoType>(*iter())
559  );
560  }
561  }
562  chemisTree().clear();
563 
564  // Pointers to chemPoint are not valid anymore, clear the list
565  MRUList_.clear();
566 
567  // Construct the tree without giving a reference to attach to it
568  // since the structure has been completely discarded
569  chemPointISAT<ThermoType>* nulPhi = 0;
570  forAll(tempList, i)
571  {
572  chemisTree().insertNewLeaf
573  (
574  tempList[i]->phi(),
575  tempList[i]->Rphi(),
576  tempList[i]->A(),
577  scaleFactor(),
578  this->tolerance(),
579  scaleFactor_.size(),
580  nulPhi
581  );
582  deleteDemandDrivenData(tempList[i]);
583  }
584  }
585 
586  // The structure has been changed, it will force the binary tree to
587  // perform a new search and find the most appropriate point still stored
588  lastSearch_ = nullptr;
589  }
590 
591  // Compute the A matrix needed to store the chemPoint.
592  label ASize = this->chemistry_.nEqns() + 1;
593  scalarSquareMatrix A(ASize, Zero);
594  computeA(A, Rphiq, li, rho, deltaT);
595 
596  chemisTree().insertNewLeaf
597  (
598  phiq,
599  Rphiq,
600  A,
601  scaleFactor(),
602  this->tolerance(),
603  scaleFactor_.size(),
604  lastSearch_ // lastSearch_ may be nullptr (handled by binaryTree)
605  );
606  if (lastSearch_ != nullptr)
607  {
608  addToMRU(lastSearch_);
609  }
610  nAdd_++;
611 
612  return growthOrAddFlag;
613 }
614 
615 
616 template<class ThermoType>
618 {
619  if (this->log())
620  {
621  nRetrievedFile_()
622  << runTime_.timeOutputValue() << " " << nRetrieved_ << endl;
623  nRetrieved_ = 0;
624 
625  nGrowthFile_()
626  << runTime_.timeOutputValue() << " " << nGrowth_ << endl;
627  nGrowth_ = 0;
628 
629  nAddFile_()
630  << runTime_.timeOutputValue() << " " << nAdd_ << endl;
631  nAdd_ = 0;
632 
633  sizeFile_()
634  << runTime_.timeOutputValue() << " " << this->size() << endl;
635  }
636 }
637 
638 
639 // ************************************************************************* //
Extends standardChemistryModel by adding the TDAC method.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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
const scalarSquareMatrix & A() const
dimensionedScalar log(const dimensionedScalar &ds)
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
Class to perform the LU decomposition on a symmetric matrix.
const iterator & end()
Definition: LList.H:286
Template class for non-intrusive linked lists.
Definition: LList.H:51
T & first()
Return the first entry added.
Definition: LList.H:139
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
bool checkSolution(const scalarField &phiq, const scalarField &Rphiq)
If phiq is not in the EOA, then the mapping is computed.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Leaf of the binary tree. The chemPoint stores the composition &#39;phi&#39;, the mapping of this composition ...
List< label > & completeToSimplifiedIndex()
const label & timeTag()
virtual bool retrieve(const Foam::scalarField &phiq, scalarField &Rphiq)
Find the closest stored leaf of phiQ and store the result in.
Definition: ISAT.C:431
ISAT(const dictionary &chemistryProperties, TDACChemistryModel< ThermoType > &chemistry)
Construct from dictionary.
Definition: ISAT.C:33
label size() const
Return the number of elements in matrix (m*n)
Definition: MatrixI.H:71
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects...
Definition: DynamicList.H:56
virtual label add(const scalarField &phiq, const scalarField &Rphiq, const label li, const scalar rho, const scalar deltaT)
Add information to the tabulation.
Definition: ISAT.C:512
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Definition: DynamicListI.H:296
bool grow(const scalarField &phiq)
More details about the minimum-volume ellipsoid covering an.
phi
Definition: correctPhi.H:3
static const zero Zero
Definition: zero.H:97
const scalarField & phi() const
basicChemistryModel & chemistry
iterator begin()
Definition: LList.H:281
const scalarField & Rphi() const
An abstract class for chemistry tabulation.
const Time & time() const
Return time.
Definition: IOobject.C:318
virtual label nEqns() const
Number of ODE&#39;s to solve.
void increaseNumRetrieve()
Increases the number of retrieves the chempoint has generated.
bool inEOA(const scalarField &phiq)
To RETRIEVE the mapping from the stored chemPoint phi, the query.
autoPtr< OFstream > logFile(const word &name) const
Create and return a TDAC log file of the given name.
Implementation of the ISAT (In-situ adaptive tabulation), for chemistry calculation.
Definition: ISAT.H:57
const dimensionedScalar phi0
Magnetic flux quantum: default SI units: [Wb].
void deleteDemandDrivenData(DataPtr &dataPtr)
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:844
void inv(scalarSquareMatrix &M) const
Set M to the inverse of this square matrix.