kOmegaSSTIDDES.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) 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::kOmegaSSTIDDES
29
30Group
31 grpDESTurbulence
32
33Description
34 k-omega-SST IDDES turbulence model for incompressible and compressible
35 flows
36
37 Reference:
38 \verbatim
39 Gritskevich, M.S., Garbaruk, A.V., Schuetze, J., Menter, F.R. (2011)
40 Development of DDES and IDDES Formulations for the k-omega
41 Shear Stress Transport Model, Flow, Turbulence and Combustion,
42 pp. 1-19
43 \endverbatim
44
45SourceFiles
46 kOmegaSSTIDDES.C
47
48\*---------------------------------------------------------------------------*/
49
50#ifndef kOmegaSSTIDDES_H
51#define kOmegaSSTIDDES_H
52
53#include "kOmegaSSTDES.H"
54
55// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
56
57namespace Foam
58{
59namespace LESModels
60{
61
62/*---------------------------------------------------------------------------*\
63 class kOmegaSSTIDDES Declaration
64\*---------------------------------------------------------------------------*/
65
66template<class BasicTurbulenceModel>
68:
69 public kOmegaSSTDES<BasicTurbulenceModel>
70{
71 // Private Member Functions
72
73 //- Check that the supplied delta is an IDDESDelta
74 const IDDESDelta& setDelta() const;
75
76 tmp<volScalarField> alpha() const;
77 tmp<volScalarField> ft(const volScalarField& magGradU) const;
78 tmp<volScalarField> fl(const volScalarField& magGradU) const;
79
81 (
82 const volScalarField& nur,
83 const volScalarField& magGradU
84 ) const;
85
86 //- Delay function
87 tmp<volScalarField> fdt(const volScalarField& magGradU) const;
88
89 //- No copy construct
90 kOmegaSSTIDDES(const kOmegaSSTIDDES&) = delete;
91
92 //- No copy assignment
93 void operator=(const kOmegaSSTIDDES&) = delete;
94
95
96protected:
97
98 // Protected data
99
100 // Model coefficients
106
107 // Fields
109 const IDDESDelta& IDDESDelta_;
110
111 // Protected Member Functions
112
113 //- Length scale
115 (
116 const volScalarField& magGradU,
117 const volScalarField& CDES
118 ) const;
119
120
121public:
123 typedef typename BasicTurbulenceModel::alphaField alphaField;
124 typedef typename BasicTurbulenceModel::rhoField rhoField;
125 typedef typename BasicTurbulenceModel::transportModel transportModel;
126
127
128 //- Runtime type information
129 TypeName("kOmegaSSTIDDES");
130
131
132 // Constructors
133
134 //- Construct from components
136 (
137 const alphaField& alpha,
138 const rhoField& rho,
139 const volVectorField& U,
140 const surfaceScalarField& alphaRhoPhi,
141 const surfaceScalarField& phi,
142 const transportModel& transport,
143 const word& propertiesName = turbulenceModel::propertiesName,
144 const word& type = typeName
145 );
146
147
148 //- Destructor
149 virtual ~kOmegaSSTIDDES() = default;
150
151
152 // Member Functions
153
154 //- Re-read model coefficients if they have changed
155 virtual bool read();
156};
157
158
159// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
160
161} // End namespace LESModels
162} // End namespace Foam
163
164// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
165
166#ifdef NoRepository
167 #include "kOmegaSSTIDDES.C"
168#endif
169
170// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
171
172#endif
173
174// ************************************************************************* //
surfaceScalarField & phi
k-omega-SST DES turbulence model for incompressible and compressible flows
Definition: kOmegaSSTDES.H:74
virtual tmp< volScalarField > CDES(const volScalarField &F1) const
Blending for CDES parameter.
Definition: kOmegaSSTDES.H:98
k-omega-SST IDDES turbulence model for incompressible and compressible flows
virtual tmp< volScalarField > dTilda(const volScalarField &magGradU, const volScalarField &CDES) const
Length scale.
BasicTurbulenceModel::alphaField alphaField
BasicTurbulenceModel::rhoField rhoField
virtual ~kOmegaSSTIDDES()=default
Destructor.
BasicTurbulenceModel::transportModel transportModel
TypeName("kOmegaSSTIDDES")
Runtime type information.
virtual bool read()
Re-read model coefficients if they have changed.
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
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition: typeInfo.H:73