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 "LRR.H"
00027 #include <OpenFOAM/addToRunTimeSelectionTable.H>
00028 #include <finiteVolume/wallFvPatch.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(LRR, 0);
00044 addToRunTimeSelectionTable(RASModel, LRR, dictionary);
00045
00046
00047
00048 LRR::LRR
00049 (
00050 const volVectorField& U,
00051 const surfaceScalarField& phi,
00052 transportModel& lamTransportModel
00053 )
00054 :
00055 RASModel(typeName, U, phi, lamTransportModel),
00056
00057 Cmu_
00058 (
00059 dimensioned<scalar>::lookupOrAddToDict
00060 (
00061 "Cmu",
00062 coeffDict_,
00063 0.09
00064 )
00065 ),
00066 Clrr1_
00067 (
00068 dimensioned<scalar>::lookupOrAddToDict
00069 (
00070 "Clrr1",
00071 coeffDict_,
00072 1.8
00073 )
00074 ),
00075 Clrr2_
00076 (
00077 dimensioned<scalar>::lookupOrAddToDict
00078 (
00079 "Clrr2",
00080 coeffDict_,
00081 0.6
00082 )
00083 ),
00084 C1_
00085 (
00086 dimensioned<scalar>::lookupOrAddToDict
00087 (
00088 "C1",
00089 coeffDict_,
00090 1.44
00091 )
00092 ),
00093 C2_
00094 (
00095 dimensioned<scalar>::lookupOrAddToDict
00096 (
00097 "C2",
00098 coeffDict_,
00099 1.92
00100 )
00101 ),
00102 Cs_
00103 (
00104 dimensioned<scalar>::lookupOrAddToDict
00105 (
00106 "Cs",
00107 coeffDict_,
00108 0.25
00109 )
00110 ),
00111 Ceps_
00112 (
00113 dimensioned<scalar>::lookupOrAddToDict
00114 (
00115 "Ceps",
00116 coeffDict_,
00117 0.15
00118 )
00119 ),
00120 sigmaEps_
00121 (
00122 dimensioned<scalar>::lookupOrAddToDict
00123 (
00124 "sigmaEps",
00125 coeffDict_,
00126 1.3
00127 )
00128 ),
00129 couplingFactor_
00130 (
00131 dimensioned<scalar>::lookupOrAddToDict
00132 (
00133 "couplingFactor",
00134 coeffDict_,
00135 0.0
00136 )
00137 ),
00138
00139 R_
00140 (
00141 IOobject
00142 (
00143 "R",
00144 runTime_.timeName(),
00145 mesh_,
00146 IOobject::NO_READ,
00147 IOobject::AUTO_WRITE
00148 ),
00149 autoCreateR("R", mesh_)
00150 ),
00151 k_
00152 (
00153 IOobject
00154 (
00155 "k",
00156 runTime_.timeName(),
00157 mesh_,
00158 IOobject::NO_READ,
00159 IOobject::AUTO_WRITE
00160 ),
00161 autoCreateK("k", mesh_)
00162 ),
00163 epsilon_
00164 (
00165 IOobject
00166 (
00167 "epsilon",
00168 runTime_.timeName(),
00169 mesh_,
00170 IOobject::NO_READ,
00171 IOobject::AUTO_WRITE
00172 ),
00173 autoCreateEpsilon("epsilon", mesh_)
00174 ),
00175 nut_
00176 (
00177 IOobject
00178 (
00179 "nut",
00180 runTime_.timeName(),
00181 mesh_,
00182 IOobject::NO_READ,
00183 IOobject::AUTO_WRITE
00184 ),
00185 autoCreateNut("nut", mesh_)
00186 )
00187 {
00188 nut_ = Cmu_*sqr(k_)/(epsilon_ + epsilonSmall_);
00189 nut_.correctBoundaryConditions();
00190
00191 if (couplingFactor_.value() < 0.0 || couplingFactor_.value() > 1.0)
00192 {
00193 FatalErrorIn
00194 (
00195 "LRR::LRR"
00196 "(const volVectorField& U, const surfaceScalarField& phi,"
00197 "transportModel& lamTransportModel)"
00198 ) << "couplingFactor = " << couplingFactor_
00199 << " is not in range 0 - 1" << nl
00200 << exit(FatalError);
00201 }
00202
00203 printCoeffs();
00204 }
00205
00206
00207
00208
00209 tmp<volSymmTensorField> LRR::devReff() const
00210 {
00211 return tmp<volSymmTensorField>
00212 (
00213 new volSymmTensorField
00214 (
00215 IOobject
00216 (
00217 "devRhoReff",
00218 runTime_.timeName(),
00219 mesh_,
00220 IOobject::NO_READ,
00221 IOobject::NO_WRITE
00222 ),
00223 R_ - nu()*dev(twoSymm(fvc::grad(U_)))
00224 )
00225 );
00226 }
00227
00228
00229 tmp<fvVectorMatrix> LRR::divDevReff(volVectorField& U) const
00230 {
00231 if (couplingFactor_.value() > 0.0)
00232 {
00233 return
00234 (
00235 fvc::div(R_ + couplingFactor_*nut_*fvc::grad(U), "div(R)")
00236 + fvc::laplacian
00237 (
00238 (1.0 - couplingFactor_)*nut_,
00239 U,
00240 "laplacian(nuEff,U)"
00241 )
00242 - fvm::laplacian(nuEff(), U)
00243 );
00244 }
00245 else
00246 {
00247 return
00248 (
00249 fvc::div(R_)
00250 + fvc::laplacian(nut_, U, "laplacian(nuEff,U)")
00251 - fvm::laplacian(nuEff(), U)
00252 );
00253 }
00254 }
00255
00256
00257 bool LRR::read()
00258 {
00259 if (RASModel::read())
00260 {
00261 Cmu_.readIfPresent(coeffDict());
00262 Clrr1_.readIfPresent(coeffDict());
00263 Clrr2_.readIfPresent(coeffDict());
00264 C1_.readIfPresent(coeffDict());
00265 C2_.readIfPresent(coeffDict());
00266 Cs_.readIfPresent(coeffDict());
00267 Ceps_.readIfPresent(coeffDict());
00268 sigmaEps_.readIfPresent(coeffDict());
00269
00270 couplingFactor_.readIfPresent(coeffDict());
00271
00272 if (couplingFactor_.value() < 0.0 || couplingFactor_.value() > 1.0)
00273 {
00274 FatalErrorIn("LRR::read()")
00275 << "couplingFactor = " << couplingFactor_
00276 << " is not in range 0 - 1"
00277 << exit(FatalError);
00278 }
00279
00280 return true;
00281 }
00282 else
00283 {
00284 return false;
00285 }
00286 }
00287
00288
00289 void LRR::correct()
00290 {
00291 RASModel::correct();
00292
00293 if (!turbulence_)
00294 {
00295 return;
00296 }
00297
00298 volSymmTensorField P = -twoSymm(R_ & fvc::grad(U_));
00299 volScalarField G("RASModel::G", 0.5*mag(tr(P)));
00300
00301
00302 epsilon_.boundaryField().updateCoeffs();
00303
00304
00305 tmp<fvScalarMatrix> epsEqn
00306 (
00307 fvm::ddt(epsilon_)
00308 + fvm::div(phi_, epsilon_)
00309 - fvm::Sp(fvc::div(phi_), epsilon_)
00310
00311 - fvm::laplacian(DepsilonEff(), epsilon_)
00312 ==
00313 C1_*G*epsilon_/k_
00314 - fvm::Sp(C2_*epsilon_/k_, epsilon_)
00315 );
00316
00317 epsEqn().relax();
00318
00319 epsEqn().boundaryManipulate(epsilon_.boundaryField());
00320
00321 solve(epsEqn);
00322 bound(epsilon_, epsilon0_);
00323
00324
00325
00326
00327 const fvPatchList& patches = mesh_.boundary();
00328
00329 forAll(patches, patchi)
00330 {
00331 const fvPatch& curPatch = patches[patchi];
00332
00333 if (isA<wallFvPatch>(curPatch))
00334 {
00335 forAll(curPatch, facei)
00336 {
00337 label faceCelli = curPatch.faceCells()[facei];
00338 P[faceCelli]
00339 *= min(G[faceCelli]/(0.5*mag(tr(P[faceCelli])) + SMALL), 1.0);
00340 }
00341 }
00342 }
00343
00344
00345 tmp<fvSymmTensorMatrix> REqn
00346 (
00347 fvm::ddt(R_)
00348 + fvm::div(phi_, R_)
00349 - fvm::Sp(fvc::div(phi_), R_)
00350
00351 - fvm::laplacian(DREff(), R_)
00352 + fvm::Sp(Clrr1_*epsilon_/k_, R_)
00353 ==
00354 P
00355 - (2.0/3.0*(1 - Clrr1_)*I)*epsilon_
00356 - Clrr2_*dev(P)
00357 );
00358
00359 REqn().relax();
00360 solve(REqn);
00361
00362 R_.max
00363 (
00364 dimensionedSymmTensor
00365 (
00366 "zero",
00367 R_.dimensions(),
00368 symmTensor
00369 (
00370 k0_.value(), -GREAT, -GREAT,
00371 k0_.value(), -GREAT,
00372 k0_.value()
00373 )
00374 )
00375 );
00376
00377 k_ = 0.5*tr(R_);
00378 bound(k_, k0_);
00379
00380
00381
00382 nut_ = Cmu_*sqr(k_)/epsilon_;
00383 nut_.correctBoundaryConditions();
00384
00385
00386
00387
00388 forAll(patches, patchi)
00389 {
00390 const fvPatch& curPatch = patches[patchi];
00391
00392 if (isA<wallFvPatch>(curPatch))
00393 {
00394 symmTensorField& Rw = R_.boundaryField()[patchi];
00395
00396 const scalarField& nutw = nut_.boundaryField()[patchi];
00397
00398 vectorField snGradU = U_.boundaryField()[patchi].snGrad();
00399
00400 const vectorField& faceAreas
00401 = mesh_.Sf().boundaryField()[patchi];
00402
00403 const scalarField& magFaceAreas
00404 = mesh_.magSf().boundaryField()[patchi];
00405
00406 forAll(curPatch, facei)
00407 {
00408
00409 tensor gradUw
00410 = (faceAreas[facei]/magFaceAreas[facei])*snGradU[facei];
00411
00412
00413 tensor tauw = -nutw[facei]*2*symm(gradUw);
00414
00415
00416 Rw[facei].xy() = tauw.xy();
00417 Rw[facei].xz() = tauw.xz();
00418 Rw[facei].yz() = tauw.yz();
00419 }
00420 }
00421 }
00422 }
00423
00424
00425
00426
00427 }
00428 }
00429 }
00430
00431