DAC.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2016-2017 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "DAC.H"
29 
30 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
31 
32 template<class CompType, class ThermoType>
34 (
35  const IOdictionary& dict,
37 )
38 :
40  searchInitSet_(this->coeffsDict_.subDict("initialSet").size()),
41  zprime_(0),
42  nbCLarge_(this->coeffsDict_.template lookupOrDefault<label>("nbCLarge", 3)),
43  sC_(this->nSpecie_, 0),
44  sH_(this->nSpecie_, 0),
45  sO_(this->nSpecie_, 0),
46  sN_(this->nSpecie_, 0),
47  CO2Id_(-1),
48  COId_(-1),
49  HO2Id_(-1),
50  H2OId_(-1),
51  NOId_(-1),
52  automaticSIS_
53  (
54  this->coeffsDict_.template lookupOrDefault<Switch>
55  (
56  "automaticSIS",
57  true
58  )
59  ),
60  phiTol_
61  (
62  this->coeffsDict_.template lookupOrDefault<scalar>
63  (
64  "phiTol", this->tolerance()
65  )
66  ),
67  NOxThreshold_
68  (
69  this->coeffsDict_.template lookupOrDefault<scalar>
70  (
71  "NOxThreshold",
72  1800
73  )
74  ),
75  CO2Name_(this->coeffsDict_.template lookupOrDefault<word>("CO2", "CO2")),
76  COName_(this->coeffsDict_.template lookupOrDefault<word>("CO", "CO")),
77  HO2Name_(this->coeffsDict_.template lookupOrDefault<word>("HO2", "HO2")),
78  H2OName_(this->coeffsDict_.template lookupOrDefault<word>("H2O", "H2O")),
79  NOName_(this->coeffsDict_.template lookupOrDefault<word>("NO", "NO")),
80  forceFuelInclusion_
81  (
82  this->coeffsDict_.template lookupOrDefault<Switch>
83  (
84  "forceFuelInclusion",
85  false
86  )
87  )
88 {
89  label j = 0;
90  dictionary initSet = this->coeffsDict_.subDict("initialSet");
91 
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 
100  if (j < searchInitSet_.size())
101  {
103  << searchInitSet_.size()-j
104  << " species in the intial set is not in the mechanism "
105  << initSet
106  << exit(FatalError);
107  }
108 
109 
110  const List<List<specieElement>>& specieComposition =
111  chemistry.specieComp();
112 
113  for (label i=0; i<this->nSpecie_; i++)
114  {
115  const List<specieElement>& curSpecieComposition =
116  specieComposition[i];
117 
118  // For all elements in the current species
119  forAll(curSpecieComposition, j)
120  {
121  const specieElement& curElement = curSpecieComposition[j];
122 
123  if (curElement.name() == "C")
124  {
125  sC_[i] = curElement.nAtoms();
126  }
127  else if (curElement.name() == "H")
128  {
129  sH_[i] = curElement.nAtoms();
130  }
131  else if (curElement.name() == "O")
132  {
133  sO_[i] = curElement.nAtoms();
134  }
135  else if (curElement.name() == "N")
136  {
137  sN_[i] = curElement.nAtoms();
138  }
139  else
140  {
141  Info<< " element " << curElement.name() << " not considered"
142  << endl;
143  }
144  }
145  if (this->chemistry_.Y()[i].member() == CO2Name_)
146  {
147  CO2Id_ = i;
148  }
149  else if (this->chemistry_.Y()[i].member() == COName_)
150  {
151  COId_ = i;
152  }
153  else if (this->chemistry_.Y()[i].member() == HO2Name_)
154  {
155  HO2Id_ = i;
156  }
157  else if (this->chemistry_.Y()[i].member() == H2OName_)
158  {
159  H2OId_ = i;
160  }
161  else if (this->chemistry_.Y()[i].member() == NOName_)
162  {
163  NOId_ = i;
164  }
165  }
166 
167  if ((CO2Id_==-1 || COId_==-1 || HO2Id_==-1 || H2OId_==-1) && automaticSIS_)
168  {
170  << "The name of the species used in automatic SIS are not found in "
171  << " the mechanism. You should either set the name for CO2, CO, H2O"
172  << " and HO2 properly or set automaticSIS to off "
173  << exit(FatalError);
174  }
175 
176  // To compute zprime, the fuel species should be specified.
177  // According to the given mass fraction, an equivalent O/C ratio is computed
178  if (automaticSIS_)
179  {
180  dictionary fuelDict;
181  if (this->coeffsDict_.found("fuelSpecies"))
182  {
183  fuelDict = this->coeffsDict_.subDict("fuelSpecies");
184  fuelSpecies_ = fuelDict.toc();
185  if (fuelSpecies_.size() == 0)
186  {
188  << "With automatic SIS, the fuel species should be "
189  << "specified in the fuelSpecies subDict"
190  << exit(FatalError);
191  }
192  }
193  else
194  {
196  << "With automatic SIS, the fuel species should be "
197  << "specified in the fuelSpecies subDict"
198  << exit(FatalError);
199  }
200 
201 
202  fuelSpeciesID_.setSize(fuelSpecies_.size());
203  fuelSpeciesProp_.setSize(fuelSpecies_.size());
204  scalar Mmtot(0.0);
205 
206  forAll(fuelSpecies_, i)
207  {
208  fuelSpeciesProp_[i] = fuelDict.get<scalar>(fuelSpecies_[i]);
209  for (label j=0; j<this->nSpecie_; j++)
210  {
211  if (this->chemistry_.Y()[j].member() == fuelSpecies_[i])
212  {
213  fuelSpeciesID_[i] = j;
214  break;
215  }
216  }
217  scalar curMm =
218  this->chemistry_.specieThermo()[fuelSpeciesID_[i]].W();
219  Mmtot += fuelSpeciesProp_[i]/curMm;
220  }
221 
222  Mmtot = 1.0/Mmtot;
223  scalar nbC(0.0);
224  scalar nbO(0.0);
225  forAll(fuelSpecies_, i)
226  {
227  label curID = fuelSpeciesID_[i];
228  scalar curMm = this->chemistry_.specieThermo()[curID].W();
229 
230  nbC += fuelSpeciesProp_[i]*Mmtot/curMm*sC_[curID];
231  nbO += fuelSpeciesProp_[i]*Mmtot/curMm*sO_[curID];
232  }
233  zprime_ = nbO/nbC;
234  }
235 }
236 
237 
238 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
239 
240 template<class CompType, class ThermoType>
242 {}
243 
244 
245 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
246 
247 
248 template<class CompType, class ThermoType>
250 (
251  const scalarField &c,
252  const scalar T,
253  const scalar p
254 )
255 {
256  scalarField& completeC(this->chemistry_.completeC());
257  scalarField c1(this->chemistry_.nEqns(), Zero);
258  for (label i=0; i<this->nSpecie_; i++)
259  {
260  c1[i] = c[i];
261  completeC[i] = c[i];
262  }
263 
264  c1[this->nSpecie_] = T;
265  c1[this->nSpecie_+1] = p;
266 
267  // Compute the rAB matrix
268  RectangularMatrix<scalar> rABNum(this->nSpecie_, this->nSpecie_, Zero);
269  scalarField PA(this->nSpecie_, Zero);
270  scalarField CA(this->nSpecie_, Zero);
271 
272  // Number of initialized rAB for each lines
273  Field<label> NbrABInit(this->nSpecie_, Zero);
274  // Position of the initialized rAB, -1 when not initialized
275  RectangularMatrix<label> rABPos(this->nSpecie_, this->nSpecie_, -1);
276  // Index of the other species involved in the rABNum
277  RectangularMatrix<label> rABOtherSpec(this->nSpecie_, this->nSpecie_, -1);
278 
279  scalar pf, cf, pr, cr;
280  label lRef, rRef;
281  for (const Reaction<ThermoType>& R : this->chemistry_.reactions())
282  {
283  // for each reaction compute omegai
284  scalar omegai = this->chemistry_.omega
285  (
286  R, c1, T, p, pf, cf, lRef, pr, cr, rRef
287  );
288 
289  // Then for each pair of species composing this reaction,
290  // compute the rAB matrix (separate the numerator and
291  // denominator)
292 
293  // While computing the rAB for all the species involved in the reaction
294  // we should consider that one can write a reaction A+B=2C as A+B=C+C
295  // In this case, the following algorithm only take once the effect
296  // of the species. It stores the species encountered in the reaction but
297  // use another list to see if this species has already been used
298 
299  DynamicList<scalar> wA(R.lhs().size() + R.rhs().size());
300  DynamicList<label> wAID(R.lhs().size() + R.rhs().size());
301 
302  forAll(R.lhs(), s) // Compute rAB for all species in the left hand side
303  {
304  label ss = R.lhs()[s].index;
305  scalar sl = -R.lhs()[s].stoichCoeff; // vAi = v''-v' => here -v'
306  List<bool> deltaBi(this->nSpecie_, false);
307  FIFOStack<label> usedIndex;
308  forAll(R.lhs(), j)
309  {
310  label sj = R.lhs()[j].index;
311  usedIndex.push(sj);
312  deltaBi[sj] = true;
313  }
314  forAll(R.rhs(), j)
315  {
316  label sj = R.rhs()[j].index;
317  usedIndex.push(sj);
318  deltaBi[sj] = true;
319  }
320 
321  // Disable for self reference (by definition rAA=0)
322  deltaBi[ss] = false;
323 
324  while (!usedIndex.empty())
325  {
326  label curIndex = usedIndex.pop();
327  if (deltaBi[curIndex])
328  {
329  // Disable to avoid counting it more than once
330  deltaBi[curIndex] = false;
331  // Test if this rAB is not initialized
332  if (rABPos(ss, curIndex) == -1)
333  {
334  // It starts at rABPos(ss, sj)=0
335  rABPos(ss, curIndex) = NbrABInit[ss];
336  NbrABInit[ss]++;
337  // to avoid overflow
338  rABNum(ss, rABPos(ss, curIndex)) = sl*omegai;
339  // store the other specie involved
340  rABOtherSpec(ss, rABPos(ss, curIndex)) = curIndex;
341  }
342  else
343  {
344  rABNum(ss, rABPos(ss, curIndex)) += sl*omegai;
345  }
346  }
347  }
348 
349  bool found(false);
350  forAll(wAID, id)
351  {
352  if (ss == wAID[id])
353  {
354  wA[id] += sl*omegai;
355  found = true;
356  }
357  }
358  if (!found)
359  {
360  wA.append(sl*omegai);
361  wAID.append(ss);
362  }
363  }
364 
365  // Compute rAB for all species in the right hand side
366  forAll(R.rhs(), s)
367  {
368  label ss = R.rhs()[s].index;
369  scalar sl = R.rhs()[s].stoichCoeff; // vAi = v''-v' => here v''
370  List<bool> deltaBi(this->nSpecie_, false);
371  FIFOStack<label> usedIndex;
372  forAll(R.lhs(), j)
373  {
374  label sj = R.lhs()[j].index;
375  usedIndex.push(sj);
376  deltaBi[sj] = true;
377  }
378  forAll(R.rhs(), j)
379  {
380  label sj = R.rhs()[j].index;
381  usedIndex.push(sj);
382  deltaBi[sj] = true;
383  }
384 
385  // Disable for self reference (by definition rAA=0)
386  deltaBi[ss] = false;
387 
388  while (!usedIndex.empty())
389  {
390  label curIndex = usedIndex.pop();
391  if (deltaBi[curIndex])
392  {
393  // Disable to avoid counting it more than once
394  deltaBi[curIndex] = false;
395 
396  // Test if this rAB is not initialized
397  if (rABPos(ss, curIndex) == -1)
398  {
399  // it starts at rABPos(ss, sj)=0
400  rABPos(ss, curIndex) = NbrABInit[ss];
401  NbrABInit[ss]++;
402  rABNum(ss, rABPos(ss, curIndex)) = sl*omegai;
403  rABOtherSpec(ss, rABPos(ss, curIndex)) = curIndex;
404  }
405  else
406  {
407  rABNum(ss, rABPos(ss, curIndex)) += sl*omegai;
408  }
409  }
410  }
411 
412  bool found(false);
413  forAll(wAID, id)
414  {
415  if (ss==wAID[id])
416  {
417  wA[id] += sl*omegai;
418  found = true;
419  }
420  }
421  if (!found)
422  {
423  wA.append(sl*omegai);
424  wAID.append(ss);
425  }
426  }
427  wAID.shrink();
428 
429  // Now that every species of the reactions has been visited, we can
430  // compute the production and consumption rate. This way, it avoids
431  // getting wrong results when species are present in both lhs and rhs
432  forAll(wAID, id)
433  {
434  if (wA[id] > 0.0)
435  {
436  if (PA[wAID[id]] == 0.0)
437  {
438  PA[wAID[id]] = wA[id];
439  }
440  else
441  {
442  PA[wAID[id]] += wA[id];
443  }
444  }
445  else
446  {
447  if (CA[wAID[id]] == 0.0)
448  {
449  CA[wAID[id]] = -wA[id];
450  }
451  else
452  {
453  CA[wAID[id]] += -wA[id];
454  }
455  }
456  }
457  }
458  // rii = 0.0 by definition
459 
460  scalar phiLarge(0.0);
461  scalar phiProgress(0.0);
462  if (automaticSIS_)
463  {
464  // Compute the progress equivalence ratio
465  // and the equivalence ratio for fuel decomposition
466  label nElements = 4; // 4 main elements (C, H, O, N)
467 
468  // Total number of C, H and O (in this order)
469  scalarList Na(nElements, Zero);
470  scalarList Nal(nElements, Zero); // for large hydrocarbons
471 
472  for (label i=0; i<this->nSpecie_; i++)
473  {
474  // Complete combustion products are not considered
475  if
476  (
477  this->chemistry_.Y()[i].member() == "CO2"
478  || this->chemistry_.Y()[i].member() == "H2O"
479  )
480  {
481  continue;
482  }
483 
484  Na[0] += sC_[i]*c[i];
485  Na[1] += sH_[i]*c[i];
486  Na[2] += sO_[i]*c[i];
487  if (sC_[i]>nbCLarge_ || this->chemistry_.Y()[i].member() == "O2")
488  {
489  Nal[0] += sC_[i]*c[i];
490  Nal[1] += sH_[i]*c[i];
491  Nal[2] += sO_[i]*c[i];
492  }
493  }
494 
495  // 2C(-CO2) + H(-H2O)/2 - z'C(-CO2)
496  // Progress equivalence ratio = ----------------------------------
497  // O(-CO2-H2O) - z' C(-CO2)
498  // where minus means that this species is not considered for the number
499  // of atoms and z' is the ratio of the number of O and C in the fuel(s)
500  phiProgress = (2*Na[0]+Na[1]/2-zprime_*Na[0])/(Na[2]-zprime_*Na[0]);
501 
502  // 2Cl + Hl/2
503  // Equivalence ratio for fuel decomposition = ----------
504  // Ol(+O2)
505  phiLarge = (2*Nal[0] + Nal[1]/2)/Nal[2];
506  }
507 
508  // Using the rAB matrix (numerator and denominator separated)
509  // compute the R value according to the search initiating set
510  scalarField Rvalue(this->nSpecie_, Zero);
511  label speciesNumber = 0;
512 
513  // Set all species to inactive and activate them according
514  // to rAB and initial set
515  for (label i=0; i<this->nSpecie_; ++i)
516  {
517  this->activeSpecies_[i] = false;
518  }
519 
520  // Initialize the FIFOStack for search set
522 
523  const labelList& SIS(searchInitSet_);
524 
525  // If automaticSIS is on, the search initiating set is selected according to
526  // phiProgress and phiLarge
527  if (automaticSIS_)
528  {
529  if (phiLarge >= phiTol_ && phiProgress >= phiTol_)
530  {
531  // When phiLarge and phiProgress >= phiTol then
532  // CO, HO2 and fuel are in the SIS
533  Q.push(COId_);
534  ++speciesNumber;
535  this->activeSpecies_[COId_] = true;
536  Rvalue[COId_] = 1.0;
537  Q.push(HO2Id_);
538  ++speciesNumber;
539  this->activeSpecies_[HO2Id_] = true;
540  Rvalue[HO2Id_] = 1.0;
541  for (const label fuelId : fuelSpeciesID_)
542  {
543  Q.push(fuelId);
544  ++speciesNumber;
545  this->activeSpecies_[fuelId] = true;
546  Rvalue[fuelId] = 1.0;
547  }
548 
549  }
550  else if (phiLarge < phiTol_ && phiProgress >= phiTol_)
551  {
552  // When phiLarge < phiTol and phiProgress >= phiTol then
553  // CO, HO2 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 
563  if (forceFuelInclusion_)
564  {
565  for (const label fuelId : fuelSpeciesID_)
566  {
567  Q.push(fuelId);
568  ++speciesNumber;
569  this->activeSpecies_[fuelId] = true;
570  Rvalue[fuelId] = 1.0;
571  }
572  }
573  }
574  else
575  {
576  // When phiLarge and phiProgress< phiTol then
577  // CO2, H2O are in the SIS
578  Q.push(CO2Id_);
579  speciesNumber++;
580  this->activeSpecies_[CO2Id_] = true;
581  Rvalue[CO2Id_] = 1.0;
582 
583  Q.push(H2OId_);
584  speciesNumber++;
585  this->activeSpecies_[H2OId_] = true;
586  Rvalue[H2OId_] = 1.0;
587  if (forceFuelInclusion_)
588  {
589  for (const label fuelId : fuelSpeciesID_)
590  {
591  Q.push(fuelId);
592  ++speciesNumber;
593  this->activeSpecies_[fuelId] = true;
594  Rvalue[fuelId] = 1.0;
595  }
596  }
597  }
598 
599  if (T > NOxThreshold_ && NOId_ != -1)
600  {
601  Q.push(NOId_);
602  ++speciesNumber;
603  this->activeSpecies_[NOId_] = true;
604  Rvalue[NOId_] = 1.0;
605  }
606  }
607  else // No automaticSIS => all species of the SIS are added
608  {
609  for (label i=0; i<SIS.size(); i++)
610  {
611  label q = SIS[i];
612  this->activeSpecies_[q] = true;
613  ++speciesNumber;
614  Q.push(q);
615  Rvalue[q] = 1.0;
616  }
617  }
618 
619  // Execute the main loop for R-value
620  while (!Q.empty())
621  {
622  label u = Q.pop();
623  scalar Den = max(PA[u], CA[u]);
624  if (Den != 0)
625  {
626  for (label v=0; v<NbrABInit[u]; ++v)
627  {
628  label otherSpec = rABOtherSpec(u, v);
629  scalar rAB = mag(rABNum(u, v))/Den;
630  if (rAB > 1)
631  {
632  rAB = 1;
633  }
634 
635  // The direct link is weaker than the user-defined tolerance
636  if (rAB >= this->tolerance())
637  {
638  scalar Rtemp = Rvalue[u]*rAB;
639  // a link analysed previously is stronger
640  // the (composed) link is stronger than the user-defined
641  // tolerance
642  if ((Rvalue[otherSpec]<Rtemp) && (Rtemp>=this->tolerance()))
643  {
644  Q.push(otherSpec);
645  Rvalue[otherSpec] = Rtemp;
646  if (!this->activeSpecies_[otherSpec])
647  {
648  this->activeSpecies_[otherSpec] = true;
649  ++speciesNumber;
650  }
651  }
652  }
653  }
654  }
655  }
656 
657  // Put a flag on the reactions containing at least one removed species
658  forAll(this->chemistry_.reactions(), i)
659  {
660  const Reaction<ThermoType>& R = this->chemistry_.reactions()[i];
661  this->chemistry_.reactionsDisabled()[i] = false;
662 
663  forAll(R.lhs(), s)
664  {
665  label ss = R.lhs()[s].index;
666  if (!this->activeSpecies_[ss])
667  {
668  // Flag the reaction to disable it
669  this->chemistry_.reactionsDisabled()[i] = true;
670  break;
671  }
672  }
673  if (!this->chemistry_.reactionsDisabled()[i])
674  {
675  forAll(R.rhs(), s)
676  {
677  label ss = R.rhs()[s].index;
678  if (!this->activeSpecies_[ss])
679  {
680  // Flag the reaction to disable it
681  this->chemistry_.reactionsDisabled()[i] = true;
682  break;
683  }
684  }
685  }
686  }
687 
688  this->NsSimp_ = speciesNumber;
689  scalarField& simplifiedC(this->chemistry_.simplifiedC());
690  simplifiedC.setSize(this->NsSimp_ + 2);
691  DynamicList<label>& s2c(this->chemistry_.simplifiedToCompleteIndex());
692  s2c.setSize(this->NsSimp_);
693  Field<label>& c2s(this->chemistry_.completeToSimplifiedIndex());
694 
695  label j = 0;
696  for (label i=0; i<this->nSpecie_; i++)
697  {
698  if (this->activeSpecies_[i])
699  {
700  s2c[j] = i;
701  simplifiedC[j] = c[i];
702  c2s[i] = j++;
703  if (!this->chemistry_.active(i))
704  {
705  this->chemistry_.setActive(i);
706  }
707  }
708  else
709  {
710  c2s[i] = -1;
711  }
712  }
713 
714  simplifiedC[this->NsSimp_] = T;
715  simplifiedC[this->NsSimp_+1] = p;
716  this->chemistry_.setNsDAC(this->NsSimp_);
717 
718  // Change temporary Ns in chemistryModel
719  // to make the function nEqns working
720  this->chemistry_.setNSpecie(this->NsSimp_);
721 }
722 
723 
724 // ************************************************************************* //
Foam::IOdictionary
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:54
p
volScalarField & p
Definition: createFieldRefs.H:8
s
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.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))
Definition: gmvOutputSpray.H:25
Foam::Zero
static constexpr const zero Zero
Global zero.
Definition: zero.H:128
Foam::DynamicList< scalar >
Foam::dictionary::found
bool found(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Search for an entry (const access) with the given keyword.
Definition: dictionary.C:359
chemistry
BasicChemistryModel< psiReactionThermo > & chemistry
Definition: createFieldRefs.H:1
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:337
Foam::dictionary::get
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:81
Foam::specieElement::name
const word & name() const
Return the name of the element.
Definition: specieElementI.H:52
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
Foam::chemistryReductionMethods::DAC::reduceMechanism
virtual void reduceMechanism(const scalarField &c, const scalar T, const scalar p)
Reduce the mechanism.
Definition: DAC.C:250
Foam::FIFOStack::push
void push(const T &element)
Push an element onto the back of the stack.
Definition: FIFOStack.H:78
R
#define R(A, B, C, D, E, F, K, M)
Foam::constant::physicoChemical::c1
const dimensionedScalar c1
First radiation constant: default SI units: [W/m2].
Foam::specieElement::nAtoms
label nAtoms() const
Return the number of atoms of this element in the specie.
Definition: specieElementI.H:64
Foam::label
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:62
Foam::Field< scalar >
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
Foam::chemistryReductionMethods::DAC::~DAC
virtual ~DAC()
Destructor.
Definition: DAC.C:241
Foam::dictionary::subDict
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
Definition: dictionary.C:523
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
Foam::RectangularMatrix< scalar >
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
T
const volScalarField & T
Definition: createFieldRefs.H:2
Foam::DynamicList::setSize
void setSize(const label nElem)
Alter addressable list size.
Definition: DynamicListI.H:282
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
found
bool found
Definition: TABSMDCalcMethod2.H:32
Foam::FIFOStack::pop
T pop()
Pop the bottom element off the stack.
Definition: FIFOStack.H:90
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:355
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:102
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::specieElement
Definition: specieElement.H:54
Foam::constant::universal::c
const dimensionedScalar c
Speed of light in a vacuum.
Foam::dictionary::toc
wordList toc() const
Return the table of contents.
Definition: dictionary.C:665
DAC.H
Foam::TDACChemistryModel< CompType, ThermoType >
Foam::FIFOStack
A FIFO stack based on a singly-linked list.
Definition: FIFOStack.H:51
Foam::chemistryReductionMethod
An abstract class for methods of chemical mechanism reduction.
Definition: chemistryReductionMethod.H:56
Foam::Reaction
Simple extension of ReactionThermo to handle reaction kinetics in addition to the equilibrium thermod...
Definition: Reaction.H:56
Foam::chemistryReductionMethods::DAC::DAC
DAC(const IOdictionary &dict, TDACChemistryModel< CompType, ThermoType > &chemistry)
Construct from components.
Definition: DAC.C:34