Forward compatibility: flex
[foam-extend-3.2.git] / src / turbulenceModels / incompressible / RAS / coupledKOmegaSST / coupledKOmegaSST.H
blob14e36bb840774d857cc7dac13c245adee5f14fbf
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::incompressible::RASModels::coupledKOmegaSST
27 Description
28     Implementation of the k-omega-SST turbulence model for incompressible
29     flows using a block-coupled solution.
31     Turbulence model described in:
32     @verbatim
33         Menter, F., Esch, T.
34         "Elements of Industrial Heat Transfer Prediction"
35         16th Brazilian Congress of Mechanical Engineering (COBEM),
36         Nov. 2001
37     @endverbatim
39     Note that this implementation is written in terms of alpha diffusion
40     coefficients rather than the more traditional sigma (alpha = 1/sigma) so
41     that the blending can be applied to all coefficuients in a consistent
42     manner.  The paper suggests that sigma is blended but this would not be
43     consistent with the blending of the k-epsilon and k-omega models.
45     Also note that the error in the last term of equation (2) relating to
46     sigma has been corrected.
48     Wall-functions are applied in this implementation by using equations (14)
49     to specify the near-wall omega as appropriate.
51     The blending functions (15) and (16) are not currently used because of the
52     uncertainty in their origin, range of applicability and that is y+ becomes
53     sufficiently small blending u_tau in this manner clearly becomes nonsense.
55     The default model coefficients correspond to the following:
56     @verbatim
57         coupledKOmegaSSTCoeffs
58         {
59             alphaK1     0.85034;
60             alphaK2     1.0;
61             alphaOmega1 0.5;
62             alphaOmega2 0.85616;
63             beta1       0.075;
64             beta2       0.0828;
65             betaStar    0.09;
66             gamma1      0.5532;
67             gamma2      0.4403;
68             a1          0.31;
69             c1          10.0;
70         }
71     @endverbatim
73 Author
74     Hrvoje Jasak, Wikki Ltd.  All rights reserved.
76 SourceFiles
77     coupledKOmegaSST.C
79 \*---------------------------------------------------------------------------*/
81 #ifndef coupledKOmegaSST_H
82 #define coupledKOmegaSST_H
84 #include "RASModel.H"
85 #include "wallDist.H"
87 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
89 namespace Foam
91 namespace incompressible
93 namespace RASModels
96 /*---------------------------------------------------------------------------*\
97                            Class kOmega Declaration
98 \*---------------------------------------------------------------------------*/
100 class coupledKOmegaSST
102     public RASModel
104     // Private data
106         // Model coefficients
107             dimensionedScalar alphaK1_;
108             dimensionedScalar alphaK2_;
110             dimensionedScalar alphaOmega1_;
111             dimensionedScalar alphaOmega2_;
113             dimensionedScalar gamma1_;
114             dimensionedScalar gamma2_;
116             dimensionedScalar beta1_;
117             dimensionedScalar beta2_;
119             dimensionedScalar betaStar_;
121             dimensionedScalar a1_;
122             dimensionedScalar c1_;
125         //- Wall distance field
126         //  Note: different to wall distance in parent RASModel
127         wallDist y_;
130         // Fields
132             volScalarField k_;
133             volScalarField omega_;
134             volScalarField nut_;
136        // Block vector field for k-omega
137        volVector2Field kOmega_;
140     // Private member functions
142         tmp<volScalarField> F1(const volScalarField& CDkOmega) const;
143         tmp<volScalarField> F2() const;
145         tmp<volScalarField> blend
146         (
147             const volScalarField& F1,
148             const dimensionedScalar& psi1,
149             const dimensionedScalar& psi2
150         ) const
151         {
152             return F1*(psi1 - psi2) + psi2;
153         }
155         tmp<volScalarField> alphaK
156         (
157             const volScalarField& F1
158         ) const
159         {
160             return blend(F1, alphaK1_, alphaK2_);
161         }
163         tmp<volScalarField> alphaOmega
164         (
165             const volScalarField& F1
166         ) const
167         {
168             return blend(F1, alphaOmega1_, alphaOmega2_);
169         }
171         tmp<volScalarField> beta
172         (
173             const volScalarField& F1
174         ) const
175         {
176             return blend(F1, beta1_, beta2_);
177         }
179         tmp<volScalarField> gamma
180         (
181             const volScalarField& F1
182         ) const
183         {
184             return blend(F1, gamma1_, gamma2_);
185         }
188 public:
190     //- Runtime type information
191     TypeName("coupledKOmegaSST");
194     // Constructors
196         //- Construct from components
197         coupledKOmegaSST
198         (
199             const volVectorField& U,
200             const surfaceScalarField& phi,
201             transportModel& transport
202         );
205     //- Destructor
207         virtual ~coupledKOmegaSST()
208         {}
211     // Member Functions
213         //- Return the turbulence viscosity
214         virtual tmp<volScalarField> nut() const
215         {
216             return nut_;
217         }
219         //- Return the effective diffusivity for k
220         tmp<volScalarField> DkEff(const volScalarField& F1) const
221         {
222             return tmp<volScalarField>
223             (
224                 new volScalarField("DkEff", alphaK(F1)*nut_ + nu())
225             );
226         }
228         //- Return the effective diffusivity for omega
229         tmp<volScalarField> DomegaEff(const volScalarField& F1) const
230         {
231             return tmp<volScalarField>
232             (
233                 new volScalarField("DomegaEff", alphaOmega(F1)*nut_ + nu())
234             );
235         }
237         //- Return the turbulence kinetic energy
238         virtual tmp<volScalarField> k() const
239         {
240             return k_;
241         }
243         //- Return the turbulence specific dissipation rate
244         virtual tmp<volScalarField> omega() const
245         {
246             return omega_;
247         }
249         //- Return the turbulence kinetic energy dissipation rate
250         virtual tmp<volScalarField> epsilon() const
251         {
252             return tmp<volScalarField>
253             (
254                 new volScalarField
255                 (
256                     IOobject
257                     (
258                         "epsilon",
259                         mesh_.time().timeName(),
260                         mesh_
261                     ),
262                     betaStar_*k_*omega_,
263                     omega_.boundaryField().types()
264                 )
265             );
266         }
268         //- Return the Reynolds stress tensor
269         virtual tmp<volSymmTensorField> R() const;
271         //- Return the effective stress tensor including the laminar stress
272         virtual tmp<volSymmTensorField> devReff() const;
274         //- Return the source term for the momentum equation
275         virtual tmp<fvVectorMatrix> divDevReff(volVectorField& U) const;
277         //- Solve the turbulence equations and correct the turbulence viscosity
278         virtual void correct();
280         //- Read RASProperties dictionary
281         virtual bool read();
285 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
287 } // End namespace RASModels
288 } // namespace incompressible
289 } // End namespace Foam
291 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
293 #endif
295 // ************************************************************************* //