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