30 template<
class ThermoType>
40 const wordHashSet initSet(this->coeffsDict_.lookup(
"initialSet"));
43 searchInitSet_.append(chemistry.
mixture().species()[iter.key()]);
50 template<
class ThermoType>
57 template<
class ThermoType>
68 for(
label i=0; i<this->nSpecie_; i++)
73 c1[this->nSpecie_] =
T;
74 c1[this->nSpecie_+1] =
p;
89 forAll(this->chemistry_.reactions(), i)
94 scalar omegaf, omegar;
95 const scalar omegai = R.
omega(p, T, c1, li, omegaf, omegar);
106 scalar sl = -R.
lhs()[
s].stoichCoeff;
119 wA.append(sl*omegai);
126 scalar sl = R.
rhs()[
s].stoichCoeff;
139 wA.append(sl*omegai);
150 label curID = wAID[id];
153 scalar curwA = ((wA[id]>=0) ? wA[
id] : -wA[
id]);
171 deltaBi[curID] =
false;
172 while(!usedIndex.empty())
176 if (deltaBi[curIndex])
179 deltaBi[curIndex] =
false;
182 if (rABPos(curID, curIndex)==-1)
184 rABPos(curID, curIndex) = NbrABInit[curID];
186 rABNum(curID, rABPos(curID, curIndex)) = curwA;
187 rABOtherSpec(curID, rABPos(curID, curIndex)) = curIndex;
191 rABNum(curID, rABPos(curID, curIndex)) += curwA;
195 if (rABDen[curID] == 0.0)
197 rABDen[curID] = curwA;
201 rABDen[curID] +=curwA;
207 label speciesNumber = 0;
211 for (
label i=0; i<this->nSpecie_; i++)
213 this->activeSpecies_[i] =
false;
220 for (
label i=0; i<searchInitSet_.size(); i++)
222 label q = searchInitSet_[i];
223 this->activeSpecies_[q] =
true;
232 scalar Den = rABDen[u];
236 for (
label v=0; v<NbrABInit[u]; v++)
238 label otherSpec = rABOtherSpec(u, v);
239 scalar rAB = rABNum(u, v)/Den;
250 rAB >= this->tolerance()
251 && !this->activeSpecies_[otherSpec]
255 this->activeSpecies_[otherSpec] =
true;
263 forAll(this->chemistry_.reactions(), i)
266 this->chemistry_.reactionsDisabled()[i] =
false;
273 if (!this->activeSpecies_[ss])
276 this->chemistry_.reactionsDisabled()[i] =
true;
282 if (!this->chemistry_.reactionsDisabled()[i])
287 if (!this->activeSpecies_[ss])
289 this->chemistry_.reactionsDisabled()[i] =
true;
296 this->NsSimp_ = speciesNumber;
297 this->chemistry_.simplifiedC().setSize(this->NsSimp_+2);
298 this->chemistry_.simplifiedToCompleteIndex().setSize(this->NsSimp_);
301 for (
label i=0; i<this->nSpecie_; i++)
303 if (this->activeSpecies_[i])
305 this->chemistry_.simplifiedToCompleteIndex()[j] = i;
306 this->chemistry_.simplifiedC()[j] = c[i];
307 this->chemistry_.completeToSimplifiedIndex()[i] = j++;
308 if (!this->chemistry_.active(i))
310 this->chemistry_.setActive(i);
315 this->chemistry_.completeToSimplifiedIndex()[i] = -1;
319 this->chemistry_.simplifiedC()[this->NsSimp_] =
T;
320 this->chemistry_.simplifiedC()[this->NsSimp_+1] =
p;
321 this->chemistry_.setNsDAC(this->NsSimp_);
325 this->chemistry_.setNSpecie(this->NsSimp_);
A FIFO stack based on a singly-linked list.
Extends standardChemistryModel by adding the TDAC method.
A HashTable with keys but without contents.
#define forAll(list, i)
Loop across all elements in list.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
const List< specieCoeffs > & lhs() const
Return the components of the left hand side.
virtual void reduceMechanism(const scalar p, const scalar T, const scalarField &c, const label li)
Reduce the mechanism.
Simple extension of ReactionThermo to handle reaction kinetics in addition to the equilibrium thermod...
void omega(const scalar p, const scalar T, const scalarField &c, const label li, scalarField &dcdt) const
Net reaction rate for individual species.
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
const dimensionedScalar c1
First radiation constant: default SI units: [W/m^2].
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))
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
basicChemistryModel & chemistry
const List< specieCoeffs > & rhs() const
Return the components of the right hand side.
DRG(const IOdictionary &dict, TDACChemistryModel< ThermoType > &chemistry)
Construct from components.
#define R(A, B, C, D, E, F, K, M)
void push(const T &a)
Push an element onto the stack.
const multiComponentMixture< ThermoType > & mixture() const
Return reference to the mixture.
T pop()
Pop the bottom element off the stack.
An abstract class for methods of chemical mechanism reduction.