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-2024 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "DRGEP.H"
27 #include "SortableListDRGEP.H"
28 
29 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
30 
31 template<class ThermoType>
33 (
34  const dictionary& dict,
36 )
37 :
39  searchInitSet_(),
40  sC_(this->nSpecie(), 0),
41  sH_(this->nSpecie(), 0),
42  sO_(this->nSpecie(), 0),
43  sN_(this->nSpecie(), 0),
44  NGroupBased_(50)
45 {
47 
48  const wordHashSet initSet(this->coeffDict(dict).lookup("initialSet"));
49  forAllConstIter(wordHashSet, initSet, iter)
50  {
51  searchInitSet_.append(chemistry.thermo().species()[iter.key()]);
52  }
53 
54  if (this->coeffDict(dict).found("NGroupBased"))
55  {
56  NGroupBased_ =
57  this->coeffDict(dict).template lookup<label>("NGroupBased");
58  }
59 
60  for (label i=0; i<this->nSpecie(); i++)
61  {
62  const List<specieElement>& curSpecieComposition =
63  chemistry.mixture().specieComposition(i);
64 
65  // for all elements in the current species
66  forAll(curSpecieComposition, j)
67  {
68  const specieElement& curElement =
69  curSpecieComposition[j];
70  if (curElement.name() == "C")
71  {
72  sC_[i] = curElement.nAtoms();
73  }
74  else if (curElement.name() == "H")
75  {
76  sH_[i] = curElement.nAtoms();
77  }
78  else if (curElement.name() == "O")
79  {
80  sO_[i] = curElement.nAtoms();
81  }
82  else if (curElement.name() == "N")
83  {
84  sN_[i] = curElement.nAtoms();
85  }
86  else
87  {
88  Info<< "element not considered"<< endl;
89  }
90  }
91  }
92 }
93 
94 
95 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
96 
97 template<class ThermoType>
99 {}
100 
101 
102 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
103 
104 template<class ThermoType>
106 (
107  const scalar p,
108  const scalar T,
109  const scalarField& c,
110  List<label>& ctos,
111  DynamicList<label>& stoc,
112  const label li
113 )
114 {
115  scalarField c1(this->chemistry_.nEqns(), 0.0);
116 
117  for (label i=0; i<this->nSpecie(); i++)
118  {
119  c1[i] = c[i];
120  }
121 
122  c1[this->nSpecie()] = T;
123  c1[this->nSpecie()+1] = p;
124 
125  // Compute the rAB matrix
126  RectangularMatrix<scalar> rABNum(this->nSpecie(),this->nSpecie(),0.0);
127  scalarField PA(this->nSpecie(),0.0);
128  scalarField CA(this->nSpecie(),0.0);
129 
130  // Number of initialised rAB for each lines
131  Field<label> NbrABInit(this->nSpecie(),0);
132  // Position of the initialised rAB, -1 when not initialised
133  RectangularMatrix<label> rABPos(this->nSpecie(), this->nSpecie(), -1);
134  // Index of the other species involved in the rABNum
135  RectangularMatrix<label> rABOtherSpec(this->nSpecie(), this->nSpecie(), -1);
136 
137  scalarField omegaV(this->chemistry_.reactions().size());
138  forAll(this->chemistry_.reactions(), i)
139  {
140  const Reaction<ThermoType>& R = this->chemistry_.reactions()[i];
141 
142  // for each reaction compute omegai
143  scalar omegaf, omegar;
144  const scalar omegai = R.omega(p, T, c1, li, omegaf, omegar);
145  omegaV[i] = omegai;
146 
147  // then for each pair of species composing this reaction,
148  // compute the rAB matrix (separate the numerator and
149  // denominator)
150 
151  DynamicList<scalar> wA(R.lhs().size()+R.rhs().size());
152  DynamicList<label> wAID(R.lhs().size()+R.rhs().size());
153  forAll(R.lhs(), s)// compute rAB for all species in the left hand side
154  {
155  label ss = R.lhs()[s].index;
156  scalar sl = -R.lhs()[s].stoichCoeff; // vAi = v''-v' => here -v'
157  List<bool> deltaBi(this->nSpecie(), false);
158  FIFOStack<label> usedIndex;
159  forAll(R.lhs(), j)
160  {
161  label sj = R.lhs()[j].index;
162  usedIndex.push(sj);
163  deltaBi[sj] = true;
164  }
165  forAll(R.rhs(), j)
166  {
167  label sj = R.rhs()[j].index;
168  usedIndex.push(sj);
169  deltaBi[sj] = true;
170  }
171 
172  // Disable for self reference (by definition rAA=0)
173  deltaBi[ss] = false;
174 
175  while(!usedIndex.empty())
176  {
177  label curIndex = usedIndex.pop();
178  if (deltaBi[curIndex])
179  {
180  // disable to avoid counting it more than once
181  deltaBi[curIndex] = false;
182  // test if this rAB is not initialised
183  if (rABPos(ss, curIndex)==-1)
184  {
185  rABPos(ss, curIndex) = NbrABInit[ss];
186  NbrABInit[ss]++;
187  rABNum(ss, rABPos(ss, curIndex)) = sl*omegai;
188  rABOtherSpec(ss, rABPos(ss, curIndex)) = curIndex;
189  }
190  else
191  {
192  rABNum(ss, rABPos(ss, curIndex)) += sl*omegai;
193  }
194  }
195  }
196  bool found(false);
197  forAll(wAID, id)
198  {
199  if (ss==wAID[id])
200  {
201  wA[id] += sl*omegai;
202  found = true;
203  }
204  }
205  if (!found)
206  {
207  wA.append(sl*omegai);
208  wAID.append(ss);
209  }
210  }
211 
212  // Compute rAB for all species in the right hand side
213  forAll(R.rhs(), s)
214  {
215  label ss = R.rhs()[s].index;
216  scalar sl = R.rhs()[s].stoichCoeff; // vAi = v''-v' => here v''
217  List<bool> deltaBi(this->nSpecie(), false);
218  FIFOStack<label> usedIndex;
219  forAll(R.lhs(), j)
220  {
221  label sj = R.lhs()[j].index;
222  usedIndex.push(sj);
223  deltaBi[sj] = true;
224  }
225  forAll(R.rhs(), j)
226  {
227  label sj = R.rhs()[j].index;
228  usedIndex.push(sj);
229  deltaBi[sj] = true;
230  }
231 
232  // Disable for self reference (by definition rAA=0)
233  deltaBi[ss] = false;
234 
235  while(!usedIndex.empty())
236  {
237  label curIndex = usedIndex.pop();
238  if (deltaBi[curIndex])
239  {
240  // disable to avoid counting it more than once
241  deltaBi[curIndex] = false;
242  // test if this rAB is not initialised
243  if (rABPos(ss, curIndex)==-1)
244  {
245  rABPos(ss, curIndex) = NbrABInit[ss];
246  NbrABInit[ss]++;
247  rABNum(ss, rABPos(ss, curIndex)) = sl*omegai;
248  rABOtherSpec(ss, rABPos(ss, curIndex)) = curIndex;
249  }
250  else
251  {
252  rABNum(ss, rABPos(ss, curIndex)) += sl*omegai;
253  }
254  }
255  }
256 
257  bool found(false);
258  forAll(wAID, id)
259  {
260  if (ss==wAID[id])
261  {
262  wA[id] += sl*omegai;
263  found = true;
264  }
265  }
266  if (!found)
267  {
268  wA.append(sl*omegai);
269  wAID.append(ss);
270  }
271  }
272 
273  wAID.shrink();
274  // Now that every species of the reactions has been visited, we can
275  // compute the production and consumption rate. This way, it avoids
276  // getting wrong results when species are present in both lhs and rhs
277  forAll(wAID, id)
278  {
279  if (wA[id] > 0.0)
280  {
281  if (PA[wAID[id]] == 0.0)
282  {
283  PA[wAID[id]] = wA[id];
284  }
285  else
286  {
287  PA[wAID[id]] += wA[id];
288  }
289  }
290  else
291  {
292  if (CA[wAID[id]] == 0.0)
293  {
294  CA[wAID[id]] = -wA[id];
295  }
296  else
297  {
298  CA[wAID[id]] += -wA[id];
299  }
300  }
301  }
302  }
303  // rii = 0.0 by definition
304 
305  // Compute the production rate of each element Pa
306  label nElements = 4; // 4 main elements (C, H, O, N)
307  scalarList Pa(nElements,0.0);
308  scalarList Ca(nElements,0.0);
309 
310  // for (label q=0; q<SIS.size(); q++)
311  for (label i=0; i<this->nSpecie(); i++)
312  {
313  Pa[0] += sC_[i]*max(0.0,(PA[i]-CA[i]));
314  Ca[0] += sC_[i]*max(0.0,-(PA[i]-CA[i]));
315  Pa[1] += sH_[i]*max(0.0,(PA[i]-CA[i]));
316  Ca[1] += sH_[i]*max(0.0,-(PA[i]-CA[i]));
317  Pa[2] += sO_[i]*max(0.0,(PA[i]-CA[i]));
318  Ca[2] += sO_[i]*max(0.0,-(PA[i]-CA[i]));
319  Pa[3] += sN_[i]*max(0.0,(PA[i]-CA[i]));
320  Ca[3] += sN_[i]*max(0.0,-(PA[i]-CA[i]));
321  }
322 
323  // Using the rAB matrix (numerator and denominator separated)
324  // compute the R value according to the search initiating set
325  scalarField Rvalue(this->nSpecie(),0.0);
326  List<bool> disabledSpecies(this->nSpecie(),false);
327 
328  // set all species to inactive and activate them according
329  // to rAB and initial set
330  label NActiveSpecies = 0;
331  for (label i=0; i<this->nSpecie(); i++)
332  {
333  this->activeSpecies_[i] = false;
334  }
335 
336  // Initialise the FIFOStack for search set
338  const labelList& SIS(this->searchInitSet_);
339  DynamicList<label> QStart(SIS.size());
340  DynamicList<scalar> alphaQ(SIS.size());
341 
342  // Compute the alpha coefficient and initialise the R value of the species
343  // in the SIS
344  for (label i=0; i<SIS.size(); i++)
345  {
346  label q = SIS[i];
347  // compute alpha
348  scalar alphaA(0.0);
349  // C
350  if (Pa[0] > vSmall)
351  {
352  scalar alphaTmp = (sC_[q]*mag(PA[q]-CA[q])/Pa[0]);
353  if (alphaTmp > alphaA)
354  {
355  alphaA = alphaTmp;
356  }
357  }
358  // H
359  if (Pa[1] > vSmall)
360  {
361  scalar alphaTmp = (sH_[q]*mag(PA[q]-CA[q])/Pa[1]);
362  if (alphaTmp > alphaA)
363  {
364  alphaA = alphaTmp;
365  }
366  }
367  // O
368  if (Pa[2] > vSmall)
369  {
370  scalar alphaTmp = (sO_[q]*mag(PA[q]-CA[q])/Pa[2]);
371  if (alphaTmp > alphaA)
372  {
373  alphaA = alphaTmp;
374  }
375  }
376  // N
377  if (Pa[3] > vSmall)
378  {
379  scalar alphaTmp = (sN_[q]*mag(PA[q]-CA[q])/Pa[3]);
380  if (alphaTmp > alphaA)
381  {
382  alphaA = alphaTmp;
383  }
384  }
385  if (alphaA > this->tolerance())
386  {
387  this->activeSpecies_[q] = true;
388  NActiveSpecies ++;
389  Q.push(q);
390  QStart.append(q);
391  alphaQ.append(1.0);
392  Rvalue[q] = 1.0;
393  }
394  else
395  {
396  Rvalue[q] = alphaA;
397  }
398  }
399 
400  // if all species from the SIS has been removed
401  // force the use of the species with maximum Rvalue
402  if (Q.empty())
403  {
404  scalar Rmax=0.0;
405  label specID=-1;
406  forAll(SIS, i)
407  {
408  if (Rvalue[SIS[i]] > Rmax)
409  {
410  Rmax = Rvalue[SIS[i]];
411  specID=SIS[i];
412  }
413  }
414  Q.push(specID);
415  QStart.append(specID);
416  alphaQ.append(1.0);
417  NActiveSpecies ++;
418  Rvalue[specID] = 1.0;
419  this->activeSpecies_[specID] = true;
420  }
421 
422  // Execute the main loop for R-value
423  while (!Q.empty())
424  {
425  label u = Q.pop();
426  scalar Den = max(PA[u],CA[u]);
427  if (Den > vSmall)
428  {
429  for (label v=0; v<NbrABInit[u]; v++)
430  {
431  label otherSpec = rABOtherSpec(u, v);
432  scalar rAB = mag(rABNum(u, v))/Den;
433  if (rAB > 1)
434  {
435  rAB = 1;
436  }
437 
438  scalar Rtemp = Rvalue[u]*rAB;
439  // a link analysed previously is stronger
440  if (Rvalue[otherSpec] < Rtemp)
441  {
442  Rvalue[otherSpec] = Rtemp;
443  // the (composed) link is stronger than the tolerance
444  if (Rtemp >= this->tolerance())
445  {
446  Q.push(otherSpec);
447  if (!this->activeSpecies_[otherSpec])
448  {
449  this->activeSpecies_[otherSpec] = true;
450  NActiveSpecies ++;
451  }
452  }
453  }
454  }
455  }
456  }
457 
458  // Group-based reduction
459  // number of species disabled in the first step
460  label NDisabledSpecies(this->nSpecie() - NActiveSpecies);
461 
462  // while the number of removed species is greater than NGroupBased, the rAB
463  // are reevaluated according to the group based definition for each loop the
464  // temporary disabled species (in the first reduction) are sorted to disable
465  // definitely the NGroupBased species with lower R then these R value a
466  // reevaluated taking into account these disabled species
467  while(NDisabledSpecies > NGroupBased_)
468  {
469  // keep track of disabled species using sortablelist to extract only
470  // NGroupBased lower R value
471  SortableListDRGEP<scalar> Rdisabled(NDisabledSpecies);
472  labelList Rindex(NDisabledSpecies);
473  label nD = 0;
474  forAll(disabledSpecies, i)
475  {
476  // if just disabled and not in a previous loop
477  if (!this->activeSpecies_[i] && !disabledSpecies[i])
478  {
479  // Note: non-reached species will be removed first (Rvalue=0)
480  Rdisabled[nD] = Rvalue[i];
481  Rindex[nD++] = i;
482  }
483  }
484  // sort the Rvalue to obtain the NGroupBased lower R value
485  Rdisabled.partialSort(NGroupBased_);
486  labelList tmpIndex(Rdisabled.indices());
487 
488  // disable definitely NGroupBased species in this loop
489  for (label i=0; i<NGroupBased_; i++)
490  {
491  disabledSpecies[Rindex[tmpIndex[i]]] = true;
492  }
493  NDisabledSpecies -= NGroupBased_;
494 
495  // reevaluate the rAB according to the group-based definition rAB{S} [1]
496  // only update the numerator
497  forAll(NbrABInit, i)
498  {
499  for (label v=0; v<NbrABInit[i]; v++)
500  {
501  rABNum(i, v) = 0.0;
502  }
503  }
504  forAll(this->chemistry_.reactions(), i)
505  {
506  const Reaction<ThermoType>& R = this->chemistry_.reactions()[i];
507 
508  scalar omegai = omegaV[i];
509 
510  forAll(R.lhs(), s)
511  {
512  label ss = R.lhs()[s].index;
513  scalar sl = -R.lhs()[s].stoichCoeff; // vAi = v''-v' => here -v'
514  List<bool> deltaBi(this->nSpecie(), false);
515  bool alreadyDisabled(false);
516  FIFOStack<label> usedIndex;
517  forAll(R.lhs(), j)
518  {
519  label sj = R.lhs()[j].index;
520  usedIndex.push(sj);
521  deltaBi[sj] = true;
522  if (disabledSpecies[sj])
523  {
524  alreadyDisabled=true;
525  }
526  }
527  forAll(R.rhs(), j)
528  {
529  label sj = R.rhs()[j].index;
530  usedIndex.push(sj);
531  deltaBi[sj] = true;
532  if (disabledSpecies[sj])
533  {
534  alreadyDisabled=true;
535  }
536  }
537 
538  deltaBi[ss] = false;
539 
540  if (alreadyDisabled)
541  {
542  // if one of the species in this reaction is disabled, all
543  // species connected to species ss are modified
544  for (label v=0; v<NbrABInit[ss]; v++)
545  {
546  rABNum(ss, v) += sl*omegai;
547  }
548  }
549  else
550  {
551  while(!usedIndex.empty())
552  {
553  label curIndex = usedIndex.pop();
554  if (deltaBi[curIndex])
555  {
556  // disable to avoid counting it more than once
557  deltaBi[curIndex] = false;
558  rABNum(ss, rABPos(ss, curIndex)) += sl*omegai;
559  }
560  }
561  }
562  }
563 
564  forAll(R.rhs(), s)
565  {
566  label ss = R.rhs()[s].index;
567  scalar sl = R.rhs()[s].stoichCoeff; // vAi = v''-v' => here v''
568  List<bool> deltaBi(this->nSpecie(), false);
569  bool alreadyDisabled(false);
570  FIFOStack<label> usedIndex;
571  forAll(R.lhs(), j)
572  {
573  label sj = R.lhs()[j].index;
574  usedIndex.push(sj);
575  deltaBi[sj] = true;
576  if (disabledSpecies[sj])
577  {
578  alreadyDisabled=true;
579  }
580  }
581  forAll(R.rhs(), j)
582  {
583  label sj = R.rhs()[j].index;
584  usedIndex.push(sj);
585  deltaBi[sj] = true;
586  if (disabledSpecies[sj])
587  {
588  alreadyDisabled=true;
589  }
590  }
591 
592  deltaBi[ss] = false;
593 
594  if (alreadyDisabled)
595  {
596  // if one of the species in this reaction is disabled, all
597  // species connected to species ss are modified
598  for (label v=0; v<NbrABInit[ss]; v++)
599  {
600  rABNum(ss, v) += sl*omegai;
601  }
602  }
603  else
604  {
605  while(!usedIndex.empty())
606  {
607  label curIndex = usedIndex.pop();
608  if (deltaBi[curIndex])
609  {
610  deltaBi[curIndex] = false;
611  rABNum(ss, rABPos(ss, curIndex)) += sl*omegai;
612  }
613  }
614  }
615  }
616  }
617 
618  forAll(QStart, qi)
619  {
620  label u = QStart[qi];
621  Q.push(u);
622  }
623 
624  while (!Q.empty())
625  {
626  label u = Q.pop();
627  scalar Den = max(PA[u],CA[u]);
628  if (Den!=0.0)
629  {
630  for (label v=0; v<NbrABInit[u]; v++)
631  {
632  label otherSpec = rABOtherSpec(u, v);
633  if (!disabledSpecies[otherSpec])
634  {
635  scalar rAB = mag(rABNum(u, v))/Den;
636  if (rAB > 1)
637  {
638  rAB = 1;
639  }
640 
641  scalar Rtemp = Rvalue[u]*rAB;
642  // a link analysed previously is stronger
643  if (Rvalue[otherSpec] < Rtemp)
644  {
645  Rvalue[otherSpec] = Rtemp;
646  if (Rtemp >= this->tolerance())
647  {
648  Q.push(otherSpec);
649  if (!this->activeSpecies_[otherSpec])
650  {
651  this->activeSpecies_[otherSpec] = true;
652  NDisabledSpecies--;
653  }
654  }
655  }
656  }
657  }
658  }
659  }
660  }
661 
663 }
664 
665 
666 // ************************************************************************* //
bool found
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:476
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Definition: DynamicListI.H:296
DynamicList< T, SizeInc, SizeMult, SizeDiv > & shrink()
Shrink the allocated space to the number of elements used.
Definition: DynamicListI.H:252
A FIFO stack based on a singly-linked list.
Definition: FIFOStack.H:54
void push(const T &a)
Push an element onto the stack.
Definition: FIFOStack.H:84
T pop()
Pop the bottom element off the stack.
Definition: FIFOStack.H:90
A HashTable with keys but without contents.
Definition: HashSet.H:62
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: List.H:91
void append(const T &)
Append an element at the end of the list.
Definition: ListI.H:178
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
Simple extension of ThermoType to handle reaction kinetics in addition to the equilibrium thermodynam...
Definition: Reaction.H:72
A list that is sorted upon construction or when explicitly requested with the sort() method.
const labelList & indices() const
Return the list of sorted indices. Updated every sort.
void partialSort(int M)
Partial sort the list (if changed after construction time)
Extends base chemistry model by adding a thermo package, and ODE functions. Introduces chemistry equa...
An abstract class for methods of chemical mechanism reduction.
const dictionary & coeffDict(const dictionary &dict) const
Protected Member Functions.
void initReduceMechanism()
Initialise reduction of the mechanism.
label nSpecie()
Return the number of species.
void endReduceMechanism(List< label > &ctos, DynamicList< label > &stoc)
End reduction of the mechanism.
virtual ~DRGEP()
Destructor.
Definition: DRGEP.C:98
DRGEP(const dictionary &dict, chemistryModel< ThermoType > &chemistry)
Construct from components.
Definition: DRGEP.C:33
virtual void reduceMechanism(const scalar p, const scalar T, const scalarField &c, List< label > &ctos, DynamicList< label > &stoc, const label li)
Reduce the mechanism.
Definition: DRGEP.C:106
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
label nAtoms() const
Return the number of atoms of this element in the specie.
const word & name() const
Return the name of the element.
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(lagrangian::Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(lagrangian::Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(lagrangian::Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.name(), lagrangian::cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
const dimensionedScalar c1
First radiation constant: default SI units: [W/m^2].
const dimensionedScalar c
Speed of light in a vacuum.
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
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:258
void T(LagrangianPatchField< Type > &f, const LagrangianPatchField< Type > &f1)
messageStream Info
void mag(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
static scalar R(const scalar a, const scalar x)
Definition: invIncGamma.C:102
const label nSpecie
dictionary dict
volScalarField & p
basicChemistryModel & chemistry