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 {
46  const wordHashSet initSet(this->coeffsDict_.lookup("initialSet"));
47  forAllConstIter(wordHashSet, initSet, iter)
48  {
49  searchInitSet_.append(chemistry.mixture().species()[iter.key()]);
50  }
51 
52  if (this->coeffsDict_.found("NGroupBased"))
53  {
54  NGroupBased_ = this->coeffsDict_.template lookup<label>("NGroupBased");
55  }
56 
57  for (label i=0; i<this->nSpecie_; i++)
58  {
59  const List<specieElement>& curSpecieComposition =
60  chemistry.mixture().specieComposition(i);
61 
62  // for all elements in the current species
63  forAll(curSpecieComposition, j)
64  {
65  const specieElement& curElement =
66  curSpecieComposition[j];
67  if (curElement.name() == "C")
68  {
69  sC_[i] = curElement.nAtoms();
70  }
71  else if (curElement.name() == "H")
72  {
73  sH_[i] = curElement.nAtoms();
74  }
75  else if (curElement.name() == "O")
76  {
77  sO_[i] = curElement.nAtoms();
78  }
79  else if (curElement.name() == "N")
80  {
81  sN_[i] = curElement.nAtoms();
82  }
83  else
84  {
85  Info<< "element not considered"<< endl;
86  }
87  }
88  }
89 }
90 
91 
92 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
93 
94 template<class ThermoType>
96 {}
97 
98 
99 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
100 
101 template<class ThermoType>
103 (
104  const scalar p,
105  const scalar T,
106  const scalarField& c,
107  const label li
108 )
109 {
110  scalarField& completeC(this->chemistry_.completeC());
111  scalarField c1(this->chemistry_.nEqns(), 0.0);
112 
113  for (label i=0; i<this->nSpecie_; i++)
114  {
115  c1[i] = c[i];
116  completeC[i] = c[i];
117  }
118 
119  c1[this->nSpecie_] = T;
120  c1[this->nSpecie_+1] = p;
121 
122  // Compute the rAB matrix
123  RectangularMatrix<scalar> rABNum(this->nSpecie_,this->nSpecie_,0.0);
124  scalarField PA(this->nSpecie_,0.0);
125  scalarField CA(this->nSpecie_,0.0);
126 
127  // Number of initialised rAB for each lines
128  Field<label> NbrABInit(this->nSpecie_,0);
129  // Position of the initialised rAB, -1 when not initialised
130  RectangularMatrix<label> rABPos(this->nSpecie_, this->nSpecie_, -1);
131  // Index of the other species involved in the rABNum
132  RectangularMatrix<label> rABOtherSpec(this->nSpecie_, this->nSpecie_, -1);
133 
134  scalarField omegaV(this->chemistry_.reactions().size());
135  forAll(this->chemistry_.reactions(), i)
136  {
137  const Reaction<ThermoType>& R = this->chemistry_.reactions()[i];
138 
139  // for each reaction compute omegai
140  scalar omegaf, omegar;
141  const scalar omegai = R.omega(p, T, c1, li, omegaf, omegar);
142  omegaV[i] = omegai;
143 
144  // then for each pair of species composing this reaction,
145  // compute the rAB matrix (separate the numerator and
146  // denominator)
147 
148  DynamicList<scalar> wA(R.lhs().size()+R.rhs().size());
149  DynamicList<label> wAID(R.lhs().size()+R.rhs().size());
150  forAll(R.lhs(), s)// compute rAB for all species in the left hand side
151  {
152  label ss = R.lhs()[s].index;
153  scalar sl = -R.lhs()[s].stoichCoeff; // vAi = v''-v' => here -v'
154  List<bool> deltaBi(this->nSpecie_, false);
155  FIFOStack<label> usedIndex;
156  forAll(R.lhs(), j)
157  {
158  label sj = R.lhs()[j].index;
159  usedIndex.push(sj);
160  deltaBi[sj] = true;
161  }
162  forAll(R.rhs(), j)
163  {
164  label sj = R.rhs()[j].index;
165  usedIndex.push(sj);
166  deltaBi[sj] = true;
167  }
168 
169  // Disable for self reference (by definition rAA=0)
170  deltaBi[ss] = false;
171 
172  while(!usedIndex.empty())
173  {
174  label curIndex = usedIndex.pop();
175  if (deltaBi[curIndex])
176  {
177  // disable to avoid counting it more than once
178  deltaBi[curIndex] = false;
179  // test if this rAB is not initialised
180  if (rABPos(ss, curIndex)==-1)
181  {
182  rABPos(ss, curIndex) = NbrABInit[ss];
183  NbrABInit[ss]++;
184  rABNum(ss, rABPos(ss, curIndex)) = sl*omegai;
185  rABOtherSpec(ss, rABPos(ss, curIndex)) = curIndex;
186  }
187  else
188  {
189  rABNum(ss, rABPos(ss, curIndex)) += sl*omegai;
190  }
191  }
192  }
193  bool found(false);
194  forAll(wAID, id)
195  {
196  if (ss==wAID[id])
197  {
198  wA[id] += sl*omegai;
199  found = true;
200  }
201  }
202  if (!found)
203  {
204  wA.append(sl*omegai);
205  wAID.append(ss);
206  }
207  }
208 
209  // Compute rAB for all species in the right hand side
210  forAll(R.rhs(), s)
211  {
212  label ss = R.rhs()[s].index;
213  scalar sl = R.rhs()[s].stoichCoeff; // vAi = v''-v' => here v''
214  List<bool> deltaBi(this->nSpecie_, false);
215  FIFOStack<label> usedIndex;
216  forAll(R.lhs(), j)
217  {
218  label sj = R.lhs()[j].index;
219  usedIndex.push(sj);
220  deltaBi[sj] = true;
221  }
222  forAll(R.rhs(), j)
223  {
224  label sj = R.rhs()[j].index;
225  usedIndex.push(sj);
226  deltaBi[sj] = true;
227  }
228 
229  // Disable for self reference (by definition rAA=0)
230  deltaBi[ss] = false;
231 
232  while(!usedIndex.empty())
233  {
234  label curIndex = usedIndex.pop();
235  if (deltaBi[curIndex])
236  {
237  // disable to avoid counting it more than once
238  deltaBi[curIndex] = false;
239  // test if this rAB is not initialised
240  if (rABPos(ss, curIndex)==-1)
241  {
242  rABPos(ss, curIndex) = NbrABInit[ss];
243  NbrABInit[ss]++;
244  rABNum(ss, rABPos(ss, curIndex)) = sl*omegai;
245  rABOtherSpec(ss, rABPos(ss, curIndex)) = curIndex;
246  }
247  else
248  {
249  rABNum(ss, rABPos(ss, curIndex)) += sl*omegai;
250  }
251  }
252  }
253 
254  bool found(false);
255  forAll(wAID, id)
256  {
257  if (ss==wAID[id])
258  {
259  wA[id] += sl*omegai;
260  found = true;
261  }
262  }
263  if (!found)
264  {
265  wA.append(sl*omegai);
266  wAID.append(ss);
267  }
268  }
269 
270  wAID.shrink();
271  // Now that every species of the reactions has been visited, we can
272  // compute the production and consumption rate. This way, it avoids
273  // getting wrong results when species are present in both lhs and rhs
274  forAll(wAID, id)
275  {
276  if (wA[id] > 0.0)
277  {
278  if (PA[wAID[id]] == 0.0)
279  {
280  PA[wAID[id]] = wA[id];
281  }
282  else
283  {
284  PA[wAID[id]] += wA[id];
285  }
286  }
287  else
288  {
289  if (CA[wAID[id]] == 0.0)
290  {
291  CA[wAID[id]] = -wA[id];
292  }
293  else
294  {
295  CA[wAID[id]] += -wA[id];
296  }
297  }
298  }
299  }
300  // rii = 0.0 by definition
301 
302  // Compute the production rate of each element Pa
303  label nElements = 4; // 4 main elements (C, H, O, N)
304  scalarList Pa(nElements,0.0);
305  scalarList Ca(nElements,0.0);
306 
307  // for (label q=0; q<SIS.size(); q++)
308  for (label i=0; i<this->nSpecie_; i++)
309  {
310  Pa[0] += sC_[i]*max(0.0,(PA[i]-CA[i]));
311  Ca[0] += sC_[i]*max(0.0,-(PA[i]-CA[i]));
312  Pa[1] += sH_[i]*max(0.0,(PA[i]-CA[i]));
313  Ca[1] += sH_[i]*max(0.0,-(PA[i]-CA[i]));
314  Pa[2] += sO_[i]*max(0.0,(PA[i]-CA[i]));
315  Ca[2] += sO_[i]*max(0.0,-(PA[i]-CA[i]));
316  Pa[3] += sN_[i]*max(0.0,(PA[i]-CA[i]));
317  Ca[3] += sN_[i]*max(0.0,-(PA[i]-CA[i]));
318  }
319 
320  // Using the rAB matrix (numerator and denominator separated)
321  // compute the R value according to the search initiating set
322  scalarField Rvalue(this->nSpecie_,0.0);
323  label speciesNumber = 0;
324  List<bool> disabledSpecies(this->nSpecie_,false);
325 
326  // set all species to inactive and activate them according
327  // to rAB and initial set
328  for (label i=0; i<this->nSpecie_; i++)
329  {
330  this->activeSpecies_[i] = false;
331  }
332  // Initialise the FIFOStack for search set
334  const labelList& SIS(this->searchInitSet_);
335  DynamicList<label> QStart(SIS.size());
336  DynamicList<scalar> alphaQ(SIS.size());
337 
338  // Compute the alpha coefficient and initialise the R value of the species
339  // in the SIS
340  for (label i=0; i<SIS.size(); i++)
341  {
342  label q = SIS[i];
343  // compute alpha
344  scalar alphaA(0.0);
345  // C
346  if (Pa[0] > vSmall)
347  {
348  scalar alphaTmp = (sC_[q]*mag(PA[q]-CA[q])/Pa[0]);
349  if (alphaTmp > alphaA)
350  {
351  alphaA = alphaTmp;
352  }
353  }
354  // H
355  if (Pa[1] > vSmall)
356  {
357  scalar alphaTmp = (sH_[q]*mag(PA[q]-CA[q])/Pa[1]);
358  if (alphaTmp > alphaA)
359  {
360  alphaA = alphaTmp;
361  }
362  }
363  // O
364  if (Pa[2] > vSmall)
365  {
366  scalar alphaTmp = (sO_[q]*mag(PA[q]-CA[q])/Pa[2]);
367  if (alphaTmp > alphaA)
368  {
369  alphaA = alphaTmp;
370  }
371  }
372  // N
373  if (Pa[3] > vSmall)
374  {
375  scalar alphaTmp = (sN_[q]*mag(PA[q]-CA[q])/Pa[3]);
376  if (alphaTmp > alphaA)
377  {
378  alphaA = alphaTmp;
379  }
380  }
381  if (alphaA > this->tolerance())
382  {
383  this->activeSpecies_[q] = true;
384  speciesNumber++;
385  Q.push(q);
386  QStart.append(q);
387  alphaQ.append(1.0);
388  Rvalue[q] = 1.0;
389  }
390  else
391  {
392  Rvalue[q] = alphaA;
393  }
394  }
395 
396  // if all species from the SIS has been removed
397  // force the use of the species with maximum Rvalue
398  if (Q.empty())
399  {
400  scalar Rmax=0.0;
401  label specID=-1;
402  forAll(SIS, i)
403  {
404  if (Rvalue[SIS[i]] > Rmax)
405  {
406  Rmax = Rvalue[SIS[i]];
407  specID=SIS[i];
408  }
409  }
410  Q.push(specID);
411  QStart.append(specID);
412  alphaQ.append(1.0);
413  speciesNumber++;
414  Rvalue[specID] = 1.0;
415  this->activeSpecies_[specID] = true;
416  }
417 
418  // Execute the main loop for R-value
419  while (!Q.empty())
420  {
421  label u = Q.pop();
422  scalar Den = max(PA[u],CA[u]);
423  if (Den > vSmall)
424  {
425  for (label v=0; v<NbrABInit[u]; v++)
426  {
427  label otherSpec = rABOtherSpec(u, v);
428  scalar rAB = mag(rABNum(u, v))/Den;
429  if (rAB > 1)
430  {
431  rAB = 1;
432  }
433 
434  scalar Rtemp = Rvalue[u]*rAB;
435  // a link analysed previously is stronger
436  if (Rvalue[otherSpec] < Rtemp)
437  {
438  Rvalue[otherSpec] = Rtemp;
439  // the (composed) link is stronger than the tolerance
440  if (Rtemp >= this->tolerance())
441  {
442  Q.push(otherSpec);
443  if (!this->activeSpecies_[otherSpec])
444  {
445  this->activeSpecies_[otherSpec] = true;
446  speciesNumber++;
447  }
448  }
449  }
450  }
451  }
452  }
453 
454  // Group-based reduction
455  // number of species disabled in the first step
456  label NDisabledSpecies(this->nSpecie_-speciesNumber);
457 
458  // while the number of removed species is greater than NGroupBased, the rAB
459  // are reevaluated according to the group based definition for each loop the
460  // temporary disabled species (in the first reduction) are sorted to disable
461  // definitely the NGroupBased species with lower R then these R value a
462  // reevaluated taking into account these disabled species
463  while(NDisabledSpecies > NGroupBased_)
464  {
465  // keep track of disabled species using sortablelist to extract only
466  // NGroupBased lower R value
467  SortableListDRGEP<scalar> Rdisabled(NDisabledSpecies);
468  labelList Rindex(NDisabledSpecies);
469  label nD = 0;
470  forAll(disabledSpecies, i)
471  {
472  // if just disabled and not in a previous loop
473  if (!this->activeSpecies_[i] && !disabledSpecies[i])
474  {
475  // Note: non-reached species will be removed first (Rvalue=0)
476  Rdisabled[nD] = Rvalue[i];
477  Rindex[nD++] = i;
478  }
479  }
480  // sort the Rvalue to obtain the NGroupBased lower R value
481  Rdisabled.partialSort(NGroupBased_);
482  labelList tmpIndex(Rdisabled.indices());
483 
484  // disable definitely NGroupBased species in this loop
485  for (label i=0; i<NGroupBased_; i++)
486  {
487  disabledSpecies[Rindex[tmpIndex[i]]] = true;
488  }
489  NDisabledSpecies -= NGroupBased_;
490 
491  // reevaluate the rAB according to the group-based definition rAB{S} [1]
492  // only update the numerator
493  forAll(NbrABInit, i)
494  {
495  for (label v=0; v<NbrABInit[i]; v++)
496  {
497  rABNum(i, v) = 0.0;
498  }
499  }
500  forAll(this->chemistry_.reactions(), i)
501  {
502  const Reaction<ThermoType>& R = this->chemistry_.reactions()[i];
503 
504  scalar omegai = omegaV[i];
505 
506  forAll(R.lhs(), s)
507  {
508  label ss = R.lhs()[s].index;
509  scalar sl = -R.lhs()[s].stoichCoeff; // vAi = v''-v' => here -v'
510  List<bool> deltaBi(this->nSpecie_, false);
511  bool alreadyDisabled(false);
512  FIFOStack<label> usedIndex;
513  forAll(R.lhs(), j)
514  {
515  label sj = R.lhs()[j].index;
516  usedIndex.push(sj);
517  deltaBi[sj] = true;
518  if (disabledSpecies[sj])
519  {
520  alreadyDisabled=true;
521  }
522  }
523  forAll(R.rhs(), j)
524  {
525  label sj = R.rhs()[j].index;
526  usedIndex.push(sj);
527  deltaBi[sj] = true;
528  if (disabledSpecies[sj])
529  {
530  alreadyDisabled=true;
531  }
532  }
533 
534  deltaBi[ss] = false;
535 
536  if (alreadyDisabled)
537  {
538  // if one of the species in this reaction is disabled, all
539  // species connected to species ss are modified
540  for (label v=0; v<NbrABInit[ss]; v++)
541  {
542  rABNum(ss, v) += sl*omegai;
543  }
544  }
545  else
546  {
547  while(!usedIndex.empty())
548  {
549  label curIndex = usedIndex.pop();
550  if (deltaBi[curIndex])
551  {
552  // disable to avoid counting it more than once
553  deltaBi[curIndex] = false;
554  rABNum(ss, rABPos(ss, curIndex)) += sl*omegai;
555  }
556  }
557  }
558  }
559 
560  forAll(R.rhs(), s)
561  {
562  label ss = R.rhs()[s].index;
563  scalar sl = R.rhs()[s].stoichCoeff; // vAi = v''-v' => here v''
564  List<bool> deltaBi(this->nSpecie_, false);
565  bool alreadyDisabled(false);
566  FIFOStack<label> usedIndex;
567  forAll(R.lhs(), j)
568  {
569  label sj = R.lhs()[j].index;
570  usedIndex.push(sj);
571  deltaBi[sj] = true;
572  if (disabledSpecies[sj])
573  {
574  alreadyDisabled=true;
575  }
576  }
577  forAll(R.rhs(), j)
578  {
579  label sj = R.rhs()[j].index;
580  usedIndex.push(sj);
581  deltaBi[sj] = true;
582  if (disabledSpecies[sj])
583  {
584  alreadyDisabled=true;
585  }
586  }
587 
588  deltaBi[ss] = false;
589 
590  if (alreadyDisabled)
591  {
592  // if one of the species in this reaction is disabled, all
593  // species connected to species ss are modified
594  for (label v=0; v<NbrABInit[ss]; v++)
595  {
596  rABNum(ss, v) += sl*omegai;
597  }
598  }
599  else
600  {
601  while(!usedIndex.empty())
602  {
603  label curIndex = usedIndex.pop();
604  if (deltaBi[curIndex])
605  {
606  deltaBi[curIndex] = false;
607  rABNum(ss, rABPos(ss, curIndex)) += sl*omegai;
608  }
609  }
610  }
611  }
612  }
613 
614  forAll(QStart, qi)
615  {
616  label u = QStart[qi];
617  Q.push(u);
618  }
619 
620  while (!Q.empty())
621  {
622  label u = Q.pop();
623  scalar Den = max(PA[u],CA[u]);
624  if (Den!=0.0)
625  {
626  for (label v=0; v<NbrABInit[u]; v++)
627  {
628  label otherSpec = rABOtherSpec(u, v);
629  if (!disabledSpecies[otherSpec])
630  {
631  scalar rAB = mag(rABNum(u, v))/Den;
632  if (rAB > 1)
633  {
634  rAB = 1;
635  }
636 
637  scalar Rtemp = Rvalue[u]*rAB;
638  // a link analysed previously is stronger
639  if (Rvalue[otherSpec] < Rtemp)
640  {
641  Rvalue[otherSpec] = Rtemp;
642  if (Rtemp >= this->tolerance())
643  {
644  Q.push(otherSpec);
645  if (!this->activeSpecies_[otherSpec])
646  {
647  this->activeSpecies_[otherSpec] = true;
648  speciesNumber++;
649  NDisabledSpecies--;
650  }
651  }
652  }
653  }
654  }
655  }
656  }
657  }
658 
659  // End of group-based reduction
660 
661  // Put a flag on the reactions containing at least one removed species
662  forAll(this->chemistry_.reactions(), i)
663  {
664  const Reaction<ThermoType>& R = this->chemistry_.reactions()[i];
665  this->chemistry_.reactionsDisabled()[i] = false;
666  forAll(R.lhs(), s)
667  {
668  label ss = R.lhs()[s].index;
669  if (!this->activeSpecies_[ss])
670  {
671  this->chemistry_.reactionsDisabled()[i] = true;
672  break;
673  }
674  }
675  if (!this->chemistry_.reactionsDisabled()[i])
676  {
677  forAll(R.rhs(), s)
678  {
679  label ss = R.rhs()[s].index;
680  if (!this->activeSpecies_[ss])
681  {
682  this->chemistry_.reactionsDisabled()[i] = true;
683  break;
684  }
685  }
686  }
687  }
688 
689  this->NsSimp_ = speciesNumber;
690  scalarField& simplifiedC(this->chemistry_.simplifiedC());
691  simplifiedC.setSize(this->NsSimp_+2);
692  DynamicList<label>& s2c(this->chemistry_.simplifiedToCompleteIndex());
693  s2c.setSize(this->NsSimp_);
694  Field<label>& c2s(this->chemistry_.completeToSimplifiedIndex());
695 
696  label j = 0;
697  for (label i=0; i<this->nSpecie_; i++)
698  {
699  if (this->activeSpecies_[i])
700  {
701  s2c[j] = i;
702  simplifiedC[j] = c[i];
703  c2s[i] = j++;
704  if (!this->chemistry_.active(i))
705  {
706  this->chemistry_.setActive(i);
707  }
708  }
709  else
710  {
711  c2s[i] = -1;
712  }
713  }
714  simplifiedC[this->NsSimp_] = T;
715  simplifiedC[this->NsSimp_+1] = p;
716  this->chemistry_.setNsDAC(this->NsSimp_);
717  // change temporary Ns in chemistryModel
718  // to make the function nEqns working
719  this->chemistry_.setNSpecie(this->NsSimp_);
720 }
721 
722 
723 // ************************************************************************* //
A FIFO stack based on a singly-linked list.
Definition: FIFOStack.H:51
Extends standardChemistryModel by adding the TDAC method.
A HashTable with keys but without contents.
Definition: HashSet.H:59
dictionary dict
#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
virtual ~DRGEP()
Destructor.
Definition: DRGEP.C:95
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
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 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
void omega(const scalar p, const scalar T, const scalarField &c, const label li, scalarField &dcdt) const
Net reaction rate for individual species.
Definition: Reaction.C:290
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))
void setSize(const label)
Alter the addressed list size.
Definition: DynamicListI.H:175
const word & name() const
Return the name of the element.
label nAtoms() const
Return the number of atoms of this element in the specie.
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
basicChemistryModel & chemistry
const volScalarField & T
void setSize(const label)
Reset size of List.
Definition: List.C:281
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)
virtual void reduceMechanism(const scalar p, const scalar T, const scalarField &c, const label li)
Reduce the mechanism.
Definition: DRGEP.C:103
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.
const multiComponentMixture< ThermoType > & mixture() const
Return reference to the mixture.
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
DRGEP(const IOdictionary &dict, TDACChemistryModel< ThermoType > &chemistry)
Construct from components.
Definition: DRGEP.C:33