Loading...
Searching...
No Matches
kOmegaSSTSAS.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-2016 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::RASModels::kOmegaSSTSAS
29
30Group
31 grpLESTurbulence
32
33Description
34 Scale-adaptive URAS model based on the k-omega-SST RAS model.
35
36 References:
37 \verbatim
38 Egorov, Y., & Menter F.R. (2008).
39 Development and Application of SST-SAS Model in the DESIDER Project.
40 Advances in Hybrid RANS-LES Modelling,
41 Notes on Num. Fluid Mech. And Multidisciplinary Design,
42 Volume 97, 261-270.
43 \endverbatim
44
45 The model coefficients are
46 \verbatim
47 kOmegaSSTSASCoeffs
48 {
49 // Default SST coefficients
50 alphaK1 0.85;
51 alphaK2 1.0;
52 alphaOmega1 0.5;
53 alphaOmega2 0.856;
54 beta1 0.075;
55 beta2 0.0828;
56 betaStar 0.09;
57 gamma1 5/9;
58 gamma2 0.44;
59 a1 0.31;
60 b1 1.0;
61 c1 10.0;
62 F3 no;
63
64 // Default SAS coefficients
65 Cs 0.11;
66 kappa 0.41;
67 zeta2 3.51;
68 sigmaPhi 2.0/3.0;
69 C 2;
70
71 // Delta must be specified for SAS e.g.
72 delta cubeRootVol;
73
74 cubeRootVolCoeffs
75 {}
76 }
77 \endverbatim
78
79SourceFiles
80 kOmegaSSTSAS.C
81
82\*---------------------------------------------------------------------------*/
83
84#ifndef kOmegaSSTSAS_H
85#define kOmegaSSTSAS_H
86
87#include "kOmegaSST.H"
88#include "LESdelta.H"
89
90// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
91
92namespace Foam
93{
94namespace RASModels
95{
96
97/*---------------------------------------------------------------------------*\
98 Class kOmegaSSTSAS Declaration
99\*---------------------------------------------------------------------------*/
100
101template<class BasicTurbulenceModel>
102class kOmegaSSTSAS
103:
104 public kOmegaSST<BasicTurbulenceModel>
105{
106 // Private Member Functions
107
108 //- No copy construct
109 kOmegaSSTSAS(const kOmegaSSTSAS&) = delete;
110
111 //- No copy assignment
112 void operator=(const kOmegaSSTSAS&) = delete;
113
114
115protected:
116
117 // Protected data
118
119 // Model constants
126
127
128 // Fields
129
130 //- Run-time selectable delta model
133
134 // Protected Member Functions
135
136 //- SAS omega source
138 (
142 ) const;
143
144
145public:
146
147 typedef typename BasicTurbulenceModel::alphaField alphaField;
148 typedef typename BasicTurbulenceModel::rhoField rhoField;
149 typedef typename BasicTurbulenceModel::transportModel transportModel;
152 //- Runtime type information
153 TypeName("kOmegaSSTSAS");
154
155
156 // Constructors
157
158 //- Construct from components
159 kOmegaSSTSAS
160 (
161 const alphaField& alpha,
162 const rhoField& rho,
163 const volVectorField& U,
164 const surfaceScalarField& alphaRhoPhi,
165 const surfaceScalarField& phi,
166 const transportModel& transport,
167 const word& propertiesName = turbulenceModel::propertiesName,
168 const word& type = typeName
169 );
170
171
172 //- Destructor
173 virtual ~kOmegaSSTSAS() = default;
174
175
176 // Member Functions
177
178 //- Re-read model coefficients if they have changed
179 virtual bool read();
180
181 //- Access function to filter width
182 inline const volScalarField& delta() const
183 {
184 return *delta_;
185 }
186};
187
188
189// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
191} // End namespace RASModels
192} // End namespace Foam
193
194// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
196#ifdef NoRepository
197 #include "kOmegaSSTSAS.C"
198#endif
199
200// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
201
202#endif
203
204// ************************************************************************* //
DimensionedField< scalar, volMesh > Internal
Foam::EddyDiffusivity< Foam::fluidThermoCompressibleTurbulenceModel >::alphaField alphaField
Foam::EddyDiffusivity< Foam::fluidThermoCompressibleTurbulenceModel >::rhoField rhoField
kOmegaSSTSAS(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 ~kOmegaSSTSAS()=default
Destructor.
TypeName("kOmegaSSTSAS")
Runtime type information.
virtual tmp< fvScalarMatrix > Qsas(const volScalarField::Internal &S2, const volScalarField::Internal &gamma, const volScalarField::Internal &beta) const
SAS omega source.
Foam::EddyDiffusivity< Foam::fluidThermoCompressibleTurbulenceModel >::transportModel transportModel
const volScalarField & delta() const
Access function to filter width.
virtual bool read()
Re-read model coefficients if they have changed.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
tmp< volScalarField::Internal > gamma(const volScalarField::Internal &F1) const
tmp< volScalarField::Internal > beta(const volScalarField::Internal &F1) const
virtual tmp< volScalarField > S2(const volTensorField &gradU) const
Return square of strain rate.
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 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