Go to the documentation of this file.
33 namespace laminarFlameSpeedModels
56 coeffsDict_(
dict.subDict(typeName +
"Coeffs").subDict(fuel_)),
57 pPoints_(coeffsDict_.lookup(
"pPoints")),
58 EqRPoints_(coeffsDict_.lookup(
"EqRPoints")),
59 alpha_(coeffsDict_.lookup(
"alpha")),
60 beta_(coeffsDict_.lookup(
"beta")),
63 checkPointsMonotonicity(
"equivalenceRatio", EqRPoints_);
64 checkPointsMonotonicity(
"pressure", pPoints_);
65 checkCoefficientArrayShape(
"alpha", alpha_);
66 checkCoefficientArrayShape(
"beta", beta_);
84 for (
label i = 1; i <
x.size(); i ++)
91 ) <<
"Data points for the " <<
name
92 <<
" do not increase monotonically" <<
endl
107 ok &=
x.size() == EqRPoints_.size() - 1;
111 ok &=
x[i].size() == pPoints_.size();
115 ok &=
x[i][j].size() ==
x[i][0].size();
124 ) <<
"Inconsistent size of " <<
name <<
" coefficients array" <<
endl
139 if (
x < xPoints.first())
143 xLim = xPoints.first();
147 else if (
x > xPoints.last())
149 xIndex = xPoints.
size() - 2;
151 xLim = xPoints.last();
155 for (xIndex = 0;
x > xPoints[xIndex+1]; xIndex ++)
160 xXi = (
x - xPoints[xIndex])/(xPoints[xIndex+1] - xPoints[xIndex]);
192 for (
label i = 1; i < coeffs.
size(); i++)
194 y += i*coeffs[i]*xPow;
203 const label EqRIndex,
212 polynomial(beta_[EqRIndex][pIndex],EqR)
220 const label EqRIndex,
228 polynomial(alpha_[EqRIndex][pIndex],EqR)
229 *THatPowB(EqRIndex, pIndex, EqR, Tu);
236 const label EqRIndex,
243 scalar
A = polynomial(alpha_[EqRIndex][pIndex], EqRLim);
244 scalar dA = dPolynomial(alpha_[EqRIndex][pIndex], EqRLim);
245 scalar dB = dPolynomial(beta_[EqRIndex][pIndex], EqRLim);
246 scalar TB = THatPowB(EqRIndex, pIndex, EqRLim, Tu);
249 return max(TB*(
A + (dA +
A*
log(Tu/TRef_)*dB)*(EqR - EqRLim)), 0.0);
262 label EqRIndex, pIndex;
267 EqRInRange = interval(EqRPoints_, EqR, EqRIndex, EqRXi, EqRLim);
269 interval(pPoints_,
p, pIndex, pXi, pLim);
271 for (
label pI = 0; pI < 2; pI ++)
275 s = correlationInRange(EqRIndex, pIndex + pI, EqR, Tu);
279 s = correlationOutOfRange(EqRIndex, pIndex + pI, EqR, EqRLim, Tu);
313 if (psiuReactionThermo_.composition().contains(
"ft"))
315 const volScalarField& ft = psiuReactionThermo_.composition().Y(
"ft");
320 psiuReactionThermo_.lookup(
"stoichiometricAirFuelMassRatio")
321 )*ft/
max(1 - ft, SMALL);
325 EqR = equivalenceRatio_;
350 Su0[cellI] = speed(EqR[cellI],
p[cellI], Tu[cellI]);
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
RaviPetersen(const RaviPetersen &)
Construct as copy (not implemented)
A class for handling words, derived from string.
defineTypeNameAndDebug(constant, 0)
#define forAll(list, i)
Loop across all elements in list.
A class for managing temporary objects.
const dimensionSet dimVelocity
Foam::psiuReactionThermo.
virtual ~RaviPetersen()
Destructor.
scalar speed(const scalar EqR, const scalar p, const scalar Tu) const
Return the laminar flame speed [m/s].
scalar THatPowB(const label EqRIndex, const label pIndex, const scalar EqR, const scalar Tu) const
Calculate normalised temperature to the power of the B polynomial.
Ostream & endl(Ostream &os)
Add newline and flush stream.
simpleMatrix< scalar > A(Nc)
scalar dPolynomial(const List< scalar > &coeffs, const scalar x) const
Evaluate a polynomial differential.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
tmp< GeometricField< Type, fvPatchField, volMesh > > Su(const GeometricField< Type, fvPatchField, volMesh > &su, const GeometricField< Type, fvPatchField, volMesh > &vf)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
tmp< volScalarField > operator()() const
Return the laminar flame speed [m/s].
void checkPointsMonotonicity(const word &name, const List< scalar > &x) const
Check that input points are ordered.
A list of keyword definitions, which are a keyword followed by any number of values (e....
dimensionedScalar log(const dimensionedScalar &ds)
Macros for easy insertion into run-time selection tables.
Generic dimensioned Type class.
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().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))
addToRunTimeSelectionTable(laminarFlameSpeed, constant, dictionary)
errorManipArg< error, int > exit(error &err, const int errNo=1)
bool interval(const List< scalar > &xPoints, const scalar x, label &xIndex, scalar &xXi, scalar &xLim) const
Find and interpolate a value in the data point arrays.
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
scalar correlationInRange(const label EqRIndex, const label pIndex, const scalar EqR, const scalar Tu) const
Return the flame speed within the correlation range.
Abstract class for laminar flame speed.
void checkCoefficientArrayShape(const word &name, const List< List< List< scalar > > > &x) const
Check that the coefficient arrays are of the correct shape.
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
scalar correlationOutOfRange(const label EqRIndex, const label pIndex, const scalar EqR, const scalar EqRLim, const scalar Tu) const
Extrapolate the flame speed correlation outside its range.
void size(const label)
Override size to be inconsistent with allocated storage.
Generic GeometricField class.
scalar polynomial(const List< scalar > &coeffs, const scalar x) const
Evaluate a polynomial.
word name(const complex &)
Return a string representation of a complex.