ISAT.C
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  Copyright (C) 2019 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "ISAT.H"
30 #include "LUscalarMatrix.H"
31 
32 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
33 
34 template<class CompType, class ThermoType>
36 (
37  const dictionary& chemistryProperties,
39 )
40 :
42  (
43  chemistryProperties,
44  chemistry
45  ),
46  chemisTree_(chemistry, this->coeffsDict_),
47  scaleFactor_(chemistry.nEqns() + ((this->variableTimeStep()) ? 1 : 0), 1),
48  runTime_(chemistry.time()),
49  chPMaxLifeTime_
50  (
51  this->coeffsDict_.lookupOrDefault("chPMaxLifeTime", INT_MAX)
52  ),
53  maxGrowth_(this->coeffsDict_.lookupOrDefault("maxGrowth", INT_MAX)),
54  checkEntireTreeInterval_
55  (
56  this->coeffsDict_.lookupOrDefault("checkEntireTreeInterval", INT_MAX)
57  ),
58  maxDepthFactor_
59  (
60  this->coeffsDict_.lookupOrDefault
61  (
62  "maxDepthFactor",
63  (chemisTree_.maxNLeafs() - 1)
64  /(log(scalar(chemisTree_.maxNLeafs()))/log(2.0))
65  )
66  ),
67  minBalanceThreshold_
68  (
69  this->coeffsDict_.lookupOrDefault
70  (
71  "minBalanceThreshold",0.1*chemisTree_.maxNLeafs()
72  )
73  ),
74  MRURetrieve_(this->coeffsDict_.lookupOrDefault("MRURetrieve", false)),
75  maxMRUSize_(this->coeffsDict_.lookupOrDefault("maxMRUSize", 0)),
76  lastSearch_(nullptr),
77  growPoints_(this->coeffsDict_.lookupOrDefault("growPoints", true)),
78  nRetrieved_(0),
79  nGrowth_(0),
80  nAdd_(0),
81  cleaningRequired_(false)
82 {
83  if (this->active_)
84  {
85  dictionary scaleDict(this->coeffsDict_.subDict("scaleFactor"));
86  const label Ysize = this->chemistry_.Y().size();
87  const scalar otherScaleFactor = scaleDict.get<scalar>("otherSpecies");
88  for (label i=0; i<Ysize; ++i)
89  {
90  if (!scaleDict.found(this->chemistry_.Y()[i].member()))
91  {
92  scaleFactor_[i] = otherScaleFactor;
93  }
94  else
95  {
96  scaleFactor_[i] =
97  scaleDict.get<scalar>(this->chemistry_.Y()[i].member());
98  }
99  }
100 
101  scaleDict.readEntry("Temperature", scaleFactor_[Ysize]);
102  scaleDict.readEntry("Pressure", scaleFactor_[Ysize + 1]);
103 
104  if (this->variableTimeStep())
105  {
106  scaleDict.readEntry("deltaT", scaleFactor_[Ysize + 2]);
107  }
108  }
109 
110  if (this->variableTimeStep())
111  {
112  nAdditionalEqns_ = 3;
113  }
114  else
115  {
116  nAdditionalEqns_ = 2;
117  }
118 
119  if (this->log())
120  {
121  nRetrievedFile_ = chemistry.logFile("found_isat.out");
122  nGrowthFile_ = chemistry.logFile("growth_isat.out");
123  nAddFile_ = chemistry.logFile("add_isat.out");
124  sizeFile_ = chemistry.logFile("size_isat.out");
125  }
126 }
127 
128 
129 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
130 
131 template<class CompType, class ThermoType>
133 {}
134 
135 
136 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
137 
138 template<class CompType, class ThermoType>
140 (
142 )
143 {
144  if (maxMRUSize_ > 0 && MRURetrieve_)
145  {
146  auto iter = MRUList_.begin();
147 
148  // First search if the chemPoint is already in the list
149  bool isInList = false;
150  for (; iter.good(); ++iter)
151  {
152  if (*iter == phi0)
153  {
154  isInList = true;
155  break;
156  }
157  }
158 
159  if (isInList)
160  {
161  // If it is in the list, then move it to front
162  if (iter != MRUList_.begin())
163  {
164  MRUList_.remove(iter);
165  MRUList_.insert(phi0);
166  }
167  }
168  else
169  {
170  // chemPoint not yet in the list, iter is last
171  if (MRUList_.size() == maxMRUSize_)
172  {
173  if (iter == MRUList_.end())
174  {
175  MRUList_.remove(iter);
176  MRUList_.insert(phi0);
177  }
178  else
179  {
181  << "Error in MRUList construction"
182  << exit(FatalError);
183  }
184  }
185  else
186  {
187  MRUList_.insert(phi0);
188  }
189  }
190  }
191 }
192 
193 
194 template<class CompType, class ThermoType>
196 (
197  chemPointISAT<CompType, ThermoType>* phi0,
198  const scalarField& phiq,
199  scalarField& Rphiq
200 )
201 {
202  label nEqns = this->chemistry_.nEqns(); // Species, T, p
203  bool mechRedActive = this->chemistry_.mechRed()->active();
204  Rphiq = phi0->Rphi();
205  scalarField dphi(phiq-phi0->phi());
206  const scalarSquareMatrix& gradientsMatrix = phi0->A();
207  List<label>& completeToSimplified(phi0->completeToSimplifiedIndex());
208 
209  // Rphiq[i]=Rphi0[i]+A(i, j)dphi[j]
210  // where Aij is dRi/dphi_j
211  for (label i=0; i<nEqns-nAdditionalEqns_; ++i)
212  {
213  if (mechRedActive)
214  {
215  label si = completeToSimplified[i];
216  // The species is active
217  if (si != -1)
218  {
219  for (label j=0; j<nEqns-2; j++)
220  {
221  label sj = completeToSimplified[j];
222  if (sj != -1)
223  {
224  Rphiq[i] += gradientsMatrix(si, sj)*dphi[j];
225  }
226  }
227  Rphiq[i] +=
228  gradientsMatrix(si, phi0->nActiveSpecies())*dphi[nEqns - 2];
229  Rphiq[i] +=
230  gradientsMatrix(si, phi0->nActiveSpecies() + 1)
231  *dphi[nEqns - 1];
232 
233  if (this->variableTimeStep())
234  {
235  Rphiq[i] +=
236  gradientsMatrix(si, phi0->nActiveSpecies() + 2)
237  *dphi[nEqns];
238  }
239 
240  // As we use an approximation of A, Rphiq should be checked for
241  // negative values
242  Rphiq[i] = max(0.0, Rphiq[i]);
243  }
244  // The species is not active A(i, j) = I(i, j)
245  else
246  {
247  Rphiq[i] += dphi[i];
248  Rphiq[i] = max(0.0, Rphiq[i]);
249  }
250  }
251  else // Mechanism reduction is not active
252  {
253  for (label j=0; j<nEqns; ++j)
254  {
255  Rphiq[i] += gradientsMatrix(i, j)*dphi[j];
256  }
257  // As we use a first order gradient matrix, Rphiq should be checked
258  // for negative values
259  Rphiq[i] = max(0.0, Rphiq[i]);
260  }
261  }
262 }
263 
264 
265 template<class CompType, class ThermoType>
267 (
268  chemPointISAT<CompType, ThermoType>* phi0,
269  const scalarField& phiq,
270  const scalarField& Rphiq
271 )
272 {
273  // If the pointer to the chemPoint is nullptr, the function stops
274  if (phi0 == nullptr)
275  {
276  return false;
277  }
278 
279  // Raise a flag when the chemPoint used has been grown more than the
280  // allowed number of time
281  if (phi0->nGrowth() > maxGrowth_)
282  {
283  cleaningRequired_ = true;
284  phi0->toRemove() = true;
285  return false;
286  }
287 
288  // If the solution RphiQ is still within the tolerance we try to grow it
289  // in some cases this might result in a failure and the grow function of
290  // the chemPoint returns false
291  if (phi0->checkSolution(phiq, Rphiq))
292  {
293  return phi0->grow(phiq);
294  }
295 
296  // The actual solution and the approximation given by ISAT are too different
297  return false;
298 }
299 
300 
301 template<class CompType, class ThermoType>
302 bool
304 {
305  bool treeModified(false);
306 
307  // Check all chemPoints to see if we need to delete some of the chemPoints
308  // according to the elapsed time and number of growths
309  chemPointISAT<CompType, ThermoType>* x = chemisTree_.treeMin();
310  while (x != nullptr)
311  {
312  chemPointISAT<CompType, ThermoType>* xtmp =
313  chemisTree_.treeSuccessor(x);
314 
315  scalar elapsedTimeSteps = this->chemistry_.timeSteps() - x->timeTag();
316 
317  if ((elapsedTimeSteps > chPMaxLifeTime_) || (x->nGrowth() > maxGrowth_))
318  {
319  chemisTree_.deleteLeaf(x);
320  treeModified = true;
321  }
322  x = xtmp;
323  }
324  // Check if the tree should be balanced according to criterion:
325  // - the depth of the tree bigger than a*log2(size), log2(size) being the
326  // ideal depth (e.g. 4 leafs can be stored in a tree of depth 2)
327  if
328  (
329  chemisTree_.size() > minBalanceThreshold_
330  && chemisTree_.depth() >
331  maxDepthFactor_*log(scalar(chemisTree_.size()))/log(2.0)
332  )
333  {
334  chemisTree_.balance();
335  MRUList_.clear();
336  treeModified = true;
337  }
338 
339  // Return a bool to specify if the tree structure has been modified and is
340  // now below the user specified limit (true if not full)
341  return (treeModified && !chemisTree_.isFull());
342 }
343 
344 
345 template<class CompType, class ThermoType>
347 (
349  const scalarField& Rphiq,
350  const scalar rhoi,
351  const scalar dt
352 )
353 {
354  bool mechRedActive = this->chemistry_.mechRed()->active();
355  label speciesNumber = this->chemistry_.nSpecie();
356  scalarField Rcq(this->chemistry_.nEqns() + nAdditionalEqns_ - 2);
357  for (label i=0; i<speciesNumber; ++i)
358  {
359  label s2c = i;
360  if (mechRedActive)
361  {
362  s2c = this->chemistry_.simplifiedToCompleteIndex()[i];
363  }
364  Rcq[i] = rhoi*Rphiq[s2c]/this->chemistry_.specieThermo()[s2c].W();
365  }
366  Rcq[speciesNumber] = Rphiq[Rphiq.size() - nAdditionalEqns_];
367  Rcq[speciesNumber+1] = Rphiq[Rphiq.size() - nAdditionalEqns_ + 1];
368  if (this->variableTimeStep())
369  {
370  Rcq[speciesNumber + 2] = Rphiq[Rphiq.size() - nAdditionalEqns_ + 2];
371  }
372 
373  // Aaa is computed implicitly,
374  // A is given by A = C(psi0, t0+dt), where C is obtained through solving
375  // d/dt C(psi0,t) = J(psi(t))C(psi0,t)
376  // If we solve it implicitly:
377  // (C(psi0, t0+dt) - C(psi0,t0))/dt = J(psi(t0+dt))C(psi0,t0+dt)
378  // The Jacobian is thus computed according to the mapping
379  // C(psi0,t0+dt)*(I-dt*J(psi(t0+dt))) = C(psi0, t0)
380  // A = C(psi0,t0)/(I-dt*J(psi(t0+dt)))
381  // where C(psi0,t0) = I
382 
383  this->chemistry_.jacobian(runTime_.value(), Rcq, A);
384 
385  // The jacobian is computed according to the molar concentration
386  // the following conversion allows the code to use A with mass fraction
387  for (label i=0; i<speciesNumber; ++i)
388  {
389  label si = i;
390 
391  if (mechRedActive)
392  {
393  si = this->chemistry_.simplifiedToCompleteIndex()[i];
394  }
395 
396  for (label j=0; j<speciesNumber; ++j)
397  {
398  label sj = j;
399  if (mechRedActive)
400  {
401  sj = this->chemistry_.simplifiedToCompleteIndex()[j];
402  }
403  A(i, j) *=
404  -dt*this->chemistry_.specieThermo()[si].W()
405  /this->chemistry_.specieThermo()[sj].W();
406  }
407 
408  A(i, i) += 1;
409  // Columns for pressure and temperature
410  A(i, speciesNumber) *=
411  -dt*this->chemistry_.specieThermo()[si].W()/rhoi;
412  A(i, speciesNumber+1) *=
413  -dt*this->chemistry_.specieThermo()[si].W()/rhoi;
414  }
415 
416  // For temperature and pressure, only unity on the diagonal
417  A(speciesNumber, speciesNumber) = 1;
418  A(speciesNumber + 1, speciesNumber + 1) = 1;
419  if (this->variableTimeStep())
420  {
421  A[speciesNumber + 2][speciesNumber + 2] = 1;
422  }
423 
424  // Inverse of (I-dt*J(psi(t0+dt)))
425  LUscalarMatrix LUA(A);
426  LUA.inv(A);
427 }
428 
429 
430 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
431 
432 template<class CompType, class ThermoType>
434 (
435  const Foam::scalarField& phiq,
436  scalarField& Rphiq
437 )
438 {
439  bool retrieved(false);
441 
442  // If the tree is not empty
443  if (chemisTree_.size())
444  {
445  chemisTree_.binaryTreeSearch(phiq, chemisTree_.root(), phi0);
446 
447  // lastSearch keeps track of the chemPoint we obtain by the regular
448  // binary tree search
449  lastSearch_ = phi0;
450  if (phi0->inEOA(phiq))
451  {
452  retrieved = true;
453  }
454  // After a successful secondarySearch, phi0 store a pointer to the
455  // found chemPoint
456  else if (chemisTree_.secondaryBTSearch(phiq, phi0))
457  {
458  retrieved = true;
459  }
460  else if (MRURetrieve_)
461  {
462  forAllConstIters(MRUList_, iter)
463  {
464  phi0 = iter();
465  if (phi0->inEOA(phiq))
466  {
467  retrieved = true;
468  break;
469  }
470  }
471  }
472  }
473  // The tree is empty, retrieved is still false
474  else
475  {
476  // There is no chempoints that we can try to grow
477  lastSearch_ = nullptr;
478  }
479 
480  if (retrieved)
481  {
482  phi0->increaseNumRetrieve();
483  scalar elapsedTimeSteps =
484  this->chemistry_.timeSteps() - phi0->timeTag();
485 
486  // Raise a flag when the chemPoint has been used more than the allowed
487  // number of time steps
488  if (elapsedTimeSteps > chPMaxLifeTime_ && !phi0->toRemove())
489  {
490  cleaningRequired_ = true;
491  phi0->toRemove() = true;
492  }
493  lastSearch_->lastTimeUsed() = this->chemistry_.timeSteps();
494  addToMRU(phi0);
495  calcNewC(phi0,phiq, Rphiq);
496  ++nRetrieved_;
497  return true;
498  }
499 
500 
501  // This point is reached when every retrieve trials have failed
502  // or if the tree is empty
503  return false;
504 }
505 
506 
507 template<class CompType, class ThermoType>
509 (
510  const scalarField& phiq,
511  const scalarField& Rphiq,
512  const scalar rho,
513  const scalar deltaT
514 )
515 {
516  label growthOrAddFlag = 1;
517 
518  // If lastSearch_ holds a valid pointer to a chemPoint AND the growPoints_
519  // option is on, the code first tries to grow the point hold by lastSearch_
520  if (lastSearch_ && growPoints_)
521  {
522  if (grow(lastSearch_,phiq, Rphiq))
523  {
524  ++nGrowth_;
525  growthOrAddFlag = 0;
526 
527  // The structure of the tree is not modified, return false
528  return growthOrAddFlag;
529  }
530  }
531 
532  // If the code reach this point, it is either because lastSearch_ is not
533  // valid, OR because growPoints_ is not on, OR because the grow operation
534  // has failed. In the three cases, a new point is added to the tree.
535  if (chemisTree().isFull())
536  {
537  // If cleanAndBalance operation do not result in a reduction of the tree
538  // size, the last possibility is to delete completely the tree.
539  // It can be partially rebuilt with the MRU list if this is used.
540  if (!cleanAndBalance())
541  {
543  if (maxMRUSize_>0)
544  {
545  // Create a copy of each chemPointISAT of the MRUList_ before
546  // they are deleted
547  forAllConstIters(MRUList_, iter)
548  {
549  tempList.append
550  (
552  );
553  }
554  }
555  chemisTree().clear();
556 
557  // Pointers to chemPoint are not valid anymore, clear the list
558  MRUList_.clear();
559 
560  // Construct the tree without giving a reference to attach to it
561  // since the structure has been completely discarded
562  chemPointISAT<CompType, ThermoType>* nulPhi = nullptr;
563  for (auto& t : tempList)
564  {
565  chemisTree().insertNewLeaf
566  (
567  t->phi(),
568  t->Rphi(),
569  t->A(),
570  scaleFactor(),
571  this->tolerance(),
572  scaleFactor_.size(),
573  nulPhi
574  );
576  }
577  }
578 
579  // The structure has been changed, it will force the binary tree to
580  // perform a new search and find the most appropriate point still stored
581  lastSearch_ = nullptr;
582  }
583 
584  // Compute the A matrix needed to store the chemPoint.
585  label ASize = this->chemistry_.nEqns() + nAdditionalEqns_ - 2;
586  scalarSquareMatrix A(ASize, Zero);
587  computeA(A, Rphiq, rho, deltaT);
588 
589  chemisTree().insertNewLeaf
590  (
591  phiq,
592  Rphiq,
593  A,
594  scaleFactor(),
595  this->tolerance(),
596  scaleFactor_.size(),
597  lastSearch_ // lastSearch_ may be nullptr (handled by binaryTree)
598  );
599 
600  ++nAdd_;
601 
602  return growthOrAddFlag;
603 }
604 
605 
606 template<class CompType, class ThermoType>
607 void
609 {
610  if (this->log())
611  {
612  nRetrievedFile_()
613  << runTime_.timeOutputValue() << " " << nRetrieved_ << endl;
614  nRetrieved_ = 0;
615 
616  nGrowthFile_()
617  << runTime_.timeOutputValue() << " " << nGrowth_ << endl;
618  nGrowth_ = 0;
619 
620  nAddFile_()
621  << runTime_.timeOutputValue() << " " << nAdd_ << endl;
622  nAdd_ = 0;
623 
624  sizeFile_()
625  << runTime_.timeOutputValue() << " " << this->size() << endl;
626  }
627 }
628 
629 
630 // ************************************************************************* //
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
Foam::chemistryTabulationMethods::ISAT
Implementation of the ISAT (In-situ adaptive tabulation), for chemistry calculation.
Definition: ISAT.H:59
Foam::Zero
static constexpr const zero Zero
Global zero.
Definition: zero.H:128
Foam::DynamicList
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:57
A
static const Foam::dimensionedScalar A("", Foam::dimPressure, 611.21)
chemistry
BasicChemistryModel< psiReactionThermo > & chemistry
Definition: createFieldRefs.H:1
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:337
LUscalarMatrix.H
Foam::dictionary::get
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:81
Foam::chemistryTabulationMethods::ISAT::add
virtual label add(const scalarField &phiq, const scalarField &Rphiq, const scalar rho, const scalar deltaT)
Add information to the tabulation.
Definition: ISAT.C:509
rho
rho
Definition: readInitialConditions.H:96
Foam::deleteDemandDrivenData
void deleteDemandDrivenData(DataPtr &dataPtr)
Definition: demandDrivenData.H:42
Foam::label
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:62
Foam::Field< scalar >
Foam::chemPointISAT
Leaf of the binary tree. The chemPoint stores the composition 'phi', the mapping of this composition ...
Definition: chemPointISAT.H:142
Foam::DynamicList::append
DynamicList< T, SizeMin > & append(const T &val)
Append an element to the end of this list.
Definition: DynamicListI.H:472
Foam::constant::electromagnetic::phi0
const dimensionedScalar phi0
Magnetic flux quantum: default SI units: [Wb].
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
Foam::chemistryTabulationMethods::ISAT::writePerformance
virtual void writePerformance()
Definition: ISAT.C:608
Foam::scalarSquareMatrix
SquareMatrix< scalar > scalarSquareMatrix
Definition: scalarMatrices.H:57
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
Foam::log
dimensionedScalar log(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:262
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::SquareMatrix< scalar >
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:355
forAllConstIters
forAllConstIters(mixture.phases(), phase)
Definition: pEqn.H:28
Foam::chemistryTabulationMethods::ISAT::retrieve
virtual bool retrieve(const Foam::scalarField &phiq, scalarField &Rphiq)
Definition: ISAT.C:434
x
x
Definition: LISASMDCalcMethod2.H:52
Foam::chemistryTabulationMethods::ISAT::~ISAT
virtual ~ISAT()
Definition: ISAT.C:132
Foam::TDACChemistryModel< CompType, ThermoType >
ISAT.H
Foam::chemistryTabulationMethod::active
bool active()
Definition: chemistryTabulationMethod.H:124
Foam::chemistryTabulationMethod
An abstract class for chemistry tabulation.
Definition: chemistryTabulationMethod.H:59