ENH: autoLayerDriver: better layering information message
[OpenFOAM-2.0.x.git] / src / turbulenceModels / incompressible / LES / kOmegaSSTSAS / kOmegaSSTSAS.H
blob2ca98380afad07e6d762abd6ba539e23e5aa5a9d
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::incompressible::LESModels::kOmegaSSTSAS
27 Description
28     kOmegaSSTSAS LES turbulence model for incompressible flows
29     based on:
31     "Evaluation of the SST-SAS model: channel flow, asymmetric diffuser
32     and axi-symmetric hill".
33     European Conference on Computational Fluid Dynamics ECCOMAS CFD 2006.
34     Lars Davidson
37     The first term of the Qsas expression is corrected following:
39     DESider A European Effort on Hybrid RANS-LES Modelling:
40     Results of the European-Union Funded Project, 2004 - 2007
41     (Notes on Numerical Fluid Mechanics and Multidisciplinary Design).
42     Chapter 2, section 8 Formulation of the Scale-Adaptive Simulation (SAS)
43     Model during the DESIDER Project. Published in Springer-Verlag Berlin
44     Heidelberg 2009.
45     F. R. Menter and Y. Egorov.
48 SourceFiles
49     kOmegaSSTSAS.C
51 \*---------------------------------------------------------------------------*/
53 #ifndef kOmegaSSTSAS_H
54 #define kOmegaSSTSAS_H
56 #include "LESModel.H"
57 #include "volFields.H"
58 #include "wallDist.H"
60 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
62 namespace Foam
64 namespace incompressible
66 namespace LESModels
69 /*---------------------------------------------------------------------------*\
70                         Class kOmegaSSTSAS Declaration
71 \*---------------------------------------------------------------------------*/
73 class kOmegaSSTSAS
75     public LESModel
77     // Private Member Functions
79         //- Update sub-grid scale fields
80         void updateSubGridScaleFields(const volScalarField& D);
82         // Disallow default bitwise copy construct and assignment
83         kOmegaSSTSAS(const kOmegaSSTSAS&);
84         kOmegaSSTSAS& operator=(const kOmegaSSTSAS&);
87 protected:
89     // Protected data
91         // Model constants
93             dimensionedScalar alphaK1_;
94             dimensionedScalar alphaK2_;
96             dimensionedScalar alphaOmega1_;
97             dimensionedScalar alphaOmega2_;
99             dimensionedScalar gamma1_;
100             dimensionedScalar gamma2_;
102             dimensionedScalar beta1_;
103             dimensionedScalar beta2_;
105             dimensionedScalar betaStar_;
107             dimensionedScalar a1_;
108             dimensionedScalar c1_;
109             dimensionedScalar Cs_;
111             dimensionedScalar alphaPhi_;
112             dimensionedScalar zetaTilda2_;
113             dimensionedScalar FSAS_;
115             dimensionedScalar omegaMin_;
117             wallDist y_;
118             dimensionedScalar Cmu_;
119             dimensionedScalar kappa_;
122         // Fields
124             volScalarField k_;
125             volScalarField omega_;
126             volScalarField nuSgs_;
129     // Protected Member Functions
131         tmp<volScalarField> Lvk2
132         (
133             const volScalarField& S2
134         ) const;
136         tmp<volScalarField> F1(const volScalarField& CDkOmega) const;
137         tmp<volScalarField> F2() const;
139         tmp<volScalarField> blend
140         (
141             const volScalarField& F1,
142             const dimensionedScalar& psi1,
143             const dimensionedScalar& psi2
144         ) const
145         {
146             return F1*(psi1 - psi2) + psi2;
147         }
149         tmp<volScalarField> alphaK
150         (
151             const volScalarField& F1
152         ) const
153         {
154             return blend(F1, alphaK1_, alphaK2_);
155         }
157         tmp<volScalarField> alphaOmega
158         (
159              const volScalarField& F1
160         ) const
161         {
162             return blend(F1, alphaOmega1_, alphaOmega2_);
163         }
165         tmp<volScalarField> beta
166         (
167             const volScalarField& F1
168         ) const
169         {
170             return blend(F1, beta1_, beta2_);
171         }
173         tmp<volScalarField> gamma
174         (
175             const volScalarField& F1
176         ) const
177         {
178             return blend(F1, gamma1_, gamma2_);
179         }
182 public:
184     //- Runtime type information
185     TypeName("kOmegaSSTSAS");
188     // Constructors
190         //- Construct from components
191         kOmegaSSTSAS
192         (
193             const volVectorField& U,
194             const surfaceScalarField& phi,
195             transportModel& transport,
196             const word& turbulenceModelName = turbulenceModel::typeName,
197             const word& modelName = typeName
198         );
201     //- Destructor
202     virtual ~kOmegaSSTSAS()
203     {}
206     // Member Functions
208         //- Return SGS kinetic energy
209         virtual tmp<volScalarField> k() const
210         {
211             return k_;
212         }
214         //- Return omega
215         virtual tmp<volScalarField> omega() const
216         {
217             return omega_;
218         }
220         //- Return the effective diffusivity for k
221         tmp<volScalarField> DkEff(const volScalarField& F1) const
222         {
223             return tmp<volScalarField>
224             (
225                  new volScalarField("DkEff", alphaK(F1)*nuSgs_ + nu())
226             );
227         }
229         //- Return the effective diffusivity for omega
230         tmp<volScalarField> DomegaEff(const volScalarField& F1) const
231         {
232             return tmp<volScalarField>
233             (
234                 new volScalarField("DomegaEff", alphaOmega(F1)*nuSgs_ + nu())
235             );
236         }
238         //- Return sub-grid disipation rate
239         virtual tmp<volScalarField> epsilon() const;
241         //- Return SGS viscosity
242         virtual tmp<volScalarField> nuSgs() const
243         {
244             return nuSgs_;
245         }
247         //- Return the sub-grid stress tensor.
248         virtual tmp<volSymmTensorField> B() const;
250         //- Return the effective sub-grid turbulence stress tensor
251         //  including the laminar stress
252         virtual tmp<volSymmTensorField> devBeff() const;
254         //- Return the deviatoric part of the divergence of Beff
255         //  i.e. the additional term in the filtered NSE.
256         virtual tmp<fvVectorMatrix> divDevBeff(volVectorField& U) const;
258         //- Solve the turbulence equations (k-w) and correct the turbulence
259         //  viscosity
260         virtual void correct(const tmp<volTensorField>& gradU);
262         //- Read LESProperties dictionary
263         virtual bool read();
267 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
269 } // End namespace LESModels
270 } // End namespace incompressible
271 } // End namespace Foam
273 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
275 #endif
277 // ************************************************************************* //