33 namespace laminarFlameSpeedModels
58 W_(coeffDict.lookup<scalar>(
"W")),
59 eta_(coeffDict.lookup<scalar>(
"eta")),
60 xi_(coeffDict.lookup<scalar>(
"xi")),
61 f_(coeffDict.lookup<scalar>(
"f")),
62 alpha_(coeffDict.lookup<scalar>(
"alpha")),
63 beta_(coeffDict.lookup<scalar>(
"beta"))
75 inline Foam::scalar Foam::laminarFlameSpeedModels::Gulder::SuRef
82 return W_*
pow(phi, eta_)*
exp(-xi_*
sqr(phi - 1.075));
91 inline Foam::scalar Foam::laminarFlameSpeedModels::Gulder::Su0pTphi
99 static const scalar Tref = 300.0;
100 static const scalar pRef = 1.013e5;
102 return SuRef(phi)*
pow((Tu/Tref), alpha_)*
pow((
p/pRef), beta_)*(1 - f_*Yres);
113 tmp<volScalarField> tSu0
127 Su0[celli] = Su0pTphi(
p[celli], Tu[celli], phi, 0);
139 p.boundaryField()[
patchi][facei],
140 Tu.boundaryField()[
patchi][facei],
158 tmp<volScalarField> tSu0
172 Su0[celli] = Su0pTphi(
p[celli], Tu[celli], phi[celli], 0);
184 p.boundaryField()[
patchi][facei],
185 Tu.boundaryField()[
patchi][facei],
186 phi.boundaryField()[
patchi][facei],
197 Foam::laminarFlameSpeedModels::Gulder::Su0pTphi
205 tmp<volScalarField> tSu0
219 Su0[celli] = Su0pTphi(
p[celli], Tu[celli], phi[celli], egr[celli]);
231 p.boundaryField()[
patchi][facei],
232 Tu.boundaryField()[
patchi][facei],
233 phi.boundaryField()[
patchi][facei],
234 egr.boundaryField()[
patchi][facei]
248 psiuMulticomponentThermo_.containsSpecie(
"ft")
249 && psiuMulticomponentThermo_.containsSpecie(
"egr")
257 psiuMulticomponentThermo_.p(),
258 psiuMulticomponentThermo_.Tu(),
261 "stoichiometricAirFuelMassRatio",
263 psiuMulticomponentThermo_.properties()
264 )*ft/
max(1 - ft - egr, small),
268 else if (psiuMulticomponentThermo_.containsSpecie(
"ft"))
274 psiuMulticomponentThermo_.p(),
275 psiuMulticomponentThermo_.Tu(),
278 "stoichiometricAirFuelMassRatio",
280 psiuMulticomponentThermo_.properties()
281 )*ft/
max(1 - ft, small)
288 psiuMulticomponentThermo_.p(),
289 psiuMulticomponentThermo_.Tu(),
#define forAll(list, i)
Loop across all elements in list.
Macros for easy insertion into run-time selection tables.
Generic GeometricField class.
GeometricBoundaryField< Type, GeoMesh, PrimitiveField > Boundary
Type of the boundary field.
static tmp< GeometricField< Type, GeoMesh, PrimitiveField > > New(const word &name, const Internal &, const PtrList< Patch > &, const HashPtrTable< Source > &=HashPtrTable< Source >())
Return a temporary field constructed from name,.
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Laminar flame speed obtained from Gulder's correlation.
Gulder(const dictionary &dict, const dictionary &coeffDict, const psiuMulticomponentThermo &)
Construct from dictionary and psiuMulticomponentThermo.
tmp< volScalarField > operator()() const
Return the laminar flame speed [m/s].
virtual ~Gulder()
Destructor.
Abstract class for laminar flame speed.
Base-class for combustion fluid thermodynamic properties based on compressibility.
A class for managing temporary objects.
addToRunTimeSelectionTable(laminarFlameSpeed, constant, dictionary)
defineTypeNameAndDebug(constant, 0)
dimensionedScalar exp(const dimensionedScalar &ds)
const dimensionSet dimless
void pow(LagrangianPatchField< typename powProduct< Type, r >::type > &f, const LagrangianPatchField< Type > &f1)
VolField< scalar > volScalarField
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
void sqr(LagrangianPatchField< typename outerProduct< Type, Type >::type > &f, const LagrangianPatchField< Type > &f1)
const dimensionSet dimVelocity
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.