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>
109 const scalar
nu =
mu[i]/(
rho[i] + SMALL);
111 const scalar Da =
clamp 114 scalar(1
e-10), scalar(10)
118 const scalar CtauI =
min(C1_/(Da*
sqrt(ReT + 1)), 2.1377);
120 const scalar CgammaI =
121 clamp(C2_*
sqrt(Da*(ReT + 1)), scalar(0.4082), scalar(5));
123 const scalar gammaL =
139 pow(gammaL, exp1_)/(1 -
pow(gammaL, exp2_)),
148 kappa_.correctBoundaryConditions();
154 const scalar
nu =
mu[i]/(
rho[i] + SMALL);
155 const scalar gammaL =
170 pow(gammaL, exp1_)/(1 -
pow(gammaL, exp2_)),
179 kappa_.correctBoundaryConditions();
182 Info<<
"Chemistry time solved max/min : " 185 this->chemistryPtr_->solve(tauStar);
190 template<
class ReactionThermo>
198 template<
class ReactionThermo>
208 this->
thermo().phasePropertyName(typeName +
":Qdot"),
222 tQdot.ref() = kappa_*this->chemistryPtr_->Qdot();
229 template<
class ReactionThermo>
243 C1_ = this->coeffs().getOrDefault(
"C1", 0.05774);
244 C2_ = this->coeffs().getOrDefault(
"C2", 0.5);
245 Cgamma_ = this->coeffs().getOrDefault(
"Cgamma", 2.1377);
246 Ctau_ = this->coeffs().getOrDefault(
"Ctau", 0.4083);
247 exp1_ = this->coeffs().getOrDefault(
"exp1",
EDCexp1[
int(version_)]);
248 exp2_ = this->coeffs().getOrDefault(
"exp2",
EDCexp2[
int(version_)]);
const EDCversions EDCdefaultVersion
virtual tmp< fvScalarMatrix > R(volScalarField &Y) const
Fuel consumption rate matrix.
Type gMin(const FieldField< Field, Type > &f)
const Enum< EDCversions > EDCversionNames
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
virtual bool read()
Update properties from given dictionary.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
compressible::turbulenceModel & turb
dimensionedScalar sqrt(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
dimensionedScalar pow025(const dimensionedScalar &ds)
tmp< volScalarField > trho
label k
Boltzmann constant.
Ignore writing from objectRegistry::writeObject()
const dimensionSet dimless
Dimensionless.
#define forAll(list, i)
Loop across all elements in list.
psiReactionThermo & thermo
const dimensionSet dimVolume(pow3(dimLength))
Laminar combustion model.
const dimensionedScalar e
Elementary charge.
A class for handling words, derived from Foam::string.
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
Abstract base class for turbulence models (RAS, LES and laminar).
Type gMax(const FieldField< Field, Type > &f)
Info<< "Reading field U\"<< 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\"<< endl;autoPtr< compressible::turbulenceModel > turbulence(compressible::turbulenceModel::New(rho, U, phi, thermo))
Eddy Dissipation Concept (EDC) turbulent combustion model.
const dimensionSet dimEnergy
const dimensionedScalar mu
Atomic mass unit.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
PtrList< volScalarField > & Y
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
virtual tmp< volScalarField > Qdot() const
Heat release rate [kg/m/s3].
messageStream Info
Information stream (stdout output on master, null elsewhere)
virtual ~EDC()
Destructor.
virtual void correct()
Correct combustion rate.
A class for managing temporary objects.
Defines the attributes of an object for which implicit objectRegistry management is supported...
virtual tmp< fvScalarMatrix > R(volScalarField &Y) const
Fuel consumption rate matrix.
Do not request registration (bool: false)
virtual bool read()
Update properties from given dictionary.
dimensionSet clamp(const dimensionSet &a, const dimensionSet &range)
static constexpr const zero Zero
Global zero (0)