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 "kOmegaSST.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(kOmegaSST, 0);
00043 addToRunTimeSelectionTable(RASModel, kOmegaSST, dictionary);
00044
00045
00046
00047 tmp<volScalarField> kOmegaSST::F1(const volScalarField& CDkOmega) const
00048 {
00049 volScalarField CDkOmegaPlus = max
00050 (
00051 CDkOmega,
00052 dimensionedScalar("1.0e-10", dimless/sqr(dimTime), 1.0e-10)
00053 );
00054
00055 volScalarField arg1 = min
00056 (
00057 min
00058 (
00059 max
00060 (
00061 (scalar(1)/betaStar_)*sqrt(k_)/(omega_*y_),
00062 scalar(500)*nu()/(sqr(y_)*omega_)
00063 ),
00064 (4*alphaOmega2_)*k_/(CDkOmegaPlus*sqr(y_))
00065 ),
00066 scalar(10)
00067 );
00068
00069 return tanh(pow4(arg1));
00070 }
00071
00072 tmp<volScalarField> kOmegaSST::F2() const
00073 {
00074 volScalarField arg2 = min
00075 (
00076 max
00077 (
00078 (scalar(2)/betaStar_)*sqrt(k_)/(omega_*y_),
00079 scalar(500)*nu()/(sqr(y_)*omega_)
00080 ),
00081 scalar(100)
00082 );
00083
00084 return tanh(sqr(arg2));
00085 }
00086
00087
00088
00089
00090 kOmegaSST::kOmegaSST
00091 (
00092 const volVectorField& U,
00093 const surfaceScalarField& phi,
00094 transportModel& lamTransportModel
00095 )
00096 :
00097 RASModel(typeName, U, phi, lamTransportModel),
00098
00099 alphaK1_
00100 (
00101 dimensioned<scalar>::lookupOrAddToDict
00102 (
00103 "alphaK1",
00104 coeffDict_,
00105 0.85034
00106 )
00107 ),
00108 alphaK2_
00109 (
00110 dimensioned<scalar>::lookupOrAddToDict
00111 (
00112 "alphaK2",
00113 coeffDict_,
00114 1.0
00115 )
00116 ),
00117 alphaOmega1_
00118 (
00119 dimensioned<scalar>::lookupOrAddToDict
00120 (
00121 "alphaOmega1",
00122 coeffDict_,
00123 0.5
00124 )
00125 ),
00126 alphaOmega2_
00127 (
00128 dimensioned<scalar>::lookupOrAddToDict
00129 (
00130 "alphaOmega2",
00131 coeffDict_,
00132 0.85616
00133 )
00134 ),
00135 gamma1_
00136 (
00137 dimensioned<scalar>::lookupOrAddToDict
00138 (
00139 "gamma1",
00140 coeffDict_,
00141 0.5532
00142 )
00143 ),
00144 gamma2_
00145 (
00146 dimensioned<scalar>::lookupOrAddToDict
00147 (
00148 "gamma2",
00149 coeffDict_,
00150 0.4403
00151 )
00152 ),
00153 beta1_
00154 (
00155 dimensioned<scalar>::lookupOrAddToDict
00156 (
00157 "beta1",
00158 coeffDict_,
00159 0.075
00160 )
00161 ),
00162 beta2_
00163 (
00164 dimensioned<scalar>::lookupOrAddToDict
00165 (
00166 "beta2",
00167 coeffDict_,
00168 0.0828
00169 )
00170 ),
00171 betaStar_
00172 (
00173 dimensioned<scalar>::lookupOrAddToDict
00174 (
00175 "betaStar",
00176 coeffDict_,
00177 0.09
00178 )
00179 ),
00180 a1_
00181 (
00182 dimensioned<scalar>::lookupOrAddToDict
00183 (
00184 "a1",
00185 coeffDict_,
00186 0.31
00187 )
00188 ),
00189 c1_
00190 (
00191 dimensioned<scalar>::lookupOrAddToDict
00192 (
00193 "c1",
00194 coeffDict_,
00195 10.0
00196 )
00197 ),
00198
00199 y_(mesh_),
00200
00201 k_
00202 (
00203 IOobject
00204 (
00205 "k",
00206 runTime_.timeName(),
00207 mesh_,
00208 IOobject::NO_READ,
00209 IOobject::AUTO_WRITE
00210 ),
00211 autoCreateK("k", mesh_)
00212 ),
00213 omega_
00214 (
00215 IOobject
00216 (
00217 "omega",
00218 runTime_.timeName(),
00219 mesh_,
00220 IOobject::NO_READ,
00221 IOobject::AUTO_WRITE
00222 ),
00223 autoCreateOmega("omega", mesh_)
00224 ),
00225 nut_
00226 (
00227 IOobject
00228 (
00229 "nut",
00230 runTime_.timeName(),
00231 mesh_,
00232 IOobject::NO_READ,
00233 IOobject::AUTO_WRITE
00234 ),
00235 autoCreateNut("nut", mesh_)
00236 )
00237 {
00238 bound(omega_, omega0_);
00239
00240 nut_ = a1_*k_/max(a1_*omega_, F2()*sqrt(2.0)*mag(symm(fvc::grad(U_))));
00241 nut_.correctBoundaryConditions();
00242
00243 printCoeffs();
00244 }
00245
00246
00247
00248
00249 tmp<volSymmTensorField> kOmegaSST::R() const
00250 {
00251 return tmp<volSymmTensorField>
00252 (
00253 new volSymmTensorField
00254 (
00255 IOobject
00256 (
00257 "R",
00258 runTime_.timeName(),
00259 mesh_,
00260 IOobject::NO_READ,
00261 IOobject::NO_WRITE
00262 ),
00263 ((2.0/3.0)*I)*k_ - nut_*twoSymm(fvc::grad(U_)),
00264 k_.boundaryField().types()
00265 )
00266 );
00267 }
00268
00269
00270 tmp<volSymmTensorField> kOmegaSST::devReff() const
00271 {
00272 return tmp<volSymmTensorField>
00273 (
00274 new volSymmTensorField
00275 (
00276 IOobject
00277 (
00278 "devRhoReff",
00279 runTime_.timeName(),
00280 mesh_,
00281 IOobject::NO_READ,
00282 IOobject::NO_WRITE
00283 ),
00284 -nuEff()*dev(twoSymm(fvc::grad(U_)))
00285 )
00286 );
00287 }
00288
00289
00290 tmp<fvVectorMatrix> kOmegaSST::divDevReff(volVectorField& U) const
00291 {
00292 return
00293 (
00294 - fvm::laplacian(nuEff(), U)
00295 - fvc::div(nuEff()*dev(fvc::grad(U)().T()))
00296 );
00297 }
00298
00299
00300 bool kOmegaSST::read()
00301 {
00302 if (RASModel::read())
00303 {
00304 alphaK1_.readIfPresent(coeffDict());
00305 alphaK2_.readIfPresent(coeffDict());
00306 alphaOmega1_.readIfPresent(coeffDict());
00307 alphaOmega2_.readIfPresent(coeffDict());
00308 gamma1_.readIfPresent(coeffDict());
00309 gamma2_.readIfPresent(coeffDict());
00310 beta1_.readIfPresent(coeffDict());
00311 beta2_.readIfPresent(coeffDict());
00312 betaStar_.readIfPresent(coeffDict());
00313 a1_.readIfPresent(coeffDict());
00314 c1_.readIfPresent(coeffDict());
00315
00316 return true;
00317 }
00318 else
00319 {
00320 return false;
00321 }
00322 }
00323
00324
00325 void kOmegaSST::correct()
00326 {
00327 RASModel::correct();
00328
00329 if (!turbulence_)
00330 {
00331 return;
00332 }
00333
00334 if (mesh_.changing())
00335 {
00336 y_.correct();
00337 }
00338
00339 volScalarField S2 = magSqr(symm(fvc::grad(U_)));
00340 volScalarField G("RASModel::G", nut_*2*S2);
00341
00342
00343 omega_.boundaryField().updateCoeffs();
00344
00345 volScalarField CDkOmega =
00346 (2*alphaOmega2_)*(fvc::grad(k_) & fvc::grad(omega_))/omega_;
00347
00348 volScalarField F1 = this->F1(CDkOmega);
00349
00350
00351 tmp<fvScalarMatrix> omegaEqn
00352 (
00353 fvm::ddt(omega_)
00354 + fvm::div(phi_, omega_)
00355 - fvm::Sp(fvc::div(phi_), omega_)
00356 - fvm::laplacian(DomegaEff(F1), omega_)
00357 ==
00358 gamma(F1)*2*S2
00359 - fvm::Sp(beta(F1)*omega_, omega_)
00360 - fvm::SuSp
00361 (
00362 (F1 - scalar(1))*CDkOmega/omega_,
00363 omega_
00364 )
00365 );
00366
00367 omegaEqn().relax();
00368
00369 omegaEqn().boundaryManipulate(omega_.boundaryField());
00370
00371 solve(omegaEqn);
00372 bound(omega_, omega0_);
00373
00374
00375 tmp<fvScalarMatrix> kEqn
00376 (
00377 fvm::ddt(k_)
00378 + fvm::div(phi_, k_)
00379 - fvm::Sp(fvc::div(phi_), k_)
00380 - fvm::laplacian(DkEff(F1), k_)
00381 ==
00382 min(G, c1_*betaStar_*k_*omega_)
00383 - fvm::Sp(betaStar_*omega_, k_)
00384 );
00385
00386 kEqn().relax();
00387 solve(kEqn);
00388 bound(k_, k0_);
00389
00390
00391
00392 nut_ = a1_*k_/max(a1_*omega_, F2()*sqrt(2*S2));
00393 nut_.correctBoundaryConditions();
00394 }
00395
00396
00397
00398
00399 }
00400 }
00401 }
00402
00403