FreeFOAM The Cross-Platform CFD Toolkit
Hosted by SourceForge:
Get FreeFOAM at SourceForge.net.
            Fast, secure and Free Open Source software downloads

oneEqEddy.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-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::compressible::LESModels::oneEqEddy
00026 
00027 Description
00028     One Equation Eddy Viscosity Model for incompressible flows
00029 
00030     Eddy viscosity SGS model using a modeled balance equation to simulate the
00031     behaviour of k, hence,
00032     @verbatim
00033         d/dt(rho*k) + div(rho*U*k) - div(muEff*grad(k))
00034         =
00035         -rho*D:B - ce*rho*k^3/2/delta
00036 
00037     and
00038 
00039         B = 2/3*k*I - 2*nuSgs*dev(D)
00040 
00041     where
00042 
00043         D = symm(grad(U));
00044         muSgs = ck*rho*sqrt(k)*delta
00045     @endverbatim
00046 
00047 
00048 SourceFiles
00049     oneEqEddy.C
00050 
00051 \*---------------------------------------------------------------------------*/
00052 
00053 #ifndef compressibleOneEqEddy_H
00054 #define compressibleOneEqEddy_H
00055 
00056 #include <compressibleLESModels/GenEddyVisc.H>
00057 
00058 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00059 
00060 namespace Foam
00061 {
00062 namespace compressible
00063 {
00064 namespace LESModels
00065 {
00066 
00067 /*---------------------------------------------------------------------------*\
00068                            Class oneEqEddy Declaration
00069 \*---------------------------------------------------------------------------*/
00070 
00071 class oneEqEddy
00072 :
00073     public GenEddyVisc
00074 {
00075     // Private data
00076 
00077         dimensionedScalar ck_;
00078 
00079 
00080     // Private Member Functions
00081 
00082         //- Update sub-grid scale fields
00083         void updateSubGridScaleFields();
00084 
00085         // Disallow default bitwise copy construct and assignment
00086         oneEqEddy(const oneEqEddy&);
00087         oneEqEddy& operator=(const oneEqEddy&);
00088 
00089 
00090 public:
00091 
00092     //- Runtime type information
00093     TypeName("oneEqEddy");
00094 
00095 
00096     // Constructors
00097 
00098         //- Constructor from components
00099         oneEqEddy
00100         (
00101             const volScalarField& rho,
00102             const volVectorField& U,
00103             const surfaceScalarField& phi,
00104             const basicThermo& thermoPhysicalModel
00105         );
00106 
00107 
00108     // Destructor
00109     virtual ~oneEqEddy()
00110     {}
00111 
00112 
00113     // Member Functions
00114 
00115         //- Return the effective diffusivity for k
00116         tmp<volScalarField> DkEff() const
00117         {
00118             return tmp<volScalarField>
00119             (
00120                 new volScalarField("DkEff", muSgs_ + mu())
00121             );
00122         }
00123 
00124         //- Correct Eddy-Viscosity and related properties
00125         virtual void correct(const tmp<volTensorField>& gradU);
00126 
00127         //- Read LESProperties dictionary
00128         virtual bool read();
00129 };
00130 
00131 
00132 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00133 
00134 } // End namespace LESModels
00135 } // End namespace compressible
00136 } // End namespace Foam
00137 
00138 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00139 
00140 #endif
00141 
00142 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines