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