Go to the documentation of this file.
35 template<
class ReactionThermo>
38 const word& modelType,
41 const word& combustionProperties
51 integrateReactionRate_
53 this->coeffs().getOrDefault(
"integrateReactionRate", true)
56 if (integrateReactionRate_)
58 Info<<
" using integrated reaction rate" <<
endl;
62 Info<<
" using instantaneous reaction rate" <<
endl;
69 template<
class ReactionThermo>
76 template<
class ReactionThermo>
80 return this->chemistryPtr_->tc();
84 template<
class ReactionThermo>
89 if (integrateReactionRate_)
97 if (this->coeffs().readIfPresent(
"maxIntegrationTime", maxTime))
99 this->chemistryPtr_->solve
101 min(1.0/rDeltaT, maxTime)()
106 this->chemistryPtr_->solve((1.0/rDeltaT)());
111 this->chemistryPtr_->solve(this->
mesh().time().deltaTValue());
116 this->chemistryPtr_->calculate();
122 template<
class ReactionThermo>
132 const label specieI =
133 this->
thermo().composition().species()[Y.member()];
135 Su += this->chemistryPtr_->RR(specieI);
142 template<
class ReactionThermo>
152 this->
thermo().phasePropertyName(typeName +
":Qdot"),
166 tQdot.
ref() = this->chemistryPtr_->Qdot();
173 template<
class ReactionThermo>
178 integrateReactionRate_ =
179 this->coeffs().getOrDefault(
"integrateReactionRate",
true);
tmp< volScalarField > tc() const
Defines the attributes of an object for which implicit objectRegistry management is supported,...
static const volScalarField & localRDeltaT(const fvMesh &mesh)
A class for handling words, derived from Foam::string.
virtual tmp< fvScalarMatrix > R(volScalarField &Y) const
A class for managing temporary objects.
static constexpr const zero Zero
const dimensionSet dimEnergy
psiReactionThermo & thermo
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
Ostream & endl(Ostream &os)
label min(const labelHashSet &set, label minValue=labelMax)
Laminar combustion model.
virtual tmp< volScalarField > Qdot() const
static bool enabled(const fvMesh &mesh)
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
fvMatrix< scalar > fvScalarMatrix
Generic templated field type.
Chemistry model wrapper for combustion models.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
PtrList< volScalarField > & Y
Calculate the matrix for implicit and explicit sources.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
const dimensionSet dimVolume(pow3(dimLength))
Generic GeometricField class.
compressible::turbulenceModel & turb
Abstract base class for turbulence models (RAS, LES and laminar).