BUG: UListIO: byteSize overflowing on really big faceLists
[OpenFOAM-2.0.x.git] / src / turbulenceModels / compressible / LES / SpalartAllmaras / SpalartAllmaras.H
blob5caa292642e91ba3cd6a0e9409ee10605af254c3
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2011 OpenFOAM Foundation
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
13     the Free Software Foundation, either version 3 of the License, or
14     (at your 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, see <http://www.gnu.org/licenses/>.
24 Class
25     Foam::compressible::LESModels::SpalartAllmaras
27 Description
28     SpalartAllmaras for compressible flows
30 SourceFiles
31     SpalartAllmaras.C
33 \*---------------------------------------------------------------------------*/
35 #ifndef compressibleSpalartAllmaras_H
36 #define compressibleSpalartAllmaras_H
38 #include "LESModel.H"
39 #include "volFields.H"
41 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
43 namespace Foam
45 namespace compressible
47 namespace LESModels
50 /*---------------------------------------------------------------------------*\
51                        Class SpalartAllmaras Declaration
52 \*---------------------------------------------------------------------------*/
54 class SpalartAllmaras
56     public LESModel
58     // Private data
60         // Model coefficients
62             dimensionedScalar sigmaNut_;
63             dimensionedScalar Prt_;
65             dimensionedScalar Cb1_;
66             dimensionedScalar Cb2_;
67             dimensionedScalar Cv1_;
68             dimensionedScalar Cv2_;
69             dimensionedScalar CDES_;
70             dimensionedScalar ck_;
71             dimensionedScalar kappa_;
72             dimensionedScalar Cw1_;
73             dimensionedScalar Cw2_;
74             dimensionedScalar Cw3_;
77     // Fields
79         volScalarField nuTilda_;
80         volScalarField dTilda_;
81         volScalarField muSgs_;
82         volScalarField alphaSgs_;
85     // Private Member Functions
87         //- Update sub-grid scale fields
88         void updateSubGridScaleFields();
90         tmp<volScalarField> fv1() const;
91         tmp<volScalarField> fv2() const;
92         tmp<volScalarField> fv3() const;
93         tmp<volScalarField> fw(const volScalarField& Stilda) const;
95         // Disallow default bitwise copy construct and assignment
96         SpalartAllmaras(const SpalartAllmaras&);
97         SpalartAllmaras& operator=(const SpalartAllmaras&);
100 public:
102     //- Runtime type information
103     TypeName("SpalartAllmaras");
106     // Constructors
108         //- Constructor from components
109         SpalartAllmaras
110         (
111             const volScalarField& rho,
112             const volVectorField& U,
113             const surfaceScalarField& phi,
114             const basicThermo& thermoPhysicalModel,
115             const word& turbulenceModelName = turbulenceModel::typeName,
116             const word& modelName = typeName
117         );
120     //- Destructor
121     virtual ~SpalartAllmaras()
122     {}
125     // Member Functions
127         tmp<volScalarField> nuTilda() const
128         {
129             return nuTilda_;
130         }
132         //- Return SGS kinetic energy
133         virtual tmp<volScalarField> k() const
134         {
135             return sqr(muSgs()/rho()/ck_/dTilda_);
136         }
138         //- Return sub-grid disipation rate
139         virtual tmp<volScalarField> epsilon() const;
141         //- Return SGS viscosity
142         virtual tmp<volScalarField> muSgs() const
143         {
144             return muSgs_;
145         }
147         //- Return SGS thermal diffusivity
148         virtual tmp<volScalarField> alphaSgs() const
149         {
150             return alphaSgs_;
151         }
153         //- Return the sub-grid stress tensor.
154         virtual tmp<volSymmTensorField> B() const;
156         //- Return the deviatoric part of the effective sub-grid
157         //  turbulence stress tensor including the laminar stress
158         virtual tmp<volSymmTensorField> devRhoBeff() const;
160         //- Returns div(rho*dev(B)).
161         // This is the additional term due to the filtering of the NSE.
162         virtual tmp<fvVectorMatrix> divDevRhoBeff(volVectorField& U) const;
164         //- Correct nuTilda and related properties
165         virtual void correct(const tmp<volTensorField>& gradU);
167         //- Read LESProperties dictionary
168         virtual bool read();
172 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
174 } // End namespace LESModels
175 } // End namespace compressible
176 } // End namespace Foam
178 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
180 #endif
182 // ************************************************************************* //