Loading...
Searching...
No Matches
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 Foam_kOmegaSSTDES_H
55#define Foam_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 //- Switch to activate grey-area enhanced sigma-(D)DES
92 // Model coefficients
93
96 //- DES coefficients
99
101 // Protected Member Functions
102
103 //- Blending for CDES parameter
104 virtual tmp<volScalarField> CDES(const volScalarField& F1) const
105 {
106 return this->blend(F1, CDESkom_, CDESkeps_);
107 }
108
109 virtual void correctNut(const volScalarField& S2);
110 virtual void correctNut();
111
113 (
114 const volScalarField& nur,
115 const volScalarField& magGradU
116 ) const;
118 //- Return square of strain rate
119 virtual tmp<volScalarField> S2
120 (
121 const volTensorField& gradU
122 ) const;
123
124 //- Return length scale
127 const volScalarField& magGradU,
128 const volScalarField& CDES
129 ) const;
130
131 //- Return epsilon/k
133 (
135 const volTensorField& gradU
136 ) const;
137
138 //- Return (G/nu)_0
140 (
141 const volTensorField& gradU,
142 const volScalarField& S2
143 ) const;
144
145 //- Return G/nu
147 (
151 ) const;
153
154public:
155
156 typedef typename BasicTurbulenceModel::alphaField alphaField;
157 typedef typename BasicTurbulenceModel::rhoField rhoField;
158 typedef typename BasicTurbulenceModel::transportModel transportModel;
159
160
161 //- Runtime type information
162 TypeName("kOmegaSSTDES");
163
164
165 // Constructors
166
167 //- Construct from components
168 kOmegaSSTDES
169 (
170 const alphaField& alpha,
171 const rhoField& rho,
173 const surfaceScalarField& alphaRhoPhi,
174 const surfaceScalarField& phi,
175 const transportModel& transport,
176 const word& propertiesName = turbulenceModel::propertiesName,
177 const word& type = typeName
178 );
180
181 //- Destructor
182 virtual ~kOmegaSSTDES() = default;
183
184
185 // Member Functions
186
187 //- Re-read model coefficients if they have changed
188 virtual bool read();
189
190 //- RAS length scale
191 virtual tmp<volScalarField> lengthScaleRAS() const;
192
193 //- LES length scale
195 (
196 const volScalarField& CDES
197 ) const;
198
199 //- Return the LES field indicator
200 virtual tmp<volScalarField> LESRegion() const;
201};
202
204// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
205
206} // End namespace LESModels
207} // End namespace Foam
208
209// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
210#ifdef NoRepository
211 #include "kOmegaSSTDES.C"
212#endif
213
214// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
215#endif
217// ************************************************************************* //
DimensionedField< scalar, volMesh > Internal
virtual tmp< volScalarField > dTilda(const volScalarField &magGradU, const volScalarField &CDES) const
Return length scale.
Foam::EddyDiffusivity< Foam::fluidThermoCompressibleTurbulenceModel >::alphaField alphaField
TypeName("kOmegaSSTDES")
Runtime type information.
Foam::EddyDiffusivity< Foam::fluidThermoCompressibleTurbulenceModel >::rhoField rhoField
virtual tmp< volScalarField > lengthScaleLES(const volScalarField &CDES) const
LES length scale.
virtual void correctNut(const volScalarField &S2)
virtual tmp< volScalarField > lengthScaleRAS() const
RAS length scale.
kOmegaSSTDES(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 ~kOmegaSSTDES()=default
Destructor.
virtual tmp< volScalarField > LESRegion() const
Return the LES field indicator.
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.
virtual tmp< volScalarField > CDES(const volScalarField &F1) const
Blending for CDES parameter.
Foam::EddyDiffusivity< Foam::fluidThermoCompressibleTurbulenceModel >::transportModel transportModel
virtual tmp< volScalarField::Internal > epsilonByk(const volScalarField &F1, const volTensorField &gradU) const
Return epsilon/k.
tmp< volScalarField > r(const volScalarField &nur, const volScalarField &magGradU) const
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
virtual tmp< volScalarField > F1(const volScalarField &CDkOmega) const
virtual tmp< volScalarField > F2() const
tmp< volScalarField > blend(const volScalarField &F1, const dimensionedScalar &psi1, const dimensionedScalar &psi2) 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