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