kOmegaSSTDDES.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::kOmegaSSTDDES
29
30Group
31 grpDESTurbulence
32
33Description
34 k-omega-SST DDES turbulence model for incompressible and compressible flows
35
36 Reference:
37 \verbatim
38 Gritskevich, M.S., Garbaruk, A.V., Schuetze, J., Menter, F.R. (2011)
39 Development of DDES and IDDES Formulations for the k-omega
40 Shear Stress Transport Model, Flow, Turbulence and Combustion,
41 pp. 1-19
42 \endverbatim
43
44SourceFiles
45 kOmegaSSTDDES.C
46
47\*---------------------------------------------------------------------------*/
48
49#ifndef kOmegaSSTDDES_H
50#define kOmegaSSTDDES_H
51
52#include "kOmegaSSTDES.H"
53
54// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
55
56namespace Foam
57{
58namespace LESModels
59{
60
61/*---------------------------------------------------------------------------*\
62 Class kOmegaSSTDDES Declaration
63\*---------------------------------------------------------------------------*/
64
65template<class BasicTurbulenceModel>
66class kOmegaSSTDDES
67:
68 public kOmegaSSTDES<BasicTurbulenceModel>
69{
70 // Private Member Functions
71
72 tmp<volScalarField> fd(const volScalarField& magGradU) const;
73
74 tmp<volScalarField> rd(const volScalarField& magGradU) const;
75
76 //- No copy construct
77 kOmegaSSTDDES(const kOmegaSSTDDES&) = delete;
78
79 //- No copy assignment
80 void operator=(const kOmegaSSTDDES&) = delete;
81
82
83protected:
84
85 // Protected data
86
87 // Model coefficients
91
92
93 // Protected Member Functions
94
95 //- Length scale
97 (
98 const volScalarField& magGradU,
99 const volScalarField& CDES
100 ) const;
101
102
103public:
105 typedef typename BasicTurbulenceModel::alphaField alphaField;
106 typedef typename BasicTurbulenceModel::rhoField rhoField;
107 typedef typename BasicTurbulenceModel::transportModel transportModel;
108
109
110 //- Runtime type information
111 TypeName("kOmegaSSTDDES");
112
113
114 // Constructors
115
116 //- Construct from components
118 (
119 const alphaField& alpha,
120 const rhoField& rho,
121 const volVectorField& U,
122 const surfaceScalarField& alphaRhoPhi,
123 const surfaceScalarField& phi,
124 const transportModel& transport,
125 const word& propertiesName = turbulenceModel::propertiesName,
126 const word& type = typeName
127 );
128
129
130 //- Destructor
131 virtual ~kOmegaSSTDDES() = default;
132
133
134 // Member Functions
135
136 //- Re-read model coefficients if they have changed
137 virtual bool read();
138};
139
140
141// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
142
143} // End namespace LESModels
144} // End namespace Foam
145
146// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
147
148#ifdef NoRepository
149# include "kOmegaSSTDDES.C"
150#endif
151
152// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
153
154#endif
155
156// ************************************************************************* //
surfaceScalarField & phi
k-omega-SST DDES turbulence model for incompressible and compressible flows
Definition: kOmegaSSTDDES.H:68
virtual tmp< volScalarField > dTilda(const volScalarField &magGradU, const volScalarField &CDES) const
Length scale.
Definition: kOmegaSSTDDES.C:82
BasicTurbulenceModel::alphaField alphaField
BasicTurbulenceModel::rhoField rhoField
virtual ~kOmegaSSTDDES()=default
Destructor.
BasicTurbulenceModel::transportModel transportModel
TypeName("kOmegaSSTDDES")
Runtime type information.
virtual bool read()
Re-read model coefficients if they have changed.
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
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