Loading...
Searching...
No Matches
SpalartAllmarasDDES.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) 2011-2016 OpenFOAM Foundation
9 Copyright (C) 2022, 2024 Upstream CFD GmbH
10 Copyright (C) 2019-2024 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::SpalartAllmarasDDES
30
31Group
32 grpDESTurbulence
33
34Description
35 SpalartAllmaras DDES turbulence model for incompressible and compressible
36 flows.
37
38 Reference:
39 \verbatim
40 Spalart, P. R., Deck, S., Shur, M. L., Squires,
41 K. D., Strelets, M. K., & Travin, A. (2006).
42 A new version of detached-eddy simulation,
43 resistant to ambiguous grid densities.
44 Theoretical and computational fluid dynamics, 20(3), 181-195.
45 DOI:10.1007/s00162-006-0015-0
46 \endverbatim
47
48 Reference for enhanced shielding function formulation:
49 \verbatim
50 Deck, S., Renard, N. (2020).
51 Towards an enhanced protection of attached boundary layers in hybrid
52 RANS/LES methods.
53 Journal of Computational Physics, 400, 108970.
54 DOI:10.1016/j.jcp.2019.108970
55 \endverbatim
56
57SourceFiles
58 SpalartAllmarasDDES.C
59
60\*---------------------------------------------------------------------------*/
61
62#ifndef Foam_SpalartAllmarasDDES_H
63#define Foam_SpalartAllmarasDDES_H
64
65#include "SpalartAllmarasBase.H"
66#include "Enum.H"
67
68// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
69
70namespace Foam
71{
72namespace LESModels
73{
74
75/*---------------------------------------------------------------------------*\
76 Class SpalartAllmarasDDES Declaration
77\*---------------------------------------------------------------------------*/
78
79template<class BasicTurbulenceModel>
80class SpalartAllmarasDDES
81:
82 public SpalartAllmarasDES<BasicTurbulenceModel>
83{
84public:
85
86 // Public enumerations
87
88 //- Shielding modes
89 enum class shieldingMode
90 {
93 };
94
95 //- Shielding mode names
97
98
99private:
100
101 // Private Member Functions
102
103 //- Return the shielding function
104 tmp<volScalarField> fd(const volScalarField& magGradU) const;
105
106 //- No copy construct
107 SpalartAllmarasDDES(const SpalartAllmarasDDES&) = delete;
108
109 //- No copy assignment
110 void operator=(const SpalartAllmarasDDES&) = delete;
111
112
113protected:
114
115 // Protected Data
116
117 //- Shielding mode
119
120 // Model coefficients
121
128
130 // Protected Member Functions
132 //- Return the production term
135 const volScalarField& chi,
136 const volScalarField& fv1,
137 const volTensorField& gradU,
139 ) const;
140
141 //- Return the length scale
143 (
144 const volScalarField& chi,
145 const volScalarField& fv1,
146 const volTensorField& gradU
147 ) const;
148
149
150public:
151
152 typedef typename BasicTurbulenceModel::alphaField alphaField;
153 typedef typename BasicTurbulenceModel::rhoField rhoField;
154 typedef typename BasicTurbulenceModel::transportModel transportModel;
155
156
157 //- Runtime type information
158 TypeName("SpalartAllmarasDDES");
159
160
161 // Constructors
162
163 //- Construct from components
164 SpalartAllmarasDDES
166 const alphaField& alpha,
167 const rhoField& rho,
168 const volVectorField& U,
169 const surfaceScalarField& alphaRhoPhi,
170 const surfaceScalarField& phi,
171 const transportModel& transport,
172 const word& propertiesName = turbulenceModel::propertiesName,
173 const word& type = typeName
174 );
175
176
177 //- Destructor
178 virtual ~SpalartAllmarasDDES() = default;
180
181 // Member Functions
182
183 //- Read from dictionary
184 virtual bool read();
185
186 //- Return the shielding function
187 virtual tmp<volScalarField> fd() const;
188};
189
190
191// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
192
193} // End namespace LESModels
194} // End namespace Foam
196// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
197
198#ifdef NoRepository
199 #include "SpalartAllmarasDDES.C"
200#endif
201
202// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
204#endif
205
206// ************************************************************************* //
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition Enum.H:57
Foam::EddyDiffusivity< Foam::fluidThermoCompressibleTurbulenceModel >::alphaField alphaField
Foam::EddyDiffusivity< Foam::fluidThermoCompressibleTurbulenceModel >::rhoField rhoField
virtual tmp< volScalarField > fd() const
Return the shielding function.
SpalartAllmarasDDES(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 ~SpalartAllmarasDDES()=default
Destructor.
Foam::EddyDiffusivity< Foam::fluidThermoCompressibleTurbulenceModel >::transportModel transportModel
TypeName("SpalartAllmarasDDES")
Runtime type information.
virtual tmp< volScalarField > Stilda(const volScalarField &chi, const volScalarField &fv1, const volTensorField &gradU, const volScalarField &dTilda) const
Return the production term.
virtual tmp< volScalarField > dTilda(const volScalarField &chi, const volScalarField &fv1, const volTensorField &gradU) const
virtual bool read()
Read from dictionary.
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")
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