initial commit for version 1.5.x patch release
[OpenFOAM-1.5.x.git] / src / turbulenceModels / LES / incompressible / SpalartAllmaras / SpalartAllmaras.H
blob523da4611c665213e53fb11d44d8781f07919055
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2008 OpenCFD Ltd.
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
9     This file is part of OpenFOAM.
11     OpenFOAM is free software; you can redistribute it and/or modify it
12     under the terms of the GNU General Public License as published by the
13     Free Software Foundation; either version 2 of the License, or (at your
14     option) any later version.
16     OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17     ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18     FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
19     for more details.
21     You should have received a copy of the GNU General Public License
22     along with OpenFOAM; if not, write to the Free Software Foundation,
23     Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 Class
26     Foam::incompressible::LESModels::SpalartAllmaras
28 Description
29     SpalartAllmaras for incompressible flows
31 SourceFiles
32     SpalartAllmaras.C
34 \*---------------------------------------------------------------------------*/
36 #ifndef SpalartAllmaras_H
37 #define SpalartAllmaras_H
39 #include "LESModel.H"
40 #include "volFields.H"
42 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
44 namespace Foam
46 namespace incompressible
48 namespace LESModels
51 /*---------------------------------------------------------------------------*\
52                            Class SpalartAllmaras Declaration
53 \*---------------------------------------------------------------------------*/
55 class SpalartAllmaras
57     public LESModel
59     // Private data
61         dimensionedScalar alphaNut_;
63         dimensionedScalar Cb1_;
64         dimensionedScalar Cb2_;
65         dimensionedScalar Cv1_;
66         dimensionedScalar Cv2_;
67         dimensionedScalar CDES_;
68         dimensionedScalar ck_;
69         dimensionedScalar kappa_;
70         dimensionedScalar Cw1_;
71         dimensionedScalar Cw2_;
72         dimensionedScalar Cw3_;
75     // Private member functions
77         tmp<volScalarField> fv1() const;
78         tmp<volScalarField> fv2() const;
79         tmp<volScalarField> fv3() const;
80         tmp<volScalarField> fw(const volScalarField& Stilda) const;
82         // Disallow default bitwise copy construct and assignment
83         SpalartAllmaras(const SpalartAllmaras&);
84         SpalartAllmaras& operator=(const SpalartAllmaras&);
86         volScalarField nuTilda_;
87         volScalarField dTilda_;
88         volScalarField nuSgs_;
91 public:
93     //- Runtime type information
94     TypeName("SpalartAllmaras");
97     // Constructors
99         //- Constructor from components
100         SpalartAllmaras
101         (
102             const volVectorField& U,
103             const surfaceScalarField& phi,
104             transportModel& transport
105         );
108     // Destructor
110         ~SpalartAllmaras()
111         {}
114     // Member Functions
116         //- Return SGS kinetic energy
117         tmp<volScalarField> k() const
118         {
119             return sqr(nuSgs()/ck_/dTilda_);
120         }
122         //- Return sub-grid disipation rate
123         tmp<volScalarField> epsilon() const;
125         tmp<volScalarField> nuTilda() const
126         {
127             return nuTilda_;
128         }
130         //- Return SGS viscosity
131         tmp<volScalarField> nuSgs() const
132         {
133             return nuSgs_;
134         }
136         //- Return the sub-grid stress tensor.
137         tmp<volSymmTensorField> B() const;
139         //- Return the effective sub-grid turbulence stress tensor
140         //  including the laminar stress
141         tmp<volSymmTensorField> devBeff() const;
143         //- Return the deviatoric part of the divergence of Beff
144         //  i.e. the additional term in the filtered NSE.
145         tmp<fvVectorMatrix> divDevBeff(volVectorField& U) const;
147         //- Correct nuTilda and related properties
148         void correct(const tmp<volTensorField>& gradU);
150         //- Read turbulenceProperties dictionary
151         bool read();
155 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
157 } // End namespace LESModels
158 } // End namespace incompressible
159 } // End namespace Foam
161 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
163 #endif
165 // ************************************************************************* //