kOmegaSSTDES.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2016-2018 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::kOmegaSSTDES
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 See also
46  Foam::kOmegaSST
47 
48 SourceFiles
49  kOmegaSST.C
50 
51 \*---------------------------------------------------------------------------*/
52 
53 #ifndef kOmegaSSTDES_H
54 #define kOmegaSSTDES_H
55 
56 #include "kOmegaSSTBase.H"
57 #include "LESModel.H"
58 #include "LESeddyViscosity.H"
59 
60 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
61 
62 namespace Foam
63 {
64 namespace LESModels
65 {
66 
67 /*---------------------------------------------------------------------------*\
68  Class kOmegaSST Declaration
69 \*---------------------------------------------------------------------------*/
70 
71 template<class BasicTurbulenceModel>
72 class kOmegaSSTDES
73 :
74  public Foam::kOmegaSST
75  <
76  LESeddyViscosity<BasicTurbulenceModel>,
77  BasicTurbulenceModel
78  >
79 {
80 
81 protected:
82 
83  // Protected data
84 
85  // Model constants
86 
87  //- DES coefficient
89 
90  //- Zonal filter choice
91  //
92  // - 0: no filtering
93  // - 1: (1 - F1)
94  // - 2: (1 - F2)
96 
97 
98  // Protected Member Functions
99 
100  //- Return the turbulent length-scale
102 
103  //- The DES dissipation-rate multiplier with options zonal filtering
104  // based on either F1 or F2
106  (
109  ) const;
110 
111  //- Return epsilon/k which for standard RAS is betaStar*omega
113  (
116  ) const;
117 
118 
119 public:
121  typedef typename BasicTurbulenceModel::alphaField alphaField;
122  typedef typename BasicTurbulenceModel::rhoField rhoField;
123  typedef typename BasicTurbulenceModel::transportModel transportModel;
124 
125 
126  //- Runtime type information
127  TypeName("kOmegaSSTDES");
128 
129 
130  // Constructors
131 
132  //- Construct from components
134  (
135  const alphaField& alpha,
136  const rhoField& rho,
137  const volVectorField& U,
138  const surfaceScalarField& alphaRhoPhi,
139  const surfaceScalarField& phi,
140  const transportModel& transport,
141  const word& propertiesName = turbulenceModel::propertiesName,
142  const word& type = typeName
143  );
144 
145 
146  //- Destructor
147  virtual ~kOmegaSSTDES()
148  {}
149 
150 
151  // Member Functions
152 
153  //- Read model coefficients if they have changed
154  virtual bool read();
155 };
156 
157 
158 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
159 
160 } // End namespace LESModels
161 } // End namespace Foam
162 
163 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
164 #ifdef NoRepository
165  #include "kOmegaSSTDES.C"
166 #endif
167 
168 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
169 #endif
170 
171 // ************************************************************************* //
virtual tmp< volScalarField::Internal > epsilonByk(const volScalarField::Internal &F1, const volScalarField::Internal &F2) const
Return epsilon/k which for standard RAS is betaStar*omega.
Definition: kOmegaSSTDES.C:70
surfaceScalarField & phi
TypeName("kOmegaSSTDES")
Runtime type information.
uint8_t direction
Definition: direction.H:45
BasicTurbulenceModel::rhoField rhoField
Definition: kOmegaSSTDES.H:121
DimensionedField< scalar, volMesh > Internal
Type of the internal field from which this GeometricField is derived.
virtual tmp< volScalarField > F1(const volScalarField &CDkOmega) const
direction FSST_
Zonal filter choice.
Definition: kOmegaSSTDES.H:94
BasicTurbulenceModel::transportModel transportModel
Definition: kOmegaSSTDES.H:122
Implementation of the k-omega-SST turbulence model for incompressible and compressible flows...
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from string.
Definition: word.H:59
virtual ~kOmegaSSTDES()
Destructor.
Definition: kOmegaSSTDES.H:146
Implementation of the k-omega-SST-DES turbulence model for incompressible and compressible flows...
Definition: kOmegaSSTDES.H:71
virtual tmp< volScalarField::Internal > FDES(const volScalarField::Internal &F1, const volScalarField::Internal &F2) const
The DES dissipation-rate multiplier with options zonal filtering.
Definition: kOmegaSSTDES.C:46
U
Definition: pEqn.H:72
tmp< volScalarField::Internal > Lt() const
Return the turbulent length-scale.
Definition: kOmegaSSTDES.C:38
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
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
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
A class for managing temporary objects.
Definition: PtrList.H:53
dimensionedScalar CDES_
DES coefficient.
Definition: kOmegaSSTDES.H:87
Namespace for OpenFOAM.
virtual bool read()
Read model coefficients if they have changed.
Definition: kOmegaSSTDES.C:131
BasicTurbulenceModel::alphaField alphaField
Definition: kOmegaSSTDES.H:120