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 "NonlinearKEShih.H"
00027 #include <OpenFOAM/addToRunTimeSelectionTable.H>
00028 #include <finiteVolume/wallFvPatch.H>
00029
00030
00031
00032 namespace Foam
00033 {
00034 namespace incompressible
00035 {
00036 namespace RASModels
00037 {
00038
00039
00040
00041 defineTypeNameAndDebug(NonlinearKEShih, 0);
00042 addToRunTimeSelectionTable(RASModel, NonlinearKEShih, dictionary);
00043
00044
00045
00046 NonlinearKEShih::NonlinearKEShih
00047 (
00048 const volVectorField& U,
00049 const surfaceScalarField& phi,
00050 transportModel& lamTransportModel
00051 )
00052 :
00053 RASModel(typeName, U, phi, lamTransportModel),
00054
00055 C1_
00056 (
00057 dimensioned<scalar>::lookupOrAddToDict
00058 (
00059 "C1",
00060 coeffDict_,
00061 1.44
00062 )
00063 ),
00064 C2_
00065 (
00066 dimensioned<scalar>::lookupOrAddToDict
00067 (
00068 "C2",
00069 coeffDict_,
00070 1.92
00071 )
00072 ),
00073 sigmak_
00074 (
00075 dimensioned<scalar>::lookupOrAddToDict
00076 (
00077 "sigmak",
00078 coeffDict_,
00079 1.0
00080 )
00081 ),
00082 sigmaEps_
00083 (
00084 dimensioned<scalar>::lookupOrAddToDict
00085 (
00086 "sigmaEps",
00087 coeffDict_,
00088 1.3
00089 )
00090 ),
00091 A1_
00092 (
00093 dimensioned<scalar>::lookupOrAddToDict
00094 (
00095 "A1",
00096 coeffDict_,
00097 1.25
00098 )
00099 ),
00100 A2_
00101 (
00102 dimensioned<scalar>::lookupOrAddToDict
00103 (
00104 "A2",
00105 coeffDict_,
00106 1000.0
00107 )
00108 ),
00109 Ctau1_
00110 (
00111 dimensioned<scalar>::lookupOrAddToDict
00112 (
00113 "Ctau1",
00114 coeffDict_,
00115 -4.0
00116 )
00117 ),
00118 Ctau2_
00119 (
00120 dimensioned<scalar>::lookupOrAddToDict
00121 (
00122 "Ctau2",
00123 coeffDict_,
00124 13.0
00125 )
00126 ),
00127 Ctau3_
00128 (
00129 dimensioned<scalar>::lookupOrAddToDict
00130 (
00131 "Ctau3",
00132 coeffDict_,
00133 -2.0
00134 )
00135 ),
00136 alphaKsi_
00137 (
00138 dimensioned<scalar>::lookupOrAddToDict
00139 (
00140 "alphaKsi",
00141 coeffDict_,
00142 0.9
00143 )
00144 ),
00145
00146 kappa_
00147 (
00148 dimensioned<scalar>::lookupOrAddToDict
00149 (
00150 "kappa_",
00151 coeffDict_,
00152 0.41
00153 )
00154 ),
00155 E_
00156 (
00157 dimensioned<scalar>::lookupOrAddToDict
00158 (
00159 "E",
00160 coeffDict_,
00161 9.8
00162 )
00163 ),
00164
00165 k_
00166 (
00167 IOobject
00168 (
00169 "k",
00170 runTime_.timeName(),
00171 mesh_,
00172 IOobject::MUST_READ,
00173 IOobject::AUTO_WRITE
00174 ),
00175 mesh_
00176 ),
00177
00178 epsilon_
00179 (
00180 IOobject
00181 (
00182 "epsilon",
00183 runTime_.timeName(),
00184 mesh_,
00185 IOobject::MUST_READ,
00186 IOobject::AUTO_WRITE
00187 ),
00188 mesh_
00189 ),
00190
00191 gradU_(fvc::grad(U)),
00192 eta_(k_/epsilon_*sqrt(2.0*magSqr(0.5*(gradU_ + gradU_.T())))),
00193 ksi_(k_/epsilon_*sqrt(2.0*magSqr(0.5*(gradU_ - gradU_.T())))),
00194 Cmu_(2.0/(3.0*(A1_ + eta_ + alphaKsi_*ksi_))),
00195 fEta_(A2_ + pow(eta_, 3.0)),
00196
00197 nut_("nut", Cmu_*sqr(k_)/(epsilon_ + epsilonSmall_)),
00198
00199 nonlinearStress_
00200 (
00201 "nonlinearStress",
00202 symm
00203 (
00204 pow(k_, 3.0)/sqr(epsilon_)
00205 *(
00206 Ctau1_/fEta_
00207 *(
00208 (gradU_ & gradU_)
00209 + (gradU_ & gradU_)().T()
00210 )
00211 + Ctau2_/fEta_*(gradU_ & gradU_.T())
00212 + Ctau3_/fEta_*(gradU_.T() & gradU_)
00213 )
00214 )
00215 )
00216 {
00217 #include <incompressibleRASModels/wallNonlinearViscosityI.H>
00218
00219 printCoeffs();
00220 }
00221
00222
00223
00224
00225 tmp<volSymmTensorField> NonlinearKEShih::R() const
00226 {
00227 volTensorField gradU_ = fvc::grad(U_);
00228
00229 return tmp<volSymmTensorField>
00230 (
00231 new volSymmTensorField
00232 (
00233 IOobject
00234 (
00235 "R",
00236 runTime_.timeName(),
00237 mesh_,
00238 IOobject::NO_READ,
00239 IOobject::NO_WRITE
00240 ),
00241 ((2.0/3.0)*I)*k_ - nut_*twoSymm(gradU_) + nonlinearStress_,
00242 k_.boundaryField().types()
00243 )
00244 );
00245 }
00246
00247
00248 tmp<volSymmTensorField> NonlinearKEShih::devReff() const
00249 {
00250 return tmp<volSymmTensorField>
00251 (
00252 new volSymmTensorField
00253 (
00254 IOobject
00255 (
00256 "devRhoReff",
00257 runTime_.timeName(),
00258 mesh_,
00259 IOobject::NO_READ,
00260 IOobject::NO_WRITE
00261 ),
00262 -nuEff()*dev(twoSymm(fvc::grad(U_))) + nonlinearStress_
00263 )
00264 );
00265 }
00266
00267
00268 tmp<fvVectorMatrix> NonlinearKEShih::divDevReff(volVectorField& U) const
00269 {
00270 return
00271 (
00272 fvc::div(nonlinearStress_)
00273 - fvm::laplacian(nuEff(), U)
00274 - fvc::div(nuEff()*dev(fvc::grad(U)().T()))
00275 );
00276 }
00277
00278
00279 bool NonlinearKEShih::read()
00280 {
00281 if (RASModel::read())
00282 {
00283 C1_.readIfPresent(coeffDict());
00284 C2_.readIfPresent(coeffDict());
00285 sigmak_.readIfPresent(coeffDict());
00286 sigmaEps_.readIfPresent(coeffDict());
00287 A1_.readIfPresent(coeffDict());
00288 A2_.readIfPresent(coeffDict());
00289 Ctau1_.readIfPresent(coeffDict());
00290 Ctau2_.readIfPresent(coeffDict());
00291 Ctau3_.readIfPresent(coeffDict());
00292 alphaKsi_.readIfPresent(coeffDict());
00293
00294 kappa_.readIfPresent(coeffDict());
00295 E_.readIfPresent(coeffDict());
00296
00297 return true;
00298 }
00299 else
00300 {
00301 return false;
00302 }
00303 }
00304
00305
00306 void NonlinearKEShih::correct()
00307 {
00308 RASModel::correct();
00309
00310 if (!turbulence_)
00311 {
00312 return;
00313 }
00314
00315 gradU_ = fvc::grad(U_);
00316
00317
00318 volScalarField S2 = symm(gradU_) && gradU_;
00319
00320 volScalarField G
00321 (
00322 "RASModel::G",
00323 Cmu_*sqr(k_)/epsilon_*S2
00324 - (nonlinearStress_ && gradU_)
00325 );
00326
00327 #include <incompressibleRASModels/nonLinearWallFunctionsI.H>
00328
00329
00330 tmp<fvScalarMatrix> epsEqn
00331 (
00332 fvm::ddt(epsilon_)
00333 + fvm::div(phi_, epsilon_)
00334 - fvm::laplacian(DepsilonEff(), epsilon_)
00335 ==
00336 C1_*G*epsilon_/k_
00337 - fvm::Sp(C2_*epsilon_/k_, epsilon_)
00338 );
00339
00340 epsEqn().relax();
00341
00342 #include <incompressibleRASModels/wallDissipationI.H>
00343
00344 solve(epsEqn);
00345 bound(epsilon_, epsilon0_);
00346
00347
00348
00349
00350 tmp<fvScalarMatrix> kEqn
00351 (
00352 fvm::ddt(k_)
00353 + fvm::div(phi_, k_)
00354 - fvm::laplacian(DkEff(), k_)
00355 ==
00356 G
00357 - fvm::Sp(epsilon_/k_, k_)
00358 );
00359
00360 kEqn().relax();
00361 solve(kEqn);
00362 bound(k_, k0_);
00363
00364
00365
00366
00367 eta_ = k_/epsilon_*sqrt(2.0*magSqr(0.5*(gradU_ + gradU_.T())));
00368 ksi_ = k_/epsilon_*sqrt(2.0*magSqr(0.5*(gradU_ - gradU_.T())));
00369 Cmu_ = 2.0/(3.0*(A1_ + eta_ + alphaKsi_*ksi_));
00370 fEta_ = A2_ + pow(eta_, 3.0);
00371
00372 nut_ = Cmu_*sqr(k_)/epsilon_;
00373
00374 #include <incompressibleRASModels/wallNonlinearViscosityI.H>
00375
00376 nonlinearStress_ = symm
00377 (
00378 pow(k_, 3.0)/sqr(epsilon_)
00379 *(
00380 Ctau1_/fEta_
00381 *(
00382 (gradU_ & gradU_)
00383 + (gradU_ & gradU_)().T()
00384 )
00385 + Ctau2_/fEta_*(gradU_ & gradU_.T())
00386 + Ctau3_/fEta_*(gradU_.T() & gradU_)
00387 )
00388 );
00389 }
00390
00391
00392
00393
00394 }
00395 }
00396 }
00397
00398