33 namespace laminarFlameSpeedModels
57 pPoints_(coeffDict.lookup(
"pPoints")),
58 EqRPoints_(coeffDict.lookup(
"EqRPoints")),
59 alpha_(coeffDict.lookup(
"alpha")),
60 beta_(coeffDict.lookup(
"beta")),
61 TRef_(coeffDict.lookup<scalar>(
"TRef"))
63 checkPointsMonotonicity(coeffDict,
"equivalenceRatio", EqRPoints_);
64 checkPointsMonotonicity(coeffDict,
"pressure", pPoints_);
65 checkCoefficientArrayShape(coeffDict,
"alpha", alpha_);
66 checkCoefficientArrayShape(coeffDict,
"beta", beta_);
78 void Foam::laminarFlameSpeedModels::RaviPetersen::checkPointsMonotonicity
85 for (
label i = 1; i <
x.size(); i ++)
92 ) <<
"Data points for the " <<
name
93 <<
" do not increase monotonically" <<
endl
100 void Foam::laminarFlameSpeedModels::RaviPetersen::checkCoefficientArrayShape
102 const dictionary& coeffDict,
104 const List<List<List<scalar>>>&
x
109 ok &=
x.size() == EqRPoints_.size() - 1;
113 ok &=
x[i].size() == pPoints_.size();
117 ok &=
x[i][j].size() ==
x[i][0].size();
126 ) <<
"Inconsistent size of " <<
name <<
" coefficients array" <<
endl
132 inline bool Foam::laminarFlameSpeedModels::RaviPetersen::interval
134 const List<scalar>& xPoints,
141 if (
x < xPoints.first())
145 xLim = xPoints.first();
149 else if (
x > xPoints.last())
151 xIndex = xPoints.size() - 2;
153 xLim = xPoints.last();
157 for (xIndex = 0;
x > xPoints[xIndex+1]; xIndex ++)
162 xXi = (
x - xPoints[xIndex])/(xPoints[xIndex+1] - xPoints[xIndex]);
169 inline Foam::scalar Foam::laminarFlameSpeedModels::RaviPetersen::polynomial
171 const List<scalar>& coeffs,
186 inline Foam::scalar Foam::laminarFlameSpeedModels::RaviPetersen::dPolynomial
188 const List<scalar>& coeffs,
194 for (
label i = 1; i < coeffs.size(); i++)
196 y += i*coeffs[i]*xPow;
203 inline Foam::scalar Foam::laminarFlameSpeedModels::RaviPetersen::THatPowB
205 const label EqRIndex,
214 polynomial(beta_[EqRIndex][pIndex],EqR)
220 Foam::laminarFlameSpeedModels::RaviPetersen::correlationInRange
222 const label EqRIndex,
230 polynomial(alpha_[EqRIndex][pIndex],EqR)
231 *THatPowB(EqRIndex, pIndex, EqR, Tu);
236 Foam::laminarFlameSpeedModels::RaviPetersen::correlationOutOfRange
238 const label EqRIndex,
245 scalar
A = polynomial(alpha_[EqRIndex][pIndex], EqRLim);
246 scalar dA = dPolynomial(alpha_[EqRIndex][pIndex], EqRLim);
247 scalar dB = dPolynomial(beta_[EqRIndex][pIndex], EqRLim);
248 scalar TB = THatPowB(EqRIndex, pIndex, EqRLim, Tu);
251 return max(TB*(
A + (dA +
A*
log(Tu/TRef_)*dB)*(EqR - EqRLim)), 0.0);
255 inline Foam::scalar Foam::laminarFlameSpeedModels::RaviPetersen::speed
264 label EqRIndex, pIndex;
269 EqRInRange = interval(EqRPoints_, EqR, EqRIndex, EqRXi, EqRLim);
271 interval(pPoints_,
p, pIndex, pXi, pLim);
273 for (
label pI = 0; pI < 2; pI ++)
277 s = correlationInRange(EqRIndex, pIndex + pI, EqR, Tu);
281 s = correlationOutOfRange(EqRIndex, pIndex + pI, EqR, EqRLim, Tu);
315 if (psiuMulticomponentThermo_.containsSpecie(
"egr"))
318 <<
"The " <<
type() <<
" model does not support EGR"
321 else if (psiuMulticomponentThermo_.containsSpecie(
"ft"))
328 "stoichiometricAirFuelMassRatio",
330 psiuMulticomponentThermo_.properties()
331 )*ft/
max(1 - ft, small);
335 EqR = equivalenceRatio_;
352 Su0[celli] = speed(EqR[celli],
p[celli], Tu[celli]);
#define forAll(list, i)
Loop across all elements in list.
Macros for easy insertion into run-time selection tables.
Generic GeometricField class.
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,.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Laminar flame speed obtained from Ravi and Petersen's correlation.
virtual ~RaviPetersen()
Destructor.
RaviPetersen(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.
T & ref() const
Return non-const reference or generate a fatal error.
A class for handling words, derived from string.
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(lagrangian::Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(lagrangian::Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(lagrangian::Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.name(), lagrangian::cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
tmp< VolField< Type > > Su(const VolField< Type > &su, const VolField< Type > &vf)
addToRunTimeSelectionTable(laminarFlameSpeed, constant, dictionary)
defineTypeNameAndDebug(constant, 0)
static const coefficient A("A", dimPressure, 611.21)
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
dimensionedScalar log(const dimensionedScalar &ds)
void pow(LagrangianPatchField< typename powProduct< Type, r >::type > &f, const LagrangianPatchField< Type > &f1)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
const dimensionSet dimVelocity
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.