![]() |
|
00001 /*---------------------------------------------------------------------------*\ 00002 ========= | 00003 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox 00004 \\ / O peration | 00005 \\ / A nd | Copyright (C) 1991-2005 OpenCFD Ltd. 00006 \\/ M anipulation | 00007 ------------------------------------------------------------------------------- 00008 License 00009 This file is part of OpenFOAM. 00010 00011 OpenFOAM is free software; you can redistribute it and/or modify it 00012 under the terms of the GNU General Public License as published by the 00013 Free Software Foundation; either version 2 of the License, or (at your 00014 option) any later version. 00015 00016 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT 00017 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or 00018 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License 00019 for more details. 00020 00021 You should have received a copy of the GNU General Public License 00022 along with OpenFOAM; if not, write to the Free Software Foundation, 00023 Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA 00024 00025 Class 00026 dynSmagorinsky 00027 00028 Description 00029 <pre> 00030 The Isochoric dynamic Smagorinsky Model 00031 ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ 00032 Algebraic eddy viscosity SGS model founded on the assumption that 00033 local equilibrium prevails, hence 00034 00035 B = 2/3*k*I - 2*nuEff*dev(D) 00036 00037 where 00038 00039 k = cI*delta^2*||D||^2 00040 nuEff = cD*delta^2*||D|| + nu 00041 00042 In the dynamic version of the choric Smagorinsky model 00043 the coefficients cI and cD are calculated during the simulation, 00044 00045 cI=<K*m>/<m*m> 00046 00047 and 00048 00049 cD=<L.M>/<M.M>, 00050 00051 where 00052 00053 K = 0.5*(F(U.U) - F(U).F(U)) 00054 m = delta^2*(4*||F(D)||^2 - F(||D||^2)) 00055 L = rho*(F(U*U) - F(U)*F(U) - 0.33*K*I) 00056 M = rho*delta^2*(F(||D||*dev(D)) - 4*||F(D)||*F(dev(D))) 00057 </pre> 00058 00059 SourceFiles 00060 dynSmagorinsky.C 00061 00062 \*---------------------------------------------------------------------------*/ 00063 00064 #ifndef dynSmagorinsky_H 00065 #define dynSmagorinsky_H 00066 00067 #include "Smagorinsky.H" 00068 #include "LESfilter.H" 00069 00070 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // 00071 00072 namespace Foam 00073 { 00074 namespace LESmodels 00075 { 00076 00077 /*---------------------------------------------------------------------------*\ 00078 Class dynSmagorinsky Declaration 00079 \*---------------------------------------------------------------------------*/ 00080 00081 class dynSmagorinsky 00082 : 00083 public GenEddyVisc 00084 { 00085 // Private data 00086 00087 volScalarField k_; 00088 00089 autoPtr<LESfilter> filterPtr_; 00090 LESfilter& filter_; 00091 00092 00093 // Private Member Functions 00094 00095 //- Calculate coefficients cD, cI from filtering velocity field 00096 dimensionedScalar cD(const volTensorField& D) const; 00097 dimensionedScalar cI(const volTensorField& D) const; 00098 00099 // Disallow default bitwise copy construct and assignment 00100 dynSmagorinsky(const dynSmagorinsky&); 00101 dynSmagorinsky& operator=(const dynSmagorinsky&); 00102 00103 00104 public: 00105 00106 //- Runtime type information 00107 TypeName("dynSmagorinsky"); 00108 00109 // Constructors 00110 00111 //- Constructor from components 00112 dynSmagorinsky 00113 ( 00114 const volVectorField& U, 00115 const surfaceScalarField& phi, 00116 transportModel& transport 00117 ); 00118 00119 00120 // Destructor 00121 00122 ~dynSmagorinsky(); 00123 00124 00125 // Member Functions 00126 00127 //- Return SGS kinetic energy 00128 tmp<volScalarField> k() const 00129 { 00130 return k_; 00131 } 00132 00133 //- Correct Eddy-Viscosity and related properties 00134 void correct(const tmp<volTensorField>& gradU); 00135 00136 //- Read turbulenceProperties dictionary 00137 bool read(); 00138 }; 00139 00140 00141 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // 00142 00143 } // End namespace LESmodels 00144 } // End namespace Foam 00145 00146 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // 00147 00148 #endif 00149 00150 // ************************************************************************* //