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