32 template<
class ReactionThermo>
35 const word& modelType,
36 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>
86 if (integrateReactionRate_)
88 if (fv::localEulerDdt::enabled(this->
mesh()))
91 fv::localEulerDdt::localRDeltaT(this->
mesh());
93 if (this->coeffs().
found(
"maxIntegrationTime"))
95 scalar maxIntegrationTime
100 this->chemistryPtr_->solve
102 min(1.0/rDeltaT, maxIntegrationTime)()
107 this->chemistryPtr_->solve((1.0/rDeltaT)());
112 this->chemistryPtr_->solve(this->
mesh().time().deltaTValue());
117 this->chemistryPtr_->calculate();
123 template<
class ReactionThermo>
133 const label specieI =
136 Su += this->chemistryPtr_->RR(specieI);
143 template<
class ReactionThermo>
153 this->
thermo().phasePropertyName(typeName +
":Qdot"),
167 tQdot.
ref() = this->chemistryPtr_->Qdot();
174 template<
class ReactionThermo>
179 integrateReactionRate_ =
180 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/s3].
static word member(const word &name)
Return member (name without the extension)
const dimensionSet dimVolume(pow3(dimLength))
Laminar combustion model.
stressControl lookup("compactNormalStress") >> compactNormalStress
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.
Abstract base class for turbulence models (RAS, LES and laminar).
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
const dimensionSet dimEnergy
virtual ~laminar()
Destructor.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
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.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
virtual tmp< fvScalarMatrix > R(volScalarField &Y) const
Fuel consumption rate matrix.
Calculate the matrix for implicit and explicit sources.
virtual bool read()
Update properties from given dictionary.