35namespace laminarFlameSpeedModels
58 coeffsDict_(
dict.optionalSubDict(typeName +
"Coeffs").subDict(fuel_)),
59 W_(coeffsDict_.get<scalar>(
"W")),
60 eta_(coeffsDict_.get<scalar>(
"eta")),
61 xi_(coeffsDict_.get<scalar>(
"xi")),
62 f_(coeffsDict_.get<scalar>(
"f")),
63 alpha_(coeffsDict_.get<scalar>(
"alpha")),
64 beta_(coeffsDict_.get<scalar>(
"beta"))
76inline Foam::scalar Foam::laminarFlameSpeedModels::GuldersEGR::SuRef
92inline Foam::scalar Foam::laminarFlameSpeedModels::GuldersEGR::Su0pTphi
100 static const scalar Tref = 300.0;
101 static const scalar pRef = 1.013e5;
103 return SuRef(
phi)*
pow((Tu/Tref), alpha_)*
pow((
p/pRef), beta_)*(1 - f_*Yres);
108Foam::laminarFlameSpeedModels::GuldersEGR::Su0pTphi
115 tmp<volScalarField> tSu0
137 Su0[celli] = Su0pTphi(
p[celli], Tu[celli],
phi, 0.0);
144 forAll(Su0Bf[patchi], facei)
146 Su0Bf[patchi][facei] =
149 p.boundaryField()[patchi][facei],
150 Tu.boundaryField()[patchi][facei],
162Foam::laminarFlameSpeedModels::GuldersEGR::Su0pTphi
170 tmp<volScalarField> tSu0
192 Su0[celli] = Su0pTphi(
p[celli], Tu[celli],
phi[celli], egr[celli]);
199 forAll(Su0Bf[patchi], facei)
201 Su0Bf[patchi][facei] =
204 p.boundaryField()[patchi][facei],
205 Tu.boundaryField()[patchi][facei],
206 phi.boundaryField()[patchi][facei],
207 egr.boundaryField()[patchi][facei]
221 psiuReactionThermo_.composition().contains(
"ft")
222 && psiuReactionThermo_.composition().contains(
"egr")
227 psiuReactionThermo_.p(),
228 psiuReactionThermo_.Tu(),
231 "stoichiometricAirFuelMassRatio",
dimless, psiuReactionThermo_
234 scalar(1)/psiuReactionThermo_.composition().Y(
"ft")
237 psiuReactionThermo_.composition().Y(
"egr")
244 psiuReactionThermo_.p(),
245 psiuReactionThermo_.Tu(),
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
Internal & ref(const bool updateAccessTime=true)
Return a reference to the dimensioned internal field.
GeometricBoundaryField< scalar, fvPatchField, volMesh > Boundary
Type of boundary fields.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Laminar flame speed obtained from Gulder's correlation with EGR modelling.
virtual ~GuldersEGR()
Destructor.
tmp< volScalarField > operator()() const
Return the laminar flame speed [m/s].
Abstract class for laminar flame speed.
Foam::psiuReactionThermo.
A class for managing temporary objects.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
dimensionedScalar exp(const dimensionedScalar &ds)
const dimensionSet dimless
Dimensionless.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
const dimensionSet dimVelocity
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
static constexpr const zero Zero
Global zero (0)
#define forAll(list, i)
Loop across all elements in list.