DAC.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 "DAC.H"
27 
28 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
29 
30 template<class CompType, class ThermoType>
32 (
33  const IOdictionary& dict,
35 )
36 :
38  searchInitSet_(this->coeffsDict_.subDict("initialSet").size()),
39  zprime_(0),
40  nbCLarge_(3),
41  sC_(this->nSpecie_,0),
42  sH_(this->nSpecie_,0),
43  sO_(this->nSpecie_,0),
44  sN_(this->nSpecie_,0),
45  CO2Id_(-1),
46  COId_(-1),
47  HO2Id_(-1),
48  H2OId_(-1),
49  NOId_(-1),
50  automaticSIS_(true),
51  phiTol_(this->tolerance()),
52  NOxThreshold_(1800),
53  CO2Name_
54  (
55  dict.subDict("reduction").lookupOrDefault<word>
56  (
57  "CO2Name","CO2"
58  )
59  ),
60  COName_
61  (
62  dict.subDict("reduction").lookupOrDefault<word>
63  (
64  "COName","CO"
65  )
66  ),
67  HO2Name_
68  (
69  dict.subDict("reduction").lookupOrDefault<word>
70  (
71  "HO2Name","HO2"
72  )
73  ),
74  H2OName_
75  (
76  dict.subDict("reduction").lookupOrDefault<word>
77  (
78  "H2OName","H2O"
79  )
80  ),
81  NOName_
82  (
83  dict.subDict("reduction").lookupOrDefault<word>
84  (
85  "NOName","NO"
86  )
87  ),
88  forceFuelInclusion_(false)
89 {
90  label j=0;
91  dictionary initSet = this->coeffsDict_.subDict("initialSet");
92  for (label i=0; i<chemistry.nSpecie(); i++)
93  {
94  if (initSet.found(chemistry.Y()[i].member()))
95  {
96  searchInitSet_[j++] = i;
97  }
98  }
99  if (j<searchInitSet_.size())
100  {
102  << searchInitSet_.size()-j
103  << " species in the initial set is not in the mechanism "
104  << initSet
105  << exit(FatalError);
106  }
107 
108  if (this->coeffsDict_.found("automaticSIS"))
109  {
110  automaticSIS_.readIfPresent("automaticSIS", this->coeffsDict_);
111  }
112 
113  if (this->coeffsDict_.found("forceFuelInclusion"))
114  {
115  forceFuelInclusion_.readIfPresent
116  (
117  "forceFuelInclusion",this->coeffsDict_
118  );
119  }
120 
121  if (this->coeffsDict_.found("phiTol"))
122  {
123  phiTol_ = readScalar(this->coeffsDict_.lookup("phiTol"));
124  }
125 
126  if (this->coeffsDict_.found("NOxThreshold"))
127  {
128  NOxThreshold_ = readScalar(this->coeffsDict_.lookup("NOxThreshold"));
129  }
130  const List<List<specieElement>>& specieComposition =
131  chemistry.specieComp();
132 
133  for (label i=0; i<this->nSpecie_; i++)
134  {
135  const List<specieElement>& curSpecieComposition =
136  specieComposition[i];
137  // For all elements in the current species
138  forAll(curSpecieComposition, j)
139  {
140  const specieElement& curElement =
141  curSpecieComposition[j];
142  if (curElement.name() == "C")
143  {
144  sC_[i] = curElement.nAtoms();
145  }
146  else if (curElement.name() == "H")
147  {
148  sH_[i] = curElement.nAtoms();
149  }
150  else if (curElement.name() == "O")
151  {
152  sO_[i] = curElement.nAtoms();
153  }
154  else if (curElement.name() == "N")
155  {
156  sN_[i] = curElement.nAtoms();
157  }
158  else
159  {
160  Info<< "element not considered"<<endl;
161  }
162  }
163  if (this->chemistry_.Y()[i].member() == CO2Name_)
164  {
165  CO2Id_ = i;
166  }
167  else if (this->chemistry_.Y()[i].member() == COName_)
168  {
169  COId_ = i;
170  }
171  else if (this->chemistry_.Y()[i].member() == HO2Name_)
172  {
173  HO2Id_ = i;
174  }
175  else if (this->chemistry_.Y()[i].member() == H2OName_)
176  {
177  H2OId_ = i;
178  }
179  else if (this->chemistry_.Y()[i].member() == NOName_)
180  {
181  NOId_ = i;
182  }
183  }
184 
185  if ((CO2Id_==-1 || COId_==-1 || HO2Id_==-1 || H2OId_==-1) && automaticSIS_)
186  {
188  << "The name of the species used in automatic SIS are not found in "
189  << " the mechanism. You should either set the name for CO2, CO, H2O"
190  << " and HO2 properly or set automaticSIS to off "
191  << exit(FatalError);
192  }
193 
194  // To compute zprime, the fuel species should be specified.
195  // According to the given mass fraction, an equivalent O/C ratio is computed
196  if (automaticSIS_)
197  {
198  dictionary fuelDict;
199  if (this->coeffsDict_.found("fuelSpecies"))
200  {
201  fuelDict = this->coeffsDict_.subDict("fuelSpecies");
202  fuelSpecies_ = fuelDict.toc();
203  if (fuelSpecies_.size() == 0)
204  {
206  << "With automatic SIS, the fuel species should be "
207  << "specified in the fuelSpecies subDict"
208  << exit(FatalError);
209  }
210  }
211  else
212  {
214  << "With automatic SIS, the fuel species should be "
215  << "specified in the fuelSpecies subDict"
216  << exit(FatalError);
217  }
218 
219  if (this->coeffsDict_.found("nbCLarge"))
220  {
221  nbCLarge_ = readLabel(fuelDict.lookup("nbCLarge"));
222  }
223 
224  fuelSpeciesID_.setSize(fuelSpecies_.size());
225  fuelSpeciesProp_.setSize(fuelSpecies_.size());
226  scalar Mmtot(0.0);
227 
228  forAll(fuelSpecies_, i)
229  {
230  fuelSpeciesProp_[i] = readScalar(fuelDict.lookup(fuelSpecies_[i]));
231  for (label j=0; j<this->nSpecie_; j++)
232  {
233  if (this->chemistry_.Y()[j].member() == fuelSpecies_[i])
234  {
235  fuelSpeciesID_[i] = j;
236  break;
237  }
238  }
239  scalar curMm =
240  this->chemistry_.specieThermo()[fuelSpeciesID_[i]].W();
241  Mmtot += fuelSpeciesProp_[i]/curMm;
242  }
243 
244  Mmtot = 1.0/Mmtot;
245  scalar nbC(0.0);
246  scalar nbO(0.0);
247  forAll(fuelSpecies_, i)
248  {
249  label curID = fuelSpeciesID_[i];
250  scalar curMm = this->chemistry_.specieThermo()[curID].W();
251 
252  nbC += fuelSpeciesProp_[i]*Mmtot/curMm*sC_[curID];
253  nbO += fuelSpeciesProp_[i]*Mmtot/curMm*sO_[curID];
254  }
255  zprime_ = nbO/nbC;
256  }
257 }
258 
259 
260 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
261 
262 template<class CompType, class ThermoType>
264 {}
265 
266 
267 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
268 
269 
270 template<class CompType, class ThermoType>
272 (
273  const scalarField &c,
274  const scalar T,
275  const scalar p
276 )
277 {
278  scalarField& completeC(this->chemistry_.completeC());
279  scalarField c1(this->chemistry_.nEqns(), 0.0);
280  for(label i=0; i<this->nSpecie_; i++)
281  {
282  c1[i] = c[i];
283  completeC[i] = c[i];
284  }
285 
286  c1[this->nSpecie_] = T;
287  c1[this->nSpecie_+1] = p;
288 
289  // Compute the rAB matrix
290  RectangularMatrix<scalar> rABNum(this->nSpecie_,this->nSpecie_,0.0);
291  scalarField PA(this->nSpecie_,0.0);
292  scalarField CA(this->nSpecie_,0.0);
293 
294  // Number of initialized rAB for each lines
295  Field<label> NbrABInit(this->nSpecie_,0);
296  // Position of the initialized rAB, -1 when not initialized
297  RectangularMatrix<label> rABPos(this->nSpecie_, this->nSpecie_, -1);
298  // Index of the other species involved in the rABNum
299  RectangularMatrix<label> rABOtherSpec(this->nSpecie_, this->nSpecie_, -1);
300 
301  scalar pf, cf, pr, cr;
302  label lRef, rRef;
303  forAll(this->chemistry_.reactions(), i)
304  {
305  const Reaction<ThermoType>& R = this->chemistry_.reactions()[i];
306  // for each reaction compute omegai
307  scalar omegai = this->chemistry_.omega
308  (
309  R, c1, T, p, pf, cf, lRef, pr, cr, rRef
310  );
311 
312  // Then for each pair of species composing this reaction,
313  // compute the rAB matrix (separate the numerator and
314  // denominator)
315 
316  // While computing the rAB for all the species involved in the reaction
317  // we should consider that one can write a reaction A+B=2C as A+B=C+C
318  // In this case, the following algorithm only take once the effect
319  // of the species. It stores the species encountered in the reaction but
320  // use another list to see if this species has already been used
321 
322  DynamicList<scalar> wA(R.lhs().size()+R.rhs().size());
323  DynamicList<label> wAID(R.lhs().size()+R.rhs().size());
324 
325  forAll(R.lhs(), s) // Compute rAB for all species in the left hand side
326  {
327  label ss = R.lhs()[s].index;
328  scalar sl = -R.lhs()[s].stoichCoeff; // vAi = v''-v' => here -v'
329  List<bool> deltaBi(this->nSpecie_, false);
330  FIFOStack<label> usedIndex;
331  forAll(R.lhs(), j)
332  {
333  label sj = R.lhs()[j].index;
334  usedIndex.push(sj);
335  deltaBi[sj] = true;
336  }
337  forAll(R.rhs(), j)
338  {
339  label sj = R.rhs()[j].index;
340  usedIndex.push(sj);
341  deltaBi[sj] = true;
342  }
343 
344  // Disable for self reference (by definition rAA=0)
345  deltaBi[ss] = false;
346 
347  while(!usedIndex.empty())
348  {
349  label curIndex = usedIndex.pop();
350  if (deltaBi[curIndex])
351  {
352  // Disable to avoid counting it more than once
353  deltaBi[curIndex] = false;
354  // Test if this rAB is not initialized
355  if (rABPos(ss, curIndex)==-1)
356  {
357  // It starts at rABPos(ss, sj)=0
358  rABPos(ss, curIndex) = NbrABInit[ss];
359  NbrABInit[ss]++;
360  // to avoid overflow
361  rABNum(ss, rABPos(ss, curIndex)) = sl*omegai;
362  // store the other specie involved
363  rABOtherSpec(ss, rABPos(ss, curIndex)) = curIndex;
364  }
365  else
366  {
367  rABNum(ss, rABPos(ss, curIndex)) += sl*omegai;
368  }
369  }
370  }
371 
372  bool found(false);
373  forAll(wAID, id)
374  {
375  if (ss==wAID[id])
376  {
377  wA[id] += sl*omegai;
378  found = true;
379  }
380  }
381  if (!found)
382  {
383  wA.append(sl*omegai);
384  wAID.append(ss);
385  }
386  }
387 
388  forAll(R.rhs(), s) // Compute rAB for all species in the right hand side
389  {
390  label ss = R.rhs()[s].index;
391  scalar sl = R.rhs()[s].stoichCoeff; // vAi = v''-v' => here v''
392  List<bool> deltaBi(this->nSpecie_, false);
393  FIFOStack<label> usedIndex;
394  forAll(R.lhs(), j)
395  {
396  label sj = R.lhs()[j].index;
397  usedIndex.push(sj);
398  deltaBi[sj] = true;
399  }
400  forAll(R.rhs(), j)
401  {
402  label sj = R.rhs()[j].index;
403  usedIndex.push(sj);
404  deltaBi[sj] = true;
405  }
406 
407  // Disable for self reference (by definition rAA=0)
408  deltaBi[ss] = false;
409 
410  while(!usedIndex.empty())
411  {
412  label curIndex = usedIndex.pop();
413  if (deltaBi[curIndex])
414  {
415  // Disable to avoid counting it more than once
416  deltaBi[curIndex] = false;
417 
418  // Test if this rAB is not initialized
419  if (rABPos(ss, curIndex) == -1)
420  {
421  // it starts at rABPos(ss, sj)=0
422  rABPos(ss, curIndex) = NbrABInit[ss];
423  NbrABInit[ss]++;
424  rABNum(ss, rABPos(ss, curIndex)) = sl*omegai;
425  rABOtherSpec(ss, rABPos(ss, curIndex)) = curIndex;
426  }
427  else
428  {
429  rABNum(ss, rABPos(ss, curIndex)) += sl*omegai;
430  }
431  }
432  }
433 
434  bool found(false);
435  forAll(wAID, id)
436  {
437  if (ss==wAID[id])
438  {
439  wA[id] += sl*omegai;
440  found = true;
441  }
442  }
443  if (!found)
444  {
445  wA.append(sl*omegai);
446  wAID.append(ss);
447  }
448  }
449  wAID.shrink();
450 
451  // Now that every species of the reactions has been visited, we can
452  // compute the production and consumption rate. This way, it avoids
453  // getting wrong results when species are present in both lhs and rhs
454  forAll(wAID, id)
455  {
456  if (wA[id] > 0.0)
457  {
458  if (PA[wAID[id]] == 0.0)
459  {
460  PA[wAID[id]] = wA[id];
461  }
462  else
463  {
464  PA[wAID[id]] += wA[id];
465  }
466  }
467  else
468  {
469  if (CA[wAID[id]] == 0.0)
470  {
471  CA[wAID[id]] = -wA[id];
472  }
473  else
474  {
475  CA[wAID[id]] += -wA[id];
476  }
477  }
478  }
479  }
480  // rii = 0.0 by definition
481 
482  scalar phiLarge(0.0);
483  scalar phiProgress(0.0);
484  if (automaticSIS_)
485  {
486  // Compute the progress equivalence ratio
487  // and the equivalence ratio for fuel decomposition
488  label nElements = 4; // 4 main elements (C, H, O, N)
489 
490  // Total number of C, H and O (in this order)
491  scalarList Na(nElements,0.0);
492  scalarList Nal(nElements,0.0); // for large hydrocarbons
493 
494  for (label i=0; i<this->nSpecie_; i++)
495  {
496  // Complete combustion products are not considered
497  if
498  (
499  this->chemistry_.Y()[i].member() == "CO2"
500  || this->chemistry_.Y()[i].member() == "H2O"
501  )
502  {
503  continue;
504  }
505  Na[0] += sC_[i]*c[i];
506  Na[1] += sH_[i]*c[i];
507  Na[2] += sO_[i]*c[i];
508  if (sC_[i]>nbCLarge_ || this->chemistry_.Y()[i].member() == "O2")
509  {
510  Nal[0] += sC_[i]*c[i];
511  Nal[1] += sH_[i]*c[i];
512  Nal[2] += sO_[i]*c[i];
513  }
514  }
515 
516  // 2C(-CO2) + H(-H2O)/2 - z'C(-CO2)
517  // Progress equivalence ratio = ----------------------------------
518  // O(-CO2-H2O) - z' C(-CO2)
519  // where minus means that this species is not considered for the number
520  // of atoms and z' is the ratio of the number of O and C in the fuel(s)
521  phiProgress = (2*Na[0]+Na[1]/2-zprime_*Na[0])/(Na[2]-zprime_*Na[0]);
522 
523  // 2Cl + Hl/2
524  // Equivalence ratio for fuel decomposition = ----------
525  // Ol(+O2)
526  phiLarge = (2*Nal[0]+Nal[1]/2)/Nal[2];
527  }
528 
529  // Using the rAB matrix (numerator and denominator separated)
530  // compute the R value according to the search initiating set
531  scalarField Rvalue(this->nSpecie_,0.0);
532  label speciesNumber = 0;
533 
534  // Set all species to inactive and activate them according
535  // to rAB and initial set
536  for (label i=0; i<this->nSpecie_; i++)
537  {
538  this->activeSpecies_[i] = false;
539  }
540 
541  // Initialize the FIFOStack for search set
543 
544  const labelList& SIS(searchInitSet_);
545 
546  // If automaticSIS is on, the search initiating set is selected according to
547  // phiProgress and phiLarge
548  if (automaticSIS_)
549  {
550  if (phiLarge >= phiTol_ && phiProgress >= phiTol_)
551  {
552  // When phiLarge and phiProgress >= phiTol then
553  // CO, HO2 and fuel are in the SIS
554  Q.push(COId_);
555  speciesNumber++;
556  this->activeSpecies_[COId_] = true;
557  Rvalue[COId_] = 1.0;
558  Q.push(HO2Id_);
559  speciesNumber++;
560  this->activeSpecies_[HO2Id_] = true;
561  Rvalue[HO2Id_] = 1.0;
562  forAll(fuelSpeciesID_,i)
563  {
564  Q.push(fuelSpeciesID_[i]);
565  speciesNumber++;
566  this->activeSpecies_[fuelSpeciesID_[i]] = true;
567  Rvalue[fuelSpeciesID_[i]] = 1.0;
568  }
569 
570  }
571  else if (phiLarge < phiTol_ && phiProgress >= phiTol_)
572  {
573  // When phiLarge < phiTol and phiProgress >= phiTol then
574  // CO, HO2 are in the SIS
575  Q.push(COId_);
576  speciesNumber++;
577  this->activeSpecies_[COId_] = true;
578  Rvalue[COId_] = 1.0;
579  Q.push(HO2Id_);
580  speciesNumber++;
581  this->activeSpecies_[HO2Id_] = true;
582  Rvalue[HO2Id_] = 1.0;
583 
584  if (forceFuelInclusion_)
585  {
586  forAll(fuelSpeciesID_,i)
587  {
588  Q.push(fuelSpeciesID_[i]);
589  speciesNumber++;
590  this->activeSpecies_[fuelSpeciesID_[i]] = true;
591  Rvalue[fuelSpeciesID_[i]] = 1.0;
592  }
593  }
594  }
595  else
596  {
597  // When phiLarge and phiProgress< phiTol then
598  // CO2, H2O are in the SIS
599  Q.push(CO2Id_);
600  speciesNumber++;
601  this->activeSpecies_[CO2Id_] = true;
602  Rvalue[CO2Id_] = 1.0;
603 
604  Q.push(H2OId_);
605  speciesNumber++;
606  this->activeSpecies_[H2OId_] = true;
607  Rvalue[H2OId_] = 1.0;
608  if (forceFuelInclusion_)
609  {
610  forAll(fuelSpeciesID_,i)
611  {
612  Q.push(fuelSpeciesID_[i]);
613  speciesNumber++;
614  this->activeSpecies_[fuelSpeciesID_[i]] = true;
615  Rvalue[fuelSpeciesID_[i]] = 1.0;
616  }
617  }
618  }
619 
620  if (T>NOxThreshold_ && NOId_!=-1)
621  {
622  Q.push(NOId_);
623  speciesNumber++;
624  this->activeSpecies_[NOId_] = true;
625  Rvalue[NOId_] = 1.0;
626  }
627  }
628  else // No automaticSIS => all species of the SIS are added
629  {
630  for (label i=0; i<SIS.size(); i++)
631  {
632  label q = SIS[i];
633  this->activeSpecies_[q] = true;
634  speciesNumber++;
635  Q.push(q);
636  Rvalue[q] = 1.0;
637  }
638  }
639 
640  // Execute the main loop for R-value
641  while (!Q.empty())
642  {
643  label u = Q.pop();
644  scalar Den = max(PA[u],CA[u]);
645  if (Den != 0)
646  {
647  for (label v=0; v<NbrABInit[u]; v++)
648  {
649  label otherSpec = rABOtherSpec(u, v);
650  scalar rAB = mag(rABNum(u, v))/Den;
651  if (rAB > 1)
652  {
653  rAB = 1;
654  }
655  // The direct link is weaker than the user-defined tolerance
656  if (rAB >= this->tolerance())
657  {
658  scalar Rtemp = Rvalue[u]*rAB;
659  // a link analysed previously is stronger
660  // the (composed) link is stronger than the user-defined
661  // tolerance
662  if ((Rvalue[otherSpec]<Rtemp) && (Rtemp>=this->tolerance()))
663  {
664  Q.push(otherSpec);
665  Rvalue[otherSpec] = Rtemp;
666  if (!this->activeSpecies_[otherSpec])
667  {
668  this->activeSpecies_[otherSpec] = true;
669  speciesNumber++;
670  }
671  }
672  }
673  }
674  }
675  }
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 
683  forAll(R.lhs(), s)
684  {
685  label ss = R.lhs()[s].index;
686  if (!this->activeSpecies_[ss])
687  {
688  // Flag the reaction to disable it
689  this->chemistry_.reactionsDisabled()[i] = true;
690  break;
691  }
692  }
693  if (!this->chemistry_.reactionsDisabled()[i])
694  {
695  forAll(R.rhs(), s)
696  {
697  label ss = R.rhs()[s].index;
698  if (!this->activeSpecies_[ss])
699  {
700  // Flag the reaction to disable it
701  this->chemistry_.reactionsDisabled()[i] = true;
702  break;
703  }
704  }
705  }
706  }
707 
708  this->NsSimp_ = speciesNumber;
709  scalarField& simplifiedC(this->chemistry_.simplifiedC());
710  simplifiedC.setSize(this->NsSimp_+2);
711  DynamicList<label>& s2c(this->chemistry_.simplifiedToCompleteIndex());
712  s2c.setSize(this->NsSimp_);
713  Field<label>& c2s(this->chemistry_.completeToSimplifiedIndex());
714 
715  label j = 0;
716  for (label i=0; i<this->nSpecie_; i++)
717  {
718  if (this->activeSpecies_[i])
719  {
720  s2c[j] = i;
721  simplifiedC[j] = c[i];
722  c2s[i] = j++;
723  if (!this->chemistry_.active(i))
724  {
725  this->chemistry_.setActive(i);
726  }
727  }
728  else
729  {
730  c2s[i] = -1;
731  }
732  }
733 
734  simplifiedC[this->NsSimp_] = T;
735  simplifiedC[this->NsSimp_+1] = p;
736  this->chemistry_.setNsDAC(this->NsSimp_);
737 
738  // Change temporary Ns in chemistryModel
739  // to make the function nEqns working
740  this->chemistry_.setNSpecie(this->NsSimp_);
741 }
742 
743 
744 // ************************************************************************* //
A FIFO stack based on a singly-linked list.
Definition: FIFOStack.H:51
#define readScalar
Definition: doubleScalar.C:38
virtual label nSpecie() const
The number of species.
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
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
BasicChemistryModel< rhoReactionThermo > & chemistry
Simple extension of ReactionThermo to handle reaction kinetics in addition to the equilibrium thermod...
Definition: Reaction.H:55
wordList toc() const
Return the table of contents.
Definition: dictionary.C:776
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
A class for handling words, derived from string.
Definition: word.H:59
DAC(const IOdictionary &dict, TDACChemistryModel< CompType, ThermoType > &chemistry)
Construct from components.
Definition: DAC.C:32
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
List< List< specieElement > > & specieComp()
const volScalarField & T
virtual ~DAC()
Destructor.
Definition: DAC.C:263
virtual void reduceMechanism(const scalarField &c, const scalar T, const scalar p)
Reduce the mechanism.
Definition: DAC.C:272
void setSize(const label)
Reset size of List.
Definition: List.C:281
#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
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
bool found
const List< specieCoeffs > & rhs() const
Return the components of the right hand side.
Definition: ReactionI.H:66
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:576