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 "LaunderSharmaKE.H"
00027 #include <OpenFOAM/addToRunTimeSelectionTable.H>
00028
00029 #include <incompressibleRASModels/backwardsCompatibilityWallFunctions.H>
00030
00031
00032
00033 namespace Foam
00034 {
00035 namespace incompressible
00036 {
00037 namespace RASModels
00038 {
00039
00040
00041
00042 defineTypeNameAndDebug(LaunderSharmaKE, 0);
00043 addToRunTimeSelectionTable(RASModel, LaunderSharmaKE, dictionary);
00044
00045
00046
00047 tmp<volScalarField> LaunderSharmaKE::fMu() const
00048 {
00049 return exp(-3.4/sqr(scalar(1) + sqr(k_)/(nu()*epsilonTilda_)/50.0));
00050 }
00051
00052
00053 tmp<volScalarField> LaunderSharmaKE::f2() const
00054 {
00055 return
00056 scalar(1)
00057 - 0.3*exp(-min(sqr(sqr(k_)/(nu()*epsilonTilda_)), scalar(50.0)));
00058 }
00059
00060
00061
00062
00063 LaunderSharmaKE::LaunderSharmaKE
00064 (
00065 const volVectorField& U,
00066 const surfaceScalarField& phi,
00067 transportModel& lamTransportModel
00068 )
00069 :
00070 RASModel(typeName, U, phi, lamTransportModel),
00071
00072 Cmu_
00073 (
00074 dimensioned<scalar>::lookupOrAddToDict
00075 (
00076 "Cmu",
00077 coeffDict_,
00078 0.09
00079 )
00080 ),
00081 C1_
00082 (
00083 dimensioned<scalar>::lookupOrAddToDict
00084 (
00085 "C1",
00086 coeffDict_,
00087 1.44
00088 )
00089 ),
00090 C2_
00091 (
00092 dimensioned<scalar>::lookupOrAddToDict
00093 (
00094 "C2",
00095 coeffDict_,
00096 1.92
00097 )
00098 ),
00099 sigmaEps_
00100 (
00101 dimensioned<scalar>::lookupOrAddToDict
00102 (
00103 "sigmaEps",
00104 coeffDict_,
00105 1.3
00106 )
00107 ),
00108
00109 k_
00110 (
00111 IOobject
00112 (
00113 "k",
00114 runTime_.timeName(),
00115 mesh_,
00116 IOobject::MUST_READ,
00117 IOobject::AUTO_WRITE
00118 ),
00119 mesh_
00120 ),
00121
00122 epsilonTilda_
00123 (
00124 IOobject
00125 (
00126 "epsilon",
00127 runTime_.timeName(),
00128 mesh_,
00129 IOobject::MUST_READ,
00130 IOobject::AUTO_WRITE
00131 ),
00132 mesh_
00133 ),
00134
00135 nut_
00136 (
00137 IOobject
00138 (
00139 "nut",
00140 runTime_.timeName(),
00141 mesh_,
00142 IOobject::NO_READ,
00143 IOobject::AUTO_WRITE
00144 ),
00145 autoCreateLowReNut("nut", mesh_)
00146 )
00147 {
00148 nut_ = Cmu_*fMu()*sqr(k_)/(epsilonTilda_ + epsilonSmall_);
00149 nut_.correctBoundaryConditions();
00150
00151 printCoeffs();
00152 }
00153
00154
00155
00156
00157 tmp<volSymmTensorField> LaunderSharmaKE::R() const
00158 {
00159 return tmp<volSymmTensorField>
00160 (
00161 new volSymmTensorField
00162 (
00163 IOobject
00164 (
00165 "R",
00166 runTime_.timeName(),
00167 mesh_,
00168 IOobject::NO_READ,
00169 IOobject::NO_WRITE
00170 ),
00171 ((2.0/3.0)*I)*k_ - nut_*twoSymm(fvc::grad(U_)),
00172 k_.boundaryField().types()
00173 )
00174 );
00175 }
00176
00177
00178 tmp<volSymmTensorField> LaunderSharmaKE::devReff() const
00179 {
00180 return tmp<volSymmTensorField>
00181 (
00182 new volSymmTensorField
00183 (
00184 IOobject
00185 (
00186 "devRhoReff",
00187 runTime_.timeName(),
00188 mesh_,
00189 IOobject::NO_READ,
00190 IOobject::NO_WRITE
00191 ),
00192 -nuEff()*dev(twoSymm(fvc::grad(U_)))
00193 )
00194 );
00195 }
00196
00197
00198 tmp<fvVectorMatrix> LaunderSharmaKE::divDevReff(volVectorField& U) const
00199 {
00200 return
00201 (
00202 - fvm::laplacian(nuEff(), U)
00203 - fvc::div(nuEff()*dev(fvc::grad(U)().T()))
00204 );
00205 }
00206
00207
00208 bool LaunderSharmaKE::read()
00209 {
00210 if (RASModel::read())
00211 {
00212 Cmu_.readIfPresent(coeffDict());
00213 C1_.readIfPresent(coeffDict());
00214 C2_.readIfPresent(coeffDict());
00215 sigmaEps_.readIfPresent(coeffDict());
00216
00217 return true;
00218 }
00219 else
00220 {
00221 return false;
00222 }
00223 }
00224
00225
00226 void LaunderSharmaKE::correct()
00227 {
00228 RASModel::correct();
00229
00230 if (!turbulence_)
00231 {
00232 return;
00233 }
00234
00235 volScalarField S2 = 2*magSqr(symm(fvc::grad(U_)));
00236
00237 volScalarField G("RASModel::G", nut_*S2);
00238
00239 volScalarField E = 2.0*nu()*nut_*fvc::magSqrGradGrad(U_);
00240 volScalarField D = 2.0*nu()*magSqr(fvc::grad(sqrt(k_)));
00241
00242
00243
00244
00245 tmp<fvScalarMatrix> epsEqn
00246 (
00247 fvm::ddt(epsilonTilda_)
00248 + fvm::div(phi_, epsilonTilda_)
00249 - fvm::laplacian(DepsilonEff(), epsilonTilda_)
00250 ==
00251 C1_*G*epsilonTilda_/k_
00252 - fvm::Sp(C2_*f2()*epsilonTilda_/k_, epsilonTilda_)
00253 + E
00254 );
00255
00256 epsEqn().relax();
00257 solve(epsEqn);
00258 bound(epsilonTilda_, epsilon0_);
00259
00260
00261
00262
00263 tmp<fvScalarMatrix> kEqn
00264 (
00265 fvm::ddt(k_)
00266 + fvm::div(phi_, k_)
00267 - fvm::laplacian(DkEff(), k_)
00268 ==
00269 G - fvm::Sp((epsilonTilda_ + D)/k_, k_)
00270 );
00271
00272 kEqn().relax();
00273 solve(kEqn);
00274 bound(k_, k0_);
00275
00276
00277
00278 nut_ == Cmu_*fMu()*sqr(k_)/epsilonTilda_;
00279 }
00280
00281
00282
00283
00284 }
00285 }
00286 }
00287
00288