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

dynOneEqEddy.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::incompressible::LESModels::dynOneEqEddy
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
00031     the behaviour of k.
00032 
00033     Thus
00034     @verbatim
00035         d/dt(k) + div(U*k) - div(nuSgs*grad(k))
00036         =
00037         -B*L - ce*k^3/2/delta
00038 
00039     and
00040 
00041         B = 2/3*k*I - 2*nuSgs*dev(D)
00042         Beff = 2/3*k*I - 2*nuEff*dev(D)
00043 
00044     where
00045 
00046         D = symm(grad(U));
00047         nuSgs = ck*sqrt(k)*delta
00048         nuEff = nuSgs + nu
00049     @endverbatim
00050 
00051 SourceFiles
00052     dynOneEqEddy.C
00053 
00054 \*---------------------------------------------------------------------------*/
00055 
00056 #ifndef dynOneEqEddy_H
00057 #define dynOneEqEddy_H
00058 
00059 #include <incompressibleLESModels/GenEddyVisc.H>
00060 #include <LESfilters/LESfilter.H>
00061 
00062 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00063 
00064 namespace Foam
00065 {
00066 namespace incompressible
00067 {
00068 namespace LESModels
00069 {
00070 
00071 /*---------------------------------------------------------------------------*\
00072                            Class dynOneEqEddy Declaration
00073 \*---------------------------------------------------------------------------*/
00074 
00075 class dynOneEqEddy
00076 :
00077     public GenEddyVisc
00078 {
00079     // Private data
00080 
00081         volScalarField k_;
00082 
00083         autoPtr<LESfilter> filterPtr_;
00084         LESfilter& filter_;
00085 
00086 
00087     // Private Member Functions
00088 
00089         //- Update sub-grid scale fields
00090         void updateSubGridScaleFields(const volSymmTensorField& D);
00091 
00092         //- Calculate ck, ce by filtering the velocity field U.
00093         dimensionedScalar ck(const volSymmTensorField& D) const;
00094         dimensionedScalar ce(const volSymmTensorField& D) const;
00095 
00096         // Disallow default bitwise copy construct and assignment
00097         dynOneEqEddy(const dynOneEqEddy&);
00098         dynOneEqEddy& operator=(const dynOneEqEddy&);
00099 
00100 
00101 public:
00102 
00103     //- Runtime type information
00104     TypeName("dynOneEqEddy");
00105 
00106     // Constructors
00107 
00108         //- Construct from components
00109         dynOneEqEddy
00110         (
00111             const volVectorField& U,
00112             const surfaceScalarField& phi,
00113             transportModel& transport
00114         );
00115 
00116 
00117     //- Destructor
00118     virtual ~dynOneEqEddy()
00119     {}
00120 
00121 
00122     // Member Functions
00123 
00124         //- Return SGS kinetic energy
00125         virtual tmp<volScalarField> k() const
00126         {
00127             return k_;
00128         }
00129 
00130         //- Return the effective diffusivity for k
00131         tmp<volScalarField> DkEff() const
00132         {
00133             return tmp<volScalarField>
00134             (
00135                 new volScalarField("DkEff", nuSgs_ + nu())
00136             );
00137         }
00138 
00139         //- Correct Eddy-Viscosity and related properties
00140         virtual void correct(const tmp<volTensorField>& gradU);
00141 
00142         //- Read LESProperties dictionary
00143         virtual bool read();
00144 };
00145 
00146 
00147 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00148 
00149 } // End namespace LESModels
00150 } // End namespace incompressible
00151 } // End namespace Foam
00152 
00153 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00154 
00155 #endif
00156 
00157 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines