32 template<
class ReactionThermo>
35 const word& modelType,
36 const ReactionThermo& thermo,
38 const word& combustionProperties
48 integrateReactionRate_
50 this->coeffs().lookupOrDefault(
"integrateReactionRate",
true)
53 if (integrateReactionRate_)
55 Info<<
" using integrated reaction rate" <<
endl;
59 Info<<
" using instantaneous reaction rate" <<
endl;
66 template<
class ReactionThermo>
73 template<
class ReactionThermo>
77 return this->chemistryPtr_->tc();
81 template<
class ReactionThermo>
84 if (integrateReactionRate_)
86 if (fv::localEulerDdt::enabled(this->
mesh()))
89 fv::localEulerDdt::localRDeltaT(this->
mesh());
91 if (this->coeffs().
found(
"maxIntegrationTime"))
93 const scalar maxIntegrationTime
95 this->coeffs().
template lookup<scalar>(
"maxIntegrationTime")
98 this->chemistryPtr_->solve
100 min(1.0/rDeltaT, maxIntegrationTime)()
105 this->chemistryPtr_->solve((1.0/rDeltaT)());
110 this->chemistryPtr_->solve(this->
mesh().time().deltaTValue());
115 this->chemistryPtr_->calculate();
120 template<
class ReactionThermo>
128 Su += this->chemistryPtr_->RR(specieI);
134 template<
class ReactionThermo>
138 return this->chemistryPtr_->Qdot();
142 template<
class ReactionThermo>
147 integrateReactionRate_ =
148 this->coeffs().lookupOrDefault(
"integrateReactionRate",
true);
virtual void correct()
Correct combustion rate.
fvMatrix< scalar > fvScalarMatrix
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
T & ref() const
Return non-const reference or generate a fatal error.
Ostream & endl(Ostream &os)
Add newline and flush stream.
tmp< volScalarField > tc() const
Return the chemical time scale.
rhoReactionThermo & thermo
virtual tmp< volScalarField > Qdot() const
Heat release rate [kg/m/s^3].
static word member(const word &name)
Return member (name without the extension)
laminar(const word &modelType, const ReactionThermo &thermo, const compressibleMomentumTransportModel &turb, const word &combustionProperties)
Construct from components.
A class for handling words, derived from string.
A special matrix type and solver, designed for finite volume solutions of scalar equations. Face addressing is used to make all matrix assembly and solution loops vectorise.
Chemistry model wrapper for combustion models.
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
virtual ~laminar()
Destructor.
PtrList< volScalarField > & Y
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
A class for managing temporary objects.
virtual tmp< fvScalarMatrix > R(volScalarField &Y) const
Fuel consumption rate matrix.
Abstract base class for turbulence models (RAS, LES and laminar).
Calculate the matrix for implicit and explicit sources.
virtual bool read()
Update properties from given dictionary.