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

GuldersEGR.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 "GuldersEGR.H"
00027 #include <OpenFOAM/addToRunTimeSelectionTable.H>
00028 
00029 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
00030 
00031 namespace Foam
00032 {
00033 namespace laminarFlameSpeedModels
00034 {
00035     defineTypeNameAndDebug(GuldersEGR, 0);
00036 
00037     addToRunTimeSelectionTable
00038     (
00039         laminarFlameSpeed,
00040         GuldersEGR,
00041         dictionary
00042     );
00043 }
00044 }
00045 
00046 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
00047 
00048 Foam::laminarFlameSpeedModels::GuldersEGR::GuldersEGR
00049 (
00050     const dictionary& dict,
00051     const hhuCombustionThermo& ct
00052 )
00053 :
00054     laminarFlameSpeed(dict, ct),
00055 
00056     coeffsDict_(dict.subDict(typeName + "Coeffs").subDict(fuel_)),
00057     W_(readScalar(coeffsDict_.lookup("W"))),
00058     eta_(readScalar(coeffsDict_.lookup("eta"))),
00059     xi_(readScalar(coeffsDict_.lookup("xi"))),
00060     f_(readScalar(coeffsDict_.lookup("f"))),
00061     alpha_(readScalar(coeffsDict_.lookup("alpha"))),
00062     beta_(readScalar(coeffsDict_.lookup("beta")))
00063 {}
00064 
00065 
00066 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
00067 
00068 Foam::laminarFlameSpeedModels::GuldersEGR::~GuldersEGR()
00069 {}
00070 
00071 
00072 // * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * * //
00073 
00074 inline Foam::scalar Foam::laminarFlameSpeedModels::GuldersEGR::SuRef
00075 (
00076     scalar phi
00077 ) const
00078 {
00079     if (phi > SMALL)
00080     {
00081         return W_*pow(phi, eta_)*exp(-xi_*sqr(phi - 1.075));
00082     }
00083     else
00084     {
00085         return 0.0;
00086     }
00087 }
00088 
00089 
00090 inline Foam::scalar Foam::laminarFlameSpeedModels::GuldersEGR::Su0pTphi
00091 (
00092     scalar p,
00093     scalar Tu,
00094     scalar phi,
00095     scalar Yres
00096 ) const
00097 {
00098     static const scalar Tref = 300.0;
00099     static const scalar pRef = 1.013e5;
00100 
00101     return SuRef(phi)*pow((Tu/Tref), alpha_)*pow((p/pRef), beta_)*(1 - f_*Yres);
00102 }
00103 
00104 
00105 Foam::tmp<Foam::volScalarField>
00106 Foam::laminarFlameSpeedModels::GuldersEGR::Su0pTphi
00107 (
00108     const volScalarField& p,
00109     const volScalarField& Tu,
00110     scalar phi
00111 ) const
00112 {
00113     tmp<volScalarField> tSu0
00114     (
00115         new volScalarField
00116         (
00117             IOobject
00118             (
00119                 "Su0",
00120                 p.time().timeName(),
00121                 p.db(),
00122                 IOobject::NO_READ,
00123                 IOobject::NO_WRITE
00124             ),
00125             p.mesh(),
00126             dimensionedScalar("Su0", dimVelocity, 0.0)
00127         )
00128     );
00129 
00130     volScalarField& Su0 = tSu0();
00131 
00132     forAll(Su0, celli)
00133     {
00134         Su0[celli] = Su0pTphi(p[celli], Tu[celli], phi, 0.0);
00135     }
00136 
00137     forAll(Su0.boundaryField(), patchi)
00138     {
00139         forAll(Su0.boundaryField()[patchi], facei)
00140         {
00141             Su0.boundaryField()[patchi][facei] =
00142                 Su0pTphi
00143                 (
00144                     p.boundaryField()[patchi][facei],
00145                     Tu.boundaryField()[patchi][facei],
00146                     phi,
00147                     0.0
00148                 );
00149         }
00150     }
00151 
00152     return tSu0;
00153 }
00154 
00155 
00156 Foam::tmp<Foam::volScalarField>
00157 Foam::laminarFlameSpeedModels::GuldersEGR::Su0pTphi
00158 (
00159     const volScalarField& p,
00160     const volScalarField& Tu,
00161     const volScalarField& phi,
00162     const volScalarField& egr
00163 ) const
00164 {
00165     tmp<volScalarField> tSu0
00166     (
00167         new volScalarField
00168         (
00169             IOobject
00170             (
00171                 "Su0",
00172                 p.time().timeName(),
00173                 p.db(),
00174                 IOobject::NO_READ,
00175                 IOobject::NO_WRITE
00176             ),
00177             p.mesh(),
00178             dimensionedScalar("Su0", dimVelocity, 0.0)
00179         )
00180     );
00181 
00182     volScalarField& Su0 = tSu0();
00183 
00184     forAll(Su0, celli)
00185     {
00186         Su0[celli] = Su0pTphi(p[celli], Tu[celli], phi[celli], egr[celli]);
00187     }
00188 
00189     forAll(Su0.boundaryField(), patchi)
00190     {
00191         forAll(Su0.boundaryField()[patchi], facei)
00192         {
00193             Su0.boundaryField()[patchi][facei] =
00194                 Su0pTphi
00195                 (
00196                     p.boundaryField()[patchi][facei],
00197                     Tu.boundaryField()[patchi][facei],
00198                     phi.boundaryField()[patchi][facei],
00199                     egr.boundaryField()[patchi][facei]
00200                 );
00201         }
00202     }
00203 
00204     return tSu0;
00205 }
00206 
00207 
00208 Foam::tmp<Foam::volScalarField>
00209 Foam::laminarFlameSpeedModels::GuldersEGR::operator()() const
00210 {
00211     if 
00212     (
00213         hhuCombustionThermo_.composition().contains("ft")
00214      && hhuCombustionThermo_.composition().contains("egr")
00215     )
00216     {
00217         return Su0pTphi
00218         (
00219             hhuCombustionThermo_.p(),
00220             hhuCombustionThermo_.Tu(),
00221             dimensionedScalar
00222             (
00223                 hhuCombustionThermo_.lookup("stoichiometricAirFuelMassRatio")
00224             )/
00225             (
00226                 scalar(1)/hhuCombustionThermo_.composition().Y("ft")
00227               - scalar(1)
00228             ),
00229             hhuCombustionThermo_.composition().Y("egr")
00230         );
00231     }
00232     else
00233     {
00234         return Su0pTphi
00235         (
00236             hhuCombustionThermo_.p(),
00237             hhuCombustionThermo_.Tu(),
00238             equivalenceRatio_
00239         );
00240     }
00241 }
00242 
00243 
00244 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines