RaviPetersen.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2013-2024 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "RaviPetersen.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace laminarFlameSpeedModels
34 {
36 
38  (
42  );
43 }
44 }
45 
46 
47 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
48 
50 (
51  const dictionary& dict,
52  const dictionary& coeffDict,
53  const psiuMulticomponentThermo& ct
54 )
55 :
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"))
62 {
63  checkPointsMonotonicity(coeffDict, "equivalenceRatio", EqRPoints_);
64  checkPointsMonotonicity(coeffDict, "pressure", pPoints_);
65  checkCoefficientArrayShape(coeffDict, "alpha", alpha_);
66  checkCoefficientArrayShape(coeffDict, "beta", beta_);
67 }
68 
69 
70 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
71 
73 {}
74 
75 
76 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
77 
78 void Foam::laminarFlameSpeedModels::RaviPetersen::checkPointsMonotonicity
79 (
80  const dictionary& coeffDict,
81  const word& name,
82  const List<scalar>& x
83 ) const
84 {
85  for (label i = 1; i < x.size(); i ++)
86  {
87  if (x[i] <= x[i-1])
88  {
90  (
91  coeffDict
92  ) << "Data points for the " << name
93  << " do not increase monotonically" << endl
94  << exit(FatalIOError);
95  }
96  }
97 }
98 
99 
100 void Foam::laminarFlameSpeedModels::RaviPetersen::checkCoefficientArrayShape
101 (
102  const dictionary& coeffDict,
103  const word& name,
104  const List<List<List<scalar>>>& x
105 ) const
106 {
107  bool ok = true;
108 
109  ok &= x.size() == EqRPoints_.size() - 1;
110 
111  forAll(x, i)
112  {
113  ok &= x[i].size() == pPoints_.size();
114 
115  forAll(x[i], j)
116  {
117  ok &= x[i][j].size() == x[i][0].size();
118  }
119  }
120 
121  if (!ok)
122  {
124  (
125  coeffDict
126  ) << "Inconsistent size of " << name << " coefficients array" << endl
127  << exit(FatalIOError);
128  }
129 }
130 
131 
132 inline bool Foam::laminarFlameSpeedModels::RaviPetersen::interval
133 (
134  const List<scalar>& xPoints,
135  const scalar x,
136  label& xIndex,
137  scalar& xXi,
138  scalar& xLim
139 ) const
140 {
141  if (x < xPoints.first())
142  {
143  xIndex = 0;
144  xXi = 0.0;
145  xLim = xPoints.first();
146  return false;
147  }
148 
149  else if (x > xPoints.last())
150  {
151  xIndex = xPoints.size() - 2;
152  xXi = 1.0;
153  xLim = xPoints.last();
154  return false;
155  }
156 
157  for (xIndex = 0; x > xPoints[xIndex+1]; xIndex ++)
158  {
159  // increment xIndex until xPoints[xIndex] < x < xPoints[xIndex+1]
160  }
161 
162  xXi = (x - xPoints[xIndex])/(xPoints[xIndex+1] - xPoints[xIndex]);
163  xLim = x;
164 
165  return true;
166 }
167 
168 
169 inline Foam::scalar Foam::laminarFlameSpeedModels::RaviPetersen::polynomial
170 (
171  const List<scalar>& coeffs,
172  const scalar x
173 ) const
174 {
175  scalar xPow = 1.0;
176  scalar y = 0.0;
177  forAll(coeffs, i)
178  {
179  y += coeffs[i]*xPow;
180  xPow *= x;
181  }
182  return y;
183 }
184 
185 
186 inline Foam::scalar Foam::laminarFlameSpeedModels::RaviPetersen::dPolynomial
187 (
188  const List<scalar>& coeffs,
189  const scalar x
190 ) const
191 {
192  scalar xPow = 1.0;
193  scalar y = 0.0;
194  for (label i = 1; i < coeffs.size(); i++)
195  {
196  y += i*coeffs[i]*xPow;
197  xPow *= x;
198  }
199  return y;
200 }
201 
202 
203 inline Foam::scalar Foam::laminarFlameSpeedModels::RaviPetersen::THatPowB
204 (
205  const label EqRIndex,
206  const label pIndex,
207  const scalar EqR,
208  const scalar Tu
209 ) const
210 {
211  return pow
212  (
213  Tu/TRef_,
214  polynomial(beta_[EqRIndex][pIndex],EqR)
215  );
216 }
217 
218 
219 inline Foam::scalar
220 Foam::laminarFlameSpeedModels::RaviPetersen::correlationInRange
221 (
222  const label EqRIndex,
223  const label pIndex,
224  const scalar EqR,
225  const scalar Tu
226 ) const
227 {
228  // standard correlation
229  return
230  polynomial(alpha_[EqRIndex][pIndex],EqR)
231  *THatPowB(EqRIndex, pIndex, EqR, Tu);
232 }
233 
234 
235 inline Foam::scalar
236 Foam::laminarFlameSpeedModels::RaviPetersen::correlationOutOfRange
237 (
238  const label EqRIndex,
239  const label pIndex,
240  const scalar EqR,
241  const scalar EqRLim,
242  const scalar Tu
243 ) const
244 {
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);
249 
250  // linear extrapolation from the bounds of the correlation
251  return max(TB*(A + (dA + A*log(Tu/TRef_)*dB)*(EqR - EqRLim)), 0.0);
252 }
253 
254 
255 inline Foam::scalar Foam::laminarFlameSpeedModels::RaviPetersen::speed
256 (
257  const scalar EqR,
258  const scalar p,
259  const scalar Tu
260 ) const
261 {
262  scalar Su = 0, s;
263 
264  label EqRIndex, pIndex;
265  scalar EqRXi, pXi;
266  scalar EqRLim, pLim;
267  bool EqRInRange;
268 
269  EqRInRange = interval(EqRPoints_, EqR, EqRIndex, EqRXi, EqRLim);
270 
271  interval(pPoints_, p, pIndex, pXi, pLim);
272 
273  for (label pI = 0; pI < 2; pI ++)
274  {
275  if (EqRInRange)
276  {
277  s = correlationInRange(EqRIndex, pIndex + pI, EqR, Tu);
278  }
279  else
280  {
281  s = correlationOutOfRange(EqRIndex, pIndex + pI, EqR, EqRLim, Tu);
282  }
283 
284  Su += (1 - pXi)*s;
285  pXi = 1 - pXi;
286  }
287 
288  return Su;
289 }
290 
291 
292 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
293 
296 {
297  const volScalarField& p = psiuMulticomponentThermo_.p();
298  const volScalarField& Tu = psiuMulticomponentThermo_.Tu();
299 
300  volScalarField EqR
301  (
302  IOobject
303  (
304  "EqR",
305  p.time().name(),
306  p.db(),
309  false
310  ),
311  p.mesh(),
313  );
314 
315  if (psiuMulticomponentThermo_.containsSpecie("egr"))
316  {
318  << "The " << type() << " model does not support EGR"
319  << exit(FatalError);
320  }
321  else if (psiuMulticomponentThermo_.containsSpecie("ft"))
322  {
323  const volScalarField& ft = psiuMulticomponentThermo_.Y("ft");
324 
325  EqR =
327  (
328  "stoichiometricAirFuelMassRatio",
329  dimless,
330  psiuMulticomponentThermo_.properties()
331  )*ft/max(1 - ft, small);
332  }
333  else
334  {
335  EqR = equivalenceRatio_;
336  }
337 
339  (
341  (
342  "Su0",
343  p.mesh(),
345  )
346  );
347 
348  volScalarField& Su0 = tSu0.ref();
349 
350  forAll(Su0, celli)
351  {
352  Su0[celli] = speed(EqR[celli], p[celli], Tu[celli]);
353  }
354 
355  return tSu0;
356 }
357 
358 
359 // ************************************************************************* //
scalar y
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
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...
Definition: IOobject.H:99
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
Laminar flame speed obtained from Ravi and Petersen's correlation.
Definition: RaviPetersen.H:64
RaviPetersen(const dictionary &dict, const dictionary &coeffDict, const psiuMulticomponentThermo &)
Construct from dictionary and psiuMulticomponentThermo.
Definition: RaviPetersen.C:50
tmp< volScalarField > operator()() const
Return the laminar flame speed [m/s].
Definition: RaviPetersen.C:295
Abstract class for laminar flame speed.
Base-class for combustion fluid thermodynamic properties based on compressibility.
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:197
A class for handling words, derived from string.
Definition: word.H:62
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:346
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
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)
Definition: fvcSup.C:44
addToRunTimeSelectionTable(laminarFlameSpeed, constant, dictionary)
static const coefficient A("A", dimPressure, 611.21)
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:258
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)
IOerror FatalIOError
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
const dimensionSet dimVelocity
error FatalError
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.
Definition: POSIX.C:488
dictionary dict
volScalarField & p