RaviPetersen.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2013-2015 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 {
35  defineTypeNameAndDebug(RaviPetersen, 0);
36 
38  (
39  laminarFlameSpeed,
40  RaviPetersen,
41  dictionary
42  );
43 }
44 }
45 
46 
47 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
48 
49 Foam::laminarFlameSpeedModels::RaviPetersen::RaviPetersen
50 (
51  const dictionary& dict,
52  const psiuReactionThermo& ct
53 )
54 :
55  laminarFlameSpeed(dict, ct),
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")),
61  TRef_(readScalar(coeffsDict_.lookup("TRef")))
62 {
63  checkPointsMonotonicity("equivalenceRatio", EqRPoints_);
64  checkPointsMonotonicity("pressure", pPoints_);
65  checkCoefficientArrayShape("alpha", alpha_);
66  checkCoefficientArrayShape("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 word& name,
81  const List<scalar>& x
82 ) const
83 {
84  for (label i = 1; i < x.size(); i ++)
85  {
86  if (x[i] <= x[i-1])
87  {
89  (
90  "laminarFlameSpeedModels::RaviPetersen::checkPointsMonotonicity"
91  "("
92  "const word&, "
93  "const List<scalar>&"
94  ") const",
95  coeffsDict_
96  ) << "Data points for the " << name
97  << " do not increase monotonically" << endl
98  << exit(FatalIOError);
99  }
100  }
101 }
102 
103 
104 void Foam::laminarFlameSpeedModels::RaviPetersen::checkCoefficientArrayShape
105 (
106  const word& name,
107  const List<List<List<scalar> > >& x
108 ) const
109 {
110  bool ok = true;
111 
112  ok &= x.size() == EqRPoints_.size() - 1;
113 
114  forAll(x, i)
115  {
116  ok &= x[i].size() == pPoints_.size();
117 
118  forAll(x[i], j)
119  {
120  ok &= x[i][j].size() == x[i][0].size();
121  }
122  }
123 
124  if (!ok)
125  {
127  (
128  "laminarFlameSpeedModels::RaviPetersen::checkCoefficientArrayShape"
129  "("
130  "const word&, "
131  "const List<List<List<scalar> > >&"
132  ") const",
133  coeffsDict_
134  ) << "Inconsistent size of " << name << " coefficients array" << endl
135  << exit(FatalIOError);
136  }
137 }
138 
139 
140 inline bool Foam::laminarFlameSpeedModels::RaviPetersen::interval
141 (
142  const List<scalar>& xPoints,
143  const scalar x,
144  label& xIndex,
145  scalar& xXi,
146  scalar& xLim
147 ) const
148 {
149  if (x < xPoints.first())
150  {
151  xIndex = 0;
152  xXi = 0.0;
153  xLim = xPoints.first();
154  return false;
155  }
156 
157  else if (x > xPoints.last())
158  {
159  xIndex = xPoints.size() - 2;
160  xXi = 1.0;
161  xLim = xPoints.last();
162  return false;
163  }
164 
165  for (xIndex = 0; x > xPoints[xIndex+1]; xIndex ++)
166  {
167  // increment xIndex until xPoints[xIndex] < x < xPoints[xIndex+1]
168  }
169 
170  xXi = (x - xPoints[xIndex])/(xPoints[xIndex+1] - xPoints[xIndex]);
171  xLim = x;
172 
173  return true;
174 }
175 
176 
177 inline Foam::scalar Foam::laminarFlameSpeedModels::RaviPetersen::polynomial
178 (
179  const List<scalar>& coeffs,
180  const scalar x
181 ) const
182 {
183  scalar xPow = 1.0;
184  scalar y = 0.0;
185  forAll(coeffs, i)
186  {
187  y += coeffs[i]*xPow;
188  xPow *= x;
189  }
190  return y;
191 }
192 
193 
194 inline Foam::scalar Foam::laminarFlameSpeedModels::RaviPetersen::dPolynomial
195 (
196  const List<scalar>& coeffs,
197  const scalar x
198 ) const
199 {
200  scalar xPow = 1.0;
201  scalar y = 0.0;
202  for (label i = 1; i < coeffs.size(); i++)
203  {
204  y += i*coeffs[i]*xPow;
205  xPow *= x;
206  }
207  return y;
208 }
209 
210 
211 inline Foam::scalar Foam::laminarFlameSpeedModels::RaviPetersen::THatPowB
212 (
213  const label EqRIndex,
214  const label pIndex,
215  const scalar EqR,
216  const scalar Tu
217 ) const
218 {
219  return pow
220  (
221  Tu/TRef_,
222  polynomial(beta_[EqRIndex][pIndex],EqR)
223  );
224 }
225 
226 
227 inline Foam::scalar
228 Foam::laminarFlameSpeedModels::RaviPetersen::correlationInRange
229 (
230  const label EqRIndex,
231  const label pIndex,
232  const scalar EqR,
233  const scalar Tu
234 ) const
235 {
236  // standard correlation
237  return
238  polynomial(alpha_[EqRIndex][pIndex],EqR)
239  *THatPowB(EqRIndex, pIndex, EqR, Tu);
240 }
241 
242 
243 inline Foam::scalar
244 Foam::laminarFlameSpeedModels::RaviPetersen::correlationOutOfRange
245 (
246  const label EqRIndex,
247  const label pIndex,
248  const scalar EqR,
249  const scalar EqRLim,
250  const scalar Tu
251 ) const
252 {
253  scalar A = polynomial(alpha_[EqRIndex][pIndex], EqRLim);
254  scalar dA = dPolynomial(alpha_[EqRIndex][pIndex], EqRLim);
255  scalar dB = dPolynomial(beta_[EqRIndex][pIndex], EqRLim);
256  scalar TB = THatPowB(EqRIndex, pIndex, EqRLim, Tu);
257 
258  // linear extrapolation from the bounds of the correlation
259  return max(TB*(A + (dA + A*log(Tu/TRef_)*dB)*(EqR - EqRLim)), 0.0);
260 }
261 
262 
263 inline Foam::scalar Foam::laminarFlameSpeedModels::RaviPetersen::speed
264 (
265  const scalar EqR,
266  const scalar p,
267  const scalar Tu
268 ) const
269 {
270  scalar Su = 0, s;
271 
272  label EqRIndex, pIndex;
273  scalar EqRXi, pXi;
274  scalar EqRLim, pLim;
275  bool EqRInRange;
276 
277  EqRInRange = interval(EqRPoints_, EqR, EqRIndex, EqRXi, EqRLim);
278 
279  interval(pPoints_, p, pIndex, pXi, pLim);
280 
281  for (label pI = 0; pI < 2; pI ++)
282  {
283  if (EqRInRange)
284  {
285  s = correlationInRange(EqRIndex, pIndex + pI, EqR, Tu);
286  }
287  else
288  {
289  s = correlationOutOfRange(EqRIndex, pIndex + pI, EqR, EqRLim, Tu);
290  }
291 
292  Su += (1 - pXi)*s;
293  pXi = 1 - pXi;
294  }
295 
296  return Su;
297 }
298 
299 
300 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
301 
304 {
305  const volScalarField& p = psiuReactionThermo_.p();
306  const volScalarField& Tu = psiuReactionThermo_.Tu();
307 
308  volScalarField EqR
309  (
310  IOobject
311  (
312  "EqR",
313  p.time().timeName(),
314  p.db(),
317  false
318  ),
319  p.mesh(),
320  dimensionedScalar("EqR", dimless, 0.0)
321  );
322 
323  if (psiuReactionThermo_.composition().contains("ft"))
324  {
325  const volScalarField& ft = psiuReactionThermo_.composition().Y("ft");
326 
327  EqR =
329  (
330  psiuReactionThermo_.lookup("stoichiometricAirFuelMassRatio")
331  )*ft/max(1 - ft, SMALL);
332  }
333  else
334  {
335  EqR = equivalenceRatio_;
336  }
337 
339  (
340  new volScalarField
341  (
342  IOobject
343  (
344  "Su0",
345  p.time().timeName(),
346  p.db(),
349  false
350  ),
351  p.mesh(),
352  dimensionedScalar("Su0", dimVelocity, 0.0)
353  )
354  );
355 
356  volScalarField& Su0 = tSu0();
357 
358  forAll(Su0, cellI)
359  {
360  Su0[cellI] = speed(EqR[cellI], p[cellI], Tu[cellI]);
361  }
362 
363  return tSu0;
364 }
365 
366 
367 // ************************************************************************* //
Foam::psiuReactionThermo.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
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 ))
T & last()
Return the last element of the list.
Definition: UListI.H:131
A class for handling words, derived from string.
Definition: word.H:59
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
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:741
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:76
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Abstract class for laminar flame speed.
const Mesh & mesh() const
Return mesh.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:638
T & first()
Return the first element of the list.
Definition: UListI.H:117
Namespace for OpenFOAM.
tmp< GeometricField< Type, fvPatchField, volMesh > > Su(const GeometricField< Type, fvPatchField, volMesh > &su, const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvcSup.C:44
dimensionedScalar log(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
IOerror FatalIOError
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define forAll(list, i)
Definition: UList.H:421
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if succesful.
Definition: doubleScalar.H:63
Macros for easy insertion into run-time selection tables.
addToRunTimeSelectionTable(laminarFlameSpeed, constant, dictionary)
scalar y
const Time & time() const
Return time.
Definition: IOobject.C:251
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
const objectRegistry & db() const
Return the local objectRegistry.
Definition: IOobject.C:245
const dimensionSet dimVelocity
#define FatalIOErrorIn(functionName, ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:325
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
A class for managing temporary objects.
Definition: PtrList.H:118
tmp< volScalarField > operator()() const
Return the laminar flame speed [m/s].
Definition: RaviPetersen.C:303