35namespace laminarFlameSpeedModels
59 coeffsDict_(
dict.optionalSubDict(typeName +
"Coeffs").subDict(fuel_)),
60 W_(coeffsDict_.get<scalar>(
"W")),
61 eta_(coeffsDict_.get<scalar>(
"eta")),
62 xi_(coeffsDict_.get<scalar>(
"xi")),
63 f_(coeffsDict_.get<scalar>(
"f")),
64 alpha_(coeffsDict_.get<scalar>(
"alpha")),
65 beta_(coeffsDict_.get<scalar>(
"beta"))
77inline Foam::scalar Foam::laminarFlameSpeedModels::Gulders::SuRef
93inline Foam::scalar Foam::laminarFlameSpeedModels::Gulders::Su0pTphi
101 static const scalar Tref = 300.0;
102 static const scalar pRef = 1.013e5;
104 return SuRef(
phi)*
pow((Tu/Tref), alpha_)*
pow((
p/pRef), beta_)*(1 - f_*Yres);
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],
168 tmp<volScalarField> tSu0
190 Su0[celli] = Su0pTphi(
p[celli], Tu[celli],
phi[celli], 0.0);
197 forAll(Su0Bf[patchi], facei)
199 Su0Bf[patchi][facei] =
202 p.boundaryField()[patchi][facei],
203 Tu.boundaryField()[patchi][facei],
204 phi.boundaryField()[patchi][facei],
217 if (psiuReactionThermo_.composition().contains(
"ft"))
219 const volScalarField& ft = psiuReactionThermo_.composition().Y(
"ft");
223 psiuReactionThermo_.p(),
224 psiuReactionThermo_.Tu(),
227 "stoichiometricAirFuelMassRatio",
dimless, psiuReactionThermo_
228 )*ft/
max(1 - ft, SMALL)
235 psiuReactionThermo_.p(),
236 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.
virtual ~Gulders()
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.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
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.