31 template<
class ThermoType>
56 NGroupBased_ = this->
coeffsDict_.template lookup<label>(
"NGroupBased");
65 forAll(curSpecieComposition, j)
68 curSpecieComposition[j];
69 if (curElement.
name() ==
"C")
71 sC_[i] = curElement.
nAtoms();
73 else if (curElement.
name() ==
"H")
75 sH_[i] = curElement.
nAtoms();
77 else if (curElement.
name() ==
"O")
79 sO_[i] = curElement.
nAtoms();
81 else if (curElement.
name() ==
"N")
83 sN_[i] = curElement.
nAtoms();
87 Info<<
"element not considered"<<
endl;
96 template<
class ThermoType>
103 template<
class ThermoType>
136 scalarField omegaV(this->chemistry_.reactions().size());
137 forAll(this->chemistry_.reactions(), i)
142 scalar omegaf, omegar;
143 const scalar omegai =
R.omega(
p,
T,
c1, li, omegaf, omegar);
155 scalar sl = -
R.lhs()[
s].stoichCoeff;
160 label sj =
R.lhs()[j].index;
166 label sj =
R.rhs()[j].index;
174 while(!usedIndex.empty())
177 if (deltaBi[curIndex])
180 deltaBi[curIndex] =
false;
182 if (rABPos(ss, curIndex)==-1)
184 rABPos(ss, curIndex) = NbrABInit[ss];
186 rABNum(ss, rABPos(ss, curIndex)) = sl*omegai;
187 rABOtherSpec(ss, rABPos(ss, curIndex)) = curIndex;
191 rABNum(ss, rABPos(ss, curIndex)) += sl*omegai;
215 scalar sl =
R.rhs()[
s].stoichCoeff;
220 label sj =
R.lhs()[j].index;
226 label sj =
R.rhs()[j].index;
234 while(!usedIndex.empty())
237 if (deltaBi[curIndex])
240 deltaBi[curIndex] =
false;
242 if (rABPos(ss, curIndex)==-1)
244 rABPos(ss, curIndex) = NbrABInit[ss];
246 rABNum(ss, rABPos(ss, curIndex)) = sl*omegai;
247 rABOtherSpec(ss, rABPos(ss, curIndex)) = curIndex;
251 rABNum(ss, rABPos(ss, curIndex)) += sl*omegai;
280 if (PA[wAID[
id]] == 0.0)
282 PA[wAID[id]] = wA[id];
286 PA[wAID[id]] += wA[id];
291 if (CA[wAID[
id]] == 0.0)
293 CA[wAID[id]] = -wA[id];
297 CA[wAID[id]] += -wA[id];
312 Pa[0] += sC_[i]*
max(0.0,(PA[i]-CA[i]));
313 Ca[0] += sC_[i]*
max(0.0,-(PA[i]-CA[i]));
314 Pa[1] += sH_[i]*
max(0.0,(PA[i]-CA[i]));
315 Ca[1] += sH_[i]*
max(0.0,-(PA[i]-CA[i]));
316 Pa[2] += sO_[i]*
max(0.0,(PA[i]-CA[i]));
317 Ca[2] += sO_[i]*
max(0.0,-(PA[i]-CA[i]));
318 Pa[3] += sN_[i]*
max(0.0,(PA[i]-CA[i]));
319 Ca[3] += sN_[i]*
max(0.0,-(PA[i]-CA[i]));
329 label NActiveSpecies = 0;
332 this->activeSpecies_[i] =
false;
337 const labelList& SIS(this->searchInitSet_);
351 scalar alphaTmp = (sC_[q]*
mag(PA[q]-CA[q])/Pa[0]);
352 if (alphaTmp > alphaA)
360 scalar alphaTmp = (sH_[q]*
mag(PA[q]-CA[q])/Pa[1]);
361 if (alphaTmp > alphaA)
369 scalar alphaTmp = (sO_[q]*
mag(PA[q]-CA[q])/Pa[2]);
370 if (alphaTmp > alphaA)
378 scalar alphaTmp = (sN_[q]*
mag(PA[q]-CA[q])/Pa[3]);
379 if (alphaTmp > alphaA)
384 if (alphaA > this->tolerance())
386 this->activeSpecies_[q] =
true;
407 if (Rvalue[SIS[i]] > Rmax)
409 Rmax = Rvalue[SIS[i]];
417 Rvalue[specID] = 1.0;
418 this->activeSpecies_[specID] =
true;
425 scalar Den =
max(PA[u],CA[u]);
428 for (
label v=0; v<NbrABInit[u]; v++)
430 label otherSpec = rABOtherSpec(u, v);
431 scalar rAB =
mag(rABNum(u, v))/Den;
437 scalar Rtemp = Rvalue[u]*rAB;
439 if (Rvalue[otherSpec] < Rtemp)
441 Rvalue[otherSpec] = Rtemp;
443 if (Rtemp >= this->tolerance())
446 if (!this->activeSpecies_[otherSpec])
448 this->activeSpecies_[otherSpec] =
true;
459 label NDisabledSpecies(this->
nSpecie() - NActiveSpecies);
466 while(NDisabledSpecies > NGroupBased_)
473 forAll(disabledSpecies, i)
476 if (!this->activeSpecies_[i] && !disabledSpecies[i])
479 Rdisabled[nD] = Rvalue[i];
488 for (
label i=0; i<NGroupBased_; i++)
490 disabledSpecies[Rindex[tmpIndex[i]]] =
true;
492 NDisabledSpecies -= NGroupBased_;
498 for (
label v=0; v<NbrABInit[i]; v++)
503 forAll(this->chemistry_.reactions(), i)
507 scalar omegai = omegaV[i];
512 scalar sl = -
R.lhs()[
s].stoichCoeff;
514 bool alreadyDisabled(
false);
518 label sj =
R.lhs()[j].index;
521 if (disabledSpecies[sj])
523 alreadyDisabled=
true;
528 label sj =
R.rhs()[j].index;
531 if (disabledSpecies[sj])
533 alreadyDisabled=
true;
543 for (
label v=0; v<NbrABInit[ss]; v++)
545 rABNum(ss, v) += sl*omegai;
550 while(!usedIndex.empty())
553 if (deltaBi[curIndex])
556 deltaBi[curIndex] =
false;
557 rABNum(ss, rABPos(ss, curIndex)) += sl*omegai;
566 scalar sl =
R.rhs()[
s].stoichCoeff;
568 bool alreadyDisabled(
false);
572 label sj =
R.lhs()[j].index;
575 if (disabledSpecies[sj])
577 alreadyDisabled=
true;
582 label sj =
R.rhs()[j].index;
585 if (disabledSpecies[sj])
587 alreadyDisabled=
true;
597 for (
label v=0; v<NbrABInit[ss]; v++)
599 rABNum(ss, v) += sl*omegai;
604 while(!usedIndex.empty())
607 if (deltaBi[curIndex])
609 deltaBi[curIndex] =
false;
610 rABNum(ss, rABPos(ss, curIndex)) += sl*omegai;
619 label u = QStart[qi];
626 scalar Den =
max(PA[u],CA[u]);
629 for (
label v=0; v<NbrABInit[u]; v++)
631 label otherSpec = rABOtherSpec(u, v);
632 if (!disabledSpecies[otherSpec])
634 scalar rAB =
mag(rABNum(u, v))/Den;
640 scalar Rtemp = Rvalue[u]*rAB;
642 if (Rvalue[otherSpec] < Rtemp)
644 Rvalue[otherSpec] = Rtemp;
645 if (Rtemp >= this->tolerance())
648 if (!this->activeSpecies_[otherSpec])
650 this->activeSpecies_[otherSpec] =
true;
#define forAll(list, i)
Loop across all elements in list.
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
DynamicList< T, SizeInc, SizeMult, SizeDiv > & shrink()
Shrink the allocated space to the number of elements used.
A FIFO stack based on a singly-linked list.
void push(const T &a)
Push an element onto the stack.
T pop()
Pop the bottom element off the stack.
A HashTable with keys but without contents.
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
void append(const T &)
Append an element at the end of the list.
void size(const label)
Override size to be inconsistent with allocated storage.
Simple extension of ThermoType to handle reaction kinetics in addition to the equilibrium thermodynam...
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.
void initReduceMechanism()
Protected Member Functions.
const dictionary coeffsDict_
Dictionary that store the algorithm data.
label nSpecie()
Return the number of species.
void endReduceMechanism(List< label > &ctos, DynamicList< label > &stoc)
End reduction of the mechanism.
virtual ~DRGEP()
Destructor.
virtual void reduceMechanism(const scalar p, const scalar T, const scalarField &c, List< label > &ctos, DynamicList< label > &stoc, const label li)
Reduce the mechanism.
DRGEP(const IOdictionary &dict, chemistryModel< ThermoType > &chemistry)
Construct from components.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
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(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.name(), 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.
Ostream & endl(Ostream &os)
Add newline and flush stream.
dimensioned< scalar > mag(const dimensioned< Type > &)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
static scalar R(const scalar a, const scalar x)
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
basicChemistryModel & chemistry