33template<
class ReactionThermo>
34Foam::combustionModels::EDC<ReactionThermo>::EDC
36 const word& modelType,
39 const word& combustionProperties
42 laminar<ReactionThermo>(modelType,
thermo,
turb, combustionProperties),
77template<
class ReactionThermo>
84template<
class ReactionThermo>
90 const auto&
epsilon = tepsilon();
93 const auto&
mu = tmu();
106 const auto& tc = ttc();
110 const scalar
nu =
mu[i]/(
rho[i] + SMALL);
112 const scalar Da =
clamp
115 scalar(1
e-10), scalar(10)
119 const scalar CtauI =
min(C1_/(Da*
sqrt(ReT + 1)), 2.1377);
121 const scalar CgammaI =
122 clamp(C2_*
sqrt(Da*(ReT + 1)), scalar(0.4082), scalar(5));
124 const scalar gammaL =
140 pow(gammaL, exp1_)/(1 -
pow(gammaL, exp2_)),
149 kappa_.correctBoundaryConditions();
155 const scalar
nu =
mu[i]/(
rho[i] + SMALL);
156 const scalar gammaL =
171 pow(gammaL, exp1_)/(1 -
pow(gammaL, exp2_)),
180 kappa_.correctBoundaryConditions();
184 Info<<
"Chemistry time solved min/max : "
192template<
class ReactionThermo>
200template<
class ReactionThermo>
214 tQdot.ref() = kappa_*this->chemistryPtr_->Qdot();
221template<
class ReactionThermo>
235 C1_ = this->coeffs().getOrDefault(
"C1", 0.05774);
236 C2_ = this->coeffs().getOrDefault(
"C2", 0.5);
237 Cgamma_ = this->coeffs().getOrDefault(
"Cgamma", 2.1377);
238 Ctau_ = this->coeffs().getOrDefault(
"Ctau", 0.4083);
239 exp1_ = this->coeffs().getOrDefault(
"exp1",
EDCexp1[
int(version_)]);
240 exp2_ = this->coeffs().getOrDefault(
"exp2",
EDCexp2[
int(version_)]);
compressible::turbulenceModel & turb
autoPtr< BasicChemistryModel< ReactionThermo > > chemistryPtr_
Pointer to chemistry model.
virtual ReactionThermo & thermo()
Return access to the thermo package.
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
@ NO_REGISTER
Do not request registration (bool: false).
@ REGISTER
Request registration (bool: true).
@ NO_READ
Nothing to be read.
@ AUTO_WRITE
Automatically write from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
const Time & time() const noexcept
Return Time associated with the objectRegistry.
const dictionary & coeffs() const
Return const dictionary of the model.
const volScalarField & rho() const
Return const access to rho.
Switch active() const noexcept
Is combustion active?
const fvMesh & mesh() const
Return const access to the mesh database.
virtual void correct()
Correct combustion rate.
virtual tmp< fvScalarMatrix > R(volScalarField &Y) const
Fuel consumption rate matrix.
virtual ~EDC()
Destructor.
virtual tmp< volScalarField > Qdot() const
Heat release rate [kg/m/s3].
virtual bool read()
Update properties from given dictionary.
tmp< volScalarField > tc() const
Return the chemical time scale.
virtual tmp< fvScalarMatrix > R(volScalarField &Y) const
Fuel consumption rate matrix.
virtual bool read()
Update properties from given dictionary.
Abstract base class for turbulence models (RAS, LES and laminar).
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T, or return the given default value. FatalIOError if it is found and the number of...
A class for managing temporary objects.
A class for handling words, derived from Foam::string.
PtrList< volScalarField > & Y
const Enum< EDCversions > EDCversionNames
const EDCversions EDCdefaultVersion
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
const dimensionSet dimless
Dimensionless.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimEnergy
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
messageStream Info
Information stream (stdout output on master, null elsewhere).
dimensionSet clamp(const dimensionSet &a, const dimensionSet &range)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
Ostream & endl(Ostream &os)
Add newline and flush stream.
dimensionedScalar sqrt(const dimensionedScalar &ds)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
MinMax< Type > gMinMax(const FieldField< Field, Type > &f)
static constexpr const zero Zero
Global zero (0).
const dimensionSet dimVolume(pow3(dimLength))
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
dimensionedScalar pow025(const dimensionedScalar &ds)
tmp< volScalarField > trho
psiReactionThermo & thermo
#define forAll(list, i)
Loop across all elements in list.