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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2013-2017 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "RaviPetersen.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace laminarFlameSpeedModels
36 {
37  defineTypeNameAndDebug(RaviPetersen, 0);
38 
40  (
41  laminarFlameSpeed,
42  RaviPetersen,
43  dictionary
44  );
45 }
46 }
47 
48 
49 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
50 
51 Foam::laminarFlameSpeedModels::RaviPetersen::RaviPetersen
52 (
53  const dictionary& dict,
54  const psiuReactionThermo& ct
55 )
56 :
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"))
64 {
65  checkPointsMonotonicity("equivalenceRatio", EqRPoints_);
66  checkPointsMonotonicity("pressure", pPoints_);
67  checkCoefficientArrayShape("alpha", alpha_);
68  checkCoefficientArrayShape("beta", beta_);
69 }
70 
71 
72 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
73 
75 {}
76 
77 
78 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
79 
80 void Foam::laminarFlameSpeedModels::RaviPetersen::checkPointsMonotonicity
81 (
82  const word& name,
83  const List<scalar>& x
84 ) const
85 {
86  for (label i = 1; i < x.size(); i ++)
87  {
88  if (x[i] <= x[i-1])
89  {
90  FatalIOErrorInFunction(coeffsDict_)
91  << "Data points for the " << name
92  << " do not increase monotonically" << nl
93  << exit(FatalIOError);
94  }
95  }
96 }
97 
98 
99 void Foam::laminarFlameSpeedModels::RaviPetersen::checkCoefficientArrayShape
100 (
101  const word& name,
102  const List<List<List<scalar>>>& x
103 ) const
104 {
105  bool ok = true;
106 
107  ok &= x.size() == EqRPoints_.size() - 1;
108 
109  forAll(x, i)
110  {
111  ok &= x[i].size() == pPoints_.size();
112 
113  forAll(x[i], j)
114  {
115  ok &= x[i][j].size() == x[i][0].size();
116  }
117  }
118 
119  if (!ok)
120  {
121  FatalIOErrorInFunction(coeffsDict_)
122  << "Inconsistent size of " << name << " coefficients array" << nl
123  << exit(FatalIOError);
124  }
125 }
126 
127 
128 inline bool Foam::laminarFlameSpeedModels::RaviPetersen::interval
129 (
130  const List<scalar>& xPoints,
131  const scalar x,
132  label& xIndex,
133  scalar& xXi,
134  scalar& xLim
135 ) const
136 {
137  if (x < xPoints.first())
138  {
139  xIndex = 0;
140  xXi = 0.0;
141  xLim = xPoints.first();
142  return false;
143  }
144 
145  else if (x > xPoints.last())
146  {
147  xIndex = xPoints.size() - 2;
148  xXi = 1.0;
149  xLim = xPoints.last();
150  return false;
151  }
152 
153  for (xIndex = 0; x > xPoints[xIndex+1]; xIndex ++)
154  {
155  // increment xIndex until xPoints[xIndex] < x < xPoints[xIndex+1]
156  }
157 
158  xXi = (x - xPoints[xIndex])/(xPoints[xIndex+1] - xPoints[xIndex]);
159  xLim = x;
160 
161  return true;
162 }
163 
164 
165 inline Foam::scalar Foam::laminarFlameSpeedModels::RaviPetersen::polynomial
166 (
167  const List<scalar>& coeffs,
168  const scalar x
169 ) const
170 {
171  scalar xPow = 1.0;
172  scalar y = 0.0;
173  forAll(coeffs, i)
174  {
175  y += coeffs[i]*xPow;
176  xPow *= x;
177  }
178  return y;
179 }
180 
181 
182 inline Foam::scalar Foam::laminarFlameSpeedModels::RaviPetersen::dPolynomial
183 (
184  const List<scalar>& coeffs,
185  const scalar x
186 ) const
187 {
188  scalar xPow = 1.0;
189  scalar y = 0.0;
190  for (label i = 1; i < coeffs.size(); i++)
191  {
192  y += i*coeffs[i]*xPow;
193  xPow *= x;
194  }
195  return y;
196 }
197 
198 
199 inline Foam::scalar Foam::laminarFlameSpeedModels::RaviPetersen::THatPowB
200 (
201  const label EqRIndex,
202  const label pIndex,
203  const scalar EqR,
204  const scalar Tu
205 ) const
206 {
207  return pow
208  (
209  Tu/TRef_,
210  polynomial(beta_[EqRIndex][pIndex],EqR)
211  );
212 }
213 
214 
215 inline Foam::scalar
216 Foam::laminarFlameSpeedModels::RaviPetersen::correlationInRange
217 (
218  const label EqRIndex,
219  const label pIndex,
220  const scalar EqR,
221  const scalar Tu
222 ) const
223 {
224  // standard correlation
225  return
226  polynomial(alpha_[EqRIndex][pIndex],EqR)
227  *THatPowB(EqRIndex, pIndex, EqR, Tu);
228 }
229 
230 
231 inline Foam::scalar
232 Foam::laminarFlameSpeedModels::RaviPetersen::correlationOutOfRange
233 (
234  const label EqRIndex,
235  const label pIndex,
236  const scalar EqR,
237  const scalar EqRLim,
238  const scalar Tu
239 ) const
240 {
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);
245 
246  // linear extrapolation from the bounds of the correlation
247  return max(TB*(A + (dA + A*log(Tu/TRef_)*dB)*(EqR - EqRLim)), 0.0);
248 }
249 
250 
251 inline Foam::scalar Foam::laminarFlameSpeedModels::RaviPetersen::speed
252 (
253  const scalar EqR,
254  const scalar p,
255  const scalar Tu
256 ) const
257 {
258  scalar Su = 0, s;
259 
260  label EqRIndex, pIndex;
261  scalar EqRXi, pXi;
262  scalar EqRLim, pLim;
263  bool EqRInRange;
264 
265  EqRInRange = interval(EqRPoints_, EqR, EqRIndex, EqRXi, EqRLim);
266 
267  interval(pPoints_, p, pIndex, pXi, pLim);
268 
269  for (label pI = 0; pI < 2; pI ++)
270  {
271  if (EqRInRange)
272  {
273  s = correlationInRange(EqRIndex, pIndex + pI, EqR, Tu);
274  }
275  else
276  {
277  s = correlationOutOfRange(EqRIndex, pIndex + pI, EqR, EqRLim, Tu);
278  }
279 
280  Su += (1 - pXi)*s;
281  pXi = 1 - pXi;
282  }
283 
284  return Su;
285 }
286 
287 
288 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
289 
292 {
293  const volScalarField& p = psiuReactionThermo_.p();
294  const volScalarField& Tu = psiuReactionThermo_.Tu();
295 
296  volScalarField EqR
297  (
298  IOobject
299  (
300  "EqR",
301  p.time().timeName(),
302  p.db(),
305  false
306  ),
307  p.mesh(),
309  );
310 
311  if (psiuReactionThermo_.composition().contains("ft"))
312  {
313  const volScalarField& ft = psiuReactionThermo_.composition().Y("ft");
314 
315  EqR =
317  (
318  "stoichiometricAirFuelMassRatio", dimless, psiuReactionThermo_
319  )*ft/max(1 - ft, SMALL);
320  }
321  else
322  {
323  EqR = equivalenceRatio_;
324  }
325 
327  (
328  new volScalarField
329  (
330  IOobject
331  (
332  "Su0",
333  p.time().timeName(),
334  p.db(),
337  false
338  ),
339  p.mesh(),
341  )
342  );
343 
344  volScalarField& Su0 = tSu0.ref();
345 
346  forAll(Su0, celli)
347  {
348  Su0[celli] = speed(EqR[celli], p[celli], Tu[celli]);
349  }
350 
351  return tSu0;
352 }
353 
354 
355 // ************************************************************************* //
Foam::IOobject::NO_WRITE
Definition: IOobject.H:195
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::laminarFlameSpeedModels::defineTypeNameAndDebug
defineTypeNameAndDebug(constant, 0)
s
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))
Definition: gmvOutputSpray.H:25
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::dimVelocity
const dimensionSet dimVelocity
Foam::psiuReactionThermo
Foam::psiuReactionThermo.
Definition: psiuReactionThermo.H:55
Foam::laminarFlameSpeedModels::RaviPetersen::~RaviPetersen
virtual ~RaviPetersen()
Destructor.
Definition: RaviPetersen.C:74
A
static const Foam::dimensionedScalar A("", Foam::dimPressure, 611.21)
Foam::FatalIOError
IOerror FatalIOError
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Su
zeroField Su
Definition: alphaSuSp.H:1
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:227
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:42
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:75
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
Foam::laminarFlameSpeedModels::RaviPetersen::operator()
tmp< volScalarField > operator()() const
Return the laminar flame speed [m/s].
Definition: RaviPetersen.C:291
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
Foam::log
dimensionedScalar log(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:262
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::laminarFlameSpeedModels::addToRunTimeSelectionTable
addToRunTimeSelectionTable(laminarFlameSpeed, constant, dictionary)
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Foam::laminarFlameSpeed
Abstract class for laminar flame speed.
Definition: laminarFlameSpeed.H:60
Foam::List< scalar >
RaviPetersen.H
x
x
Definition: LISASMDCalcMethod2.H:52
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
FatalIOErrorInFunction
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:473
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::IOobject::NO_READ
Definition: IOobject.H:188
Foam::dimless
const dimensionSet dimless
Dimensionless.
Definition: dimensionSets.C:189
y
scalar y
Definition: LISASMDCalcMethod1.H:14