SCOPELaminarFlameSpeed.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) 2011-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 "IFstream.H"
27 #include "SCOPELaminarFlameSpeed.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace laminarFlameSpeedModels
35 {
37 
39  (
41  SCOPE,
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(polyDict.lookup<scalar>("lowerLimit")),
57  ul(polyDict.lookup<scalar>("upperLimit")),
58  llv(polyPhi(ll, *this)),
59  ulv(polyPhi(ul, *this)),
60  lu(0)
61 {}
62 
63 
65 (
66  const dictionary& dict,
67  const dictionary& coeffDict,
68  const psiuMulticomponentThermo& ct
69 )
70 :
72 
73  LFL_(coeffDict.lookup<scalar>("lowerFlamabilityLimit")),
74  UFL_(coeffDict.lookup<scalar>("upperFlamabilityLimit")),
75  SuPolyL_(coeffDict.subDict("lowerSuPolynomial")),
76  SuPolyU_(coeffDict.subDict("upperSuPolynomial")),
77  Texp_(coeffDict.lookup<scalar>("Texp")),
78  pexp_(coeffDict.lookup<scalar>("pexp")),
79  MaPolyL_(coeffDict.subDict("lowerMaPolynomial")),
80  MaPolyU_(coeffDict.subDict("upperMaPolynomial"))
81 {
82  SuPolyL_.ll = max(SuPolyL_.ll, LFL_) + small;
83  SuPolyU_.ul = min(SuPolyU_.ul, UFL_) - small;
84 
85  SuPolyL_.lu = 0.5*(SuPolyL_.ul + SuPolyU_.ll);
86  SuPolyU_.lu = SuPolyL_.lu - small;
87 
88  MaPolyL_.lu = 0.5*(MaPolyL_.ul + MaPolyU_.ll);
89  MaPolyU_.lu = MaPolyL_.lu - small;
90 
91  if (debug)
92  {
93  Info<< "phi Su (T = Tref, p = pref)" << endl;
94  label n = 200;
95  for (int i=0; i<n; i++)
96  {
97  scalar phi = (2.0*i)/n;
98  Info<< phi << token::TAB << SuRef(phi) << endl;
99  }
100  }
101 }
102 
103 
104 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
105 
107 {}
108 
109 
110 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
111 
112 inline Foam::scalar Foam::laminarFlameSpeedModels::SCOPE::polyPhi
113 (
114  scalar phi,
115  const polynomial& a
116 )
117 {
118  scalar x = phi - 1.0;
119 
120  return
121  a[0]
122  *(
123  scalar(1)
124  + x*(a[1] + x*(a[2] + x*(a[3] + x*(a[4] + x*(a[5] + x*a[6])))))
125  );
126 }
127 
128 
129 inline Foam::scalar Foam::laminarFlameSpeedModels::SCOPE::SuRef
130 (
131  scalar phi
132 ) const
133 {
134  if (phi < LFL_ || phi > UFL_)
135  {
136  // Return 0 beyond the flamibility limits
137  return scalar(0);
138  }
139  else if (phi < SuPolyL_.ll)
140  {
141  // Use linear interpolation between the low end of the
142  // lower polynomial and the lower flamibility limit
143  return SuPolyL_.llv*(phi - LFL_)/(SuPolyL_.ll - LFL_);
144  }
145  else if (phi > SuPolyU_.ul)
146  {
147  // Use linear interpolation between the upper end of the
148  // upper polynomial and the upper flamibility limit
149  return SuPolyU_.ulv*(UFL_ - phi)/(UFL_ - SuPolyU_.ul);
150  }
151  else if (phi < SuPolyL_.lu)
152  {
153  // Evaluate the lower polynomial
154  return polyPhi(phi, SuPolyL_);
155  }
156  else if (phi > SuPolyU_.lu)
157  {
158  // Evaluate the upper polynomial
159  return polyPhi(phi, SuPolyU_);
160  }
161  else
162  {
164  << "phi = " << phi
165  << " cannot be handled by SCOPE function with the "
166  "given coefficients"
167  << exit(FatalError);
168 
169  return scalar(0);
170  }
171 }
172 
173 
175 (
176  scalar phi
177 ) const
178 {
179  if (phi < MaPolyL_.ll)
180  {
181  // Beyond the lower limit assume Ma is constant
182  return MaPolyL_.llv;
183  }
184  else if (phi > MaPolyU_.ul)
185  {
186  // Beyond the upper limit assume Ma is constant
187  return MaPolyU_.ulv;
188  }
189  else if (phi < SuPolyL_.lu)
190  {
191  // Evaluate the lower polynomial
192  return polyPhi(phi, MaPolyL_);
193  }
194  else if (phi > SuPolyU_.lu)
195  {
196  // Evaluate the upper polynomial
197  return polyPhi(phi, MaPolyU_);
198  }
199  else
200  {
202  << "phi = " << phi
203  << " cannot be handled by SCOPE function with the "
204  "given coefficients"
205  << exit(FatalError);
206 
207  return scalar(0);
208  }
209 }
210 
211 
212 inline Foam::scalar Foam::laminarFlameSpeedModels::SCOPE::Su0pTphi
213 (
214  scalar p,
215  scalar Tu,
216  scalar phi
217 ) const
218 {
219  static const scalar Tref = 300.0;
220  static const scalar pRef = 1.013e5;
221 
222  return SuRef(phi)*pow((Tu/Tref), Texp_)*pow((p/pRef), pexp_);
223 }
224 
225 
226 Foam::tmp<Foam::volScalarField> Foam::laminarFlameSpeedModels::SCOPE::Su0pTphi
227 (
228  const volScalarField& p,
229  const volScalarField& Tu,
230  scalar phi
231 ) const
232 {
233  tmp<volScalarField> tSu0
234  (
236  (
237  "Su0",
238  p.mesh(),
240  )
241  );
242 
243  volScalarField& Su0 = tSu0.ref();
244 
245  forAll(Su0, celli)
246  {
247  Su0[celli] = Su0pTphi(p[celli], Tu[celli], phi);
248  }
249 
250  volScalarField::Boundary& Su0Bf = Su0.boundaryFieldRef();
251 
252  forAll(Su0Bf, patchi)
253  {
254  scalarField& Su0p = Su0Bf[patchi];
255  const scalarField& pp = p.boundaryField()[patchi];
256  const scalarField& Tup = Tu.boundaryField()[patchi];
257 
258  forAll(Su0p, facei)
259  {
260  Su0p[facei] = Su0pTphi(pp[facei], Tup[facei], phi);
261  }
262  }
263 
264  return tSu0;
265 }
266 
267 
268 Foam::tmp<Foam::volScalarField> Foam::laminarFlameSpeedModels::SCOPE::Su0pTphi
269 (
270  const volScalarField& p,
271  const volScalarField& Tu,
272  const volScalarField& phi
273 ) const
274 {
275  tmp<volScalarField> tSu0
276  (
278  (
279  "Su0",
280  p.mesh(),
282  )
283  );
284 
285  volScalarField& Su0 = tSu0.ref();
286 
287  forAll(Su0, celli)
288  {
289  Su0[celli] = Su0pTphi(p[celli], Tu[celli], phi[celli]);
290  }
291 
292  volScalarField::Boundary& Su0Bf = Su0.boundaryFieldRef();
293 
294  forAll(Su0Bf, patchi)
295  {
296  scalarField& Su0p = Su0Bf[patchi];
297  const scalarField& pp = p.boundaryField()[patchi];
298  const scalarField& Tup = Tu.boundaryField()[patchi];
299  const scalarField& phip = phi.boundaryField()[patchi];
300 
301  forAll(Su0p, facei)
302  {
303  Su0p[facei] =
304  Su0pTphi
305  (
306  pp[facei],
307  Tup[facei],
308  phip[facei]
309  );
310  }
311  }
312 
313  return tSu0;
314 }
315 
316 
318 (
319  const volScalarField& phi
320 ) const
321 {
322  tmp<volScalarField> tMa
323  (
325  (
326  "Ma",
327  phi.mesh(),
329  )
330  );
331 
332  volScalarField& ma = tMa.ref();
333 
334  forAll(ma, celli)
335  {
336  ma[celli] = Ma(phi[celli]);
337  }
338 
339  volScalarField::Boundary& maBf = ma.boundaryFieldRef();
340 
341  forAll(maBf, patchi)
342  {
343  scalarField& map = maBf[patchi];
344  const scalarField& phip = phi.boundaryField()[patchi];
345 
346  forAll(map, facei)
347  {
348  map[facei] = Ma(phip[facei]);
349  }
350  }
351 
352  return tMa;
353 }
354 
355 
358 {
359  if (psiuMulticomponentThermo_.containsSpecie("ft"))
360  {
361  const volScalarField& ft = psiuMulticomponentThermo_.Y("ft");
362 
363  return Ma
364  (
366  (
367  "stoichiometricAirFuelMassRatio",
368  dimless,
369  psiuMulticomponentThermo_.properties()
370  )*ft/(scalar(1) - ft)
371  );
372  }
373  else
374  {
375  const fvMesh& mesh = psiuMulticomponentThermo_.p().mesh();
376 
377  return tmp<volScalarField>
378  (
380  (
381  "Ma",
382  mesh,
383  dimensionedScalar(dimless, Ma(equivalenceRatio_))
384  )
385  );
386  }
387 }
388 
389 
392 {
393  if (psiuMulticomponentThermo_.containsSpecie("ft"))
394  {
395  const volScalarField& ft = psiuMulticomponentThermo_.Y("ft");
396 
397  return Su0pTphi
398  (
399  psiuMulticomponentThermo_.p(),
400  psiuMulticomponentThermo_.Tu(),
402  (
403  "stoichiometricAirFuelMassRatio",
404  dimless,
405  psiuMulticomponentThermo_.properties()
406  )*ft/(scalar(1) - ft)
407  );
408  }
409  else
410  {
411  return Su0pTphi
412  (
413  psiuMulticomponentThermo_.p(),
414  psiuMulticomponentThermo_.Tu(),
415  equivalenceRatio_
416  );
417  }
418 }
419 
420 
421 // ************************************************************************* //
label n
#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.
GeometricBoundaryField< Type, GeoMesh, PrimitiveField > Boundary
Type of the boundary field.
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,.
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:96
const polyMesh & mesh() const
Return reference to polyMesh.
Definition: fvMesh.H:443
Laminar flame speed obtained from the SCOPE correlation.
tmp< volScalarField > Ma() const
Return the Markstein number.
SCOPE(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.
Definition: tmp.H:55
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
label patchi
addToRunTimeSelectionTable(laminarFlameSpeed, constant, dictionary)
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
messageStream Info
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
void pow(LagrangianPatchField< typename powProduct< Type, r >::type > &f, const LagrangianPatchField< Type > &f1)
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:62
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
const dimensionSet dimVelocity
error FatalError
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
dictionary dict
volScalarField & p