Forward compatibility: flex
[foam-extend-3.2.git] / src / turbulenceModels / incompressible / LES / kOmegaSSTSAS / kOmegaSSTSAS.H
blob3d1314d197b3c3f97dc1514918dbc6536069c4ea
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | foam-extend: Open Source CFD
4    \\    /   O peration     | Version:     3.2
5     \\  /    A nd           | Web:         http://www.foam-extend.org
6      \\/     M anipulation  | For copyright notice see file Copyright
7 -------------------------------------------------------------------------------
8 License
9     This file is part of foam-extend.
11     foam-extend 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 3 of the License, or (at your
14     option) any later version.
16     foam-extend is distributed in the hope that it will be useful, but
17     WITHOUT ANY WARRANTY; without even the implied warranty of
18     MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19     General Public License for more details.
21     You should have received a copy of the GNU General Public License
22     along with foam-extend.  If not, see <http://www.gnu.org/licenses/>.
24 Class
25     Foam::LESmodels::kOmegaSSTSAS
27 Description
28     kOmegaSSTSAS LES turbulence model for incompressible flows
30     Reference:
31     DESider A European Effort on Hybrid RANS-LES Modelling:
32     Results of the European-Union Funded Project, 2004 - 2007
33     (Notes on Numerical Fluid Mechanics and Multidisciplinary Design).
34     Chapter 8 Formulation of the Scale-Adaptive Simulation (SAS) Model during
35     the DESIDER Project.
36     F. R. Menter and Y. Egorov.
38 SourceFiles
39     kOmegaSSTSAS.C
41 \*---------------------------------------------------------------------------*/
43 #ifndef kOmegaSSTSAS_H
44 #define kOmegaSSTSAS_H
46 #include "LESModel.H"
47 #include "volFields.H"
48 #include "wallDist.H"
50 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
52 namespace Foam
54 namespace incompressible
56 namespace LESModels
59 /*---------------------------------------------------------------------------*\
60                         Class kOmegaSSTSAS Declaration
61 \*---------------------------------------------------------------------------*/
63 class kOmegaSSTSAS
65     public LESModel
67     // Private member functions
69         //- Update sub-grid scale fields
70         void updateSubGridScaleFields(const volScalarField& D);
72         // Disallow default bitwise copy construct and assignment
73         kOmegaSSTSAS(const kOmegaSSTSAS&);
74         kOmegaSSTSAS& operator=(const kOmegaSSTSAS&);
77 protected:
79     // Protected data
81         // Model constants
83             dimensionedScalar alphaK1_;
84             dimensionedScalar alphaK2_;
86             dimensionedScalar alphaOmega1_;
87             dimensionedScalar alphaOmega2_;
89             dimensionedScalar gamma1_;
90             dimensionedScalar gamma2_;
92             dimensionedScalar beta1_;
93             dimensionedScalar beta2_;
95             dimensionedScalar betaStar_;
97             dimensionedScalar a1_;
98             dimensionedScalar c1_;
99             dimensionedScalar Cs_;
101             dimensionedScalar alphaPhi_;
102             dimensionedScalar zetaTilda2_;
103             dimensionedScalar FSAS_;
105             dimensionedScalar omega0_;
106             dimensionedScalar omegaSmall_;
108             wallDist y_;
109             dimensionedScalar Cmu_;
110             dimensionedScalar kappa_;
113         // Fields
115             volScalarField k_;
116             volScalarField omega_;
117             volScalarField nuSgs_;
120     // Protected member functions
122         tmp<volScalarField> Lvk2
123         (
124             const volScalarField& S2
125         ) const;
127         tmp<volScalarField> F1(const volScalarField& CDkOmega) const;
128         tmp<volScalarField> F2() const;
130         tmp<volScalarField> blend
131         (
132             const volScalarField& F1,
133             const dimensionedScalar& psi1,
134             const dimensionedScalar& psi2
135         ) const
136         {
137             return F1*(psi1 - psi2) + psi2;
138         }
140         tmp<volScalarField> alphaK
141         (
142             const volScalarField& F1
143         ) const
144         {
145             return blend(F1, alphaK1_, alphaK2_);
146         }
148         tmp<volScalarField> alphaOmega
149         (
150              const volScalarField& F1
151         ) const
152         {
153             return blend(F1, alphaOmega1_, alphaOmega2_);
154         }
156         tmp<volScalarField> beta
157         (
158             const volScalarField& F1
159         ) const
160         {
161             return blend(F1, beta1_, beta2_);
162         }
164         tmp<volScalarField> gamma
165         (
166             const volScalarField& F1
167         ) const
168         {
169             return blend(F1, gamma1_, gamma2_);
170         }
173 public:
175     //- Runtime type information
176     TypeName("kOmegaSSTSAS");
179     // Constructors
181         //- Construct from components
182         kOmegaSSTSAS
183         (
184             const volVectorField& U,
185             const surfaceScalarField& phi,
186             transportModel& transport,
187             const word& modelName = typeName
188         );
191     //- Destructor
192     virtual ~kOmegaSSTSAS()
193     {}
196     // Member Functions
198         //- Return SGS kinetic energy
199         virtual tmp<volScalarField> k() const
200         {
201             return k_;
202         }
204         //- Return omega
205         virtual tmp<volScalarField> omega() const
206         {
207             return omega_;
208         }
210         //- Return the effective diffusivity for k
211         tmp<volScalarField> DkEff(const volScalarField& F1) const
212         {
213             return tmp<volScalarField>
214             (
215                  new volScalarField("DkEff", alphaK(F1)*nuSgs_ + nu())
216             );
217         }
219         //- Return the effective diffusivity for omega
220         tmp<volScalarField> DomegaEff(const volScalarField& F1) const
221         {
222             return tmp<volScalarField>
223             (
224                 new volScalarField("DomegaEff", alphaOmega(F1)*nuSgs_ + nu())
225             );
226         }
228         //- Return sub-grid disipation rate
229         virtual tmp<volScalarField> epsilon() const;
231         //- Return SGS viscosity
232         virtual tmp<volScalarField> nuSgs() const
233         {
234             return nuSgs_;
235         }
237         //- Return the sub-grid stress tensor.
238         virtual tmp<volSymmTensorField> B() const;
240         //- Return the effective sub-grid turbulence stress tensor
241         //  including the laminar stress
242         virtual tmp<volSymmTensorField> devBeff() const;
244         //- Return the deviatoric part of the divergence of Beff
245         //  i.e. the additional term in the filtered NSE.
246         virtual tmp<fvVectorMatrix> divDevBeff(volVectorField& U) const;
248         //- Solve the turbulence equations (k-w) and correct the turbulence
249         //  viscosity
250         virtual void correct(const tmp<volTensorField>& gradU);
252         //- Read LESProperties dictionary
253         virtual bool read();
257 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
259 } // End namespace LESModels
260 } // End namespace incompressible
261 } // End namespace Foam
263 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
265 #endif
267 // ************************************************************************* //