DRGEP.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 "DRGEP.H"
27 #include "SortableListDRGEP.H"
28 
29 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
30 
31 template<class ThermoType>
33 (
34  const IOdictionary& dict,
36 )
37 :
39  searchInitSet_(),
40  sC_(this->nSpecie(), 0),
41  sH_(this->nSpecie(), 0),
42  sO_(this->nSpecie(), 0),
43  sN_(this->nSpecie(), 0),
44  NGroupBased_(50)
45 {
47 
48  const wordHashSet initSet(this->coeffsDict_.lookup("initialSet"));
49  forAllConstIter(wordHashSet, initSet, iter)
50  {
51  searchInitSet_.append(chemistry.mixture().species()[iter.key()]);
52  }
53 
54  if (this->coeffsDict_.found("NGroupBased"))
55  {
56  NGroupBased_ = this->coeffsDict_.template lookup<label>("NGroupBased");
57  }
58 
59  for (label i=0; i<this->nSpecie(); i++)
60  {
61  const List<specieElement>& curSpecieComposition =
62  chemistry.mixture().specieComposition(i);
63 
64  // for all elements in the current species
65  forAll(curSpecieComposition, j)
66  {
67  const specieElement& curElement =
68  curSpecieComposition[j];
69  if (curElement.name() == "C")
70  {
71  sC_[i] = curElement.nAtoms();
72  }
73  else if (curElement.name() == "H")
74  {
75  sH_[i] = curElement.nAtoms();
76  }
77  else if (curElement.name() == "O")
78  {
79  sO_[i] = curElement.nAtoms();
80  }
81  else if (curElement.name() == "N")
82  {
83  sN_[i] = curElement.nAtoms();
84  }
85  else
86  {
87  Info<< "element not considered"<< endl;
88  }
89  }
90  }
91 }
92 
93 
94 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
95 
96 template<class ThermoType>
98 {}
99 
100 
101 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
102 
103 template<class ThermoType>
105 (
106  const scalar p,
107  const scalar T,
108  const scalarField& c,
109  List<label>& ctos,
110  DynamicList<label>& stoc,
111  const label li
112 )
113 {
114  scalarField c1(this->chemistry_.nEqns(), 0.0);
115 
116  for (label i=0; i<this->nSpecie(); i++)
117  {
118  c1[i] = c[i];
119  }
120 
121  c1[this->nSpecie()] = T;
122  c1[this->nSpecie()+1] = p;
123 
124  // Compute the rAB matrix
125  RectangularMatrix<scalar> rABNum(this->nSpecie(),this->nSpecie(),0.0);
126  scalarField PA(this->nSpecie(),0.0);
127  scalarField CA(this->nSpecie(),0.0);
128 
129  // Number of initialised rAB for each lines
130  Field<label> NbrABInit(this->nSpecie(),0);
131  // Position of the initialised rAB, -1 when not initialised
132  RectangularMatrix<label> rABPos(this->nSpecie(), this->nSpecie(), -1);
133  // Index of the other species involved in the rABNum
134  RectangularMatrix<label> rABOtherSpec(this->nSpecie(), this->nSpecie(), -1);
135 
136  scalarField omegaV(this->chemistry_.reactions().size());
137  forAll(this->chemistry_.reactions(), i)
138  {
139  const Reaction<ThermoType>& R = this->chemistry_.reactions()[i];
140 
141  // for each reaction compute omegai
142  scalar omegaf, omegar;
143  const scalar omegai = R.omega(p, T, c1, li, omegaf, omegar);
144  omegaV[i] = omegai;
145 
146  // then for each pair of species composing this reaction,
147  // compute the rAB matrix (separate the numerator and
148  // denominator)
149 
150  DynamicList<scalar> wA(R.lhs().size()+R.rhs().size());
151  DynamicList<label> wAID(R.lhs().size()+R.rhs().size());
152  forAll(R.lhs(), s)// compute rAB for all species in the left hand side
153  {
154  label ss = R.lhs()[s].index;
155  scalar sl = -R.lhs()[s].stoichCoeff; // vAi = v''-v' => here -v'
156  List<bool> deltaBi(this->nSpecie(), false);
157  FIFOStack<label> usedIndex;
158  forAll(R.lhs(), j)
159  {
160  label sj = R.lhs()[j].index;
161  usedIndex.push(sj);
162  deltaBi[sj] = true;
163  }
164  forAll(R.rhs(), j)
165  {
166  label sj = R.rhs()[j].index;
167  usedIndex.push(sj);
168  deltaBi[sj] = true;
169  }
170 
171  // Disable for self reference (by definition rAA=0)
172  deltaBi[ss] = false;
173 
174  while(!usedIndex.empty())
175  {
176  label curIndex = usedIndex.pop();
177  if (deltaBi[curIndex])
178  {
179  // disable to avoid counting it more than once
180  deltaBi[curIndex] = false;
181  // test if this rAB is not initialised
182  if (rABPos(ss, curIndex)==-1)
183  {
184  rABPos(ss, curIndex) = NbrABInit[ss];
185  NbrABInit[ss]++;
186  rABNum(ss, rABPos(ss, curIndex)) = sl*omegai;
187  rABOtherSpec(ss, rABPos(ss, curIndex)) = curIndex;
188  }
189  else
190  {
191  rABNum(ss, rABPos(ss, curIndex)) += sl*omegai;
192  }
193  }
194  }
195  bool found(false);
196  forAll(wAID, id)
197  {
198  if (ss==wAID[id])
199  {
200  wA[id] += sl*omegai;
201  found = true;
202  }
203  }
204  if (!found)
205  {
206  wA.append(sl*omegai);
207  wAID.append(ss);
208  }
209  }
210 
211  // Compute rAB for all species in the right hand side
212  forAll(R.rhs(), s)
213  {
214  label ss = R.rhs()[s].index;
215  scalar sl = R.rhs()[s].stoichCoeff; // vAi = v''-v' => here v''
216  List<bool> deltaBi(this->nSpecie(), false);
217  FIFOStack<label> usedIndex;
218  forAll(R.lhs(), j)
219  {
220  label sj = R.lhs()[j].index;
221  usedIndex.push(sj);
222  deltaBi[sj] = true;
223  }
224  forAll(R.rhs(), j)
225  {
226  label sj = R.rhs()[j].index;
227  usedIndex.push(sj);
228  deltaBi[sj] = true;
229  }
230 
231  // Disable for self reference (by definition rAA=0)
232  deltaBi[ss] = false;
233 
234  while(!usedIndex.empty())
235  {
236  label curIndex = usedIndex.pop();
237  if (deltaBi[curIndex])
238  {
239  // disable to avoid counting it more than once
240  deltaBi[curIndex] = false;
241  // test if this rAB is not initialised
242  if (rABPos(ss, curIndex)==-1)
243  {
244  rABPos(ss, curIndex) = NbrABInit[ss];
245  NbrABInit[ss]++;
246  rABNum(ss, rABPos(ss, curIndex)) = sl*omegai;
247  rABOtherSpec(ss, rABPos(ss, curIndex)) = curIndex;
248  }
249  else
250  {
251  rABNum(ss, rABPos(ss, curIndex)) += sl*omegai;
252  }
253  }
254  }
255 
256  bool found(false);
257  forAll(wAID, id)
258  {
259  if (ss==wAID[id])
260  {
261  wA[id] += sl*omegai;
262  found = true;
263  }
264  }
265  if (!found)
266  {
267  wA.append(sl*omegai);
268  wAID.append(ss);
269  }
270  }
271 
272  wAID.shrink();
273  // Now that every species of the reactions has been visited, we can
274  // compute the production and consumption rate. This way, it avoids
275  // getting wrong results when species are present in both lhs and rhs
276  forAll(wAID, id)
277  {
278  if (wA[id] > 0.0)
279  {
280  if (PA[wAID[id]] == 0.0)
281  {
282  PA[wAID[id]] = wA[id];
283  }
284  else
285  {
286  PA[wAID[id]] += wA[id];
287  }
288  }
289  else
290  {
291  if (CA[wAID[id]] == 0.0)
292  {
293  CA[wAID[id]] = -wA[id];
294  }
295  else
296  {
297  CA[wAID[id]] += -wA[id];
298  }
299  }
300  }
301  }
302  // rii = 0.0 by definition
303 
304  // Compute the production rate of each element Pa
305  label nElements = 4; // 4 main elements (C, H, O, N)
306  scalarList Pa(nElements,0.0);
307  scalarList Ca(nElements,0.0);
308 
309  // for (label q=0; q<SIS.size(); q++)
310  for (label i=0; i<this->nSpecie(); i++)
311  {
312  Pa[0] += sC_[i]*max(0.0,(PA[i]-CA[i]));
313  Ca[0] += sC_[i]*max(0.0,-(PA[i]-CA[i]));
314  Pa[1] += sH_[i]*max(0.0,(PA[i]-CA[i]));
315  Ca[1] += sH_[i]*max(0.0,-(PA[i]-CA[i]));
316  Pa[2] += sO_[i]*max(0.0,(PA[i]-CA[i]));
317  Ca[2] += sO_[i]*max(0.0,-(PA[i]-CA[i]));
318  Pa[3] += sN_[i]*max(0.0,(PA[i]-CA[i]));
319  Ca[3] += sN_[i]*max(0.0,-(PA[i]-CA[i]));
320  }
321 
322  // Using the rAB matrix (numerator and denominator separated)
323  // compute the R value according to the search initiating set
324  scalarField Rvalue(this->nSpecie(),0.0);
325  List<bool> disabledSpecies(this->nSpecie(),false);
326 
327  // set all species to inactive and activate them according
328  // to rAB and initial set
329  label NActiveSpecies = 0;
330  for (label i=0; i<this->nSpecie(); i++)
331  {
332  this->activeSpecies_[i] = false;
333  }
334 
335  // Initialise the FIFOStack for search set
337  const labelList& SIS(this->searchInitSet_);
338  DynamicList<label> QStart(SIS.size());
339  DynamicList<scalar> alphaQ(SIS.size());
340 
341  // Compute the alpha coefficient and initialise the R value of the species
342  // in the SIS
343  for (label i=0; i<SIS.size(); i++)
344  {
345  label q = SIS[i];
346  // compute alpha
347  scalar alphaA(0.0);
348  // C
349  if (Pa[0] > vSmall)
350  {
351  scalar alphaTmp = (sC_[q]*mag(PA[q]-CA[q])/Pa[0]);
352  if (alphaTmp > alphaA)
353  {
354  alphaA = alphaTmp;
355  }
356  }
357  // H
358  if (Pa[1] > vSmall)
359  {
360  scalar alphaTmp = (sH_[q]*mag(PA[q]-CA[q])/Pa[1]);
361  if (alphaTmp > alphaA)
362  {
363  alphaA = alphaTmp;
364  }
365  }
366  // O
367  if (Pa[2] > vSmall)
368  {
369  scalar alphaTmp = (sO_[q]*mag(PA[q]-CA[q])/Pa[2]);
370  if (alphaTmp > alphaA)
371  {
372  alphaA = alphaTmp;
373  }
374  }
375  // N
376  if (Pa[3] > vSmall)
377  {
378  scalar alphaTmp = (sN_[q]*mag(PA[q]-CA[q])/Pa[3]);
379  if (alphaTmp > alphaA)
380  {
381  alphaA = alphaTmp;
382  }
383  }
384  if (alphaA > this->tolerance())
385  {
386  this->activeSpecies_[q] = true;
387  NActiveSpecies ++;
388  Q.push(q);
389  QStart.append(q);
390  alphaQ.append(1.0);
391  Rvalue[q] = 1.0;
392  }
393  else
394  {
395  Rvalue[q] = alphaA;
396  }
397  }
398 
399  // if all species from the SIS has been removed
400  // force the use of the species with maximum Rvalue
401  if (Q.empty())
402  {
403  scalar Rmax=0.0;
404  label specID=-1;
405  forAll(SIS, i)
406  {
407  if (Rvalue[SIS[i]] > Rmax)
408  {
409  Rmax = Rvalue[SIS[i]];
410  specID=SIS[i];
411  }
412  }
413  Q.push(specID);
414  QStart.append(specID);
415  alphaQ.append(1.0);
416  NActiveSpecies ++;
417  Rvalue[specID] = 1.0;
418  this->activeSpecies_[specID] = true;
419  }
420 
421  // Execute the main loop for R-value
422  while (!Q.empty())
423  {
424  label u = Q.pop();
425  scalar Den = max(PA[u],CA[u]);
426  if (Den > vSmall)
427  {
428  for (label v=0; v<NbrABInit[u]; v++)
429  {
430  label otherSpec = rABOtherSpec(u, v);
431  scalar rAB = mag(rABNum(u, v))/Den;
432  if (rAB > 1)
433  {
434  rAB = 1;
435  }
436 
437  scalar Rtemp = Rvalue[u]*rAB;
438  // a link analysed previously is stronger
439  if (Rvalue[otherSpec] < Rtemp)
440  {
441  Rvalue[otherSpec] = Rtemp;
442  // the (composed) link is stronger than the tolerance
443  if (Rtemp >= this->tolerance())
444  {
445  Q.push(otherSpec);
446  if (!this->activeSpecies_[otherSpec])
447  {
448  this->activeSpecies_[otherSpec] = true;
449  NActiveSpecies ++;
450  }
451  }
452  }
453  }
454  }
455  }
456 
457  // Group-based reduction
458  // number of species disabled in the first step
459  label NDisabledSpecies(this->nSpecie() - NActiveSpecies);
460 
461  // while the number of removed species is greater than NGroupBased, the rAB
462  // are reevaluated according to the group based definition for each loop the
463  // temporary disabled species (in the first reduction) are sorted to disable
464  // definitely the NGroupBased species with lower R then these R value a
465  // reevaluated taking into account these disabled species
466  while(NDisabledSpecies > NGroupBased_)
467  {
468  // keep track of disabled species using sortablelist to extract only
469  // NGroupBased lower R value
470  SortableListDRGEP<scalar> Rdisabled(NDisabledSpecies);
471  labelList Rindex(NDisabledSpecies);
472  label nD = 0;
473  forAll(disabledSpecies, i)
474  {
475  // if just disabled and not in a previous loop
476  if (!this->activeSpecies_[i] && !disabledSpecies[i])
477  {
478  // Note: non-reached species will be removed first (Rvalue=0)
479  Rdisabled[nD] = Rvalue[i];
480  Rindex[nD++] = i;
481  }
482  }
483  // sort the Rvalue to obtain the NGroupBased lower R value
484  Rdisabled.partialSort(NGroupBased_);
485  labelList tmpIndex(Rdisabled.indices());
486 
487  // disable definitely NGroupBased species in this loop
488  for (label i=0; i<NGroupBased_; i++)
489  {
490  disabledSpecies[Rindex[tmpIndex[i]]] = true;
491  }
492  NDisabledSpecies -= NGroupBased_;
493 
494  // reevaluate the rAB according to the group-based definition rAB{S} [1]
495  // only update the numerator
496  forAll(NbrABInit, i)
497  {
498  for (label v=0; v<NbrABInit[i]; v++)
499  {
500  rABNum(i, v) = 0.0;
501  }
502  }
503  forAll(this->chemistry_.reactions(), i)
504  {
505  const Reaction<ThermoType>& R = this->chemistry_.reactions()[i];
506 
507  scalar omegai = omegaV[i];
508 
509  forAll(R.lhs(), s)
510  {
511  label ss = R.lhs()[s].index;
512  scalar sl = -R.lhs()[s].stoichCoeff; // vAi = v''-v' => here -v'
513  List<bool> deltaBi(this->nSpecie(), false);
514  bool alreadyDisabled(false);
515  FIFOStack<label> usedIndex;
516  forAll(R.lhs(), j)
517  {
518  label sj = R.lhs()[j].index;
519  usedIndex.push(sj);
520  deltaBi[sj] = true;
521  if (disabledSpecies[sj])
522  {
523  alreadyDisabled=true;
524  }
525  }
526  forAll(R.rhs(), j)
527  {
528  label sj = R.rhs()[j].index;
529  usedIndex.push(sj);
530  deltaBi[sj] = true;
531  if (disabledSpecies[sj])
532  {
533  alreadyDisabled=true;
534  }
535  }
536 
537  deltaBi[ss] = false;
538 
539  if (alreadyDisabled)
540  {
541  // if one of the species in this reaction is disabled, all
542  // species connected to species ss are modified
543  for (label v=0; v<NbrABInit[ss]; v++)
544  {
545  rABNum(ss, v) += sl*omegai;
546  }
547  }
548  else
549  {
550  while(!usedIndex.empty())
551  {
552  label curIndex = usedIndex.pop();
553  if (deltaBi[curIndex])
554  {
555  // disable to avoid counting it more than once
556  deltaBi[curIndex] = false;
557  rABNum(ss, rABPos(ss, curIndex)) += sl*omegai;
558  }
559  }
560  }
561  }
562 
563  forAll(R.rhs(), s)
564  {
565  label ss = R.rhs()[s].index;
566  scalar sl = R.rhs()[s].stoichCoeff; // vAi = v''-v' => here v''
567  List<bool> deltaBi(this->nSpecie(), false);
568  bool alreadyDisabled(false);
569  FIFOStack<label> usedIndex;
570  forAll(R.lhs(), j)
571  {
572  label sj = R.lhs()[j].index;
573  usedIndex.push(sj);
574  deltaBi[sj] = true;
575  if (disabledSpecies[sj])
576  {
577  alreadyDisabled=true;
578  }
579  }
580  forAll(R.rhs(), j)
581  {
582  label sj = R.rhs()[j].index;
583  usedIndex.push(sj);
584  deltaBi[sj] = true;
585  if (disabledSpecies[sj])
586  {
587  alreadyDisabled=true;
588  }
589  }
590 
591  deltaBi[ss] = false;
592 
593  if (alreadyDisabled)
594  {
595  // if one of the species in this reaction is disabled, all
596  // species connected to species ss are modified
597  for (label v=0; v<NbrABInit[ss]; v++)
598  {
599  rABNum(ss, v) += sl*omegai;
600  }
601  }
602  else
603  {
604  while(!usedIndex.empty())
605  {
606  label curIndex = usedIndex.pop();
607  if (deltaBi[curIndex])
608  {
609  deltaBi[curIndex] = false;
610  rABNum(ss, rABPos(ss, curIndex)) += sl*omegai;
611  }
612  }
613  }
614  }
615  }
616 
617  forAll(QStart, qi)
618  {
619  label u = QStart[qi];
620  Q.push(u);
621  }
622 
623  while (!Q.empty())
624  {
625  label u = Q.pop();
626  scalar Den = max(PA[u],CA[u]);
627  if (Den!=0.0)
628  {
629  for (label v=0; v<NbrABInit[u]; v++)
630  {
631  label otherSpec = rABOtherSpec(u, v);
632  if (!disabledSpecies[otherSpec])
633  {
634  scalar rAB = mag(rABNum(u, v))/Den;
635  if (rAB > 1)
636  {
637  rAB = 1;
638  }
639 
640  scalar Rtemp = Rvalue[u]*rAB;
641  // a link analysed previously is stronger
642  if (Rvalue[otherSpec] < Rtemp)
643  {
644  Rvalue[otherSpec] = Rtemp;
645  if (Rtemp >= this->tolerance())
646  {
647  Q.push(otherSpec);
648  if (!this->activeSpecies_[otherSpec])
649  {
650  this->activeSpecies_[otherSpec] = true;
651  NDisabledSpecies--;
652  }
653  }
654  }
655  }
656  }
657  }
658  }
659  }
660 
662 }
663 
664 
665 // ************************************************************************* //
A FIFO stack based on a singly-linked list.
Definition: FIFOStack.H:51
A HashTable with keys but without contents.
Definition: HashSet.H:59
dictionary dict
virtual void reduceMechanism(const scalar p, const scalar T, const scalarField &c, List< label > &ctos, DynamicList< label > &stoc, const label li)
Reduce the mechanism.
Definition: DRGEP.C:105
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
scalar omega(const scalar p, const scalar T, const scalarField &c, const label li, scalar &omegaf, scalar &omegar) const
Net reaction rate.
Definition: Reaction.C:319
virtual ~DRGEP()
Destructor.
Definition: DRGEP.C:97
DRGEP(const IOdictionary &dict, chemistryModel< ThermoType > &chemistry)
Construct from components.
Definition: DRGEP.C:33
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
A list that is sorted upon construction or when explicitly requested with the sort() method...
const label nSpecie
const List< specieCoeffs > & lhs() const
Return the components of the left hand side.
Definition: reactionI.H:36
Simple extension of ReactionThermo to handle reaction kinetics in addition to the equilibrium thermod...
Definition: Reaction.H:55
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:53
const dimensionedScalar c1
First radiation constant: default SI units: [W/m^2].
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects...
Definition: DynamicList.H:56
const word & name() const
Return the name of the element.
label nAtoms() const
Return the number of atoms of this element in the specie.
basicChemistryModel & chemistry
const multiComponentMixture< ThermoType > & mixture() const
Return reference to the mixture.
const volScalarField & T
Extends base chemistry model by adding a thermo package, and ODE functions. Introduces chemistry equa...
const List< specieCoeffs > & rhs() const
Return the components of the right hand side.
Definition: reactionI.H:42
#define R(A, B, C, D, E, F, K, M)
void push(const T &a)
Push an element onto the stack.
Definition: FIFOStack.H:84
messageStream Info
const labelList & indices() const
Return the list of sorted indices. Updated every sort.
dimensioned< scalar > mag(const dimensioned< Type > &)
T pop()
Pop the bottom element off the stack.
Definition: FIFOStack.H:90
An abstract class for methods of chemical mechanism reduction.
volScalarField & p
void partialSort(int M)
Partial sort the list (if changed after construction time)
bool found