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

RASModel.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 Namespace
00025     Foam::compressible::RASModels
00026 
00027 Description
00028     Namespace for compressible RAS turbulence models.
00029 
00030 
00031 Class
00032     Foam::compressible::RASModel
00033 
00034 Description
00035     Abstract base class for turbulence models for compressible and combusting
00036     flows.
00037 
00038 SourceFiles
00039     RASModel.C
00040     newTurbulenceModel.C
00041 
00042 \*---------------------------------------------------------------------------*/
00043 
00044 #ifndef compressibleRASModel_H
00045 #define compressibleRASModel_H
00046 
00047 #include <compressibleTurbulenceModel/turbulenceModel.H>
00048 #include <finiteVolume/volFields.H>
00049 #include <finiteVolume/surfaceFields.H>
00050 #include <finiteVolume/nearWallDist.H>
00051 #include <finiteVolume/fvm.H>
00052 #include <finiteVolume/fvc.H>
00053 #include <finiteVolume/fvMatrices.H>
00054 #include <basicThermophysicalModels/basicThermo.H>
00055 #include <OpenFOAM/IOdictionary.H>
00056 #include <OpenFOAM/Switch.H>
00057 #include <finiteVolume/bound.H>
00058 #include <OpenFOAM/autoPtr.H>
00059 #include <OpenFOAM/runTimeSelectionTables.H>
00060 
00061 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00062 
00063 namespace Foam
00064 {
00065 namespace compressible
00066 {
00067 
00068 /*---------------------------------------------------------------------------*\
00069                            Class RASModel Declaration
00070 \*---------------------------------------------------------------------------*/
00071 
00072 class RASModel
00073 :
00074     public turbulenceModel,
00075     public IOdictionary
00076 {
00077 
00078 protected:
00079 
00080     // Protected data
00081 
00082         //- Turbulence on/off flag
00083         Switch turbulence_;
00084 
00085         //- Flag to print the model coeffs at run-time
00086         Switch printCoeffs_;
00087 
00088         //- Model coefficients dictionary
00089         dictionary coeffDict_;
00090 
00091         //- Value of y+ at the edge of the laminar sublayer
00092         scalar yPlusLam_;
00093 
00094         //- Lower limit of k
00095         dimensionedScalar k0_;
00096 
00097         //- Lower limit of epsilon
00098         dimensionedScalar epsilon0_;
00099 
00100         //- Small epsilon value used to avoid divide by zero
00101         dimensionedScalar epsilonSmall_;
00102 
00103         //- Lower limit for omega
00104         dimensionedScalar omega0_;
00105 
00106         //- Small omega value used to avoid divide by zero
00107         dimensionedScalar omegaSmall_;
00108 
00109         //- Near wall distance boundary field
00110         nearWallDist y_;
00111 
00112 
00113     // Protected member functions
00114 
00115         //- Print model coefficients
00116         virtual void printCoeffs();
00117 
00118 
00119 private:
00120 
00121     // Private Member Functions
00122 
00123         //- Disallow default bitwise copy construct
00124         RASModel(const RASModel&);
00125 
00126         //- Disallow default bitwise assignment
00127         void operator=(const RASModel&);
00128 
00129 
00130 public:
00131 
00132     //- Runtime type information
00133     TypeName("RASModel");
00134 
00135 
00136     // Declare run-time constructor selection table
00137 
00138         declareRunTimeSelectionTable
00139         (
00140             autoPtr,
00141             RASModel,
00142             dictionary,
00143             (
00144                 const volScalarField& rho,
00145                 const volVectorField& U,
00146                 const surfaceScalarField& phi,
00147                 const basicThermo& thermoPhysicalModel
00148             ),
00149             (rho, U, phi, thermoPhysicalModel)
00150         );
00151 
00152 
00153     // Constructors
00154 
00155         //- Construct from components
00156         RASModel
00157         (
00158             const word& type,
00159             const volScalarField& rho,
00160             const volVectorField& U,
00161             const surfaceScalarField& phi,
00162             const basicThermo& thermoPhysicalModel
00163         );
00164 
00165 
00166     // Selectors
00167 
00168         //- Return a reference to the selected turbulence model
00169         static autoPtr<RASModel> New
00170         (
00171             const volScalarField& rho,
00172             const volVectorField& U,
00173             const surfaceScalarField& phi,
00174             const basicThermo& thermoPhysicalModel
00175         );
00176 
00177 
00178     //- Destructor
00179     virtual ~RASModel()
00180     {}
00181 
00182 
00183     // Member Functions
00184 
00185         // Access
00186 
00187             //- Return the value of k0 which k is not allowed to be less than
00188             const dimensionedScalar& k0() const
00189             {
00190                 return k0_;
00191             }
00192 
00193             //- Return the value of epsilon0 which epsilon is not allowed to be
00194             //  less than
00195             const dimensionedScalar& epsilon0() const
00196             {
00197                 return epsilon0_;
00198             }
00199 
00200             //- Return the value of epsilonSmall which is added to epsilon when
00201             //  calculating nut
00202             const dimensionedScalar& epsilonSmall() const
00203             {
00204                 return epsilonSmall_;
00205             }
00206 
00207             //- Return the value of omega0 which epsilon is not allowed to be
00208             //  less than
00209             const dimensionedScalar& omega0() const
00210             {
00211                 return omega0_;
00212             }
00213 
00214             //- Return the value of omegaSmall which is added to epsilon when
00215             //  calculating nut
00216             const dimensionedScalar& omegaSmall() const
00217             {
00218                 return omegaSmall_;
00219             }
00220 
00221             //- Allow k0 to be changed
00222             dimensionedScalar& k0()
00223             {
00224                 return k0_;
00225             }
00226 
00227             //- Allow epsilon0 to be changed
00228             dimensionedScalar& epsilon0()
00229             {
00230                 return epsilon0_;
00231             }
00232 
00233             //- Allow epsilonSmall to be changed
00234             dimensionedScalar& epsilonSmall()
00235             {
00236                 return epsilonSmall_;
00237             }
00238 
00239             //- Allow omega0 to be changed
00240             dimensionedScalar& omega0()
00241             {
00242                 return omega0_;
00243             }
00244 
00245             //- Allow omegaSmall to be changed
00246             dimensionedScalar& omegaSmall()
00247             {
00248                 return omegaSmall_;
00249             }
00250 
00251             //- Return the near wall distances
00252             const nearWallDist& y() const
00253             {
00254                 return y_;
00255             }
00256 
00257             //- Calculate y+ at the edge of the laminar sublayer
00258             scalar yPlusLam(const scalar kappa, const scalar E) const;
00259 
00260             //- Const access to the coefficients dictionary
00261             const dictionary& coeffDict() const
00262             {
00263                 return coeffDict_;
00264             }
00265 
00266 
00267         //- Return the turbulence viscosity
00268         virtual tmp<volScalarField> mut() const = 0;
00269 
00270         //- Return the effective viscosity
00271         virtual tmp<volScalarField> muEff() const
00272         {
00273             return tmp<volScalarField>
00274             (
00275                 new volScalarField("muEff", mut() + mu())
00276             );
00277         }
00278 
00279         //- Return the effective turbulent thermal diffusivity
00280         virtual tmp<volScalarField> alphaEff() const = 0;
00281 
00282         //- Return the turbulence kinetic energy
00283         virtual tmp<volScalarField> k() const = 0;
00284 
00285         //- Return the turbulence kinetic energy dissipation rate
00286         virtual tmp<volScalarField> epsilon() const = 0;
00287 
00288         //- Return the Reynolds stress tensor
00289         virtual tmp<volSymmTensorField> R() const = 0;
00290 
00291         //- Return the effective stress tensor including the laminar stress
00292         virtual tmp<volSymmTensorField> devRhoReff() const = 0;
00293 
00294         //- Return the source term for the momentum equation
00295         virtual tmp<fvVectorMatrix> divDevRhoReff(volVectorField& U) const = 0;
00296 
00297         //- Return yPlus for the given patch
00298         virtual tmp<scalarField> yPlus
00299         (
00300             const label patchI,
00301             const scalar Cmu
00302         ) const;
00303 
00304         //- Solve the turbulence equations and correct the turbulence viscosity
00305         virtual void correct() = 0;
00306 
00307         //- Read RASProperties dictionary
00308         virtual bool read() = 0;
00309 };
00310 
00311 
00312 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00313 
00314 } // End namespace compressible
00315 } // End namespace Foam
00316 
00317 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00318 
00319 #endif
00320 
00321 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines