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 "LienCubicKELowRe.H"
00027 #include <finiteVolume/wallFvPatch.H>
00028 #include <OpenFOAM/addToRunTimeSelectionTable.H>
00029
00030 #include <incompressibleRASModels/backwardsCompatibilityWallFunctions.H>
00031
00032
00033
00034 namespace Foam
00035 {
00036 namespace incompressible
00037 {
00038 namespace RASModels
00039 {
00040
00041
00042
00043 defineTypeNameAndDebug(LienCubicKELowRe, 0);
00044 addToRunTimeSelectionTable(RASModel, LienCubicKELowRe, dictionary);
00045
00046
00047
00048 LienCubicKELowRe::LienCubicKELowRe
00049 (
00050 const volVectorField& U,
00051 const surfaceScalarField& phi,
00052 transportModel& lamTransportModel
00053 )
00054 :
00055 RASModel(typeName, U, phi, lamTransportModel),
00056
00057 C1_
00058 (
00059 dimensioned<scalar>::lookupOrAddToDict
00060 (
00061 "C1",
00062 coeffDict_,
00063 1.44
00064 )
00065 ),
00066 C2_
00067 (
00068 dimensioned<scalar>::lookupOrAddToDict
00069 (
00070 "C2",
00071 coeffDict_,
00072 1.92
00073 )
00074 ),
00075 sigmak_
00076 (
00077 dimensioned<scalar>::lookupOrAddToDict
00078 (
00079 "sigmak",
00080 coeffDict_,
00081 1.0
00082 )
00083 ),
00084 sigmaEps_
00085 (
00086 dimensioned<scalar>::lookupOrAddToDict
00087 (
00088 "sigmaEps",
00089 coeffDict_,
00090 1.3
00091 )
00092 ),
00093 A1_
00094 (
00095 dimensioned<scalar>::lookupOrAddToDict
00096 (
00097 "A1",
00098 coeffDict_,
00099 1.25
00100 )
00101 ),
00102 A2_
00103 (
00104 dimensioned<scalar>::lookupOrAddToDict
00105 (
00106 "A2",
00107 coeffDict_,
00108 1000.0
00109 )
00110 ),
00111 Ctau1_
00112 (
00113 dimensioned<scalar>::lookupOrAddToDict
00114 (
00115 "Ctau1",
00116 coeffDict_,
00117 -4.0
00118 )
00119 ),
00120 Ctau2_
00121 (
00122 dimensioned<scalar>::lookupOrAddToDict
00123 (
00124 "Ctau2",
00125 coeffDict_,
00126 13.0
00127 )
00128 ),
00129 Ctau3_
00130 (
00131 dimensioned<scalar>::lookupOrAddToDict
00132 (
00133 "Ctau3",
00134 coeffDict_,
00135 -2.0
00136 )
00137 ),
00138 alphaKsi_
00139 (
00140 dimensioned<scalar>::lookupOrAddToDict
00141 (
00142 "alphaKsi",
00143 coeffDict_,
00144 0.9
00145 )
00146 ),
00147 CmuWall_
00148 (
00149 dimensioned<scalar>::lookupOrAddToDict
00150 (
00151 "Cmu",
00152 coeffDict_,
00153 0.09
00154 )
00155 ),
00156 kappa_
00157 (
00158 dimensioned<scalar>::lookupOrAddToDict
00159 (
00160 "kappa",
00161 coeffDict_,
00162 0.41
00163 )
00164 ),
00165 Am_
00166 (
00167 dimensioned<scalar>::lookupOrAddToDict
00168 (
00169 "Am",
00170 coeffDict_,
00171 0.016
00172 )
00173 ),
00174 Aepsilon_
00175 (
00176 dimensioned<scalar>::lookupOrAddToDict
00177 (
00178 "Aepsilon",
00179 coeffDict_,
00180 0.263
00181 )
00182 ),
00183 Amu_
00184 (
00185 dimensioned<scalar>::lookupOrAddToDict
00186 (
00187 "Amu",
00188 coeffDict_,
00189 0.00222
00190 )
00191 ),
00192
00193 k_
00194 (
00195 IOobject
00196 (
00197 "k",
00198 runTime_.timeName(),
00199 mesh_,
00200 IOobject::MUST_READ,
00201 IOobject::AUTO_WRITE
00202 ),
00203 mesh_
00204 ),
00205
00206 epsilon_
00207 (
00208 IOobject
00209 (
00210 "epsilon",
00211 runTime_.timeName(),
00212 mesh_,
00213 IOobject::MUST_READ,
00214 IOobject::AUTO_WRITE
00215 ),
00216 mesh_
00217 ),
00218
00219 y_(mesh_),
00220
00221 gradU_(fvc::grad(U)),
00222 eta_(k_/epsilon_*sqrt(2.0*magSqr(0.5*(gradU_ + gradU_.T())))),
00223 ksi_(k_/epsilon_*sqrt(2.0*magSqr(0.5*(gradU_ - gradU_.T())))),
00224 Cmu_(2.0/(3.0*(A1_ + eta_ + alphaKsi_*ksi_))),
00225 fEta_(A2_ + pow(eta_, 3.0)),
00226
00227 C5viscosity_
00228 (
00229 -2.0*pow(Cmu_, 3.0)*pow(k_, 4.0)/pow(epsilon_, 3.0)
00230 *(magSqr(gradU_ + gradU_.T()) - magSqr(gradU_ - gradU_.T()))
00231 ),
00232
00233 yStar_(sqrt(k_)*y_/nu() + SMALL),
00234
00235 nut_
00236 (
00237 IOobject
00238 (
00239 "nut",
00240 runTime_.timeName(),
00241 mesh_,
00242 IOobject::NO_READ,
00243 IOobject::AUTO_WRITE
00244 ),
00245 autoCreateLowReNut("nut", mesh_)
00246 ),
00247
00248 nonlinearStress_
00249 (
00250 "nonlinearStress",
00251 symm
00252 (
00253
00254 pow(k_, 3.0)/sqr(epsilon_)
00255 *(
00256 Ctau1_/fEta_
00257 *(
00258 (gradU_ & gradU_)
00259 + (gradU_ & gradU_)().T()
00260 )
00261 + Ctau2_/fEta_*(gradU_ & gradU_.T())
00262 + Ctau3_/fEta_*(gradU_.T() & gradU_)
00263 )
00264
00265 - 20.0*pow(k_, 4.0)/pow(epsilon_, 3.0)
00266 *pow(Cmu_, 3.0)
00267 *(
00268 ((gradU_ & gradU_) & gradU_.T())
00269 + ((gradU_ & gradU_.T()) & gradU_.T())
00270 - ((gradU_.T() & gradU_) & gradU_)
00271 - ((gradU_.T() & gradU_.T()) & gradU_)
00272 )
00273
00274 + min
00275 (
00276 C5viscosity_,
00277 dimensionedScalar("0", C5viscosity_.dimensions(), 0.0)
00278 )*gradU_
00279 )
00280 )
00281 {
00282 nut_ = Cmu_
00283 *(
00284 scalar(1) - exp(-Am_*yStar_))
00285 /(scalar(1) - exp(-Aepsilon_*yStar_) + SMALL
00286 )
00287 *sqr(k_)/(epsilon_ + epsilonSmall_)
00288
00289 + max
00290 (
00291 C5viscosity_,
00292 dimensionedScalar("0", C5viscosity_.dimensions(), 0.0)
00293 );
00294
00295 nut_.correctBoundaryConditions();
00296
00297 printCoeffs();
00298 }
00299
00300
00301
00302
00303 tmp<volSymmTensorField> LienCubicKELowRe::R() const
00304 {
00305 return tmp<volSymmTensorField>
00306 (
00307 new volSymmTensorField
00308 (
00309 IOobject
00310 (
00311 "R",
00312 runTime_.timeName(),
00313 mesh_,
00314 IOobject::NO_READ,
00315 IOobject::NO_WRITE
00316 ),
00317 ((2.0/3.0)*I)*k_ - nut_*twoSymm(gradU_) + nonlinearStress_,
00318 k_.boundaryField().types()
00319 )
00320 );
00321 }
00322
00323
00324 tmp<volSymmTensorField> LienCubicKELowRe::devReff() const
00325 {
00326 return tmp<volSymmTensorField>
00327 (
00328 new volSymmTensorField
00329 (
00330 IOobject
00331 (
00332 "devRhoReff",
00333 runTime_.timeName(),
00334 mesh_,
00335 IOobject::NO_READ,
00336 IOobject::NO_WRITE
00337 ),
00338 -nuEff()*dev(twoSymm(fvc::grad(U_))) + nonlinearStress_
00339 )
00340 );
00341 }
00342
00343
00344 tmp<fvVectorMatrix> LienCubicKELowRe::divDevReff(volVectorField& U) const
00345 {
00346 return
00347 (
00348 fvc::div(nonlinearStress_)
00349 - fvm::laplacian(nuEff(), U)
00350 - fvc::div(nuEff()*dev(fvc::grad(U)().T()))
00351 );
00352 }
00353
00354
00355 bool LienCubicKELowRe::read()
00356 {
00357 if (RASModel::read())
00358 {
00359 C1_.readIfPresent(coeffDict());
00360 C2_.readIfPresent(coeffDict());
00361 sigmak_.readIfPresent(coeffDict());
00362 sigmaEps_.readIfPresent(coeffDict());
00363 A1_.readIfPresent(coeffDict());
00364 A2_.readIfPresent(coeffDict());
00365 Ctau1_.readIfPresent(coeffDict());
00366 Ctau2_.readIfPresent(coeffDict());
00367 Ctau3_.readIfPresent(coeffDict());
00368 alphaKsi_.readIfPresent(coeffDict());
00369 CmuWall_.readIfPresent(coeffDict());
00370 kappa_.readIfPresent(coeffDict());
00371 Am_.readIfPresent(coeffDict());
00372 Aepsilon_.readIfPresent(coeffDict());
00373 Amu_.readIfPresent(coeffDict());
00374
00375 return true;
00376 }
00377 else
00378 {
00379 return false;
00380 }
00381 }
00382
00383
00384 void LienCubicKELowRe::correct()
00385 {
00386 RASModel::correct();
00387
00388 if (!turbulence_)
00389 {
00390 return;
00391 }
00392
00393 if (mesh_.changing())
00394 {
00395 y_.correct();
00396 }
00397
00398 gradU_ = fvc::grad(U_);
00399
00400
00401 volScalarField S2 = symm(gradU_) && gradU_;
00402
00403 yStar_ = sqrt(k_)*y_/nu() + SMALL;
00404 volScalarField Rt = sqr(k_)/(nu()*epsilon_);
00405
00406 volScalarField fMu =
00407 (scalar(1) - exp(-Am_*yStar_))
00408 /(scalar(1) - exp(-Aepsilon_*yStar_) + SMALL);
00409
00410 volScalarField f2 = scalar(1) - 0.3*exp(-sqr(Rt));
00411
00412 volScalarField G =
00413 Cmu_*fMu*sqr(k_)/epsilon_*S2 - (nonlinearStress_ && gradU_);
00414
00415
00416 tmp<fvScalarMatrix> epsEqn
00417 (
00418 fvm::ddt(epsilon_)
00419 + fvm::div(phi_, epsilon_)
00420 - fvm::laplacian(DepsilonEff(), epsilon_)
00421 ==
00422 C1_*G*epsilon_/k_
00423
00424 + C2_*f2*pow(Cmu_, 0.75)*pow(k_, scalar(0.5))
00425 /(kappa_*y_*(scalar(1) - exp(-Aepsilon_*yStar_)))
00426 *exp(-Amu_*sqr(yStar_))*epsilon_
00427 - fvm::Sp(C2_*f2*epsilon_/k_, epsilon_)
00428 );
00429
00430 epsEqn().relax();
00431
00432 # include "LienCubicKELowReSetWallDissipation.H"
00433 # include <incompressibleRASModels/wallDissipationI.H>
00434
00435 solve(epsEqn);
00436 bound(epsilon_, epsilon0_);
00437
00438
00439
00440
00441 tmp<fvScalarMatrix> kEqn
00442 (
00443 fvm::ddt(k_)
00444 + fvm::div(phi_, k_)
00445 - fvm::laplacian(DkEff(), k_)
00446 ==
00447 G
00448 - fvm::Sp(epsilon_/k_, k_)
00449 );
00450
00451 kEqn().relax();
00452 solve(kEqn);
00453 bound(k_, k0_);
00454
00455
00456
00457
00458 eta_ = k_/epsilon_*sqrt(2.0*magSqr(0.5*(gradU_ + gradU_.T())));
00459 ksi_ = k_/epsilon_*sqrt(2.0*magSqr(0.5*(gradU_ - gradU_.T())));
00460 Cmu_ = 2.0/(3.0*(A1_ + eta_ + alphaKsi_*ksi_));
00461 fEta_ = A2_ + pow(eta_, 3.0);
00462
00463 C5viscosity_ =
00464 - 2.0*pow(Cmu_, 3.0)*pow(k_, 4.0)/pow(epsilon_, 3.0)
00465 *(magSqr(gradU_ + gradU_.T()) - magSqr(gradU_ - gradU_.T()));
00466
00467 nut_ =
00468 Cmu_*fMu*sqr(k_)/epsilon_
00469
00470 + max
00471 (
00472 C5viscosity_,
00473 dimensionedScalar("0", C5viscosity_.dimensions(), 0.0)
00474 );
00475
00476 nonlinearStress_ = symm
00477 (
00478
00479 pow(k_, 3.0)/sqr(epsilon_)
00480 *(
00481 Ctau1_/fEta_
00482 *(
00483 (gradU_ & gradU_)
00484 + (gradU_ & gradU_)().T()
00485 )
00486 + Ctau2_/fEta_*(gradU_ & gradU_.T())
00487 + Ctau3_/fEta_*(gradU_.T() & gradU_)
00488 )
00489
00490 - 20.0*pow(k_, 4.0)/pow(epsilon_, 3.0)
00491 *pow(Cmu_, 3.0)
00492 *(
00493 ((gradU_ & gradU_) & gradU_.T())
00494 + ((gradU_ & gradU_.T()) & gradU_.T())
00495 - ((gradU_.T() & gradU_) & gradU_)
00496 - ((gradU_.T() & gradU_.T()) & gradU_)
00497 )
00498
00499 + min
00500 (
00501 C5viscosity_,
00502 dimensionedScalar("0", C5viscosity_.dimensions(), 0.0)
00503 )*gradU_
00504 );
00505 }
00506
00507
00508
00509
00510 }
00511 }
00512 }
00513
00514