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