OpenFOAM logo
Open Source CFD Toolkit

lowReOneEqEddy.H

Go to the documentation of this file.
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     lowReOneEqEddy
00027 
00028 Description
00029 <pre>
00030     One Equation Eddy Viscosity Model
00031     ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
00032         d/dt(rho*k) + div(rho*U*k) - div(muEff*grad(k))
00033         =
00034         -rho*B*L - ce*rho*k^3/2/delta
00035 
00036     and
00037 
00038         B = 2/3*k*I - 2*nuSgs*dev(D)
00039 
00040     where
00041 
00042         nuSgsHiRe = ck*sqrt(k)*delta
00043         nuSgs = (nu/beta)*(1 - exp(-beta*nuSgsHiRe/nu));
00044         muEff = rho*nuSgs + mu
00045 </pre>
00046 
00047 SourceFiles
00048     lowReOneEqEddy.C
00049 
00050 \*---------------------------------------------------------------------------*/
00051 
00052 #ifndef compressibleLowReOneEqEddy_H
00053 #define compressibleLowReOneEqEddy_H
00054 
00055 #include "GenEddyVisc.H"
00056 
00057 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00058 
00059 namespace Foam
00060 {
00061 namespace compressible
00062 {
00063 namespace LESmodels
00064 {
00065 
00066 /*---------------------------------------------------------------------------*\
00067                            Class lowReOneEqEddy Declaration
00068 \*---------------------------------------------------------------------------*/
00069 
00070 class lowReOneEqEddy
00071 :
00072     public GenEddyVisc
00073 {
00074     // Private data
00075 
00076         dimensionedScalar ck_;
00077         dimensionedScalar beta_;
00078 
00079     // Private Member Functions
00080 
00081         // Disallow default bitwise copy construct and assignment
00082         lowReOneEqEddy(const lowReOneEqEddy&);
00083         lowReOneEqEddy& operator=(const lowReOneEqEddy&);
00084 
00085 
00086 public:
00087 
00088     //- Runtime type information
00089     TypeName("lowReOneEqEddy");
00090 
00091 
00092     // Constructors
00093 
00094         //- Constructor from components
00095         lowReOneEqEddy
00096         (
00097             const volScalarField& rho,
00098             const volVectorField& U,
00099             const surfaceScalarField& phi,
00100             const basicThermo& thermoPhysicalModel
00101         );
00102 
00103 
00104     // Destructor
00105 
00106         ~lowReOneEqEddy()
00107         {}
00108 
00109 
00110     // Member Functions
00111 
00112         //- Return the effective diffusivity for k
00113         tmp<volScalarField> DkEff() const
00114         {
00115             return tmp<volScalarField> 
00116             (
00117                 new volScalarField("DkEff", muSgs_ + mu())
00118             );
00119         }
00120 
00121         //- Correct Eddy-Viscosity and related properties
00122         void correct(const tmp<volTensorField>& gradU);
00123 
00124         //- Read turbulenceProperties dictionary
00125         bool read();
00126 };
00127 
00128 
00129 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00130 
00131 } // End namespace LESmodels
00132 } // End namespace compressible
00133 } // End namespace Foam
00134 
00135 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00136 
00137 #endif
00138 
00139 // ************************************************************************* //
For further information go to www.openfoam.org