FreeFOAM The Cross-Platform CFD Toolkit
Hosted by SourceForge:
Get FreeFOAM at SourceForge.net.
            Fast, secure and Free Open Source software downloads

SCOPELaminarFlameSpeed.C

Go to the documentation of this file.
00001 /*---------------------------------------------------------------------------*\
00002   =========                 |
00003   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
00004    \\    /   O peration     |
00005     \\  /    A nd           | Copyright (C) 1991-2010 OpenCFD Ltd.
00006      \\/     M anipulation  |
00007 -------------------------------------------------------------------------------
00008 License
00009     This file is part of OpenFOAM.
00010 
00011     OpenFOAM is free software: you can redistribute it and/or modify it
00012     under the terms of the GNU General Public License as published by
00013     the Free Software Foundation, either version 3 of the License, or
00014     (at your option) any later version.
00015 
00016     OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
00017     ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
00018     FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
00019     for more details.
00020 
00021     You should have received a copy of the GNU General Public License
00022     along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.
00023 
00024 \*---------------------------------------------------------------------------*/
00025 
00026 #include "SCOPELaminarFlameSpeed.H"
00027 #include <OpenFOAM/addToRunTimeSelectionTable.H>
00028 
00029 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
00030 
00031 namespace Foam
00032 {
00033 namespace laminarFlameSpeedModels
00034 {
00035     defineTypeNameAndDebug(SCOPE, 0);
00036 
00037     addToRunTimeSelectionTable
00038     (
00039         laminarFlameSpeed,
00040         SCOPE,
00041         dictionary
00042     );
00043 }
00044 }
00045 
00046 
00047 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
00048 
00049 Foam::laminarFlameSpeedModels::SCOPE::polynomial::polynomial
00050 (
00051     const dictionary& polyDict
00052 )
00053 :
00054     FixedList<scalar, 7>(polyDict.lookup("coefficients")),
00055     ll(readScalar(polyDict.lookup("lowerLimit"))),
00056     ul(readScalar(polyDict.lookup("upperLimit"))),
00057     llv(polyPhi(ll, *this)),
00058     ulv(polyPhi(ul, *this)),
00059     lu(0)
00060 {}
00061 
00062 
00063 Foam::laminarFlameSpeedModels::SCOPE::SCOPE
00064 (
00065     const dictionary& dict,
00066     const hhuCombustionThermo& ct
00067 )
00068 :
00069     laminarFlameSpeed(dict, ct),
00070 
00071     coeffsDict_(dict.subDict(typeName + "Coeffs").subDict(fuel_)),
00072     LFL_(readScalar(coeffsDict_.lookup("lowerFlamabilityLimit"))),
00073     UFL_(readScalar(coeffsDict_.lookup("upperFlamabilityLimit"))),
00074     SuPolyL_(coeffsDict_.subDict("lowerSuPolynomial")),
00075     SuPolyU_(coeffsDict_.subDict("upperSuPolynomial")),
00076     Texp_(readScalar(coeffsDict_.lookup("Texp"))),
00077     pexp_(readScalar(coeffsDict_.lookup("pexp"))),
00078     MaPolyL_(coeffsDict_.subDict("lowerMaPolynomial")),
00079     MaPolyU_(coeffsDict_.subDict("upperMaPolynomial"))
00080 {
00081     SuPolyL_.ll = max(SuPolyL_.ll, LFL_) + SMALL;
00082     SuPolyU_.ul = min(SuPolyU_.ul, UFL_) - SMALL;
00083 
00084     SuPolyL_.lu = 0.5*(SuPolyL_.ul + SuPolyU_.ll);
00085     SuPolyU_.lu = SuPolyL_.lu - SMALL;
00086 
00087     MaPolyL_.lu = 0.5*(MaPolyL_.ul + MaPolyU_.ll);
00088     MaPolyU_.lu = MaPolyL_.lu - SMALL;
00089 
00090     if (debug)
00091     {
00092         Info<< "phi     Su  (T = Tref, p = pref)" << endl;
00093         label n = 200;
00094         for (int i=0; i<n; i++)
00095         {
00096             scalar phi = (2.0*i)/n;
00097             Info<< phi << token::TAB << SuRef(phi) << endl;
00098         }
00099     }
00100 }
00101 
00102 
00103 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
00104 
00105 Foam::laminarFlameSpeedModels::SCOPE::~SCOPE()
00106 {}
00107 
00108 
00109 // * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * * //
00110 
00111 inline Foam::scalar Foam::laminarFlameSpeedModels::SCOPE::polyPhi
00112 (
00113     scalar phi,
00114     const polynomial& a
00115 )
00116 {
00117     scalar x = phi - 1.0;
00118 
00119     return 
00120         a[0]
00121        *(
00122            scalar(1)
00123          + x*(a[1] + x*(a[2] + x*(a[3] + x*(a[4] + x*(a[5] + x*a[6])))))
00124         );
00125 }
00126 
00127 
00128 inline Foam::scalar Foam::laminarFlameSpeedModels::SCOPE::SuRef
00129 (
00130     scalar phi
00131 ) const
00132 {
00133     if (phi < LFL_ || phi > UFL_)
00134     {
00135         // Return 0 beyond the flamibility limits
00136         return scalar(0);
00137     }
00138     else if (phi < SuPolyL_.ll)
00139     {
00140         // Use linear interpolation between the low end of the
00141         // lower polynomial and the lower flamibility limit
00142         return SuPolyL_.llv*(phi - LFL_)/(SuPolyL_.ll - LFL_);
00143     }
00144     else if (phi > SuPolyU_.ul)
00145     {
00146         // Use linear interpolation between the upper end of the
00147         // upper polynomial and the upper flamibility limit
00148         return SuPolyU_.ulv*(UFL_ - phi)/(UFL_ - SuPolyU_.ul);
00149     }
00150     else if (phi < SuPolyL_.lu)
00151     {
00152         // Evaluate the lower polynomial
00153         return polyPhi(phi, SuPolyL_);
00154     }
00155     else if (phi > SuPolyU_.lu)
00156     {
00157         // Evaluate the upper polynomial
00158         return polyPhi(phi, SuPolyU_);
00159     }
00160     else
00161     {
00162         FatalErrorIn("laminarFlameSpeedModels::SCOPE::SuRef(scalar phi)")
00163             << "phi = " << phi
00164             << " cannot be handled by SCOPE function with the "
00165                "given coefficients"
00166             << exit(FatalError);
00167 
00168         return scalar(0);
00169     }
00170 }
00171 
00172 
00173 inline Foam::scalar Foam::laminarFlameSpeedModels::SCOPE::Ma
00174 (
00175     scalar phi
00176 ) const
00177 {
00178     if (phi < MaPolyL_.ll)
00179     {
00180         // Beyond the lower limit assume Ma is constant
00181         return MaPolyL_.llv;
00182     }
00183     else if (phi > MaPolyU_.ul)
00184     {
00185         // Beyond the upper limit assume Ma is constant
00186         return MaPolyU_.ulv;
00187     }
00188     else if (phi < SuPolyL_.lu)
00189     {
00190         // Evaluate the lower polynomial
00191         return polyPhi(phi, MaPolyL_);
00192     }
00193     else if (phi > SuPolyU_.lu)
00194     {
00195         // Evaluate the upper polynomial
00196         return polyPhi(phi, MaPolyU_);
00197     }
00198     else
00199     {
00200         FatalErrorIn("laminarFlameSpeedModels::SCOPE::Ma(scalar phi)")
00201             << "phi = " << phi
00202             << " cannot be handled by SCOPE function with the "
00203                "given coefficients"
00204             << exit(FatalError);
00205 
00206         return scalar(0);
00207     }
00208 }
00209 
00210 
00211 inline Foam::scalar Foam::laminarFlameSpeedModels::SCOPE::Su0pTphi
00212 (
00213     scalar p,
00214     scalar Tu,
00215     scalar phi
00216 ) const
00217 {
00218     static const scalar Tref = 300.0;
00219     static const scalar pRef = 1.013e5;
00220 
00221     return SuRef(phi)*pow((Tu/Tref), Texp_)*pow((p/pRef), pexp_);
00222 }
00223 
00224 
00225 Foam::tmp<Foam::volScalarField> Foam::laminarFlameSpeedModels::SCOPE::Su0pTphi
00226 (
00227     const volScalarField& p,
00228     const volScalarField& Tu,
00229     scalar phi
00230 ) const
00231 {
00232     tmp<volScalarField> tSu0
00233     (
00234         new volScalarField
00235         (
00236             IOobject
00237             (
00238                 "Su0",
00239                 p.time().timeName(),
00240                 p.db(),
00241                 IOobject::NO_READ,
00242                 IOobject::NO_WRITE
00243             ),
00244             p.mesh(),
00245             dimensionedScalar("Su0", dimVelocity, 0.0)
00246         )
00247     );
00248 
00249     volScalarField& Su0 = tSu0();
00250 
00251     forAll(Su0, celli)
00252     {
00253         Su0[celli] = Su0pTphi(p[celli], Tu[celli], phi);
00254     }
00255 
00256     forAll(Su0.boundaryField(), patchi)
00257     {
00258         scalarField& Su0p = Su0.boundaryField()[patchi];
00259         const scalarField& pp = p.boundaryField()[patchi];
00260         const scalarField& Tup = Tu.boundaryField()[patchi];
00261 
00262         forAll(Su0p, facei)
00263         {
00264             Su0p[facei] = Su0pTphi(pp[facei], Tup[facei], phi);
00265         }
00266     }
00267 
00268     return tSu0;
00269 }
00270 
00271 
00272 Foam::tmp<Foam::volScalarField> Foam::laminarFlameSpeedModels::SCOPE::Su0pTphi
00273 (
00274     const volScalarField& p,
00275     const volScalarField& Tu,
00276     const volScalarField& phi
00277 ) const
00278 {
00279     tmp<volScalarField> tSu0
00280     (
00281         new volScalarField
00282         (
00283             IOobject
00284             (
00285                 "Su0",
00286                 p.time().timeName(),
00287                 p.db(),
00288                 IOobject::NO_READ,
00289                 IOobject::NO_WRITE
00290             ),
00291             p.mesh(),
00292             dimensionedScalar("Su0", dimVelocity, 0.0)
00293         )
00294     );
00295 
00296     volScalarField& Su0 = tSu0();
00297 
00298     forAll(Su0, celli)
00299     {
00300         Su0[celli] = Su0pTphi(p[celli], Tu[celli], phi[celli]);
00301     }
00302 
00303     forAll(Su0.boundaryField(), patchi)
00304     {
00305         scalarField& Su0p = Su0.boundaryField()[patchi];
00306         const scalarField& pp = p.boundaryField()[patchi];
00307         const scalarField& Tup = Tu.boundaryField()[patchi];
00308         const scalarField& phip = phi.boundaryField()[patchi];
00309 
00310         forAll(Su0p, facei)
00311         {
00312             Su0p[facei] =
00313                 Su0pTphi
00314                 (
00315                     pp[facei],
00316                     Tup[facei],
00317                     phip[facei]
00318                 );
00319         }
00320     }
00321 
00322     return tSu0;
00323 }
00324 
00325 
00326 Foam::tmp<Foam::volScalarField> Foam::laminarFlameSpeedModels::SCOPE::Ma
00327 (
00328     const volScalarField& phi
00329 ) const
00330 {
00331     tmp<volScalarField> tMa
00332     (
00333         new volScalarField
00334         (
00335             IOobject
00336             (
00337                 "Ma",
00338                 phi.time().timeName(),
00339                 phi.db(),
00340                 IOobject::NO_READ,
00341                 IOobject::NO_WRITE
00342             ),
00343             phi.mesh(),
00344             dimensionedScalar("Ma", dimless, 0.0)
00345         )
00346     );
00347 
00348     volScalarField& ma = tMa();
00349 
00350     forAll(ma, celli)
00351     {
00352         ma[celli] = Ma(phi[celli]);
00353     }
00354 
00355     forAll(ma.boundaryField(), patchi)
00356     {
00357         scalarField& map = ma.boundaryField()[patchi];
00358         const scalarField& phip = phi.boundaryField()[patchi];
00359 
00360         forAll(map, facei)
00361         {
00362             map[facei] = Ma(phip[facei]);
00363         }
00364     }
00365 
00366     return tMa;
00367 }
00368 
00369 
00370 Foam::tmp<Foam::volScalarField>
00371 Foam::laminarFlameSpeedModels::SCOPE::Ma() const
00372 {
00373     if (hhuCombustionThermo_.composition().contains("ft"))
00374     {
00375         const volScalarField& ft = hhuCombustionThermo_.composition().Y("ft");
00376 
00377         return Ma
00378         (
00379             dimensionedScalar
00380             (
00381                 hhuCombustionThermo_.lookup("stoichiometricAirFuelMassRatio")
00382             )*ft/(scalar(1) - ft)
00383         );
00384     }
00385     else
00386     {
00387         const fvMesh& mesh = hhuCombustionThermo_.p().mesh();
00388 
00389         return tmp<volScalarField>
00390         (
00391             new volScalarField
00392             (
00393                 IOobject
00394                 (
00395                     "Ma",
00396                     mesh.time().timeName(),
00397                     mesh,
00398                     IOobject::NO_READ,
00399                     IOobject::NO_WRITE
00400                 ),
00401                 mesh,
00402                 dimensionedScalar("Ma", dimless, Ma(equivalenceRatio_))
00403             )
00404         );
00405     }
00406 }
00407 
00408 
00409 Foam::tmp<Foam::volScalarField>
00410 Foam::laminarFlameSpeedModels::SCOPE::operator()() const
00411 {
00412     if (hhuCombustionThermo_.composition().contains("ft"))
00413     {
00414         const volScalarField& ft = hhuCombustionThermo_.composition().Y("ft");
00415 
00416         return Su0pTphi
00417         (
00418             hhuCombustionThermo_.p(),
00419             hhuCombustionThermo_.Tu(),
00420             dimensionedScalar
00421             (
00422                 hhuCombustionThermo_.lookup("stoichiometricAirFuelMassRatio")
00423             )*ft/(scalar(1) - ft)
00424         );
00425     }
00426     else
00427     {
00428         return Su0pTphi
00429         (
00430             hhuCombustionThermo_.p(),
00431             hhuCombustionThermo_.Tu(),
00432             equivalenceRatio_
00433         );
00434     }
00435 }
00436 
00437 
00438 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines