Loading...
Searching...
No Matches
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) 2022 Upstream CFD GmbH
10 Copyright (C) 2019, 2022 OpenCFD Ltd.
11-------------------------------------------------------------------------------
12License
13 This file is part of OpenFOAM.
14
15 OpenFOAM is free software: you can redistribute it and/or modify it
16 under the terms of the GNU General Public License as published by
17 the Free Software Foundation, either version 3 of the License, or
18 (at your option) any later version.
19
20 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
21 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
22 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
23 for more details.
24
25 You should have received a copy of the GNU General Public License
26 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
27
28Class
29 Foam::LESModels::kOmegaSSTDDES
30
31Group
32 grpDESTurbulence
33
34Description
35 k-omega-SST DDES turbulence model for 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 kOmegaSSTDDES.C
49
50\*---------------------------------------------------------------------------*/
51
52#ifndef Foam_kOmegaSSTDDES_H
53#define Foam_kOmegaSSTDDES_H
54
55#include "kOmegaSSTDES.H"
56
57// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
58
59namespace Foam
60{
61namespace LESModels
62{
63
64/*---------------------------------------------------------------------------*\
65 Class kOmegaSSTDDES Declaration
66\*---------------------------------------------------------------------------*/
67
68template<class BasicTurbulenceModel>
69class kOmegaSSTDDES
70:
71 public kOmegaSSTDES<BasicTurbulenceModel>
72{
73 // Private Member Functions
74
75 //- Return the shielding function field
76 tmp<volScalarField> fd(const volScalarField& magGradU) const;
77
78 //- No copy construct
79 kOmegaSSTDDES(const kOmegaSSTDDES&) = delete;
80
81 //- No copy assignment
82 void operator=(const kOmegaSSTDDES&) = delete;
83
84
85protected:
86
87 // Protected Data
88
89 // Model coefficients
90
94
95 // Protected Member Functions
96
97 //- Return square of strain rate
99 (
100 const volTensorField& gradU
101 ) const;
102
103 //- Length scale
105 (
106 const volScalarField& magGradU,
107 const volScalarField& CDES
108 ) const;
110 //- Return (G/nu)_0
112 (
113 const volTensorField& gradU,
114 const volScalarField& S2
115 ) const;
116
117 //- Return G/nu
119 (
123 ) const;
124
125public:
126
127 typedef typename BasicTurbulenceModel::alphaField alphaField;
128 typedef typename BasicTurbulenceModel::rhoField rhoField;
129 typedef typename BasicTurbulenceModel::transportModel transportModel;
130
131
132 //- Runtime type information
133 TypeName("kOmegaSSTDDES");
134
135
136 // Constructors
138 //- Construct from components
139 kOmegaSSTDDES
140 (
141 const alphaField& alpha,
142 const rhoField& rho,
143 const volVectorField& U,
144 const surfaceScalarField& alphaRhoPhi,
145 const surfaceScalarField& phi,
146 const transportModel& transport,
147 const word& propertiesName = turbulenceModel::propertiesName,
148 const word& type = typeName
149 );
150
151
152 //- Destructor
153 virtual ~kOmegaSSTDDES() = default;
154
155
156 // Member Functions
157
158 //- Re-read model coefficients if they have changed
159 virtual bool read();
160
161 //- Return the shielding function
162 virtual tmp<volScalarField> fd() const;
163};
164
165
166// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
167
168} // End namespace LESModels
169} // End namespace Foam
170
171// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
172
173#ifdef NoRepository
174# include "kOmegaSSTDDES.C"
175#endif
177// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
178
179#endif
180
181// ************************************************************************* //
DimensionedField< scalar, volMesh > Internal
virtual tmp< volScalarField > dTilda(const volScalarField &magGradU, const volScalarField &CDES) const
Length scale.
Foam::EddyDiffusivity< Foam::fluidThermoCompressibleTurbulenceModel >::alphaField alphaField
Foam::EddyDiffusivity< Foam::fluidThermoCompressibleTurbulenceModel >::rhoField rhoField
virtual tmp< volScalarField > fd() const
Return the shielding function.
virtual ~kOmegaSSTDDES()=default
Destructor.
kOmegaSSTDDES(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 tmp< volScalarField > S2(const volTensorField &gradU) const
Return square of strain rate.
virtual tmp< volScalarField::Internal > GbyNu(const volScalarField::Internal &GbyNu0, const volScalarField::Internal &F2, const volScalarField::Internal &S2) const
Return G/nu.
virtual tmp< volScalarField::Internal > GbyNu0(const volTensorField &gradU, const volScalarField &S2) const
Return (G/nu)_0.
Foam::EddyDiffusivity< Foam::fluidThermoCompressibleTurbulenceModel >::transportModel transportModel
TypeName("kOmegaSSTDDES")
Runtime type information.
virtual bool read()
Re-read model coefficients if they have changed.
virtual tmp< volScalarField > F2() const
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")
GeometricField< tensor, fvPatchField, volMesh > volTensorField
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