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 "locDynOneEqEddy.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(locDynOneEqEddy, 0);
00041 addToRunTimeSelectionTable(LESModel, locDynOneEqEddy, dictionary);
00042
00043
00044
00045 void locDynOneEqEddy::updateSubGridScaleFields
00046 (
00047 const volSymmTensorField& D,
00048 const volScalarField& KK
00049 )
00050 {
00051 nuSgs_ = ck(D, KK)*sqrt(k_)*delta();
00052 nuSgs_.correctBoundaryConditions();
00053 }
00054
00055
00056 volScalarField locDynOneEqEddy::ck
00057 (
00058 const volSymmTensorField& D,
00059 const volScalarField& KK
00060 ) const
00061 {
00062 volSymmTensorField LL =
00063 simpleFilter_(dev(filter_(sqr(U())) - (sqr(filter_(U())))));
00064
00065 volSymmTensorField MM = simpleFilter_(-2.0*delta()*pow(KK, 0.5)*filter_(D));
00066
00067 volScalarField ck =
00068 simpleFilter_(0.5*(LL && MM))
00069 /(
00070 simpleFilter_(magSqr(MM))
00071 + dimensionedScalar("small", sqr(MM.dimensions()), VSMALL)
00072 );
00073
00074 return 0.5*(mag(ck) + ck);
00075 }
00076
00077
00078 volScalarField locDynOneEqEddy::ce
00079 (
00080 const volSymmTensorField& D,
00081 const volScalarField& KK
00082 ) const
00083 {
00084 volScalarField ce =
00085 simpleFilter_(nuEff()*(filter_(magSqr(D)) - magSqr(filter_(D))))
00086 /simpleFilter_(pow(KK, 1.5)/(2.0*delta()));
00087
00088 return 0.5*(mag(ce) + ce);
00089 }
00090
00091
00092
00093
00094 locDynOneEqEddy::locDynOneEqEddy
00095 (
00096 const volVectorField& U,
00097 const surfaceScalarField& phi,
00098 transportModel& transport
00099 )
00100 :
00101 LESModel(typeName, U, phi, transport),
00102 GenEddyVisc(U, phi, transport),
00103
00104 k_
00105 (
00106 IOobject
00107 (
00108 "k",
00109 runTime_.timeName(),
00110 mesh_,
00111 IOobject::MUST_READ,
00112 IOobject::AUTO_WRITE
00113 ),
00114 mesh_
00115 ),
00116
00117 simpleFilter_(U.mesh()),
00118 filterPtr_(LESfilter::New(U.mesh(), coeffDict())),
00119 filter_(filterPtr_())
00120 {
00121 volScalarField KK = 0.5*(filter_(magSqr(U)) - magSqr(filter_(U)));
00122 updateSubGridScaleFields(symm(fvc::grad(U)), KK);
00123
00124 printCoeffs();
00125 }
00126
00127
00128
00129
00130 void locDynOneEqEddy::correct(const tmp<volTensorField>& gradU)
00131 {
00132 LESModel::correct(gradU);
00133
00134 volSymmTensorField D = symm(gradU);
00135
00136 volScalarField KK = 0.5*(filter_(magSqr(U())) - magSqr(filter_(U())));
00137 KK.max(dimensionedScalar("small", KK.dimensions(), SMALL));
00138
00139 volScalarField P = 2.0*nuSgs_*magSqr(D);
00140
00141 fvScalarMatrix kEqn
00142 (
00143 fvm::ddt(k_)
00144 + fvm::div(phi(), k_)
00145 - fvm::laplacian(DkEff(), k_)
00146 ==
00147 P
00148 - fvm::Sp(ce(D, KK)*sqrt(k_)/delta(), k_)
00149 );
00150
00151 kEqn.relax();
00152 kEqn.solve();
00153
00154 bound(k_, k0());
00155
00156 updateSubGridScaleFields(D, KK);
00157 }
00158
00159
00160 bool locDynOneEqEddy::read()
00161 {
00162 if (GenEddyVisc::read())
00163 {
00164 filter_.read(coeffDict());
00165
00166 return true;
00167 }
00168 else
00169 {
00170 return false;
00171 }
00172 }
00173
00174
00175
00176
00177 }
00178 }
00179 }
00180
00181