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