Go to the documentation of this file.
33 template<
class ReactionThermo>
36 const word& modelType,
39 const word& combustionProperties
52 C1_(this->coeffs().getOrDefault(
"C1", 0.05774)),
53 C2_(this->coeffs().getOrDefault(
"C2", 0.5)),
54 Cgamma_(this->coeffs().getOrDefault(
"Cgamma", 2.1377)),
55 Ctau_(this->coeffs().getOrDefault(
"Ctau", 0.4083)),
56 exp1_(this->coeffs().getOrDefault(
"exp1",
EDCexp1[
int(version_)])),
57 exp2_(this->coeffs().getOrDefault(
"exp2",
EDCexp2[
int(version_)])),
62 this->
thermo().phasePropertyName(typeName +
":kappa"),
76 template<
class ReactionThermo>
83 template<
class ReactionThermo>
102 if (version_ == EDCversions::v2016)
109 const scalar
nu =
mu[i]/(
rho[i] + SMALL);
115 const scalar CtauI =
min(C1_/(Da*
sqrt(ReT + 1)), 2.1377);
117 const scalar CgammaI =
120 const scalar gammaL =
136 pow(gammaL, exp1_)/(1 -
pow(gammaL, exp2_)),
148 const scalar
nu =
mu[i]/(
rho[i] + SMALL);
149 const scalar gammaL =
164 pow(gammaL, exp1_)/(1 -
pow(gammaL, exp2_)),
173 Info<<
"Chemistry time solved max/min : "
176 this->chemistryPtr_->solve(tauStar);
181 template<
class ReactionThermo>
189 template<
class ReactionThermo>
199 this->
thermo().phasePropertyName(typeName +
":Qdot"),
213 tQdot.
ref() = kappa_*this->chemistryPtr_->Qdot();
220 template<
class ReactionThermo>
234 C1_ = this->coeffs().getOrDefault(
"C1", 0.05774);
235 C2_ = this->coeffs().getOrDefault(
"C2", 0.5);
236 Cgamma_ = this->coeffs().getOrDefault(
"Cgamma", 2.1377);
237 Ctau_ = this->coeffs().getOrDefault(
"Ctau", 0.4083);
238 exp1_ = this->coeffs().getOrDefault(
"exp1",
EDCexp1[
int(version_)]);
239 exp2_ = this->coeffs().getOrDefault(
"exp2",
EDCexp2[
int(version_)]);
Defines the attributes of an object for which implicit objectRegistry management is supported,...
A class for handling words, derived from Foam::string.
Info<< "Reading field U\n"<< endl;volVectorField U(IOobject("U", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);volScalarField rho(IOobject("rho", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), thermo.rho());volVectorField rhoU(IOobject("rhoU", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *U);volScalarField rhoE(IOobject("rhoE", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *(e+0.5 *magSqr(U)));surfaceScalarField pos(IOobject("pos", runTime.timeName(), mesh), mesh, dimensionedScalar("pos", dimless, 1.0));surfaceScalarField neg(IOobject("neg", runTime.timeName(), mesh), mesh, dimensionedScalar("neg", dimless, -1.0));surfaceScalarField phi("phi", fvc::flux(rhoU));Info<< "Creating turbulence model\n"<< endl;autoPtr< compressible::turbulenceModel > turbulence(compressible::turbulenceModel::New(rho, U, phi, thermo))
const dimensionedScalar mu
Atomic mass unit.
A class for managing temporary objects.
static constexpr const zero Zero
Global zero (0)
const dimensionSet dimEnergy
psiReactionThermo & thermo
virtual tmp< volScalarField > Qdot() const
Heat release rate [kg/m/s3].
const EDCversions EDCdefaultVersion
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
Ostream & endl(Ostream &os)
Add newline and flush stream.
virtual void correct()
Correct combustion rate.
virtual ~EDC()
Destructor.
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
const Enum< EDCversions > EDCversionNames
#define forAll(list, i)
Loop across all elements in list.
Laminar combustion model.
dimensionedScalar pow025(const dimensionedScalar &ds)
tmp< volScalarField > trho
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
#define R(A, B, C, D, E, F, K, M)
messageStream Info
Information stream (stdout output on master, null elsewhere)
Eddy Dissipation Concept (EDC) turbulent combustion model.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
virtual bool read()
Update properties from given dictionary.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
PtrList< volScalarField > & Y
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar sqrt(const dimensionedScalar &ds)
label k
Boltzmann constant.
const dimensionedScalar e
Elementary charge.
Type gMin(const FieldField< Field, Type > &f)
const dimensionSet dimVolume(pow3(dimLength))
compressible::turbulenceModel & turb
virtual tmp< fvScalarMatrix > R(volScalarField &Y) const
Fuel consumption rate matrix.
Abstract base class for turbulence models (RAS, LES and laminar).
Type gMax(const FieldField< Field, Type > &f)
const dimensionSet dimless
Dimensionless.