00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026 #include "Gulders.H"
00027 #include <OpenFOAM/addToRunTimeSelectionTable.H>
00028
00029
00030
00031 namespace Foam
00032 {
00033 namespace laminarFlameSpeedModels
00034 {
00035 defineTypeNameAndDebug(Gulders, 0);
00036
00037 addToRunTimeSelectionTable
00038 (
00039 laminarFlameSpeed,
00040 Gulders,
00041 dictionary
00042 );
00043 }
00044 }
00045
00046
00047
00048
00049 Foam::laminarFlameSpeedModels::Gulders::Gulders
00050 (
00051 const dictionary& dict,
00052 const hhuCombustionThermo& ct
00053 )
00054 :
00055 laminarFlameSpeed(dict, ct),
00056
00057 coeffsDict_(dict.subDict(typeName + "Coeffs").subDict(fuel_)),
00058 W_(readScalar(coeffsDict_.lookup("W"))),
00059 eta_(readScalar(coeffsDict_.lookup("eta"))),
00060 xi_(readScalar(coeffsDict_.lookup("xi"))),
00061 f_(readScalar(coeffsDict_.lookup("f"))),
00062 alpha_(readScalar(coeffsDict_.lookup("alpha"))),
00063 beta_(readScalar(coeffsDict_.lookup("beta")))
00064 {}
00065
00066
00067
00068
00069 Foam::laminarFlameSpeedModels::Gulders::~Gulders()
00070 {}
00071
00072
00073
00074
00075 inline Foam::scalar Foam::laminarFlameSpeedModels::Gulders::SuRef
00076 (
00077 scalar phi
00078 ) const
00079 {
00080 if (phi > SMALL)
00081 {
00082 return W_*pow(phi, eta_)*exp(-xi_*sqr(phi - 1.075));
00083 }
00084 else
00085 {
00086 return 0.0;
00087 }
00088 }
00089
00090
00091 inline Foam::scalar Foam::laminarFlameSpeedModels::Gulders::Su0pTphi
00092 (
00093 scalar p,
00094 scalar Tu,
00095 scalar phi,
00096 scalar Yres
00097 ) const
00098 {
00099 static const scalar Tref = 300.0;
00100 static const scalar pRef = 1.013e5;
00101
00102 return SuRef(phi)*pow((Tu/Tref), alpha_)*pow((p/pRef), beta_)*(1 - f_*Yres);
00103 }
00104
00105
00106 Foam::tmp<Foam::volScalarField> Foam::laminarFlameSpeedModels::Gulders::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> Foam::laminarFlameSpeedModels::Gulders::Su0pTphi
00157 (
00158 const volScalarField& p,
00159 const volScalarField& Tu,
00160 const volScalarField& phi
00161 ) const
00162 {
00163 tmp<volScalarField> tSu0
00164 (
00165 new volScalarField
00166 (
00167 IOobject
00168 (
00169 "Su0",
00170 p.time().timeName(),
00171 p.db(),
00172 IOobject::NO_READ,
00173 IOobject::NO_WRITE
00174 ),
00175 p.mesh(),
00176 dimensionedScalar("Su0", dimVelocity, 0.0)
00177 )
00178 );
00179
00180 volScalarField& Su0 = tSu0();
00181
00182 forAll(Su0, celli)
00183 {
00184 Su0[celli] = Su0pTphi(p[celli], Tu[celli], phi[celli], 0.0);
00185 }
00186
00187 forAll(Su0.boundaryField(), patchi)
00188 {
00189 forAll(Su0.boundaryField()[patchi], facei)
00190 {
00191 Su0.boundaryField()[patchi][facei] =
00192 Su0pTphi
00193 (
00194 p.boundaryField()[patchi][facei],
00195 Tu.boundaryField()[patchi][facei],
00196 phi.boundaryField()[patchi][facei],
00197 0.0
00198 );
00199 }
00200 }
00201
00202 return tSu0;
00203 }
00204
00205
00206 Foam::tmp<Foam::volScalarField>
00207 Foam::laminarFlameSpeedModels::Gulders::operator()() const
00208 {
00209 if (hhuCombustionThermo_.composition().contains("ft"))
00210 {
00211 const volScalarField& ft = hhuCombustionThermo_.composition().Y("ft");
00212
00213 return Su0pTphi
00214 (
00215 hhuCombustionThermo_.p(),
00216 hhuCombustionThermo_.Tu(),
00217 dimensionedScalar
00218 (
00219 hhuCombustionThermo_.lookup("stoichiometricAirFuelMassRatio")
00220 )*ft/((1 + SMALL) - ft)
00221 );
00222 }
00223 else
00224 {
00225 return Su0pTphi
00226 (
00227 hhuCombustionThermo_.p(),
00228 hhuCombustionThermo_.Tu(),
00229 equivalenceRatio_
00230 );
00231 }
00232 }
00233
00234
00235