Go to the documentation of this file.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 "dynOneEqEddy.H"
00027 #include <OpenFOAM/addToRunTimeSelectionTable.H>
00028
00029
00030
00031 namespace Foam
00032 {
00033 namespace compressible
00034 {
00035 namespace LESModels
00036 {
00037
00038
00039
00040 defineTypeNameAndDebug(dynOneEqEddy, 0);
00041 addToRunTimeSelectionTable(LESModel, dynOneEqEddy, dictionary);
00042
00043
00044
00045 void dynOneEqEddy::updateSubGridScaleFields(const volSymmTensorField& D)
00046 {
00047 muSgs_ = ck_(D)*rho()*sqrt(k_)*delta();
00048 muSgs_.correctBoundaryConditions();
00049
00050 alphaSgs_ = muSgs_/Prt_;
00051 alphaSgs_.correctBoundaryConditions();
00052 }
00053
00054
00055 dimensionedScalar dynOneEqEddy::ck_(const volSymmTensorField& D) const
00056 {
00057 volScalarField KK = 0.5*(filter_(magSqr(U())) - magSqr(filter_(U())));
00058
00059 volSymmTensorField LL = dev(filter_(sqr(U())) - (sqr(filter_(U()))));
00060
00061 volSymmTensorField MM =
00062 delta()*(filter_(sqrt(k_)*D) - 2*sqrt(KK + filter_(k_))*filter_(D));
00063
00064 return average(LL && MM)/average(magSqr(MM));
00065 }
00066
00067
00068 dimensionedScalar dynOneEqEddy::ce_(const volSymmTensorField& D) const
00069 {
00070 volScalarField KK = 0.5*(filter_(magSqr(U())) - magSqr(filter_(U())));
00071
00072 volScalarField mm =
00073 pow(KK + filter_(k_), 1.5)/(2*delta()) - filter_(pow(k_, 1.5))/delta();
00074
00075 volScalarField ee =
00076 2*delta()*ck_(D)*
00077 (
00078 filter_(sqrt(k_)*magSqr(D))
00079 - 2*sqrt(KK + filter_(k_))*magSqr(filter_(D))
00080 );
00081
00082 return average(ee*mm)/average(mm*mm);
00083 }
00084
00085
00086
00087
00088 dynOneEqEddy::dynOneEqEddy
00089 (
00090 const volScalarField& rho,
00091 const volVectorField& U,
00092 const surfaceScalarField& phi,
00093 const basicThermo& thermoPhysicalModel
00094 )
00095 :
00096 LESModel(typeName, rho, U, phi, thermoPhysicalModel),
00097 GenEddyVisc(rho, U, phi, thermoPhysicalModel),
00098
00099 filterPtr_(LESfilter::New(U.mesh(), coeffDict())),
00100 filter_(filterPtr_())
00101 {
00102 updateSubGridScaleFields(dev(symm(fvc::grad(U))));
00103
00104 printCoeffs();
00105 }
00106
00107
00108
00109
00110 void dynOneEqEddy::correct(const tmp<volTensorField>& tgradU)
00111 {
00112 const volTensorField& gradU = tgradU();
00113
00114 GenEddyVisc::correct(gradU);
00115
00116 volSymmTensorField D = dev(symm(gradU));
00117 volScalarField divU = fvc::div(phi()/fvc::interpolate(rho()));
00118 volScalarField G = 2*muSgs_*(gradU && D);
00119
00120 solve
00121 (
00122 fvm::ddt(rho(), k_)
00123 + fvm::div(phi(), k_)
00124 - fvm::laplacian(DkEff(), k_)
00125 ==
00126 G
00127 - fvm::SuSp(2.0/3.0*rho()*divU, k_)
00128 - fvm::Sp(ce_(D)*rho()*sqrt(k_)/delta(), k_)
00129 );
00130
00131 bound(k_, dimensionedScalar("0", k_.dimensions(), 1.0e-10));
00132
00133 updateSubGridScaleFields(D);
00134 }
00135
00136
00137 bool dynOneEqEddy::read()
00138 {
00139 if (GenEddyVisc::read())
00140 {
00141 filter_.read(coeffDict());
00142
00143 return true;
00144 }
00145 else
00146 {
00147 return false;
00148 }
00149 }
00150
00151
00152
00153
00154 }
00155 }
00156 }
00157
00158