kOmegaSSTDES.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2016 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
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.
15 
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.
20 
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/>.
23 
24 Class
25  Foam::LESModels::kOmegaSST
26 
27 Description
28  Implementation of the k-omega-SST-DES turbulence model for
29  incompressible and compressible flows.
30 
31  DES model described in:
32  \verbatim
33  Menter, F. R., Kuntz, M., and Langtry, R. (2003).
34  Ten Years of Industrial Experience with the SST Turbulence Model.
35  Turbulence, Heat and Mass Transfer 4, ed: K. Hanjalic, Y. Nagano,
36  & M. Tummers, Begell House, Inc., 625 - 632.
37  \endverbatim
38 
39  Optional support for zonal filtering based on F1 or F2 is provided as
40  described in the paper.
41 
42  For further details of the implementation of the base k-omega-SST model
43  see Foam::kOmegaSST.
44 
45 Group
46  grpLESTurbulence
47 
48 See also
49  Foam::kOmegaSST
50 
51 SourceFiles
52  kOmegaSST.C
53 
54 \*---------------------------------------------------------------------------*/
55 
56 #ifndef kOmegaSSTDES_H
57 #define kOmegaSSTDES_H
58 
59 #include "kOmegaSSTBase.H"
60 #include "LESModel.H"
61 #include "LESeddyViscosity.H"
62 
63 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
64 
65 namespace Foam
66 {
67 namespace LESModels
68 {
69 
70 /*---------------------------------------------------------------------------*\
71  Class kOmegaSST Declaration
72 \*---------------------------------------------------------------------------*/
73 
74 template<class BasicTurbulenceModel>
75 class kOmegaSSTDES
76 :
77  public Foam::kOmegaSST
78  <
79  LESeddyViscosity<BasicTurbulenceModel>,
80  BasicTurbulenceModel
81  >
82 {
83 
84 protected:
85 
86  // Protected data
87 
88  // Model constants
89 
90  //- DES coefficient
92 
93  //- Zonal filter choice
94  //
95  // - 0: no filtering
96  // - 1: (1 - F1)
97  // - 2: (1 - F2)
99 
100 
101  // Protected Member Functions
102 
103  //- Return the turbulent length-scale
104  tmp<volScalarField> Lt() const;
105 
106  //- The DES dissipation-rate multiplier with options zonal filtering
107  // based on either F1 or F2
108  virtual tmp<volScalarField> FDES
109  (
110  const volScalarField& F1,
111  const volScalarField& F2
112  ) const;
113 
114  //- Return epsilon/k which for standard RAS is betaStar*omega
116  (
117  const volScalarField& F1,
118  const volScalarField& F2
119  ) const;
120 
121 
122 public:
124  typedef typename BasicTurbulenceModel::alphaField alphaField;
125  typedef typename BasicTurbulenceModel::rhoField rhoField;
126  typedef typename BasicTurbulenceModel::transportModel transportModel;
127 
128 
129  //- Runtime type information
130  TypeName("kOmegaSSTDES");
131 
132 
133  // Constructors
134 
135  //- Construct from components
137  (
138  const alphaField& alpha,
139  const rhoField& rho,
140  const volVectorField& U,
141  const surfaceScalarField& alphaRhoPhi,
142  const surfaceScalarField& phi,
143  const transportModel& transport,
144  const word& propertiesName = turbulenceModel::propertiesName,
145  const word& type = typeName
146  );
147 
148 
149  //- Destructor
150  virtual ~kOmegaSSTDES()
151  {}
152 
153 
154  // Member Functions
155 
156  //- Read model coefficients if they have changed
157  virtual bool read();
158 };
159 
160 
161 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
162 
163 } // End namespace LESModels
164 } // End namespace Foam
165 
166 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
167 #ifdef NoRepository
168  #include "kOmegaSSTDES.C"
169 #endif
170 
171 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
172 #endif
173 
174 // ************************************************************************* //
surfaceScalarField & phi
U
Definition: pEqn.H:83
uint8_t direction
Definition: direction.H:46
TypeName("kOmegaSSTDES")
Runtime type information.
BasicTurbulenceModel::rhoField rhoField
Definition: kOmegaSSTDES.H:124
direction FSST_
Zonal filter choice.
Definition: kOmegaSSTDES.H:97
BasicTurbulenceModel::transportModel transportModel
Definition: kOmegaSSTDES.H:125
virtual tmp< volScalarField > epsilonByk(const volScalarField &F1, const volScalarField &F2) const
Return epsilon/k which for standard RAS is betaStar*omega.
Definition: kOmegaSSTDES.C:70
Implementation of the k-omega-SST turbulence model for incompressible and compressible flows...
tmp< volScalarField > F1(const volScalarField &CDkOmega) const
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from string.
Definition: word.H:59
tmp< volScalarField > Lt() const
Return the turbulent length-scale.
Definition: kOmegaSSTDES.C:38
virtual tmp< volScalarField > FDES(const volScalarField &F1, const volScalarField &F2) const
The DES dissipation-rate multiplier with options zonal filtering.
Definition: kOmegaSSTDES.C:46
virtual ~kOmegaSSTDES()
Destructor.
Definition: kOmegaSSTDES.H:149
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:461
kOmegaSSTDES(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName=turbulenceModel::propertiesName, const word &type=typeName)
Construct from components.
Definition: kOmegaSSTDES.C:83
A class for managing temporary objects.
Definition: PtrList.H:54
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
dimensionedScalar CDES_
DES coefficient.
Definition: kOmegaSSTDES.H:90
Namespace for OpenFOAM.
virtual bool read()
Read model coefficients if they have changed.
Definition: kOmegaSSTDES.C:131
BasicTurbulenceModel::alphaField alphaField
Definition: kOmegaSSTDES.H:123