Loading...
Searching...
No Matches
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-2022 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
35 incompressible and compressible flows.
36
37 Reference:
38 \verbatim
39 Gritskevich, M. S., Garbaruk, A. V.,
40 Schütze, J., & Menter, F. R. (2012).
41 Development of DDES and IDDES formulations for
42 the k-ω shear stress transport model.
43 Flow, turbulence and combustion, 88(3), 431-449.
44 DOI:10.1007/s10494-011-9378-4
45 \endverbatim
46
47SourceFiles
48 kOmegaSSTIDDES.C
49
50\*---------------------------------------------------------------------------*/
51
52#ifndef Foam_kOmegaSSTIDDES_H
53#define Foam_kOmegaSSTIDDES_H
54
55#include "kOmegaSSTDES.H"
56
57// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
58
59namespace Foam
60{
61namespace LESModels
62{
63
64/*---------------------------------------------------------------------------*\
65 class kOmegaSSTIDDES Declaration
66\*---------------------------------------------------------------------------*/
67
68template<class BasicTurbulenceModel>
69class kOmegaSSTIDDES
70:
71 public kOmegaSSTDES<BasicTurbulenceModel>
72{
73 // Private Member Functions
74
75 //- Check that the supplied delta is an IDDESDelta
76 const IDDESDelta& setDelta() const;
77
79 tmp<volScalarField> ft(const volScalarField& magGradU) const;
80 tmp<volScalarField> fl(const volScalarField& magGradU) const;
81
82 //- Delay function
83 tmp<volScalarField> fdt(const volScalarField& magGradU) const;
84
85 //- No copy construct
86 kOmegaSSTIDDES(const kOmegaSSTIDDES&) = delete;
87
88 //- No copy assignment
89 void operator=(const kOmegaSSTIDDES&) = delete;
90
91
92protected:
93
94 // Protected Data
95
96 // Model coefficients
97
105 //- IDDES delta
106 const IDDESDelta& IDDESDelta_;
107
108
109 // Protected Member Functions
110
111 //- Return the length scale
113 (
114 const volScalarField& magGradU,
115 const volScalarField& CDES
116 ) const;
117
118
119public:
120
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("kOmegaSSTIDDES");
130 // Constructors
131
132 //- Construct from components
133 kOmegaSSTIDDES
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 );
145
146 //- Destructor
147 virtual ~kOmegaSSTIDDES() = default;
148
149
150 // Member Functions
151
152 //- Re-read model coefficients if they have changed
153 virtual bool read();
154
155 //- Return the shielding function
156 virtual tmp<volScalarField> fd() const;
157};
158
159
160// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
161
162} // End namespace LESModels
163} // End namespace Foam
164
165// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
166
167#ifdef NoRepository
169#endif
170
171// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
172
173#endif
174
175// ************************************************************************* //
virtual tmp< volScalarField > dTilda(const volScalarField &magGradU, const volScalarField &CDES) const
Return the length scale.
Foam::EddyDiffusivity< Foam::fluidThermoCompressibleTurbulenceModel >::alphaField alphaField
Foam::EddyDiffusivity< Foam::fluidThermoCompressibleTurbulenceModel >::rhoField rhoField
virtual tmp< volScalarField > fd() const
Return the shielding function.
kOmegaSSTIDDES(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.
virtual ~kOmegaSSTIDDES()=default
Destructor.
Foam::EddyDiffusivity< Foam::fluidThermoCompressibleTurbulenceModel >::transportModel transportModel
TypeName("kOmegaSSTIDDES")
Runtime type information.
virtual bool read()
Re-read model coefficients if they have changed.
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition Switch.H:81
A class for managing temporary objects.
Definition tmp.H:75
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from Foam::string.
Definition word.H:66
U
Definition pEqn.H:72
Namespace for LES SGS models.
Namespace for OpenFOAM.
GeometricField< vector, fvPatchField, volMesh > volVectorField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
volScalarField & alpha
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition typeInfo.H:68