SCOPELaminarFlameSpeed.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) 2011-2016 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 "IFstream.H"
27 #include "SCOPELaminarFlameSpeed.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace laminarFlameSpeedModels
35 {
36  defineTypeNameAndDebug(SCOPE, 0);
37 
39  (
40  laminarFlameSpeed,
41  SCOPE,
42  dictionary
43  );
44 }
45 }
46 
47 
48 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
49 
50 Foam::laminarFlameSpeedModels::SCOPE::polynomial::polynomial
51 (
52  const dictionary& polyDict
53 )
54 :
55  FixedList<scalar, 7>(polyDict.lookup("coefficients")),
56  ll(readScalar(polyDict.lookup("lowerLimit"))),
57  ul(readScalar(polyDict.lookup("upperLimit"))),
58  llv(polyPhi(ll, *this)),
59  ulv(polyPhi(ul, *this)),
60  lu(0)
61 {}
62 
63 
64 Foam::laminarFlameSpeedModels::SCOPE::SCOPE
65 (
66  const dictionary& dict,
67  const psiuReactionThermo& ct
68 )
69 :
70  laminarFlameSpeed(dict, ct),
71 
72  coeffsDict_
73  (
74  dictionary
75  (
76  IFstream
77  (
78  fileName
79  (
80  dict.lookup("fuelFile")
81  )
82  )()
83  ).subDict(typeName + "Coeffs")
84  ),
85  LFL_(readScalar(coeffsDict_.lookup("lowerFlamabilityLimit"))),
86  UFL_(readScalar(coeffsDict_.lookup("upperFlamabilityLimit"))),
87  SuPolyL_(coeffsDict_.subDict("lowerSuPolynomial")),
88  SuPolyU_(coeffsDict_.subDict("upperSuPolynomial")),
89  Texp_(readScalar(coeffsDict_.lookup("Texp"))),
90  pexp_(readScalar(coeffsDict_.lookup("pexp"))),
91  MaPolyL_(coeffsDict_.subDict("lowerMaPolynomial")),
92  MaPolyU_(coeffsDict_.subDict("upperMaPolynomial"))
93 {
94  SuPolyL_.ll = max(SuPolyL_.ll, LFL_) + SMALL;
95  SuPolyU_.ul = min(SuPolyU_.ul, UFL_) - SMALL;
96 
97  SuPolyL_.lu = 0.5*(SuPolyL_.ul + SuPolyU_.ll);
98  SuPolyU_.lu = SuPolyL_.lu - SMALL;
99 
100  MaPolyL_.lu = 0.5*(MaPolyL_.ul + MaPolyU_.ll);
101  MaPolyU_.lu = MaPolyL_.lu - SMALL;
102 
103  if (debug)
104  {
105  Info<< "phi Su (T = Tref, p = pref)" << endl;
106  label n = 200;
107  for (int i=0; i<n; i++)
108  {
109  scalar phi = (2.0*i)/n;
110  Info<< phi << token::TAB << SuRef(phi) << endl;
111  }
112  }
113 }
114 
115 
116 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
117 
119 {}
120 
121 
122 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
123 
124 inline Foam::scalar Foam::laminarFlameSpeedModels::SCOPE::polyPhi
125 (
126  scalar phi,
127  const polynomial& a
128 )
129 {
130  scalar x = phi - 1.0;
131 
132  return
133  a[0]
134  *(
135  scalar(1)
136  + x*(a[1] + x*(a[2] + x*(a[3] + x*(a[4] + x*(a[5] + x*a[6])))))
137  );
138 }
139 
140 
141 inline Foam::scalar Foam::laminarFlameSpeedModels::SCOPE::SuRef
142 (
143  scalar phi
144 ) const
145 {
146  if (phi < LFL_ || phi > UFL_)
147  {
148  // Return 0 beyond the flamibility limits
149  return scalar(0);
150  }
151  else if (phi < SuPolyL_.ll)
152  {
153  // Use linear interpolation between the low end of the
154  // lower polynomial and the lower flamibility limit
155  return SuPolyL_.llv*(phi - LFL_)/(SuPolyL_.ll - LFL_);
156  }
157  else if (phi > SuPolyU_.ul)
158  {
159  // Use linear interpolation between the upper end of the
160  // upper polynomial and the upper flamibility limit
161  return SuPolyU_.ulv*(UFL_ - phi)/(UFL_ - SuPolyU_.ul);
162  }
163  else if (phi < SuPolyL_.lu)
164  {
165  // Evaluate the lower polynomial
166  return polyPhi(phi, SuPolyL_);
167  }
168  else if (phi > SuPolyU_.lu)
169  {
170  // Evaluate the upper polynomial
171  return polyPhi(phi, SuPolyU_);
172  }
173  else
174  {
176  << "phi = " << phi
177  << " cannot be handled by SCOPE function with the "
178  "given coefficients"
179  << exit(FatalError);
180 
181  return scalar(0);
182  }
183 }
184 
185 
187 (
188  scalar phi
189 ) const
190 {
191  if (phi < MaPolyL_.ll)
192  {
193  // Beyond the lower limit assume Ma is constant
194  return MaPolyL_.llv;
195  }
196  else if (phi > MaPolyU_.ul)
197  {
198  // Beyond the upper limit assume Ma is constant
199  return MaPolyU_.ulv;
200  }
201  else if (phi < SuPolyL_.lu)
202  {
203  // Evaluate the lower polynomial
204  return polyPhi(phi, MaPolyL_);
205  }
206  else if (phi > SuPolyU_.lu)
207  {
208  // Evaluate the upper polynomial
209  return polyPhi(phi, MaPolyU_);
210  }
211  else
212  {
214  << "phi = " << phi
215  << " cannot be handled by SCOPE function with the "
216  "given coefficients"
217  << exit(FatalError);
218 
219  return scalar(0);
220  }
221 }
222 
223 
224 inline Foam::scalar Foam::laminarFlameSpeedModels::SCOPE::Su0pTphi
225 (
226  scalar p,
227  scalar Tu,
228  scalar phi
229 ) const
230 {
231  static const scalar Tref = 300.0;
232  static const scalar pRef = 1.013e5;
233 
234  return SuRef(phi)*pow((Tu/Tref), Texp_)*pow((p/pRef), pexp_);
235 }
236 
237 
238 Foam::tmp<Foam::volScalarField> Foam::laminarFlameSpeedModels::SCOPE::Su0pTphi
239 (
240  const volScalarField& p,
241  const volScalarField& Tu,
242  scalar phi
243 ) const
244 {
245  tmp<volScalarField> tSu0
246  (
247  new volScalarField
248  (
249  IOobject
250  (
251  "Su0",
252  p.time().timeName(),
253  p.db(),
256  ),
257  p.mesh(),
258  dimensionedScalar("Su0", dimVelocity, 0.0)
259  )
260  );
261 
262  volScalarField& Su0 = tSu0.ref();
263 
264  forAll(Su0, celli)
265  {
266  Su0[celli] = Su0pTphi(p[celli], Tu[celli], phi);
267  }
268 
269  volScalarField::Boundary& Su0Bf = Su0.boundaryFieldRef();
270 
271  forAll(Su0Bf, patchi)
272  {
273  scalarField& Su0p = Su0Bf[patchi];
274  const scalarField& pp = p.boundaryField()[patchi];
275  const scalarField& Tup = Tu.boundaryField()[patchi];
276 
277  forAll(Su0p, facei)
278  {
279  Su0p[facei] = Su0pTphi(pp[facei], Tup[facei], phi);
280  }
281  }
282 
283  return tSu0;
284 }
285 
286 
287 Foam::tmp<Foam::volScalarField> Foam::laminarFlameSpeedModels::SCOPE::Su0pTphi
288 (
289  const volScalarField& p,
290  const volScalarField& Tu,
291  const volScalarField& phi
292 ) const
293 {
294  tmp<volScalarField> tSu0
295  (
296  new volScalarField
297  (
298  IOobject
299  (
300  "Su0",
301  p.time().timeName(),
302  p.db(),
305  ),
306  p.mesh(),
307  dimensionedScalar("Su0", dimVelocity, 0.0)
308  )
309  );
310 
311  volScalarField& Su0 = tSu0.ref();
312 
313  forAll(Su0, celli)
314  {
315  Su0[celli] = Su0pTphi(p[celli], Tu[celli], phi[celli]);
316  }
317 
318  volScalarField::Boundary& Su0Bf = Su0.boundaryFieldRef();
319 
320  forAll(Su0Bf, patchi)
321  {
322  scalarField& Su0p = Su0Bf[patchi];
323  const scalarField& pp = p.boundaryField()[patchi];
324  const scalarField& Tup = Tu.boundaryField()[patchi];
325  const scalarField& phip = phi.boundaryField()[patchi];
326 
327  forAll(Su0p, facei)
328  {
329  Su0p[facei] =
330  Su0pTphi
331  (
332  pp[facei],
333  Tup[facei],
334  phip[facei]
335  );
336  }
337  }
338 
339  return tSu0;
340 }
341 
342 
344 (
345  const volScalarField& phi
346 ) const
347 {
348  tmp<volScalarField> tMa
349  (
350  new volScalarField
351  (
352  IOobject
353  (
354  "Ma",
355  phi.time().timeName(),
356  phi.db(),
359  ),
360  phi.mesh(),
361  dimensionedScalar("Ma", dimless, 0.0)
362  )
363  );
364 
365  volScalarField& ma = tMa.ref();
366 
367  forAll(ma, celli)
368  {
369  ma[celli] = Ma(phi[celli]);
370  }
371 
372  volScalarField::Boundary& maBf = ma.boundaryFieldRef();
373 
374  forAll(maBf, patchi)
375  {
376  scalarField& map = maBf[patchi];
377  const scalarField& phip = phi.boundaryField()[patchi];
378 
379  forAll(map, facei)
380  {
381  map[facei] = Ma(phip[facei]);
382  }
383  }
384 
385  return tMa;
386 }
387 
388 
391 {
393  {
394  const volScalarField& ft = psiuReactionThermo_.composition().Y("ft");
395 
396  return Ma
397  (
399  (
400  psiuReactionThermo_.lookup("stoichiometricAirFuelMassRatio")
401  )*ft/(scalar(1) - ft)
402  );
403  }
404  else
405  {
406  const fvMesh& mesh = psiuReactionThermo_.p().mesh();
407 
408  return tmp<volScalarField>
409  (
410  new volScalarField
411  (
412  IOobject
413  (
414  "Ma",
415  mesh.time().timeName(),
416  mesh,
419  ),
420  mesh,
422  )
423  );
424  }
425 }
426 
427 
430 {
432  {
433  const volScalarField& ft = psiuReactionThermo_.composition().Y("ft");
434 
435  return Su0pTphi
436  (
440  (
441  psiuReactionThermo_.lookup("stoichiometricAirFuelMassRatio")
442  )*ft/(scalar(1) - ft)
443  );
444  }
445  else
446  {
447  return Su0pTphi
448  (
452  );
453  }
454 }
455 
456 
457 // ************************************************************************* //
#define readScalar
Definition: doubleScalar.C:38
surfaceScalarField & phi
bool contains(const word &specieName) const
Does the mixture include this specie?
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
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
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
scalar equivalenceRatio_
Equivalence ratio of a homogeneous mixture.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
virtual const volScalarField & Tu() const =0
Unburnt gas temperature [K].
PtrList< volScalarField > & Y()
Return the mass-fraction fields.
Macros for easy insertion into run-time selection tables.
virtual volScalarField & p()
Pressure [Pa].
Definition: basicThermo.C:466
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
stressControl lookup("compactNormalStress") >> compactNormalStress
dynamicFvMesh & mesh
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
addToRunTimeSelectionTable(laminarFlameSpeed, constant, dictionary)
tmp< volScalarField > Ma() const
Return the Markstein number.
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if succesful.
Definition: doubleScalar.H:63
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
Internal & ref()
Return a reference to the dimensioned internal field.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
label patchi
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const Mesh & mesh() const
Return mesh.
const psiuReactionThermo & psiuReactionThermo_
messageStream Info
label n
virtual basicSpecieMixture & composition()=0
Return the composition of the multi-component mixture.
A class for managing temporary objects.
Definition: PtrList.H:54
tmp< volScalarField > operator()() const
Return the laminar flame speed [m/s].
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:451
const dimensionSet dimVelocity