Loading...
Searching...
No Matches
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 | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2017 OpenFOAM Foundation
9 Copyright (C) 2023 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "IFstream.H"
32
33// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34
35namespace Foam
36{
38{
40
42 (
43 laminarFlameSpeed,
44 SCOPE,
46 );
47}
48}
49
50
51// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
52
53Foam::laminarFlameSpeedModels::SCOPE::polynomial::polynomial
54(
55 const dictionary& polyDict
56)
57:
58 FixedList<scalar, 7>(polyDict.lookup("coefficients")),
59 ll(polyDict.get<scalar>("lowerLimit")),
60 ul(polyDict.get<scalar>("upperLimit")),
61 llv(polyPhi(ll, *this)),
62 ulv(polyPhi(ul, *this)),
63 lu(0)
64{}
65
66
67Foam::laminarFlameSpeedModels::SCOPE::SCOPE
68(
69 const dictionary& dict,
70 const psiuReactionThermo& ct
71)
72:
73 laminarFlameSpeed(dict, ct),
74
75 coeffsDict_
76 (
78 (
79 IFstream
80 (
81 dict.get<fileName>("fuelFile")
82 )()
83 ).optionalSubDict(typeName + "Coeffs")
84 ),
85 LFL_
86 (
87 coeffsDict_.getCompat<scalar>
88 (
89 "lowerFlammabilityLimit",
90 {{"lowerFlamabilityLimit", 1712}}
91 )
92 ),
93 UFL_
94 (
95 coeffsDict_.getCompat<scalar>
96 (
97 "upperFlammabilityLimit",
98 {{"upperFlamabilityLimit", 1712}}
99 )
100 ),
101 SuPolyL_(coeffsDict_.subDict("lowerSuPolynomial")),
102 SuPolyU_(coeffsDict_.subDict("upperSuPolynomial")),
103 Texp_(coeffsDict_.get<scalar>("Texp")),
104 pexp_(coeffsDict_.get<scalar>("pexp")),
105 MaPolyL_(coeffsDict_.subDict("lowerMaPolynomial")),
106 MaPolyU_(coeffsDict_.subDict("upperMaPolynomial"))
107{
108 SuPolyL_.ll = max(SuPolyL_.ll, LFL_) + SMALL;
109 SuPolyU_.ul = min(SuPolyU_.ul, UFL_) - SMALL;
110
111 SuPolyL_.lu = 0.5*(SuPolyL_.ul + SuPolyU_.ll);
112 SuPolyU_.lu = SuPolyL_.lu - SMALL;
113
114 MaPolyL_.lu = 0.5*(MaPolyL_.ul + MaPolyU_.ll);
115 MaPolyU_.lu = MaPolyL_.lu - SMALL;
116
117 if (debug)
118 {
119 Info<< "phi Su (T = Tref, p = pref)" << endl;
120 label n = 200;
121 for (int i=0; i<n; i++)
122 {
123 scalar phi = (2.0*i)/n;
124 Info<< phi << token::TAB << SuRef(phi) << endl;
125 }
126 }
127}
128
129
130// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
131
133{}
134
135
136// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
137
138inline Foam::scalar Foam::laminarFlameSpeedModels::SCOPE::polyPhi
139(
140 scalar phi,
141 const polynomial& a
142)
143{
144 scalar x = phi - 1.0;
145
146 return
147 a[0]
148 *(
149 scalar(1)
150 + x*(a[1] + x*(a[2] + x*(a[3] + x*(a[4] + x*(a[5] + x*a[6])))))
151 );
152}
153
154
155inline Foam::scalar Foam::laminarFlameSpeedModels::SCOPE::SuRef
156(
157 scalar phi
158) const
159{
160 if (phi < LFL_ || phi > UFL_)
161 {
162 // Return 0 beyond the flammability limits
163 return scalar(0);
164 }
165 else if (phi < SuPolyL_.ll)
166 {
167 // Use linear interpolation between the low end of the
168 // lower polynomial and the lower flammability limit
169 return SuPolyL_.llv*(phi - LFL_)/(SuPolyL_.ll - LFL_);
170 }
171 else if (phi > SuPolyU_.ul)
172 {
173 // Use linear interpolation between the upper end of the
174 // upper polynomial and the upper flammability limit
175 return SuPolyU_.ulv*(UFL_ - phi)/(UFL_ - SuPolyU_.ul);
176 }
177 else if (phi < SuPolyL_.lu)
178 {
179 // Evaluate the lower polynomial
180 return polyPhi(phi, SuPolyL_);
181 }
182 else if (phi > SuPolyU_.lu)
183 {
184 // Evaluate the upper polynomial
185 return polyPhi(phi, SuPolyU_);
186 }
187 else
188 {
190 << "phi = " << phi
191 << " cannot be handled by SCOPE function with the "
192 "given coefficients"
193 << exit(FatalError);
194
195 return scalar(0);
196 }
197}
198
199
201(
202 scalar phi
203) const
204{
205 if (phi < MaPolyL_.ll)
206 {
207 // Beyond the lower limit assume Ma is constant
208 return MaPolyL_.llv;
209 }
210 else if (phi > MaPolyU_.ul)
211 {
212 // Beyond the upper limit assume Ma is constant
213 return MaPolyU_.ulv;
214 }
215 else if (phi < SuPolyL_.lu)
216 {
217 // Evaluate the lower polynomial
218 return polyPhi(phi, MaPolyL_);
219 }
220 else if (phi > SuPolyU_.lu)
221 {
222 // Evaluate the upper polynomial
223 return polyPhi(phi, MaPolyU_);
224 }
225 else
226 {
228 << "phi = " << phi
229 << " cannot be handled by SCOPE function with the "
230 "given coefficients"
231 << exit(FatalError);
232
233 return scalar(0);
234 }
235}
236
237
238inline Foam::scalar Foam::laminarFlameSpeedModels::SCOPE::Su0pTphi
239(
240 scalar p,
241 scalar Tu,
242 scalar phi
243) const
244{
245 static const scalar Tref = 300.0;
246 static const scalar pRef = 1.013e5;
247
248 return SuRef(phi)*pow((Tu/Tref), Texp_)*pow((p/pRef), pexp_);
249}
250
251
252Foam::tmp<Foam::volScalarField> Foam::laminarFlameSpeedModels::SCOPE::Su0pTphi
253(
254 const volScalarField& p,
255 const volScalarField& Tu,
256 scalar phi
257) const
258{
259 auto tSu0 = volScalarField::New
260 (
261 "Su0",
263 p.mesh(),
265 );
266 auto& Su0 = tSu0.ref();
267
268 forAll(Su0, celli)
269 {
270 Su0[celli] = Su0pTphi(p[celli], Tu[celli], phi);
271 }
272
273 volScalarField::Boundary& Su0Bf = Su0.boundaryFieldRef();
274
275 forAll(Su0Bf, patchi)
276 {
277 scalarField& Su0p = Su0Bf[patchi];
278 const scalarField& pp = p.boundaryField()[patchi];
279 const scalarField& Tup = Tu.boundaryField()[patchi];
280
281 forAll(Su0p, facei)
282 {
283 Su0p[facei] = Su0pTphi(pp[facei], Tup[facei], phi);
284 }
285 }
286
287 return tSu0;
288}
289
290
291Foam::tmp<Foam::volScalarField> Foam::laminarFlameSpeedModels::SCOPE::Su0pTphi
292(
293 const volScalarField& p,
294 const volScalarField& Tu,
295 const volScalarField& phi
296) const
297{
298 auto tSu0 = volScalarField::New
299 (
300 "Su0",
302 p.mesh(),
304 );
305 auto& Su0 = tSu0.ref();
306
307 forAll(Su0, celli)
308 {
309 Su0[celli] = Su0pTphi(p[celli], Tu[celli], phi[celli]);
310 }
311
312 volScalarField::Boundary& Su0Bf = Su0.boundaryFieldRef();
313
314 forAll(Su0Bf, patchi)
315 {
316 scalarField& Su0p = Su0Bf[patchi];
317 const scalarField& pp = p.boundaryField()[patchi];
318 const scalarField& Tup = Tu.boundaryField()[patchi];
319 const scalarField& phip = phi.boundaryField()[patchi];
320
321 forAll(Su0p, facei)
322 {
323 Su0p[facei] =
324 Su0pTphi
325 (
326 pp[facei],
327 Tup[facei],
328 phip[facei]
329 );
330 }
331 }
332
333 return tSu0;
334}
335
336
337Foam::tmp<Foam::volScalarField> Foam::laminarFlameSpeedModels::SCOPE::Ma
338(
339 const volScalarField& phi
340) const
341{
342 auto tMa = volScalarField::New
343 (
344 "Ma",
346 phi.mesh(),
348 );
349 auto& ma = tMa.ref();
350
351 forAll(ma, celli)
352 {
353 ma[celli] = Ma(phi[celli]);
354 }
355
356 volScalarField::Boundary& maBf = ma.boundaryFieldRef();
357
358 forAll(maBf, patchi)
359 {
360 scalarField& map = maBf[patchi];
361 const scalarField& phip = phi.boundaryField()[patchi];
362
363 forAll(map, facei)
364 {
365 map[facei] = Ma(phip[facei]);
366 }
367 }
368
369 return tMa;
370}
371
372
373Foam::tmp<Foam::volScalarField>
375{
376 if (psiuReactionThermo_.composition().contains("ft"))
377 {
378 const volScalarField& ft = psiuReactionThermo_.composition().Y("ft");
379
380 return Ma
381 (
383 (
384 "stoichiometricAirFuelMassRatio", dimless, psiuReactionThermo_
385 )*ft/(scalar(1) - ft)
386 );
387 }
388 else
389 {
390 const fvMesh& mesh = psiuReactionThermo_.p().mesh();
391
393 (
394 "Ma",
396 mesh,
397 dimensionedScalar("Ma", dimless, Ma(equivalenceRatio_))
398 );
399 }
400}
401
402
403Foam::tmp<Foam::volScalarField>
405{
406 if (psiuReactionThermo_.composition().contains("ft"))
407 {
408 const volScalarField& ft = psiuReactionThermo_.composition().Y("ft");
409
410 return Su0pTphi
411 (
412 psiuReactionThermo_.p(),
413 psiuReactionThermo_.Tu(),
415 (
416 "stoichiometricAirFuelMassRatio", dimless, psiuReactionThermo_
417 )*ft/(scalar(1) - ft)
418 );
419 }
420 else
421 {
422 return Su0pTphi
423 (
424 psiuReactionThermo_.p(),
425 psiuReactionThermo_.Tu(),
426 equivalenceRatio_
427 );
428 }
429}
430
431
432// ************************************************************************* //
label n
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
GeometricBoundaryField< scalar, fvPatchField, volMesh > Boundary
@ NO_REGISTER
Do not request registration (bool: false).
Laminar flame speed obtained from the SCOPE correlation.
tmp< volScalarField > operator()() const
Return the laminar flame speed [m/s].
tmp< volScalarField > Ma() const
Return the Markstein number.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
volScalarField & p
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
Namespace for laminar flame speed models.
Definition constant.C:29
Namespace for OpenFOAM.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
const dimensionSet dimless
Dimensionless.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
const dimensionSet dimVelocity
messageStream Info
Information stream (stdout output on master, null elsewhere).
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299