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

SpalartAllmaras.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::LESmodels::SpalartAllmaras
00026 
00027 Description
00028     SpalartAllmaras DES (SA + LES) turbulence model for incompressible flows
00029 
00030 SourceFiles
00031     SpalartAllmaras.C
00032 
00033 \*---------------------------------------------------------------------------*/
00034 
00035 #ifndef SpalartAllmaras_H
00036 #define SpalartAllmaras_H
00037 
00038 #include <incompressibleLESModels/LESModel.H>
00039 #include <finiteVolume/volFields.H>
00040 #include <finiteVolume/wallDist.H>
00041 
00042 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00043 
00044 namespace Foam
00045 {
00046 namespace incompressible
00047 {
00048 namespace LESModels
00049 {
00050 
00051 /*---------------------------------------------------------------------------*\
00052                         Class SpalartAllmaras Declaration
00053 \*---------------------------------------------------------------------------*/
00054 
00055 class SpalartAllmaras
00056 :
00057     public LESModel
00058 {
00059     // Private member functions
00060 
00061         //- Update sub-grid scale fields
00062         void updateSubGridScaleFields();
00063 
00064         // Disallow default bitwise copy construct and assignment
00065         SpalartAllmaras(const SpalartAllmaras&);
00066         SpalartAllmaras& operator=(const SpalartAllmaras&);
00067 
00068 
00069 protected:
00070 
00071     // Protected data
00072 
00073         dimensionedScalar sigmaNut_;
00074         dimensionedScalar kappa_;
00075 
00076 
00077         // Model constants
00078 
00079             dimensionedScalar Cb1_;
00080             dimensionedScalar Cb2_;
00081             dimensionedScalar Cv1_;
00082             dimensionedScalar Cv2_;
00083             dimensionedScalar CDES_;
00084             dimensionedScalar ck_;
00085             dimensionedScalar Cw1_;
00086             dimensionedScalar Cw2_;
00087             dimensionedScalar Cw3_;
00088 
00089 
00090         // Fields
00091 
00092             wallDist y_;
00093             volScalarField nuTilda_;
00094             volScalarField nuSgs_;
00095 
00096 
00097     // Protected member functions
00098 
00099         virtual tmp<volScalarField> fv1() const;
00100         virtual tmp<volScalarField> fv2() const;
00101         virtual tmp<volScalarField> fv3() const;
00102         virtual tmp<volScalarField> S(const volTensorField& gradU) const;
00103 
00104         virtual tmp<volScalarField> STilda
00105         (
00106             const volScalarField& S,
00107             const volScalarField& dTilda
00108         ) const;
00109 
00110         virtual tmp<volScalarField> r
00111         (
00112             const volScalarField& visc,
00113             const volScalarField& S,
00114             const volScalarField& dTilda
00115         ) const;
00116 
00117         virtual tmp<volScalarField> fw
00118         (
00119             const volScalarField& S,
00120             const volScalarField& dTilda
00121         ) const;
00122 
00123         //- Length scale
00124         virtual tmp<volScalarField> dTilda(const volScalarField& S) const;
00125 
00126 
00127 public:
00128 
00129     //- Runtime type information
00130     TypeName("SpalartAllmaras");
00131 
00132 
00133     // Constructors
00134 
00135         //- Construct from components
00136         SpalartAllmaras
00137         (
00138             const volVectorField& U,
00139             const surfaceScalarField& phi,
00140             transportModel& transport,
00141             const word& modelName = typeName
00142         );
00143 
00144 
00145     //- Destructor
00146     virtual ~SpalartAllmaras()
00147     {}
00148 
00149 
00150     // Member Functions
00151 
00152         //- Return SGS kinetic energy
00153         virtual tmp<volScalarField> k() const;
00154 
00155         //- Return sub-grid disipation rate
00156         virtual tmp<volScalarField> epsilon() const;
00157 
00158         tmp<volScalarField> nuTilda() const
00159         {
00160             return nuTilda_;
00161         }
00162 
00163         //- Return SGS viscosity
00164         virtual tmp<volScalarField> nuSgs() const
00165         {
00166             return nuSgs_;
00167         }
00168 
00169         //- Return the sub-grid stress tensor.
00170         virtual tmp<volSymmTensorField> B() const;
00171 
00172         //- Return the effective sub-grid turbulence stress tensor
00173         //  including the laminar stress
00174         virtual tmp<volSymmTensorField> devBeff() const;
00175 
00176         //- Return the deviatoric part of the divergence of Beff
00177         //  i.e. the additional term in the filtered NSE.
00178         virtual tmp<fvVectorMatrix> divDevBeff(volVectorField& U) const;
00179 
00180         //- Correct nuTilda and related properties
00181         virtual void correct(const tmp<volTensorField>& gradU);
00182 
00183         //- Read LESProperties dictionary
00184         virtual bool read();
00185 };
00186 
00187 
00188 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00189 
00190 } // End namespace LESModels
00191 } // End namespace incompressible
00192 } // End namespace Foam
00193 
00194 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00195 
00196 #endif
00197 
00198 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines