00001 /*---------------------------------------------------------------------------*\ 00002 ========= | 00003 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox 00004 \\ / O peration | 00005 \\ / A nd | Copyright (C) 1991-2010 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 00013 the Free Software Foundation, either version 3 of the License, or 00014 (at your 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, see <http://www.gnu.org/licenses/>. 00023 00024 Class 00025 Foam::incompressible::LESModels::mixedSmagorinsky 00026 00027 Description 00028 The mixed Isochoric Smagorinsky Model for incompressible flows. 00029 00030 The mixed model is a linear combination of an eddy viscosity model 00031 (Smagorinsky) with a scale similarity model. Hence 00032 @verbatim 00033 B = (L + C) + R = (F(v*v) - F(v)*F(v)) + R 00034 @endverbatim 00035 00036 The algebraic eddy viscosity SGS model is founded on the assumption 00037 that local equilibrium prevails, hence 00038 @verbatim 00039 R = 2/3*k*I - 2*nuEff*dev(D) 00040 where 00041 k = cI*delta^2*||D||^2 00042 nuEff = ck*sqrt(k)*delta + nu 00043 @endverbatim 00044 00045 The Leonard and cross contributions are incorporated 00046 by adding, 00047 @verbatim 00048 + div(((filter(U*U) - filter(U)*filter(U)) - 00049 0.333*I*tr(filter(U*U) - filter(U)*filter(U)))) 00050 + div((filter(U*epsilon) - filter(U)*filter(epsilon))) 00051 @endverbatim 00052 to the rhs. of the equations. 00053 00054 SourceFiles 00055 mixedSmagorinsky.C 00056 00057 \*---------------------------------------------------------------------------*/ 00058 00059 #ifndef mixedSmagorinsky_H 00060 #define mixedSmagorinsky_H 00061 00062 #include <incompressibleLESModels/scaleSimilarity.H> 00063 #include <incompressibleLESModels/Smagorinsky.H> 00064 00065 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // 00066 00067 namespace Foam 00068 { 00069 namespace incompressible 00070 { 00071 namespace LESModels 00072 { 00073 00074 /*---------------------------------------------------------------------------*\ 00075 Class mixedSmagorinsky Declaration 00076 \*---------------------------------------------------------------------------*/ 00077 00078 class mixedSmagorinsky 00079 : 00080 public scaleSimilarity, 00081 public Smagorinsky 00082 { 00083 // Private Member Functions 00084 00085 // Disallow default bitwise copy construct and assignment 00086 mixedSmagorinsky(const mixedSmagorinsky&); 00087 mixedSmagorinsky& operator=(const mixedSmagorinsky&); 00088 00089 public: 00090 00091 //- Runtime type information 00092 TypeName("mixedSmagorinsky"); 00093 00094 00095 // Constructors 00096 00097 //- Construct from components 00098 mixedSmagorinsky 00099 ( 00100 const volVectorField& U, 00101 const surfaceScalarField& phi, 00102 transportModel& transport 00103 ); 00104 00105 00106 //- Destructor 00107 virtual ~mixedSmagorinsky() 00108 {} 00109 00110 00111 // Member Functions 00112 00113 //- Return the SGS turbulent kinetic energy. 00114 virtual tmp<volScalarField> k() const; 00115 00116 //- Return the SGS turbulent disipation rate. 00117 virtual tmp<volScalarField> epsilon() const; 00118 00119 //- Return the SGS viscosity. 00120 virtual tmp<volScalarField> nuSgs() const 00121 { 00122 return nuSgs_; 00123 } 00124 00125 //- Return the sub-grid stress tensor. 00126 virtual tmp<volSymmTensorField> B() const; 00127 00128 //- Return the effective sub-grid turbulence stress tensor 00129 // including the laminar stress 00130 virtual tmp<volSymmTensorField> devBeff() const; 00131 00132 //- Implementation of div(B). This is necessary to override 00133 // (and include) the div(B) terms from both the parent classes. 00134 virtual tmp<fvVectorMatrix> divDevBeff(volVectorField& U) const; 00135 00136 //- Correct Eddy-Viscosity and related properties 00137 virtual void correct(const tmp<volTensorField>& gradU); 00138 00139 //- Read LESProperties dictionary 00140 virtual bool read(); 00141 }; 00142 00143 00144 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // 00145 00146 } // End namespace LESModels 00147 } // End namespace incompressible 00148 } // End namespace Foam 00149 00150 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // 00151 00152 #endif 00153 00154 // ************************ vim: set sw=4 sts=4 et: ************************ //