initial commit for version 1.6.x patch release
[OpenFOAM-1.6.x.git] / src / turbulenceModels / incompressible / LES / kOmegaSSTSAS / kOmegaSSTSAS.H
blob05632f924259b7536c22277961971970cc27c114
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2008-2009 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::LESmodels::kOmegaSSTSAS
28 Description
29     kOmegaSSTSAS LES turbulence model for incompressible flows
31     Note: does not have an explicit dependency on spatial discretisation
32           i.e. LESdelta not used.
34 SourceFiles
35     kOmegaSSTSAS.C
37 \*---------------------------------------------------------------------------*/
39 #ifndef kOmegaSSTSAS_H
40 #define kOmegaSSTSAS_H
42 #include "LESModel.H"
43 #include "volFields.H"
44 #include "wallDist.H"
46 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
48 namespace Foam
50 namespace incompressible
52 namespace LESModels
55 /*---------------------------------------------------------------------------*\
56                         Class kOmegaSSTSAS Declaration
57 \*---------------------------------------------------------------------------*/
59 class kOmegaSSTSAS
61     public LESModel
63     // Private member functions
65         //- Update sub-grid scale fields
66         void updateSubGridScaleFields(const volScalarField& D);
68         // Disallow default bitwise copy construct and assignment
69         kOmegaSSTSAS(const kOmegaSSTSAS&);
70         kOmegaSSTSAS& operator=(const kOmegaSSTSAS&);
73 protected:
75     // Protected data
77         // Model constants
79             dimensionedScalar alphaK1_;
80             dimensionedScalar alphaK2_;
82             dimensionedScalar alphaOmega1_;
83             dimensionedScalar alphaOmega2_;
85             dimensionedScalar gamma1_;
86             dimensionedScalar gamma2_;
88             dimensionedScalar beta1_;
89             dimensionedScalar beta2_;
91             dimensionedScalar betaStar_;
93             dimensionedScalar a1_;
94             dimensionedScalar c1_;
96             dimensionedScalar alphaPhi_;
97             dimensionedScalar zetaTilda2_;
98             dimensionedScalar FSAS_;
100             dimensionedScalar omega0_;
101             dimensionedScalar omegaSmall_;
103             wallDist y_;
104             dimensionedScalar Cmu_;
105             dimensionedScalar kappa_;
108         // Fields
110             volScalarField k_;
111             volScalarField omega_;
112             volScalarField nuSgs_;
115     // Protected member functions
117         tmp<volScalarField> Lvk2
118         (
119             const volScalarField& S2
120         ) const;
122         tmp<volScalarField> F1(const volScalarField& CDkOmega) const;
123         tmp<volScalarField> F2() const;
125         tmp<volScalarField> blend
126         (
127             const volScalarField& F1,
128             const dimensionedScalar& psi1,
129             const dimensionedScalar& psi2
130         ) const
131         {
132             return F1*(psi1 - psi2) + psi2;
133         }
135         tmp<volScalarField> alphaK
136         (
137             const volScalarField& F1
138         ) const
139         {
140             return blend(F1, alphaK1_, alphaK2_);
141         }
143         tmp<volScalarField> alphaOmega
144         (
145              const volScalarField& F1
146         ) const
147         {
148             return blend(F1, alphaOmega1_, alphaOmega2_);
149         }
151         tmp<volScalarField> beta
152         (
153             const volScalarField& F1
154         ) const
155         {
156             return blend(F1, beta1_, beta2_);
157         }
159         tmp<volScalarField> gamma
160         (
161             const volScalarField& F1
162         ) const
163         {
164             return blend(F1, gamma1_, gamma2_);
165         }
168 public:
170     //- Runtime type information
171     TypeName("kOmegaSSTSAS");
174     // Constructors
176         //- Construct from components
177         kOmegaSSTSAS
178         (
179             const volVectorField& U,
180             const surfaceScalarField& phi,
181             transportModel& transport,
182             const word& modelName = typeName
183         );
186     //- Destructor
187     virtual ~kOmegaSSTSAS()
188     {}
191     // Member Functions
193         //- Return SGS kinetic energy
194         virtual tmp<volScalarField> k() const
195         {
196             return k_;
197         }
199         //- Return omega
200         virtual tmp<volScalarField> omega() const
201         {
202             return omega_;
203         }
205         //- Return the effective diffusivity for k
206         tmp<volScalarField> DkEff(const volScalarField& F1) const
207         {
208             return tmp<volScalarField>
209             (
210                  new volScalarField("DkEff", alphaK(F1)*nuSgs_ + nu())
211             );
212         }
214         //- Return the effective diffusivity for omega
215         tmp<volScalarField> DomegaEff(const volScalarField& F1) const
216         {
217             return tmp<volScalarField>
218             (
219                 new volScalarField("DomegaEff", alphaOmega(F1)*nuSgs_ + nu())
220             );
221         }
223         //- Return sub-grid disipation rate
224         virtual tmp<volScalarField> epsilon() const;
226         //- Return SGS viscosity
227         virtual tmp<volScalarField> nuSgs() const
228         {
229             return nuSgs_;
230         }
232         //- Return the sub-grid stress tensor.
233         virtual tmp<volSymmTensorField> B() const;
235         //- Return the effective sub-grid turbulence stress tensor
236         //  including the laminar stress
237         virtual tmp<volSymmTensorField> devBeff() const;
239         //- Return the deviatoric part of the divergence of Beff
240         //  i.e. the additional term in the filtered NSE.
241         virtual tmp<fvVectorMatrix> divDevBeff(volVectorField& U) const;
243         //- Solve the turbulence equations (k-w) and correct the turbulence
244         //  viscosity
245         virtual void correct(const tmp<volTensorField>& gradU);
247         //- Read LESProperties dictionary
248         virtual bool read();
252 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
254 } // End namespace LESModels
255 } // End namespace incompressible
256 } // End namespace Foam
258 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
260 #endif
262 // ************************************************************************* //