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

kOmegaSSTSAS.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) 2008-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::kOmegaSSTSAS
00026 
00027 Description
00028     kOmegaSSTSAS LES turbulence model for incompressible flows
00029 
00030     Reference:
00031 
00032     DESider A European Effort on Hybrid RANS-LES Modelling:
00033     Results of the European-Union Funded Project, 2004 - 2007
00034     (Notes on Numerical Fluid Mechanics and Multidisciplinary Design).
00035     Chapter 8 Formulation of the Scale-Adaptive Simulation (SAS) Model during
00036     the DESIDER Project. Published in Springer-Verlag Berlin Heidelberg 2009.
00037     F. R. Menter and Y. Egorov.
00038 
00039 SourceFiles
00040     kOmegaSSTSAS.C
00041 
00042 \*---------------------------------------------------------------------------*/
00043 
00044 #ifndef kOmegaSSTSAS_H
00045 #define kOmegaSSTSAS_H
00046 
00047 #include <incompressibleLESModels/LESModel.H>
00048 #include <finiteVolume/volFields.H>
00049 #include <finiteVolume/wallDist.H>
00050 
00051 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00052 
00053 namespace Foam
00054 {
00055 namespace incompressible
00056 {
00057 namespace LESModels
00058 {
00059 
00060 /*---------------------------------------------------------------------------*\
00061                         Class kOmegaSSTSAS Declaration
00062 \*---------------------------------------------------------------------------*/
00063 
00064 class kOmegaSSTSAS
00065 :
00066     public LESModel
00067 {
00068     // Private member functions
00069 
00070         //- Update sub-grid scale fields
00071         void updateSubGridScaleFields(const volScalarField& D);
00072 
00073         // Disallow default bitwise copy construct and assignment
00074         kOmegaSSTSAS(const kOmegaSSTSAS&);
00075         kOmegaSSTSAS& operator=(const kOmegaSSTSAS&);
00076 
00077 
00078 protected:
00079 
00080     // Protected data
00081 
00082         // Model constants
00083 
00084             dimensionedScalar alphaK1_;
00085             dimensionedScalar alphaK2_;
00086 
00087             dimensionedScalar alphaOmega1_;
00088             dimensionedScalar alphaOmega2_;
00089 
00090             dimensionedScalar gamma1_;
00091             dimensionedScalar gamma2_;
00092 
00093             dimensionedScalar beta1_;
00094             dimensionedScalar beta2_;
00095 
00096             dimensionedScalar betaStar_;
00097 
00098             dimensionedScalar a1_;
00099             dimensionedScalar c1_;
00100             dimensionedScalar Cs_;
00101 
00102             dimensionedScalar alphaPhi_;
00103             dimensionedScalar zetaTilda2_;
00104             dimensionedScalar FSAS_;
00105 
00106             dimensionedScalar omega0_;
00107             dimensionedScalar omegaSmall_;
00108 
00109             wallDist y_;
00110             dimensionedScalar Cmu_;
00111             dimensionedScalar kappa_;
00112 
00113 
00114         // Fields
00115 
00116             volScalarField k_;
00117             volScalarField omega_;
00118             volScalarField nuSgs_;
00119 
00120 
00121     // Protected member functions
00122 
00123         tmp<volScalarField> Lvk2
00124         (
00125             const volScalarField& S2
00126         ) const;
00127 
00128         tmp<volScalarField> F1(const volScalarField& CDkOmega) const;
00129         tmp<volScalarField> F2() const;
00130 
00131         tmp<volScalarField> blend
00132         (
00133             const volScalarField& F1,
00134             const dimensionedScalar& psi1,
00135             const dimensionedScalar& psi2
00136         ) const
00137         {
00138             return F1*(psi1 - psi2) + psi2;
00139         }
00140 
00141         tmp<volScalarField> alphaK
00142         (
00143             const volScalarField& F1
00144         ) const
00145         {
00146             return blend(F1, alphaK1_, alphaK2_);
00147         }
00148 
00149         tmp<volScalarField> alphaOmega
00150         (
00151              const volScalarField& F1
00152         ) const
00153         {
00154             return blend(F1, alphaOmega1_, alphaOmega2_);
00155         }
00156 
00157         tmp<volScalarField> beta
00158         (
00159             const volScalarField& F1
00160         ) const
00161         {
00162             return blend(F1, beta1_, beta2_);
00163         }
00164 
00165         tmp<volScalarField> gamma
00166         (
00167             const volScalarField& F1
00168         ) const
00169         {
00170             return blend(F1, gamma1_, gamma2_);
00171         }
00172 
00173 
00174 public:
00175 
00176     //- Runtime type information
00177     TypeName("kOmegaSSTSAS");
00178 
00179 
00180     // Constructors
00181 
00182         //- Construct from components
00183         kOmegaSSTSAS
00184         (
00185             const volVectorField& U,
00186             const surfaceScalarField& phi,
00187             transportModel& transport,
00188             const word& modelName = typeName
00189         );
00190 
00191 
00192     //- Destructor
00193     virtual ~kOmegaSSTSAS()
00194     {}
00195 
00196 
00197     // Member Functions
00198 
00199         //- Return SGS kinetic energy
00200         virtual tmp<volScalarField> k() const
00201         {
00202             return k_;
00203         }
00204 
00205         //- Return omega
00206         virtual tmp<volScalarField> omega() const
00207         {
00208             return omega_;
00209         }
00210 
00211         //- Return the effective diffusivity for k
00212         tmp<volScalarField> DkEff(const volScalarField& F1) const
00213         {
00214             return tmp<volScalarField>
00215             (
00216                  new volScalarField("DkEff", alphaK(F1)*nuSgs_ + nu())
00217             );
00218         }
00219 
00220         //- Return the effective diffusivity for omega
00221         tmp<volScalarField> DomegaEff(const volScalarField& F1) const
00222         {
00223             return tmp<volScalarField>
00224             (
00225                 new volScalarField("DomegaEff", alphaOmega(F1)*nuSgs_ + nu())
00226             );
00227         }
00228 
00229         //- Return sub-grid disipation rate
00230         virtual tmp<volScalarField> epsilon() const;
00231 
00232         //- Return SGS viscosity
00233         virtual tmp<volScalarField> nuSgs() const
00234         {
00235             return nuSgs_;
00236         }
00237 
00238         //- Return the sub-grid stress tensor.
00239         virtual tmp<volSymmTensorField> B() const;
00240 
00241         //- Return the effective sub-grid turbulence stress tensor
00242         //  including the laminar stress
00243         virtual tmp<volSymmTensorField> devBeff() const;
00244 
00245         //- Return the deviatoric part of the divergence of Beff
00246         //  i.e. the additional term in the filtered NSE.
00247         virtual tmp<fvVectorMatrix> divDevBeff(volVectorField& U) const;
00248 
00249         //- Solve the turbulence equations (k-w) and correct the turbulence
00250         //  viscosity
00251         virtual void correct(const tmp<volTensorField>& gradU);
00252 
00253         //- Read LESProperties dictionary
00254         virtual bool read();
00255 };
00256 
00257 
00258 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00259 
00260 } // End namespace LESModels
00261 } // End namespace incompressible
00262 } // End namespace Foam
00263 
00264 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00265 
00266 #endif
00267 
00268 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines