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 "kOmega.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(kOmega, 0);
00043 addToRunTimeSelectionTable(RASModel, kOmega, dictionary);
00044
00045
00046
00047 kOmega::kOmega
00048 (
00049 const volVectorField& U,
00050 const surfaceScalarField& phi,
00051 transportModel& lamTransportModel
00052 )
00053 :
00054 RASModel(typeName, U, phi, lamTransportModel),
00055
00056 Cmu_
00057 (
00058 dimensioned<scalar>::lookupOrAddToDict
00059 (
00060 "betaStar",
00061 coeffDict_,
00062 0.09
00063 )
00064 ),
00065 beta_
00066 (
00067 dimensioned<scalar>::lookupOrAddToDict
00068 (
00069 "beta",
00070 coeffDict_,
00071 0.072
00072 )
00073 ),
00074 alpha_
00075 (
00076 dimensioned<scalar>::lookupOrAddToDict
00077 (
00078 "alpha",
00079 coeffDict_,
00080 0.52
00081 )
00082 ),
00083 alphaK_
00084 (
00085 dimensioned<scalar>::lookupOrAddToDict
00086 (
00087 "alphaK",
00088 coeffDict_,
00089 0.5
00090 )
00091 ),
00092 alphaOmega_
00093 (
00094 dimensioned<scalar>::lookupOrAddToDict
00095 (
00096 "alphaOmega",
00097 coeffDict_,
00098 0.5
00099 )
00100 ),
00101
00102 k_
00103 (
00104 IOobject
00105 (
00106 "k",
00107 runTime_.timeName(),
00108 mesh_,
00109 IOobject::NO_READ,
00110 IOobject::AUTO_WRITE
00111 ),
00112 autoCreateK("k", mesh_)
00113 ),
00114 omega_
00115 (
00116 IOobject
00117 (
00118 "omega",
00119 runTime_.timeName(),
00120 mesh_,
00121 IOobject::NO_READ,
00122 IOobject::AUTO_WRITE
00123 ),
00124 autoCreateOmega("omega", mesh_)
00125 ),
00126 nut_
00127 (
00128 IOobject
00129 (
00130 "nut",
00131 runTime_.timeName(),
00132 mesh_,
00133 IOobject::NO_READ,
00134 IOobject::AUTO_WRITE
00135 ),
00136 autoCreateNut("nut", mesh_)
00137 )
00138 {
00139 nut_ = k_/(omega_ + omegaSmall_);
00140 nut_.correctBoundaryConditions();
00141
00142 printCoeffs();
00143 }
00144
00145
00146
00147
00148 tmp<volSymmTensorField> kOmega::R() const
00149 {
00150 return tmp<volSymmTensorField>
00151 (
00152 new volSymmTensorField
00153 (
00154 IOobject
00155 (
00156 "R",
00157 runTime_.timeName(),
00158 mesh_,
00159 IOobject::NO_READ,
00160 IOobject::NO_WRITE
00161 ),
00162 ((2.0/3.0)*I)*k_ - nut_*twoSymm(fvc::grad(U_)),
00163 k_.boundaryField().types()
00164 )
00165 );
00166 }
00167
00168
00169 tmp<volSymmTensorField> kOmega::devReff() const
00170 {
00171 return tmp<volSymmTensorField>
00172 (
00173 new volSymmTensorField
00174 (
00175 IOobject
00176 (
00177 "devRhoReff",
00178 runTime_.timeName(),
00179 mesh_,
00180 IOobject::NO_READ,
00181 IOobject::NO_WRITE
00182 ),
00183 -nuEff()*dev(twoSymm(fvc::grad(U_)))
00184 )
00185 );
00186 }
00187
00188
00189 tmp<fvVectorMatrix> kOmega::divDevReff(volVectorField& U) const
00190 {
00191 return
00192 (
00193 - fvm::laplacian(nuEff(), U)
00194 - fvc::div(nuEff()*dev(fvc::grad(U)().T()))
00195 );
00196 }
00197
00198
00199 bool kOmega::read()
00200 {
00201 if (RASModel::read())
00202 {
00203 Cmu_.readIfPresent(coeffDict());
00204 beta_.readIfPresent(coeffDict());
00205 alphaK_.readIfPresent(coeffDict());
00206 alphaOmega_.readIfPresent(coeffDict());
00207
00208 return true;
00209 }
00210 else
00211 {
00212 return false;
00213 }
00214 }
00215
00216
00217 void kOmega::correct()
00218 {
00219 RASModel::correct();
00220
00221 if (!turbulence_)
00222 {
00223 return;
00224 }
00225
00226 volScalarField G("RASModel::G", nut_*2*magSqr(symm(fvc::grad(U_))));
00227
00228
00229 omega_.boundaryField().updateCoeffs();
00230
00231
00232 tmp<fvScalarMatrix> omegaEqn
00233 (
00234 fvm::ddt(omega_)
00235 + fvm::div(phi_, omega_)
00236 - fvm::Sp(fvc::div(phi_), omega_)
00237 - fvm::laplacian(DomegaEff(), omega_)
00238 ==
00239 alpha_*G*omega_/k_
00240 - fvm::Sp(beta_*omega_, omega_)
00241 );
00242
00243 omegaEqn().relax();
00244
00245 omegaEqn().boundaryManipulate(omega_.boundaryField());
00246
00247 solve(omegaEqn);
00248 bound(omega_, omega0_);
00249
00250
00251
00252 tmp<fvScalarMatrix> kEqn
00253 (
00254 fvm::ddt(k_)
00255 + fvm::div(phi_, k_)
00256 - fvm::Sp(fvc::div(phi_), k_)
00257 - fvm::laplacian(DkEff(), k_)
00258 ==
00259 G
00260 - fvm::Sp(Cmu_*omega_, k_)
00261 );
00262
00263 kEqn().relax();
00264 solve(kEqn);
00265 bound(k_, k0_);
00266
00267
00268
00269 nut_ = k_/(omega_ + omegaSmall_);
00270 nut_.correctBoundaryConditions();
00271 }
00272
00273
00274
00275
00276 }
00277 }
00278 }
00279
00280