35template<
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;
69template<
class ReactionThermo>
76template<
class ReactionThermo>
80 return this->chemistryPtr_->tc();
84template<
class ReactionThermo>
89 if (integrateReactionRate_)
94 fv::localEulerDdt::localRDeltaT(this->
mesh());
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();
122template<
class ReactionThermo>
132 const label specieI =
133 this->
thermo().composition().species().find(
Y.member());
135 Su += this->chemistryPtr_->RR(specieI);
142template<
class ReactionThermo>
152 this->
thermo().phasePropertyName(typeName +
":Qdot"),
166 tQdot.
ref() = this->chemistryPtr_->Qdot();
173template<
class ReactionThermo>
178 integrateReactionRate_ =
179 this->coeffs().getOrDefault(
"integrateReactionRate",
true);
#define R(A, B, C, D, E, F, K, M)
compressible::turbulenceModel & turb
Chemistry model wrapper for combustion models.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Laminar combustion model.
tmp< volScalarField > tc() const
Return the chemical time scale.
virtual void correct()
Correct combustion rate.
virtual ~laminar()
Destructor.
virtual tmp< volScalarField > Qdot() const
Heat release rate [kg/m/s3].
virtual bool read()
Update properties from given dictionary.
Abstract base class for turbulence models (RAS, LES and laminar).
virtual bool enabled() const
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
A class for managing temporary objects.
A class for handling words, derived from Foam::string.
PtrList< volScalarField > & Y
Calculate the finiteVolume matrix for implicit and explicit sources.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const dimensionSet dimEnergy
fvMatrix< scalar > fvScalarMatrix
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
messageStream Info
Information stream (stdout output on master, null elsewhere)
Ostream & endl(Ostream &os)
Add newline and flush stream.
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
static constexpr const zero Zero
Global zero (0)
const dimensionSet dimVolume(pow3(dimLength))
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
propsDict readIfPresent("fields", acceptFields)