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-2024 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 "odeChemistryModel.H"
28 #include "LUscalarMatrix.H"
30 
31 
32 /* * * * * * * * * * * * * * * Static Member Data * * * * * * * * * * * * * */
33 
34 namespace Foam
35 {
36 namespace chemistryTabulationMethods
37 {
40 }
41 }
42 
43 
44 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
45 
46 Foam::chemistryTabulationMethods::ISAT::ISAT
47 (
48  const dictionary& chemistryProperties,
49  const dictionary& coeffDict,
50  const odeChemistryModel& chemistry
51 )
52 :
53  chemistryTabulationMethod
54  (
55  chemistryProperties,
56  chemistry
57  ),
58  chemistry_(chemistry),
59  log_(coeffDict.lookupOrDefault<Switch>("log", false)),
60  reduction_(chemistry_.reduction()),
61  chemisTree_(*this, coeffDict),
62  scaleFactor_(chemistry.nEqns() + 1, 1),
63  runTime_(chemistry.time()),
64  timeSteps_(0),
65  chPMaxLifeTime_
66  (
67  coeffDict.lookupOrDefault("chPMaxLifeTime", INT_MAX)
68  ),
69  maxGrowth_(coeffDict.lookupOrDefault("maxGrowth", INT_MAX)),
70  checkEntireTreeInterval_
71  (
72  coeffDict.lookupOrDefault("checkEntireTreeInterval", INT_MAX)
73  ),
74  maxDepthFactor_
75  (
76  coeffDict.lookupOrDefault
77  (
78  "maxDepthFactor",
79  (chemisTree_.maxNLeafs() - 1)
80  /(Foam::log(scalar(chemisTree_.maxNLeafs()))/Foam::log(2.0))
81  )
82  ),
83  minBalanceThreshold_
84  (
85  coeffDict.lookupOrDefault
86  (
87  "minBalanceThreshold", 0.1*chemisTree_.maxNLeafs()
88  )
89  ),
90  MRURetrieve_(coeffDict.lookupOrDefault("MRURetrieve", false)),
91  maxMRUSize_(coeffDict.lookupOrDefault("maxMRUSize", 0)),
92  lastSearch_(nullptr),
93  growPoints_(coeffDict.lookupOrDefault("growPoints", true)),
94  tolerance_(coeffDict.lookupOrDefault("tolerance", 1e-4)),
95  nRetrieved_(0),
96  nGrowth_(0),
97  nAdd_(0),
98  addNewLeafCpuTime_(0),
99  growCpuTime_(0),
100  searchISATCpuTime_(0),
101  tabulationResults_
102  (
103  IOobject
104  (
105  chemistry.thermo().phasePropertyName("TabulationResults"),
106  chemistry.time().name(),
107  chemistry.mesh(),
108  IOobject::NO_READ,
109  IOobject::AUTO_WRITE
110  ),
111  chemistry.mesh(),
112  scalar(0)
113  ),
114 
115  cleaningRequired_(false)
116 {
117  const dictionary& scaleDict(coeffDict.subDict("scaleFactor"));
118  label Ysize = chemistry_.Y().size();
119  scalar otherScaleFactor = scaleDict.lookup<scalar>("otherSpecies");
120  for (label i=0; i<Ysize; i++)
121  {
122  if (!scaleDict.found(chemistry_.Y()[i].member()))
123  {
124  scaleFactor_[i] = otherScaleFactor;
125  }
126  else
127  {
128  scaleFactor_[i] =
129  scaleDict.lookup<scalar>(chemistry_.Y()[i].member());
130  }
131  }
132  scaleFactor_[Ysize] = scaleDict.lookup<scalar>("Temperature");
133  scaleFactor_[Ysize + 1] = scaleDict.lookup<scalar>("Pressure");
134  scaleFactor_[Ysize + 2] = scaleDict.lookup<scalar>("deltaT");
135 
136  if (log_)
137  {
138  nRetrievedFile_ = chemistry.logFile("found_isat.out");
139  nGrowthFile_ = chemistry.logFile("growth_isat.out");
140  nAddFile_ = chemistry.logFile("add_isat.out");
141  sizeFile_ = chemistry.logFile("size_isat.out");
142 
143  cpuAddFile_ = chemistry.logFile("cpu_add.out");
144  cpuGrowFile_ = chemistry.logFile("cpu_grow.out");
145  cpuRetrieveFile_ = chemistry.logFile("cpu_retrieve.out");
146  }
147 }
148 
149 
150 Foam::chemistryTabulationMethods::ISAT::ISAT
151 (
152  const dictionary& chemistryProperties,
154 )
155 :
156  ISAT
157  (
158  chemistryProperties,
159  chemistryProperties.subDict("tabulation"),
160  chemistry
161  )
162 {}
163 
164 
165 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
166 
168 {}
169 
170 
171 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
172 
173 void Foam::chemistryTabulationMethods::ISAT::addToMRU
174 (
176 )
177 {
178  if (maxMRUSize_ > 0 && MRURetrieve_)
179  {
180  // First search if the chemPoint is already in the list
181  bool isInList = false;
182  typename SLList <chemPointISAT*>::iterator iter =
183  MRUList_.begin();
184  for ( ; iter != MRUList_.end(); ++iter)
185  {
186  if (iter() == phi0)
187  {
188  isInList = true;
189  break;
190  }
191  }
192  // If it is in the list, then move it to front
193  if (isInList)
194  {
195  if (iter() != MRUList_.first())
196  {
197  // iter hold the position of the element to move
198  MRUList_.remove(iter);
199 
200  // Insert the element in front of the list
201  MRUList_.insert(phi0);
202  }
203  }
204  else // chemPoint not yet in the list, iter is last
205  {
206  if (MRUList_.size() == maxMRUSize_)
207  {
208  MRUList_.remove(iter);
209  MRUList_.insert(phi0);
210  }
211  else
212  {
213  MRUList_.insert(phi0);
214  }
215  }
216  }
217 }
218 
219 
220 void Foam::chemistryTabulationMethods::ISAT::calcNewC
221 (
222  chemPointISAT* phi0,
223  const scalarField& phiq,
224  scalarField& Rphiq
225 )
226 {
227  const label nEqns = chemistry_.nEqns(); // Species, T, p
228  const List<label>& completeToSimplified = phi0->completeToSimplifiedIndex();
229 
230  const scalarField dphi(phiq - phi0->phi());
231  const scalarSquareMatrix& gradientsMatrix = phi0->A();
232 
233  // Linear extrapolation:
234  //
235  // Rphiq[i] = Rphi0[i] + (dR_i/dphi_j)*dphi[j]
236  // = Rphi0[i] + A(i, j)*dphi[j]
237  //
238 
239  Rphiq = phi0->Rphi();
240  for (label i=0; i<nEqns - 2; i++)
241  {
242  if (reduction_)
243  {
244  const label si = completeToSimplified[i];
245 
246  if (si != -1)
247  {
248  // If specie is active, or T or p, then extrapolate using the
249  // gradients matrix
250  for (label j=0; j<nEqns + 1; j++)
251  {
252  const label sj =
253  j < nEqns - 2
254  ? completeToSimplified[j]
255  : j - (nEqns - 2) + phi0->nActive();
256 
257  if (sj != -1)
258  {
259  Rphiq[i] += gradientsMatrix(si, sj)*dphi[j];
260  }
261  }
262  }
263  else
264  {
265  // If specie is inactive then use the tabulated value directly
266  Rphiq[i] += dphi[i];
267  }
268  }
269  else
270  {
271  // Extrapolate using the gradients matrix
272  for (label j=0; j<nEqns + 1; j++)
273  {
274  Rphiq[i] += gradientsMatrix(i, j)*dphi[j];
275  }
276  }
277 
278  // Clip
279  Rphiq[i] = max(0, Rphiq[i]);
280  }
281 }
282 
283 
284 bool Foam::chemistryTabulationMethods::ISAT::grow
285 (
286  chemPointISAT* phi0,
287  const scalarField& phiq,
288  const scalarField& Rphiq
289 )
290 {
291  // If the pointer to the chemPoint is nullptr, the function stops
292  if (!phi0)
293  {
294  return false;
295  }
296 
297  // Raise a flag when the chemPoint used has been grown more than the
298  // allowed number of time
299  if (phi0->nGrowth() > maxGrowth_)
300  {
301  cleaningRequired_ = true;
302  phi0->toRemove() = true;
303  return false;
304  }
305 
306  // If the solution RphiQ is still within the tolerance we try to grow it
307  // in some cases this might result in a failure and the grow function of
308  // the chemPoint returns false
309  if (phi0->checkSolution(phiq, Rphiq))
310  {
311  return phi0->grow(phiq);
312  }
313  // The actual solution and the approximation given by ISAT are too different
314  else
315  {
316  return false;
317  }
318 }
319 
320 
321 bool Foam::chemistryTabulationMethods::ISAT::cleanAndBalance()
322 {
323  bool treeModified(false);
324 
325  // Check all chemPoints to see if we need to delete some of the chemPoints
326  // according to the elapsed time and number of growths
327  chemPointISAT* x = chemisTree_.treeMin();
328  while(x != nullptr)
329  {
330  chemPointISAT* xtmp = chemisTree_.treeSuccessor(x);
331 
332  const scalar elapsedTimeSteps = timeSteps() - x->timeTag();
333 
334  if ((elapsedTimeSteps > chPMaxLifeTime_) || (x->nGrowth() > maxGrowth_))
335  {
336  chemisTree_.deleteLeaf(x);
337  treeModified = true;
338  }
339  x = xtmp;
340  }
341 
342  MRUList_.clear();
343 
344  // Check if the tree should be balanced according to criterion:
345  // -the depth of the tree bigger than a*log2(size), log2(size) being the
346  // ideal depth (e.g. 4 leafs can be stored in a tree of depth 2)
347  if
348  (
349  chemisTree_.size() > minBalanceThreshold_
350  && chemisTree_.depth() >
351  maxDepthFactor_*Foam::log(scalar(chemisTree_.size()))/Foam::log(2.0)
352  )
353  {
354  chemisTree_.balance();
355  treeModified = true;
356  }
357 
358  // Return a bool to specify if the tree structure has been modified and is
359  // now below the user specified limit (true if not full)
360  return (treeModified && !chemisTree_.isFull());
361 }
362 
363 
364 void Foam::chemistryTabulationMethods::ISAT::computeA
365 (
367  const scalarField& Rphiq,
368  const label li,
369  const scalar dt
370 )
371 {
372  const label nSpecie = chemistry_.nSpecie();
373 
374  scalarField Rphiqs(chemistry_.nEqns() + 1);
375  for (label i=0; i<nSpecie; i++)
376  {
377  const label si = chemistry_.sToc(i);
378  Rphiqs[i] = Rphiq[si];
379  }
380  Rphiqs[nSpecie] = Rphiq[Rphiq.size() - 3];
381  Rphiqs[nSpecie + 1] = Rphiq[Rphiq.size() - 2];
382  Rphiqs[nSpecie + 2] = Rphiq[Rphiq.size() - 1];
383 
384  // Aaa is computed implicitly,
385  // A is given by A = C(psi0, t0+dt), where C is obtained through solving
386  // d/dt C(psi0, t) = J(psi(t))C(psi0, t)
387  // If we solve it implicitly:
388  // (C(psi0, t0+dt) - C(psi0, t0))/dt = J(psi(t0+dt))C(psi0, t0+dt)
389  // The Jacobian is thus computed according to the mapping
390  // C(psi0,t0+dt)*(I-dt*J(psi(t0+dt))) = C(psi0, t0)
391  // A = C(psi0,t0)/(I-dt*J(psi(t0+dt)))
392  // where C(psi0,t0) = I
393  scalarField dYTpdt(nSpecie + 2, Zero);
394  chemistry_.jacobian(runTime_.value(), Rphiqs, li, dYTpdt, A);
395 
396  // Inverse of I - dt*J(psi(t0 + dt))
397  for (label i=0; i<nSpecie + 2; i++)
398  {
399  for (label j=0; j<nSpecie + 2; j++)
400  {
401  A(i, j) *= -dt;
402  }
403  A(i, i) += 1;
404  }
405  A(nSpecie + 2, nSpecie + 2) = 1;
406  LUscalarMatrix LUA(A);
407  LUA.inv(A);
408 
409  // After inversion, lines of p and T are set to 0 except diagonal. This
410  // avoid skewness of the ellipsoid of accuracy and potential issues in the
411  // binary tree.
412  for (label i=0; i<nSpecie; i++)
413  {
414  A(nSpecie, i) = 0;
415  A(nSpecie + 1, i) = 0;
416  }
417 }
418 
419 
420 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
421 
423 (
424  const Foam::scalarField& phiq,
425  scalarField& Rphiq
426 )
427 {
428  if (log_)
429  {
430  cpuTime_.cpuTimeIncrement();
431  }
432 
433  bool retrieved(false);
435 
436  // If the tree is not empty
437  if (chemisTree_.size())
438  {
439  chemisTree_.binaryTreeSearch(phiq, chemisTree_.root(), phi0);
440 
441  // lastSearch keeps track of the chemPoint we obtain by the regular
442  // binary tree search
443  lastSearch_ = phi0;
444  if (phi0->inEOA(phiq))
445  {
446  retrieved = true;
447  }
448  // After a successful secondarySearch, phi0 store a pointer to the
449  // found chemPoint
450  else if (chemisTree_.secondaryBTSearch(phiq, phi0))
451  {
452  retrieved = true;
453  }
454  else if (MRURetrieve_)
455  {
456  typename SLList
457  <
459  >::iterator iter = MRUList_.begin();
460 
461  for ( ; iter != MRUList_.end(); ++iter)
462  {
463  phi0 = iter();
464  if (phi0->inEOA(phiq))
465  {
466  retrieved = true;
467  break;
468  }
469  }
470  }
471  }
472  // The tree is empty, retrieved is still false
473  else
474  {
475  // There is no chempoints that we can try to grow
476  lastSearch_ = nullptr;
477  }
478 
479  if (retrieved)
480  {
481  phi0->increaseNumRetrieve();
482  const scalar elapsedTimeSteps = timeSteps() - phi0->timeTag();
483 
484  // Raise a flag when the chemPoint has been used more than the allowed
485  // number of time steps
486  if (elapsedTimeSteps > chPMaxLifeTime_ && !phi0->toRemove())
487  {
488  cleaningRequired_ = true;
489  phi0->toRemove() = true;
490  }
491  lastSearch_->lastTimeUsed() = timeSteps();
492  addToMRU(phi0);
493  calcNewC(phi0, phiq, Rphiq);
494  nRetrieved_++;
495  }
496 
497  if (log_)
498  {
499  searchISATCpuTime_ += cpuTime_.cpuTimeIncrement();
500  }
501 
502  return retrieved;
503 }
504 
505 
507 (
508  const scalarField& phiq,
509  const scalarField& Rphiq,
510  const label nActive,
511  const label li,
512  const scalar deltaT
513 )
514 {
515  if (log_)
516  {
517  cpuTime_.cpuTimeIncrement();
518  }
519 
520  label growthOrAddFlag = 1;
521 
522  // If lastSearch_ holds a valid pointer to a chemPoint AND the growPoints_
523  // option is on, the code first tries to grow the point hold by lastSearch_
524  if (lastSearch_ && growPoints_)
525  {
526  if (grow(lastSearch_, phiq, Rphiq))
527  {
528  nGrowth_++;
529  growthOrAddFlag = 0;
530  addToMRU(lastSearch_);
531 
532  tabulationResults_[li] = 1;
533 
534  if (log_)
535  {
536  growCpuTime_ += cpuTime_.cpuTimeIncrement();
537  }
538 
539  // the structure of the tree is not modified, return false
540  return growthOrAddFlag;
541  }
542  }
543 
544  // If the code reach this point, it is either because lastSearch_ is not
545  // valid, OR because growPoints_ is not on, OR because the grow operation
546  // has failed. In the three cases, a new point is added to the tree.
547  if (chemisTree().isFull())
548  {
549  // If cleanAndBalance operation do not result in a reduction of the tree
550  // size, the last possibility is to delete completely the tree.
551  // It can be partially rebuild with the MRU list if this is used.
552  if (!cleanAndBalance())
553  {
555  if (maxMRUSize_>0)
556  {
557  // Create a copy of each chemPointISAT of the MRUList_ before
558  // they are deleted
559  typename SLList
560  <
562  >::iterator iter = MRUList_.begin();
563  for ( ; iter != MRUList_.end(); ++iter)
564  {
565  tempList.append
566  (
567  new chemPointISAT(*iter())
568  );
569  }
570  }
571  chemisTree().clear();
572 
573  // Pointers to chemPoint are not valid anymore, clear the list
574  MRUList_.clear();
575 
576  // Construct the tree without giving a reference to attach to it
577  // since the structure has been completely discarded
578  chemPointISAT* nulPhi = 0;
579  forAll(tempList, i)
580  {
581  chemisTree().insertNewLeaf
582  (
583  tempList[i]->phi(),
584  tempList[i]->Rphi(),
585  tempList[i]->A(),
586  scaleFactor(),
587  tolerance_,
588  scaleFactor_.size(),
589  nActive,
590  nulPhi
591  );
592  deleteDemandDrivenData(tempList[i]);
593  }
594  }
595 
596  // The structure has been changed, it will force the binary tree to
597  // perform a new search and find the most appropriate point still stored
598  lastSearch_ = nullptr;
599  }
600 
601  // Compute the A matrix needed to store the chemPoint.
602  const label ASize = chemistry_.nEqns() + 1;
603  scalarSquareMatrix A(ASize, Zero);
604  computeA(A, Rphiq, li, deltaT);
605 
606  chemisTree().insertNewLeaf
607  (
608  phiq,
609  Rphiq,
610  A,
611  scaleFactor(),
612  tolerance_,
613  scaleFactor_.size(),
614  nActive,
615  lastSearch_ // lastSearch_ may be nullptr (handled by binaryTree)
616  );
617  if (lastSearch_ != nullptr)
618  {
619  addToMRU(lastSearch_);
620  }
621  nAdd_++;
622 
623  tabulationResults_[li] = 0;
624 
625  if (log_)
626  {
627  addNewLeafCpuTime_ += cpuTime_.cpuTimeIncrement();
628  }
629 
630  return growthOrAddFlag;
631 }
632 
633 
635 {
636  if (log_)
637  {
638  nRetrievedFile_()
639  << runTime_.userTimeValue() << " " << nRetrieved_ << endl;
640  nRetrieved_ = 0;
641 
642  nGrowthFile_()
643  << runTime_.userTimeValue() << " " << nGrowth_ << endl;
644  nGrowth_ = 0;
645 
646  nAddFile_()
647  << runTime_.userTimeValue() << " " << nAdd_ << endl;
648  nAdd_ = 0;
649 
650  sizeFile_()
651  << runTime_.userTimeValue() << " "
652  << chemisTree_.size() << endl;
653 
654  cpuRetrieveFile_()
655  << runTime_.userTimeValue()
656  << " " << searchISATCpuTime_ << endl;
657  searchISATCpuTime_ = 0;
658 
659  cpuGrowFile_()
660  << runTime_.userTimeValue()
661  << " " << growCpuTime_ << endl;
662  growCpuTime_ = 0;
663 
664  cpuAddFile_()
665  << runTime_.userTimeValue()
666  << " " << addNewLeafCpuTime_ << endl;
667  addNewLeafCpuTime_ = 0;
668  }
669 }
670 
671 
673 {
674  // Increment counter of time-step
675  timeSteps_++;
676 
677  forAll(tabulationResults_, i)
678  {
679  tabulationResults_[i] = 2;
680  }
681 }
682 
683 
685 {
686  bool updated = cleanAndBalance();
687  writePerformance();
688  return updated;
689 }
690 
691 
692 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
Macros for easy insertion into run-time selection tables.
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:78
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Definition: DynamicListI.H:296
An STL-conforming iterator.
Definition: LList.H:249
Template class for non-intrusive linked lists.
Definition: LList.H:76
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:63
virtual bool retrieve(const Foam::scalarField &phiq, scalarField &Rphiq)
Find the closest stored leaf of phiQ and store the result in.
Definition: ISAT.C:423
virtual label add(const scalarField &phiq, const scalarField &Rphiq, const label nActive, const label li, const scalar deltaT)
Add information to the tabulation.
Definition: ISAT.C:507
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
Extends base chemistry model adding an ODESystem and the reduction maps needed for tabulation.
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
addToRunTimeSelectionTable(chemistryTabulationMethod, ISAT, dictionary)
const dimensionedScalar phi0
Magnetic flux quantum: default SI units: [Wb].
const dimensionedScalar e
Elementary charge.
static const coefficient A("A", dimPressure, 611.21)
Namespace for OpenFOAM.
static const zero Zero
Definition: zero.H:97
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
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:258
void deleteDemandDrivenData(DataType *&dataPtr)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensionedScalar log(const dimensionedScalar &ds)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
SquareMatrix< scalar > scalarSquareMatrix
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
const label nSpecie
basicChemistryModel & chemistry
fluidMulticomponentThermo & thermo
Definition: createFields.H:31