32template<
class CompType,
class ThermoType>
40 searchInitSet_(this->coeffsDict_.subDict(
"initialSet").size())
44 for (label i=0; i<
chemistry.nSpecie(); i++)
48 searchInitSet_[j++] = i;
51 if (j<searchInitSet_.
size())
54 << searchInitSet_.
size()-j
55 <<
" species in the initial set is not in the mechanism "
64template<
class CompType,
class ThermoType>
71template<
class CompType,
class ThermoType>
81 for(label i=0; i<this->nSpecie_; i++)
86 c1[this->nSpecie_] =
T;
87 c1[this->nSpecie_+1] =
p;
102 scalar pf, cf, pr, cr;
104 forAll(this->chemistry_.reactions(), i)
108 scalar omegai = this->chemistry_.omega
110 R, c1,
T,
p, pf, cf, lRef, pr, cr, rRef
122 label ss =
R.lhs()[
s].index;
123 scalar sl = -
R.lhs()[
s].stoichCoeff;
142 label ss =
R.rhs()[
s].index;
143 scalar sl =
R.rhs()[
s].stoichCoeff;
167 label curID = wAID[id];
170 scalar curwA = ((wA[id]>=0) ? wA[
id] : -wA[
id]);
176 label sj =
R.lhs()[j].index;
182 label sj =
R.rhs()[j].index;
188 deltaBi[curID] =
false;
189 while(!usedIndex.empty())
191 label curIndex = usedIndex.
pop();
193 if (deltaBi[curIndex])
196 deltaBi[curIndex] =
false;
199 if (rABPos(curID, curIndex)==-1)
201 rABPos(curID, curIndex) = NbrABInit[curID];
203 rABNum(curID, rABPos(curID, curIndex)) = curwA;
204 rABOtherSpec(curID, rABPos(curID, curIndex)) = curIndex;
208 rABNum(curID, rABPos(curID, curIndex)) += curwA;
212 if (rABDen[curID] == 0.0)
214 rABDen[curID] = curwA;
218 rABDen[curID] +=curwA;
224 label speciesNumber = 0;
228 for (label i=0; i<this->nSpecie_; i++)
230 this->activeSpecies_[i] =
false;
237 for (label i=0; i<searchInitSet_.size(); i++)
239 label q = searchInitSet_[i];
240 this->activeSpecies_[q] =
true;
249 scalar Den = rABDen[u];
253 for (label v=0; v<NbrABInit[u]; v++)
255 label otherSpec = rABOtherSpec(u, v);
256 scalar rAB = rABNum(u, v)/Den;
267 rAB >= this->tolerance()
268 && !this->activeSpecies_[otherSpec]
272 this->activeSpecies_[otherSpec] =
true;
280 forAll(this->chemistry_.reactions(), i)
283 this->chemistry_.reactionsDisabled()[i] =
false;
287 label ss =
R.lhs()[
s].index;
290 if (!this->activeSpecies_[ss])
293 this->chemistry_.reactionsDisabled()[i] =
true;
299 if (!this->chemistry_.reactionsDisabled()[i])
303 label ss =
R.rhs()[
s].index;
304 if (!this->activeSpecies_[ss])
306 this->chemistry_.reactionsDisabled()[i] =
true;
313 this->NsSimp_ = speciesNumber;
314 this->chemistry_.simplifiedC().setSize(this->NsSimp_+2);
315 this->chemistry_.simplifiedToCompleteIndex().setSize(this->NsSimp_);
318 for (label i=0; i<this->nSpecie_; i++)
320 if (this->activeSpecies_[i])
322 this->chemistry_.simplifiedToCompleteIndex()[j] = i;
323 this->chemistry_.simplifiedC()[j] = c[i];
324 this->chemistry_.completeToSimplifiedIndex()[i] = j++;
325 if (!this->chemistry_.active(i))
327 this->chemistry_.setActive(i);
332 this->chemistry_.completeToSimplifiedIndex()[i] = -1;
336 this->chemistry_.simplifiedC()[this->NsSimp_] =
T;
337 this->chemistry_.simplifiedC()[this->NsSimp_+1] =
p;
338 this->chemistry_.setNsDAC(this->NsSimp_);
342 this->chemistry_.setNSpecie(this->NsSimp_);
#define R(A, B, C, D, E, F, K, M)
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
DynamicList< T, SizeMin > & shrink()
Shrink the allocated space to the number of elements used.
void append(const T &val)
Copy append an element to the end of this list.
A FIFO stack based on a singly-linked list.
void push(const T &elem)
Push an element onto the back of the stack.
T pop()
Pop the bottom element off the stack.
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...
Simple extension of ReactionThermo to handle reaction kinetics in addition to the equilibrium thermod...
A templated (M x N) rectangular matrix of objects of <Type>, containing M*N elements,...
Extends StandardChemistryModel by adding the TDAC method.
void size(const label n)
Older name for setAddressableSize.
An abstract class for methods of chemical mechanism reduction.
const dictionary coeffsDict_
Dictionary that store the algorithm data.
Implementation of the Directed Relation Graph (DRG) method.
virtual void reduceMechanism(const scalarField &c, const scalar T, const scalar p)
Reduce the mechanism.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
bool found(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Search for an entry (const access) with the given keyword.
BasicChemistryModel< psiReactionThermo > & chemistry
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
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))
static constexpr const zero Zero
Global zero (0)
errorManipArg< error, int > exit(error &err, const int errNo=1)
#define forAll(list, i)
Loop across all elements in list.