34 namespace laminarFlameSpeedModels
50 Foam::laminarFlameSpeedModels::SCOPE::polynomial::polynomial
52 const dictionary& polyDict
55 FixedList<scalar, 7>(polyDict.lookup(
"coefficients")),
56 ll(polyDict.lookup<scalar>(
"lowerLimit")),
57 ul(polyDict.lookup<scalar>(
"upperLimit")),
58 llv(polyPhi(ll, *this)),
59 ulv(polyPhi(ul, *this)),
73 LFL_(coeffDict.lookup<scalar>(
"lowerFlamabilityLimit")),
74 UFL_(coeffDict.lookup<scalar>(
"upperFlamabilityLimit")),
75 SuPolyL_(coeffDict.subDict(
"lowerSuPolynomial")),
76 SuPolyU_(coeffDict.subDict(
"upperSuPolynomial")),
77 Texp_(coeffDict.lookup<scalar>(
"Texp")),
78 pexp_(coeffDict.lookup<scalar>(
"pexp")),
79 MaPolyL_(coeffDict.subDict(
"lowerMaPolynomial")),
80 MaPolyU_(coeffDict.subDict(
"upperMaPolynomial"))
82 SuPolyL_.ll =
max(SuPolyL_.ll, LFL_) + small;
83 SuPolyU_.ul =
min(SuPolyU_.ul, UFL_) - small;
85 SuPolyL_.lu = 0.5*(SuPolyL_.ul + SuPolyU_.ll);
86 SuPolyU_.lu = SuPolyL_.lu - small;
88 MaPolyL_.lu = 0.5*(MaPolyL_.ul + MaPolyU_.ll);
89 MaPolyU_.lu = MaPolyL_.lu - small;
93 Info<<
"phi Su (T = Tref, p = pref)" <<
endl;
95 for (
int i=0; i<
n; i++)
97 scalar phi = (2.0*i)/
n;
112 inline Foam::scalar Foam::laminarFlameSpeedModels::SCOPE::polyPhi
118 scalar
x = phi - 1.0;
124 +
x*(a[1] +
x*(a[2] +
x*(a[3] +
x*(a[4] +
x*(a[5] +
x*a[6])))))
129 inline Foam::scalar Foam::laminarFlameSpeedModels::SCOPE::SuRef
134 if (phi < LFL_ || phi > UFL_)
139 else if (phi < SuPolyL_.ll)
143 return SuPolyL_.llv*(phi - LFL_)/(SuPolyL_.ll - LFL_);
145 else if (phi > SuPolyU_.ul)
149 return SuPolyU_.ulv*(UFL_ - phi)/(UFL_ - SuPolyU_.ul);
151 else if (phi < SuPolyL_.lu)
154 return polyPhi(phi, SuPolyL_);
156 else if (phi > SuPolyU_.lu)
159 return polyPhi(phi, SuPolyU_);
165 <<
" cannot be handled by SCOPE function with the "
179 if (phi < MaPolyL_.ll)
184 else if (phi > MaPolyU_.ul)
189 else if (phi < SuPolyL_.lu)
192 return polyPhi(phi, MaPolyL_);
194 else if (phi > SuPolyU_.lu)
197 return polyPhi(phi, MaPolyU_);
203 <<
" cannot be handled by SCOPE function with the "
212 inline Foam::scalar Foam::laminarFlameSpeedModels::SCOPE::Su0pTphi
219 static const scalar Tref = 300.0;
220 static const scalar pRef = 1.013e5;
222 return SuRef(phi)*
pow((Tu/Tref), Texp_)*
pow((
p/pRef), pexp_);
233 tmp<volScalarField> tSu0
247 Su0[celli] = Su0pTphi(
p[celli], Tu[celli], phi);
260 Su0p[facei] = Su0pTphi(pp[facei], Tup[facei], phi);
275 tmp<volScalarField> tSu0
289 Su0[celli] = Su0pTphi(
p[celli], Tu[celli], phi[celli]);
322 tmp<volScalarField> tMa
336 ma[celli] = Ma(phi[celli]);
348 map[facei] = Ma(phip[facei]);
359 if (psiuMulticomponentThermo_.containsSpecie(
"ft"))
367 "stoichiometricAirFuelMassRatio",
369 psiuMulticomponentThermo_.properties()
370 )*ft/(scalar(1) - ft)
393 if (psiuMulticomponentThermo_.containsSpecie(
"ft"))
399 psiuMulticomponentThermo_.p(),
400 psiuMulticomponentThermo_.Tu(),
403 "stoichiometricAirFuelMassRatio",
405 psiuMulticomponentThermo_.properties()
406 )*ft/(scalar(1) - ft)
413 psiuMulticomponentThermo_.p(),
414 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.
Mesh data needed to do the Finite Volume discretisation.
const polyMesh & mesh() const
Return reference to polyMesh.
Laminar flame speed obtained from the SCOPE correlation.
tmp< volScalarField > Ma() const
Return the Markstein number.
SCOPE(const dictionary &dict, const dictionary &coeffDict, const psiuMulticomponentThermo &)
Construct from dictionary and psiuMulticomponentThermo.
tmp< volScalarField > operator()() const
Return the laminar flame speed [m/s].
Abstract class for laminar flame speed.
Base-class for combustion fluid thermodynamic properties based on compressibility.
A class for managing temporary objects.
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
addToRunTimeSelectionTable(laminarFlameSpeed, constant, dictionary)
defineTypeNameAndDebug(constant, 0)
errorManipArg< error, int > exit(error &err, const int errNo=1)
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Ostream & endl(Ostream &os)
Add newline and flush stream.
const dimensionSet dimless
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
void pow(LagrangianPatchField< typename powProduct< Type, r >::type > &f, const LagrangianPatchField< Type > &f1)
VolField< scalar > volScalarField
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
const dimensionSet dimVelocity
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.