36template<
class ReactionThermo,
class ThermoType>
37void Foam::combustionModels::
38diffusionMulticomponent<ReactionThermo, ThermoType>::init()
41 this->coeffs().readIfPresent(
"Ci", Ci_);
42 this->coeffs().readIfPresent(
"YoxStream", YoxStream_);
43 this->coeffs().readIfPresent(
"YfStream", YfStream_);
44 this->coeffs().readIfPresent(
"sigma", sigma_);
45 this->coeffs().readIfPresent(
"ftCorr", ftCorr_);
46 this->coeffs().readIfPresent(
"alpha", alpha_);
47 this->coeffs().readIfPresent(
"laminarIgn", laminarIgn_);
54 const label nReactions = reactions_.size();
56 for (label
k=0;
k < nReactions;
k++)
70 RijPtr_[
k].storePrevIter();
75 const label fuelIndex =
species.find(fuelNames_[
k]);
76 const label oxidantIndex =
species.find(oxidantNames_[
k]);
78 const scalar Wu = specieThermo_[fuelIndex].W();
79 const scalar Wox = specieThermo_[oxidantIndex].W();
83 const label specieI = lhs[i].index;
84 specieStoichCoeffs[specieI] = -lhs[i].stoichCoeff;
86 specieThermo_[specieI].hc()*lhs[i].stoichCoeff/Wu;
91 const label specieI =
rhs[i].index;
92 specieStoichCoeffs[specieI] =
rhs[i].stoichCoeff;
94 specieThermo_[specieI].hc()*
rhs[i].stoichCoeff/Wu;
97 Info <<
"Fuel heat of combustion : " << qFuel_[
k] <<
endl;
100 (Wox*
mag(specieStoichCoeffs[oxidantIndex]))
101 / (Wu*
mag(specieStoichCoeffs[fuelIndex]));
103 Info <<
"stoichiometric oxygen-fuel ratio : " << s_[
k] <<
endl;
105 stoicRatio_[
k] = s_[
k]*YfStream_[
k]/YoxStream_[
k];
107 Info <<
"stoichiometric air-fuel ratio : " << stoicRatio_[
k] <<
endl;
109 const scalar fStoich = 1.0/(1.0 + stoicRatio_[
k]);
111 Info <<
"stoichiometric mixture fraction : " << fStoich <<
endl;
118template<
class ReactionThermo,
class ThermoType>
122 const word& modelType,
125 const word& combustionProperties
144 RijPtr_(reactions_.size()),
145 Ci_(reactions_.size(), 1.0),
146 fuelNames_(this->coeffs().
lookup(
"fuels")),
147 oxidantNames_(this->coeffs().
lookup(
"oxidants")),
148 qFuel_(reactions_.size()),
149 stoicRatio_(reactions_.size()),
150 s_(reactions_.size()),
151 YoxStream_(reactions_.size(), 0.23),
152 YfStream_(reactions_.size(), 1.0),
153 sigma_(reactions_.size(), 0.02),
154 oxidantRes_(this->coeffs().
lookup(
"oxidantRes")),
155 ftCorr_(reactions_.size(),
Zero),
165template<
class ReactionThermo,
class ThermoType>
167diffusionMulticomponent<ReactionThermo, ThermoType>::tc()
const
173template<
class ReactionThermo,
class ThermoType>
182 const label nReactions = reactions_.
size();
186 for (label
k=0;
k < nReactions;
k++)
203 const label fuelIndex =
species.find(fuelNames_[
k]);
207 Rijl.
ref() = -this->chemistryPtr_->calculateRR(
k, fuelIndex);
218 const label lIndex = lhs[l].index;
219 this->chemistryPtr_->RR(lIndex) =
Zero;
224 const label rIndex =
rhs[l].index;
225 this->chemistryPtr_->RR(rIndex) =
Zero;
230 for (label
k=0;
k < nReactions;
k++)
232 const label fuelIndex = species.
find(fuelNames_[
k]);
233 const label oxidantIndex = species.
find(oxidantNames_[
k]);
245 s_[
k]*Yfuel - (Yox - YoxStream_[
k])
248 s_[
k]*YfStream_[
k] + YoxStream_[
k]
252 const scalar
sigma = sigma_[
k];
254 const scalar fStoich = 1.0/(1.0 + stoicRatio_[
k]) + ftCorr_[
k];
259 Yox/
max(oxidantRes_[
k], 1
e-3)
265 1.0 +
sqr(OAvailScaled)
293 min(RijkDiff, topHatFilter*RijlPtr[
k]*
pos(Yox)*
pos(Yfuel));
302 if (debug && this->mesh_.time().writeTime())
309 const List<specieCoeffs>&
rhs = reactions_[
k].rhs();
310 const List<specieCoeffs>& lhs = reactions_[
k].lhs();
312 scalar fuelStoic = 1.0;
315 const label lIndex = lhs[l].index;
316 if (lIndex == fuelIndex)
318 fuelStoic = lhs[l].stoichCoeff;
323 const scalar MwFuel = specieThermo_[fuelIndex].W();
328 const label lIndex = lhs[l].index;
330 const scalar stoichCoeff = lhs[l].stoichCoeff;
332 this->chemistryPtr_->RR(lIndex) +=
333 -Rijk*stoichCoeff*specieThermo_[lIndex].W()/fuelStoic/MwFuel;
340 const label rIndex =
rhs[r].index;
342 const scalar stoichCoeff =
rhs[r].stoichCoeff;
344 this->chemistryPtr_->RR(rIndex) +=
345 Rijk*stoichCoeff*specieThermo_[rIndex].W()/fuelStoic/MwFuel;
352template<
class ReactionThermo,
class ThermoType>
360 auto&
Su = tSu.ref();
364 const label specieI =
367 Su += this->chemistryPtr_->RR(specieI);
374template<
class ReactionThermo,
class ThermoType>
391 Qdot = this->chemistryPtr_->Qdot();
398template<
class ReactionThermo,
class ThermoType>
404 this->coeffs().readIfPresent(
"Ci", Ci_);
405 this->coeffs().readIfPresent(
"sigma", sigma_);
406 this->coeffs().readIfPresent(
"oxidantRes", oxidantRes_);
407 this->coeffs().readIfPresent(
"ftCorr", ftCorr_);
408 this->coeffs().readIfPresent(
"alpha", alpha_);
409 this->coeffs().readIfPresent(
"laminarIgn", laminarIgn_);
compressible::turbulenceModel & turb
Chemistry model wrapper for combustion models.
autoPtr< BasicChemistryModel< ReactionThermo > > chemistryPtr_
Pointer to chemistry model.
virtual ReactionThermo & thermo()
Return access to the thermo package.
label size() const noexcept
The number of elements in list.
void relax(const scalar alpha)
Relax field (for steady-state solution).
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
Internal & ref(const bool updateAccessTime=true)
Same as internalFieldRef().
@ NO_REGISTER
Do not request registration (bool: false).
static word member(const word &name)
Return member (name without the extension).
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
const T * set(const label i) const
Return const pointer to element (can be nullptr), or nullptr for out-of-range access (ie,...
void size(const label n)
Older name for setAddressableSize.
const speciesTable & species() const
Return the table of species.
PtrList< volScalarField > & Y()
Return the mass-fraction fields.
const fvMesh & mesh_
Reference to the mesh database.
const compressibleTurbulenceModel & turbulence() const
Return access to turbulence.
const dictionary & coeffs() const
Return const dictionary of the model.
Switch active() const noexcept
Is combustion active?
virtual bool read()
Update properties from given dictionary.
Diffusion based turbulent combustion model for multicomponent species.
virtual void correct()
Correct combustion rate.
virtual tmp< fvScalarMatrix > R(volScalarField &Y) const
Fuel consumption rate matrix.
virtual tmp< volScalarField > Qdot() const
Heat release rate calculated from fuel consumption rate matrix.
virtual bool read()
Update properties from given dictionary.
Abstract base class for turbulence models (RAS, LES and laminar).
static const word & zeroGradientType() noexcept
The type name for zeroGradient patch fields.
label find(const word &val) const
Find index of the value (searches the hash).
virtual basicSpecieMixture & composition()=0
Return the composition of the multi-component mixture.
Lookup type of boundary radiation properties.
virtual bool write(const bool writeOnProc=true) const
Write using setting from DB.
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
A class for managing temporary objects.
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
A class for handling words, derived from Foam::string.
PtrList< volScalarField > & Y
compressible::turbulenceModel & turbulence
Calculate the gradient of the given field.
constexpr scalar pi(M_PI)
Namespace for handling debugging switches.
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
dimensionedScalar pos(const dimensionedScalar &ds)
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
dimensionedScalar exp(const dimensionedScalar &ds)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimEnergy
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
messageStream Info
Information stream (stdout output on master, null elsewhere).
hashedWordList speciesTable
A table of species as a hashedWordList.
Ostream & endl(Ostream &os)
Add newline and flush stream.
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
void rhs(fvMatrix< typename Expr::value_type > &m, const Expr &expression)
static constexpr const zero Zero
Global zero (0).
const dimensionSet dimVolume(pow3(dimLength))
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
List< scalar > scalarList
List of scalar.
void Su(fvMatrix< typename Expr::value_type > &m, const Expr &expression)
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
psiReactionThermo & thermo
dimensionedScalar sigma("sigma", dimMass/sqr(dimTime), transportProperties)
#define forAll(list, i)
Loop across all elements in list.
Hold specie index and its coefficients in the reaction rate expression.