32namespace combustionModels
37template<
class ReactionThermo,
class ThermoType>
38eddyDissipationModelBase<ReactionThermo, ThermoType>::eddyDissipationModelBase
40 const word& modelType,
43 const word& combustionProperties
53 CEDC_(this->coeffs().getScalar(
"CEDC"))
59template<
class ReactionThermo,
class ThermoType>
67template<
class ReactionThermo,
class ThermoType>
81template<
class ReactionThermo,
class ThermoType>
88 this->singleMixturePtr_->fresCorrect();
90 const label fuelI = this->singleMixturePtr_->fuelIndex();
92 const volScalarField& YFuel = this->thermo_.composition().Y()[fuelI];
96 if (this->thermo_.composition().contains(
"O2"))
102 *
min(YFuel, YO2/
s.value())
108 <<
"You selected a combustion model that requires O2 mass"
109 <<
" to be present in the mixture"
116template<
class ReactionThermo,
class ThermoType>
121 this->coeffs().readEntry(
"CEDC", CEDC_);
compressible::turbulenceModel & turb
Standard Eddy Dissipation Model based on the assumption that the reaction rates are controlled by the...
tmp< volScalarField > rtTurb() const
Return the reciprocal of the turbulent mixing time scale.
void correct()
Correct combustion rate.
virtual ~eddyDissipationModelBase()
Destructor.
virtual bool read()
Update properties.
Base class for combustion models using singleStepReactingMixture.
Abstract base class for turbulence models (RAS, LES and laminar).
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.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
compressible::turbulenceModel & turbulence
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
const dimensionSet dimVelocity
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))
errorManipArg< error, int > exit(error &err, const int errNo=1)
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)