51Foam::laminarFlameSpeedModels::RaviPetersen::RaviPetersen
58 coeffsDict_(
dict.optionalSubDict(
typeName +
"Coeffs").subDict(fuel_)),
59 pPoints_(coeffsDict_.
lookup(
"pPoints")),
60 EqRPoints_(coeffsDict_.
lookup(
"EqRPoints")),
61 alpha_(coeffsDict_.
lookup(
"alpha")),
62 beta_(coeffsDict_.
lookup(
"beta")),
63 TRef_(coeffsDict_.get<scalar>(
"TRef"))
65 checkPointsMonotonicity(
"equivalenceRatio", EqRPoints_);
66 checkPointsMonotonicity(
"pressure", pPoints_);
67 checkCoefficientArrayShape(
"alpha", alpha_);
68 checkCoefficientArrayShape(
"beta", beta_);
80void Foam::laminarFlameSpeedModels::RaviPetersen::checkPointsMonotonicity
86 for (label i = 1; i <
x.size(); i ++)
91 <<
"Data points for the " <<
name
92 <<
" do not increase monotonically" <<
nl
99void Foam::laminarFlameSpeedModels::RaviPetersen::checkCoefficientArrayShape
102 const List<List<List<scalar>>>&
x
107 ok &=
x.size() == EqRPoints_.size() - 1;
111 ok &=
x[i].size() == pPoints_.size();
115 ok &=
x[i][j].size() ==
x[i][0].size();
122 <<
"Inconsistent size of " <<
name <<
" coefficients array" <<
nl
128inline bool Foam::laminarFlameSpeedModels::RaviPetersen::interval
130 const List<scalar>& xPoints,
137 if (
x < xPoints.first())
141 xLim = xPoints.first();
145 else if (
x > xPoints.last())
147 xIndex = xPoints.size() - 2;
149 xLim = xPoints.last();
153 for (xIndex = 0;
x > xPoints[xIndex+1]; xIndex ++)
158 xXi = (
x - xPoints[xIndex])/(xPoints[xIndex+1] - xPoints[xIndex]);
165inline Foam::scalar Foam::laminarFlameSpeedModels::RaviPetersen::polynomial
167 const List<scalar>& coeffs,
182inline Foam::scalar Foam::laminarFlameSpeedModels::RaviPetersen::dPolynomial
184 const List<scalar>& coeffs,
190 for (label i = 1; i < coeffs.size(); i++)
192 y += i*coeffs[i]*xPow;
199inline Foam::scalar Foam::laminarFlameSpeedModels::RaviPetersen::THatPowB
201 const label EqRIndex,
210 polynomial(beta_[EqRIndex][pIndex],EqR)
216Foam::laminarFlameSpeedModels::RaviPetersen::correlationInRange
218 const label EqRIndex,
226 polynomial(alpha_[EqRIndex][pIndex],EqR)
227 *THatPowB(EqRIndex, pIndex, EqR, Tu);
232Foam::laminarFlameSpeedModels::RaviPetersen::correlationOutOfRange
234 const label EqRIndex,
241 scalar
A = polynomial(alpha_[EqRIndex][pIndex], EqRLim);
242 scalar dA = dPolynomial(alpha_[EqRIndex][pIndex], EqRLim);
243 scalar dB = dPolynomial(beta_[EqRIndex][pIndex], EqRLim);
244 scalar TB = THatPowB(EqRIndex, pIndex, EqRLim, Tu);
247 return max(TB*(
A + (dA +
A*
log(Tu/TRef_)*dB)*(EqR - EqRLim)), 0.0);
251inline Foam::scalar Foam::laminarFlameSpeedModels::RaviPetersen::speed
260 label EqRIndex, pIndex;
265 EqRInRange = interval(EqRPoints_, EqR, EqRIndex, EqRXi, EqRLim);
267 interval(pPoints_,
p, pIndex, pXi, pLim);
269 for (label pI = 0; pI < 2; pI ++)
273 s = correlationInRange(EqRIndex, pIndex + pI, EqR, Tu);
277 s = correlationOutOfRange(EqRIndex, pIndex + pI, EqR, EqRLim, Tu);
298 p.db().newIOobject(
"EqR"),
303 if (psiuReactionThermo_.composition().contains(
"ft"))
305 const volScalarField& ft = psiuReactionThermo_.composition().Y(
"ft");
310 "stoichiometricAirFuelMassRatio",
dimless, psiuReactionThermo_
311 )*ft/
max(1 - ft, SMALL);
315 EqR = equivalenceRatio_;
325 auto& Su0 = tSu0.ref();
329 Su0[celli] = speed(EqR[celli],
p[celli], Tu[celli]);
static const Foam::dimensionedScalar A("", Foam::dimPressure, 611.21)
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
@ NO_REGISTER
Do not request registration (bool: false).
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Laminar flame speed obtained from Ravi and Petersen's correlation.
virtual ~RaviPetersen()
Destructor.
tmp< volScalarField > operator()() const
Return the laminar flame speed [m/s].
Abstract class for laminar flame speed.
const psiuReactionThermo & psiuReactionThermo_
scalar equivalenceRatio_
Equivalence ratio of a homogeneous mixture.
Foam::psiuReactionThermo.
Lookup type of boundary radiation properties.
A class for managing temporary objects.
A class for handling words, derived from Foam::string.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.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))
Namespace for laminar flame speed models.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
const dimensionSet dimless
Dimensionless.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
const dimensionSet dimVelocity
dimensionedScalar log(const dimensionedScalar &ds)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
IOerror FatalIOError
Error stream (stdout output on all processes), with additional 'FOAM FATAL IO ERROR' header text and ...
static constexpr const zero Zero
Global zero (0).
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
errorManipArg< error, int > exit(error &err, const int errNo=1)
void Su(fvMatrix< typename Expr::value_type > &m, const Expr &expression)
constexpr char nl
The newline '\n' character (0x0a).
#define forAll(list, i)
Loop across all elements in list.