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 "kOmegaSSTSAS.H"
00027 #include <OpenFOAM/addToRunTimeSelectionTable.H>
00028 #include <finiteVolume/wallDist.H>
00029
00030
00031
00032 namespace Foam
00033 {
00034 namespace incompressible
00035 {
00036 namespace LESModels
00037 {
00038
00039
00040
00041 defineTypeNameAndDebug(kOmegaSSTSAS, 0);
00042 addToRunTimeSelectionTable(LESModel, kOmegaSSTSAS, dictionary);
00043
00044
00045
00046 void kOmegaSSTSAS::updateSubGridScaleFields(const volScalarField& S2)
00047 {
00048 nuSgs_ == a1_*k_/max(a1_*omega_, F2()*sqrt(S2));
00049 nuSgs_.correctBoundaryConditions();
00050 }
00051
00052
00053 tmp<volScalarField> kOmegaSSTSAS::F1(const volScalarField& CDkOmega) const
00054 {
00055 volScalarField CDkOmegaPlus = max
00056 (
00057 CDkOmega,
00058 dimensionedScalar("1.0e-10", dimless/sqr(dimTime), 1.0e-10)
00059 );
00060
00061 volScalarField arg1 = min
00062 (
00063 min
00064 (
00065 max
00066 (
00067 (scalar(1)/betaStar_)*sqrt(k_)/(omega_*y_),
00068 scalar(500)*nu()/(sqr(y_)*omega_)
00069 ),
00070 (4*alphaOmega2_)*k_/(CDkOmegaPlus*sqr(y_))
00071 ),
00072 scalar(10)
00073 );
00074
00075 return tanh(pow4(arg1));
00076 }
00077
00078
00079 tmp<volScalarField> kOmegaSSTSAS::F2() const
00080 {
00081 volScalarField arg2 = min
00082 (
00083 max
00084 (
00085 (scalar(2)/betaStar_)*sqrt(k_)/(omega_*y_),
00086 scalar(500)*nu()/(sqr(y_)*omega_)
00087 ),
00088 scalar(100)
00089 );
00090
00091 return tanh(sqr(arg2));
00092 }
00093
00094
00095 tmp<volScalarField> kOmegaSSTSAS::Lvk2
00096 (
00097 const volScalarField& S2
00098 ) const
00099 {
00100 return max
00101 (
00102 kappa_*sqrt(S2)
00103 /(
00104 mag(fvc::laplacian(U()))
00105 + dimensionedScalar
00106 (
00107 "ROOTVSMALL",
00108 dimensionSet(0, -1 , -1, 0, 0, 0, 0),
00109 ROOTVSMALL
00110 )
00111 ),
00112 Cs_*delta()
00113 );
00114 }
00115
00116
00117
00118
00119 kOmegaSSTSAS::kOmegaSSTSAS
00120 (
00121 const volVectorField& U,
00122 const surfaceScalarField& phi,
00123 transportModel& transport,
00124 const word& modelName
00125 )
00126 :
00127 LESModel(modelName, U, phi, transport),
00128
00129 alphaK1_
00130 (
00131 dimensioned<scalar>::lookupOrAddToDict
00132 (
00133 "alphaK1",
00134 coeffDict_,
00135 0.85034
00136 )
00137 ),
00138 alphaK2_
00139 (
00140 dimensioned<scalar>::lookupOrAddToDict
00141 (
00142 "alphaK2",
00143 coeffDict_,
00144 1.0
00145 )
00146 ),
00147 alphaOmega1_
00148 (
00149 dimensioned<scalar>::lookupOrAddToDict
00150 (
00151 "alphaOmega1",
00152 coeffDict_,
00153 0.5
00154 )
00155 ),
00156 alphaOmega2_
00157 (
00158 dimensioned<scalar>::lookupOrAddToDict
00159 (
00160 "alphaOmega2",
00161 coeffDict_,
00162 0.85616
00163 )
00164 ),
00165 gamma1_
00166 (
00167 dimensioned<scalar>::lookupOrAddToDict
00168 (
00169 "gamma1",
00170 coeffDict_,
00171 0.5532
00172 )
00173 ),
00174 gamma2_
00175 (
00176 dimensioned<scalar>::lookupOrAddToDict
00177 (
00178 "gamma2",
00179 coeffDict_,
00180 0.4403
00181 )
00182 ),
00183 beta1_
00184 (
00185 dimensioned<scalar>::lookupOrAddToDict
00186 (
00187 "beta1",
00188 coeffDict_,
00189 0.075
00190 )
00191 ),
00192 beta2_
00193 (
00194 dimensioned<scalar>::lookupOrAddToDict
00195 (
00196 "beta2",
00197 coeffDict_,
00198 0.0828
00199 )
00200 ),
00201 betaStar_
00202 (
00203 dimensioned<scalar>::lookupOrAddToDict
00204 (
00205 "betaStar",
00206 coeffDict_,
00207 0.09
00208 )
00209 ),
00210 a1_
00211 (
00212 dimensioned<scalar>::lookupOrAddToDict
00213 (
00214 "a1",
00215 coeffDict_,
00216 0.31
00217 )
00218 ),
00219 c1_
00220 (
00221 dimensioned<scalar>::lookupOrAddToDict
00222 (
00223 "c1",
00224 coeffDict_,
00225 10.0
00226 )
00227 ),
00228 Cs_
00229 (
00230 dimensioned<scalar>::lookupOrAddToDict
00231 (
00232 "Cs",
00233 coeffDict_,
00234 0.262
00235 )
00236 ),
00237 alphaPhi_
00238 (
00239 dimensioned<scalar>::lookupOrAddToDict
00240 (
00241 "alphaPhi",
00242 coeffDict_,
00243 0.666667
00244 )
00245 ),
00246 zetaTilda2_
00247 (
00248 dimensioned<scalar>::lookupOrAddToDict
00249 (
00250 "zetaTilda2",
00251 coeffDict_,
00252 1.755
00253 )
00254 ),
00255 FSAS_
00256 (
00257 dimensioned<scalar>::lookupOrAddToDict
00258 (
00259 "FSAS",
00260 coeffDict_,
00261 1.25
00262 )
00263 ),
00264
00265 omega0_("omega0", dimless/dimTime, SMALL),
00266 omegaSmall_("omegaSmall", dimless/dimTime, SMALL),
00267 y_(mesh_),
00268 Cmu_
00269 (
00270 dimensioned<scalar>::lookupOrAddToDict
00271 (
00272 "Cmu",
00273 coeffDict_,
00274 0.09
00275 )
00276 ),
00277 kappa_
00278 (
00279 dimensioned<scalar>::lookupOrAddToDict
00280 (
00281 "kappa",
00282 *this,
00283 0.41
00284 )
00285 ),
00286
00287 k_
00288 (
00289 IOobject
00290 (
00291 "k",
00292 runTime_.timeName(),
00293 mesh_,
00294 IOobject::MUST_READ,
00295 IOobject::AUTO_WRITE
00296 ),
00297 mesh_
00298 ),
00299
00300 omega_
00301 (
00302 IOobject
00303 (
00304 "omega",
00305 runTime_.timeName(),
00306 mesh_,
00307 IOobject::MUST_READ,
00308 IOobject::AUTO_WRITE
00309 ),
00310 mesh_
00311 ),
00312
00313 nuSgs_
00314 (
00315 IOobject
00316 (
00317 "nuSgs",
00318 runTime_.timeName(),
00319 mesh_,
00320 IOobject::MUST_READ,
00321 IOobject::AUTO_WRITE
00322 ),
00323 mesh_
00324 )
00325 {
00326 updateSubGridScaleFields(magSqr(2.0*symm(fvc::grad(U))));
00327
00328 printCoeffs();
00329 }
00330
00331
00332
00333
00334 void kOmegaSSTSAS::correct(const tmp<volTensorField>& gradU)
00335 {
00336 LESModel::correct(gradU);
00337
00338 if (mesh_.changing())
00339 {
00340 y_.correct();
00341 }
00342
00343 volScalarField S2 = magSqr(2.0*symm(gradU()));
00344 gradU.clear();
00345
00346 volVectorField gradK = fvc::grad(k_);
00347 volVectorField gradOmega = fvc::grad(omega_);
00348 volScalarField L = sqrt(k_)/(pow(Cmu_, 0.25)*(omega_ + omegaSmall_));
00349 volScalarField CDkOmega =
00350 (2.0*alphaOmega2_)*(gradK & gradOmega)/(omega_ + omegaSmall_);
00351 volScalarField F1 = this->F1(CDkOmega);
00352 volScalarField G = nuSgs_*0.5*S2;
00353
00354
00355 {
00356 fvScalarMatrix kEqn
00357 (
00358 fvm::ddt(k_)
00359 + fvm::div(phi(), k_)
00360 - fvm::Sp(fvc::div(phi()), k_)
00361 - fvm::laplacian(DkEff(F1), k_)
00362 ==
00363 min(G, c1_*betaStar_*k_*omega_)
00364 - fvm::Sp(betaStar_*omega_, k_)
00365 );
00366
00367 kEqn.relax();
00368 kEqn.solve();
00369 }
00370 bound(k_, k0());
00371
00372 volScalarField grad_omega_k = max
00373 (
00374 magSqr(gradOmega)/
00375 sqr(omega_ + omegaSmall_),
00376 magSqr(gradK)/
00377 sqr(k_ + k0())
00378 );
00379
00380
00381 {
00382 fvScalarMatrix omegaEqn
00383 (
00384 fvm::ddt(omega_)
00385 + fvm::div(phi(), omega_)
00386 - fvm::Sp(fvc::div(phi()), omega_)
00387 - fvm::laplacian(DomegaEff(F1), omega_)
00388 ==
00389 gamma(F1)*0.5*S2
00390 - fvm::Sp(beta(F1)*omega_, omega_)
00391 - fvm::SuSp
00392 (
00393 (F1 - scalar(1))*CDkOmega/omega_,
00394 omega_
00395 )
00396 + FSAS_
00397 *max
00398 (
00399 dimensionedScalar("zero",dimensionSet(0, 0 , -2, 0, 0),0. ),
00400 zetaTilda2_*kappa_*S2*sqr(L/Lvk2(S2))
00401 - 2.0/alphaPhi_*k_*grad_omega_k
00402 )
00403 );
00404
00405 omegaEqn.relax();
00406 omegaEqn.solve();
00407 }
00408 bound(omega_, omega0_);
00409
00410 updateSubGridScaleFields(S2);
00411 }
00412
00413
00414 tmp<volScalarField> kOmegaSSTSAS::epsilon() const
00415 {
00416 return 2.0*nuEff()*magSqr(symm(fvc::grad(U())));
00417 }
00418
00419
00420 tmp<volSymmTensorField> kOmegaSSTSAS::B() const
00421 {
00422 return ((2.0/3.0)*I)*k() - nuSgs()*twoSymm(fvc::grad(U()));
00423 }
00424
00425
00426 tmp<volSymmTensorField> kOmegaSSTSAS::devBeff() const
00427 {
00428 return -nuEff()*dev(twoSymm(fvc::grad(U())));
00429 }
00430
00431
00432 tmp<fvVectorMatrix> kOmegaSSTSAS::divDevBeff(volVectorField& U) const
00433 {
00434 return
00435 (
00436 - fvm::laplacian(nuEff(), U) - fvc::div(nuEff()*dev(fvc::grad(U)().T()))
00437 );
00438 }
00439
00440
00441 bool kOmegaSSTSAS::read()
00442 {
00443 if (LESModel::read())
00444 {
00445 alphaK1_.readIfPresent(coeffDict());
00446 alphaK2_.readIfPresent(coeffDict());
00447 alphaOmega1_.readIfPresent(coeffDict());
00448 alphaOmega2_.readIfPresent(coeffDict());
00449 gamma1_.readIfPresent(coeffDict());
00450 gamma2_.readIfPresent(coeffDict());
00451 beta1_.readIfPresent(coeffDict());
00452 beta2_.readIfPresent(coeffDict());
00453 betaStar_.readIfPresent(coeffDict());
00454 a1_.readIfPresent(coeffDict());
00455 c1_.readIfPresent(coeffDict());
00456 Cs_.readIfPresent(coeffDict());
00457 alphaPhi_.readIfPresent(coeffDict());
00458 zetaTilda2_.readIfPresent(coeffDict());
00459 FSAS_.readIfPresent(coeffDict());
00460
00461 return true;
00462 }
00463 else
00464 {
00465 return false;
00466 }
00467 }
00468
00469
00470
00471
00472 }
00473 }
00474 }
00475
00476