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 "SpalartAllmaras.H"
00027 #include <OpenFOAM/addToRunTimeSelectionTable.H>
00028
00029
00030
00031 namespace Foam
00032 {
00033 namespace incompressible
00034 {
00035 namespace LESModels
00036 {
00037
00038
00039
00040 defineTypeNameAndDebug(SpalartAllmaras, 0);
00041 addToRunTimeSelectionTable(LESModel, SpalartAllmaras, dictionary);
00042
00043
00044
00045 void SpalartAllmaras::updateSubGridScaleFields()
00046 {
00047 nuSgs_.internalField() = fv1()*nuTilda_.internalField();
00048 nuSgs_.correctBoundaryConditions();
00049 }
00050
00051
00052
00053
00054 tmp<volScalarField> SpalartAllmaras::fv1() const
00055 {
00056 volScalarField chi3 = pow3(nuTilda_/nu());
00057 return chi3/(chi3 + pow3(Cv1_));
00058 }
00059
00060
00061 tmp<volScalarField> SpalartAllmaras::fv2() const
00062 {
00063 return 1/pow3(scalar(1) + nuTilda_/(Cv2_*nu()));
00064 }
00065
00066
00067 tmp<volScalarField> SpalartAllmaras::fv3() const
00068 {
00069 volScalarField chi = nuTilda_/nu();
00070 volScalarField chiByCv2 = (1/Cv2_)*chi;
00071
00072 return
00073 (scalar(1) + chi*fv1())
00074 *(1/Cv2_)
00075 *(3*(scalar(1) + chiByCv2) + sqr(chiByCv2))
00076 /pow3(scalar(1) + chiByCv2);
00077 }
00078
00079
00080 tmp<volScalarField> SpalartAllmaras::S(const volTensorField& gradU) const
00081 {
00082 return sqrt(2.0)*mag(skew(gradU));
00083 }
00084
00085
00086 tmp<volScalarField> SpalartAllmaras::STilda
00087 (
00088 const volScalarField& S,
00089 const volScalarField& dTilda
00090 ) const
00091 {
00092 return fv3()*S + fv2()*nuTilda_/sqr(kappa_*dTilda);
00093 }
00094
00095
00096 tmp<volScalarField> SpalartAllmaras::r
00097 (
00098 const volScalarField& visc,
00099 const volScalarField& S,
00100 const volScalarField& dTilda
00101 ) const
00102 {
00103 return min
00104 (
00105 visc
00106 /(
00107 max
00108 (
00109 S,
00110 dimensionedScalar("SMALL", S.dimensions(), SMALL)
00111 )
00112 *sqr(kappa_*dTilda)
00113 + dimensionedScalar
00114 (
00115 "ROOTVSMALL",
00116 dimensionSet(0, 2 , -1, 0, 0),
00117 ROOTVSMALL
00118 )
00119 ),
00120 scalar(10)
00121 );
00122 }
00123
00124
00125 tmp<volScalarField> SpalartAllmaras::fw
00126 (
00127 const volScalarField& S,
00128 const volScalarField& dTilda
00129 ) const
00130 {
00131 volScalarField r = this->r(nuTilda_, S, dTilda);
00132
00133 volScalarField g = r + Cw2_*(pow6(r) - r);
00134
00135 return g*pow((1 + pow6(Cw3_))/(pow6(g) + pow6(Cw3_)), 1.0/6.0);
00136 }
00137
00138
00139 tmp<volScalarField> SpalartAllmaras::dTilda(const volScalarField&) const
00140 {
00141 return min(CDES_*delta(), y_);
00142 }
00143
00144
00145
00146
00147 SpalartAllmaras::SpalartAllmaras
00148 (
00149 const volVectorField& U,
00150 const surfaceScalarField& phi,
00151 transportModel& transport,
00152 const word& modelName
00153 )
00154 :
00155 LESModel(modelName, U, phi, transport),
00156
00157 sigmaNut_
00158 (
00159 dimensioned<scalar>::lookupOrAddToDict
00160 (
00161 "sigmaNut",
00162 coeffDict_,
00163 0.66666
00164 )
00165 ),
00166 kappa_
00167 (
00168 dimensioned<scalar>::lookupOrAddToDict
00169 (
00170 "kappa",
00171 coeffDict_,
00172 0.41
00173 )
00174 ),
00175 Cb1_
00176 (
00177 dimensioned<scalar>::lookupOrAddToDict
00178 (
00179 "Cb1",
00180 coeffDict_,
00181 0.1355
00182 )
00183 ),
00184 Cb2_
00185 (
00186 dimensioned<scalar>::lookupOrAddToDict
00187 (
00188 "Cb2",
00189 coeffDict_,
00190 0.622
00191 )
00192 ),
00193 Cv1_
00194 (
00195 dimensioned<scalar>::lookupOrAddToDict
00196 (
00197 "Cv1",
00198 coeffDict_,
00199 7.1
00200 )
00201 ),
00202 Cv2_
00203 (
00204 dimensioned<scalar>::lookupOrAddToDict
00205 (
00206 "Cv2",
00207 coeffDict_,
00208 5.0
00209 )
00210 ),
00211 CDES_
00212 (
00213 dimensioned<scalar>::lookupOrAddToDict
00214 (
00215 "CDES",
00216 coeffDict_,
00217 0.65
00218 )
00219 ),
00220 ck_
00221 (
00222 dimensioned<scalar>::lookupOrAddToDict
00223 (
00224 "ck",
00225 coeffDict_,
00226 0.07
00227 )
00228 ),
00229 Cw1_(Cb1_/sqr(kappa_) + (1.0 + Cb2_)/sigmaNut_),
00230 Cw2_
00231 (
00232 dimensioned<scalar>::lookupOrAddToDict
00233 (
00234 "Cw2",
00235 coeffDict_,
00236 0.3
00237 )
00238 ),
00239 Cw3_
00240 (
00241 dimensioned<scalar>::lookupOrAddToDict
00242 (
00243 "Cw3",
00244 coeffDict_,
00245 2.0
00246 )
00247 ),
00248
00249 y_(mesh_),
00250
00251 nuTilda_
00252 (
00253 IOobject
00254 (
00255 "nuTilda",
00256 runTime_.timeName(),
00257 mesh_,
00258 IOobject::MUST_READ,
00259 IOobject::AUTO_WRITE
00260 ),
00261 mesh_
00262 ),
00263
00264 nuSgs_
00265 (
00266 IOobject
00267 (
00268 "nuSgs",
00269 runTime_.timeName(),
00270 mesh_,
00271 IOobject::MUST_READ,
00272 IOobject::AUTO_WRITE
00273 ),
00274 mesh_
00275 )
00276 {
00277 updateSubGridScaleFields();
00278 }
00279
00280
00281
00282
00283 void SpalartAllmaras::correct(const tmp<volTensorField>& gradU)
00284 {
00285 LESModel::correct(gradU);
00286
00287 if (mesh_.changing())
00288 {
00289 y_.correct();
00290 y_.boundaryField() = max(y_.boundaryField(), VSMALL);
00291 }
00292
00293 const volScalarField S = this->S(gradU);
00294 const volScalarField dTilda = this->dTilda(S);
00295 const volScalarField STilda = this->STilda(S, dTilda);
00296
00297 fvScalarMatrix nuTildaEqn
00298 (
00299 fvm::ddt(nuTilda_)
00300 + fvm::div(phi(), nuTilda_)
00301 - fvm::laplacian
00302 (
00303 (nuTilda_ + nu())/sigmaNut_,
00304 nuTilda_,
00305 "laplacian(DnuTildaEff,nuTilda)"
00306 )
00307 - Cb2_/sigmaNut_*magSqr(fvc::grad(nuTilda_))
00308 ==
00309 Cb1_*STilda*nuTilda_
00310 - fvm::Sp(Cw1_*fw(STilda, dTilda)*nuTilda_/sqr(dTilda), nuTilda_)
00311 );
00312
00313 nuTildaEqn.relax();
00314 nuTildaEqn.solve();
00315
00316 bound(nuTilda_, dimensionedScalar("zero", nuTilda_.dimensions(), 0.0));
00317 nuTilda_.correctBoundaryConditions();
00318
00319 updateSubGridScaleFields();
00320 }
00321
00322
00323 tmp<volScalarField> SpalartAllmaras::k() const
00324 {
00325 return sqr(nuSgs()/ck_/dTilda(S(fvc::grad(U()))));
00326 }
00327
00328
00329 tmp<volScalarField> SpalartAllmaras::epsilon() const
00330 {
00331 return 2*nuEff()*magSqr(symm(fvc::grad(U())));
00332 }
00333
00334
00335 tmp<volSymmTensorField> SpalartAllmaras::B() const
00336 {
00337 return ((2.0/3.0)*I)*k() - nuSgs()*twoSymm(fvc::grad(U()));
00338 }
00339
00340
00341 tmp<volSymmTensorField> SpalartAllmaras::devBeff() const
00342 {
00343 return -nuEff()*dev(twoSymm(fvc::grad(U())));
00344 }
00345
00346
00347 tmp<fvVectorMatrix> SpalartAllmaras::divDevBeff(volVectorField& U) const
00348 {
00349 return
00350 (
00351 - fvm::laplacian(nuEff(), U) - fvc::div(nuEff()*dev(fvc::grad(U)().T()))
00352 );
00353 }
00354
00355
00356 bool SpalartAllmaras::read()
00357 {
00358 if (LESModel::read())
00359 {
00360 sigmaNut_.readIfPresent(coeffDict());
00361 kappa_.readIfPresent(*this);
00362 Cb1_.readIfPresent(coeffDict());
00363 Cb2_.readIfPresent(coeffDict());
00364 Cv1_.readIfPresent(coeffDict());
00365 Cv2_.readIfPresent(coeffDict());
00366 CDES_.readIfPresent(coeffDict());
00367 ck_.readIfPresent(coeffDict());
00368 Cw1_ = Cb1_/sqr(kappa_) + (1.0 + Cb2_)/sigmaNut_;
00369 Cw2_.readIfPresent(coeffDict());
00370 Cw3_.readIfPresent(coeffDict());
00371
00372 return true;
00373 }
00374 else
00375 {
00376 return false;
00377 }
00378 }
00379
00380
00381
00382
00383 }
00384 }
00385 }
00386
00387