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 | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2015 OpenFOAM Foundation
9 Copyright (C) 2016-2019 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27Class
28 Foam::LESModels::kOmegaSSTDES
29
30Group
31 grpDESTurbulence
32
33Description
34 k-omega-SST DES turbulence model for incompressible and compressible flows
35
36 Reference:
37 \verbatim
38 Strelets, M. (2001)
39 Detached Eddy Simulation of Massively Separated Flows,
40 39th AIAA Aerospace Sciences Meeting and Exhibit, Reno, NV
41 \endverbatim
42
43Note
44 The default values of the DES constants implemented are code-specific
45 values calibrated for OpenFOAM using decaying isotropic turbulence, and
46 hence deviate slightly from the values suggested in the reference
47 publication.
48
49SourceFiles
50 kOmegaSSTDES.C
51
52\*---------------------------------------------------------------------------*/
53
54#ifndef kOmegaSSTDES_H
55#define kOmegaSSTDES_H
56
57#include "DESModel.H"
58#include "kOmegaSSTBase.H"
59
60// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
61
62namespace Foam
63{
64namespace LESModels
65{
66
67/*---------------------------------------------------------------------------*\
68 class kOmegaSSTDES Declaration
69\*---------------------------------------------------------------------------*/
70
71template<class BasicTurbulenceModel>
72class kOmegaSSTDES
73:
74 public kOmegaSSTBase<DESModel<BasicTurbulenceModel>>
75{
76 // Private Member Functions
77
78 //- No copy construct
79 kOmegaSSTDES(const kOmegaSSTDES&) = delete;
80
81 //- No copy assignment
82 void operator=(const kOmegaSSTDES&) = delete;
83
84
85protected:
86
87 // Protected data
88
89 // Model coefficients
94
95
96 // Protected Member Functions
97
98 //- Blending for CDES parameter
99 virtual tmp<volScalarField> CDES(const volScalarField& F1) const
100 {
101 return this->blend(F1, CDESkom_, CDESkeps_);
102 }
103
104 virtual void correctNut(const volScalarField& S2);
105 virtual void correctNut();
106
107 //- Length scale
109 (
110 const volScalarField& magGradU,
111 const volScalarField& CDES
112 ) const;
113
114 //- Return epsilon/k
116 (
117 const volScalarField& F1,
118 const volTensorField& gradU
119 ) const;
120
121 //- Return G/nu
123 (
124 const volScalarField::Internal& GbyNu0,
127 ) const;
128
129
130public:
132 typedef typename BasicTurbulenceModel::alphaField alphaField;
133 typedef typename BasicTurbulenceModel::rhoField rhoField;
134 typedef typename BasicTurbulenceModel::transportModel transportModel;
135
136
137 //- Runtime type information
138 TypeName("kOmegaSSTDES");
139
140
141 // Constructors
142
143 //- Construct from components
145 (
146 const alphaField& alpha,
147 const rhoField& rho,
148 const volVectorField& U,
149 const surfaceScalarField& alphaRhoPhi,
150 const surfaceScalarField& phi,
151 const transportModel& transport,
152 const word& propertiesName = turbulenceModel::propertiesName,
153 const word& type = typeName
154 );
155
156
157 //- Destructor
158 virtual ~kOmegaSSTDES() = default;
159
160
161 // Member Functions
162
163 //- Re-read model coefficients if they have changed
164 virtual bool read();
165
166 //- Return the LES field indicator
167 virtual tmp<volScalarField> LESRegion() const;
168};
169
170
171// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
172
173} // End namespace LESModels
174} // End namespace Foam
175
176// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
177#ifdef NoRepository
178 #include "kOmegaSSTDES.C"
179#endif
180
181// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
182#endif
183
184// ************************************************************************* //
#define F1(B, C, D)
Definition: SHA1.C:150
surfaceScalarField & phi
k-omega-SST DES turbulence model for incompressible and compressible flows
Definition: kOmegaSSTDES.H:74
virtual tmp< volScalarField > dTilda(const volScalarField &magGradU, const volScalarField &CDES) const
Length scale.
Definition: kOmegaSSTDES.C:60
BasicTurbulenceModel::alphaField alphaField
Definition: kOmegaSSTDES.H:131
TypeName("kOmegaSSTDES")
Runtime type information.
BasicTurbulenceModel::rhoField rhoField
Definition: kOmegaSSTDES.H:132
dimensionedScalar kappa_
Definition: kOmegaSSTDES.H:90
virtual ~kOmegaSSTDES()=default
Destructor.
dimensionedScalar CDESkeps_
Definition: kOmegaSSTDES.H:92
virtual tmp< volScalarField > LESRegion() const
Return the LES field indicator.
Definition: kOmegaSSTDES.C:177
virtual tmp< volScalarField::Internal > GbyNu(const volScalarField::Internal &GbyNu0, const volScalarField::Internal &F2, const volScalarField::Internal &S2) const
Return G/nu.
Definition: kOmegaSSTDES.C:86
virtual tmp< volScalarField > CDES(const volScalarField &F1) const
Blending for CDES parameter.
Definition: kOmegaSSTDES.H:98
dimensionedScalar CDESkom_
Definition: kOmegaSSTDES.H:91
BasicTurbulenceModel::transportModel transportModel
Definition: kOmegaSSTDES.H:133
virtual tmp< volScalarField::Internal > epsilonByk(const volScalarField &F1, const volTensorField &gradU) const
Return epsilon/k.
Definition: kOmegaSSTDES.C:74
virtual bool read()
Re-read model coefficients if they have changed.
Definition: kOmegaSSTDES.C:161
Base class implementation of the k-omega-SST turbulence model for incompressible and compressible flo...
virtual tmp< volScalarField > F2() const
Definition: kOmegaSSTBase.C:72
tmp< volScalarField > blend(const volScalarField &F1, const dimensionedScalar &psi1, const dimensionedScalar &psi2) const
A class for managing temporary objects.
Definition: tmp.H:65
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from Foam::string.
Definition: word.H:68
U
Definition: pEqn.H:72
Namespace for OpenFOAM.
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:598
volScalarField & alpha
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition: typeInfo.H:73