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 | Copyright (C) 2015-2016 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 Class
25  Foam::RASModels::kOmegaSSTSAS
26 
27 Group
28  grpLESTurbulence
29 
30 Description
31  Scale-adaptive URAS model based on the k-omega-SST RAS model.
32 
33  References:
34  \verbatim
35  Egorov, Y., & Menter F.R. (2008).
36  Development and Application of SST-SAS Model in the DESIDER Project.
37  Advances in Hybrid RANS-LES Modelling,
38  Notes on Num. Fluid Mech. And Multidisciplinary Design,
39  Volume 97, 261-270.
40  \endverbatim
41 
42  The model coefficients are
43  \verbatim
44  kOmegaSSTSASCoeffs
45  {
46  // Default SST coefficients
47  alphaK1 0.85;
48  alphaK2 1.0;
49  alphaOmega1 0.5;
50  alphaOmega2 0.856;
51  beta1 0.075;
52  beta2 0.0828;
53  betaStar 0.09;
54  gamma1 5/9;
55  gamma2 0.44;
56  a1 0.31;
57  b1 1.0;
58  c1 10.0;
59  F3 no;
60 
61  // Default SAS coefficients
62  Cs 0.11;
63  kappa 0.41;
64  zeta2 3.51;
65  sigmaPhi 2.0/3.0;
66  C 2;
67 
68  // Delta must be specified for SAS e.g.
69  delta cubeRootVol;
70 
71  cubeRootVolCoeffs
72  {}
73  }
74  \endverbatim
75 
76 SourceFiles
77  kOmegaSSTSAS.C
78 
79 \*---------------------------------------------------------------------------*/
80 
81 #ifndef kOmegaSSTSAS_H
82 #define kOmegaSSTSAS_H
83 
84 #include "kOmegaSST.H"
85 #include "LESdelta.H"
86 
87 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
88 
89 namespace Foam
90 {
91 namespace RASModels
92 {
93 
94 /*---------------------------------------------------------------------------*\
95  Class kOmegaSSTSAS Declaration
96 \*---------------------------------------------------------------------------*/
97 
98 template<class BasicTurbulenceModel>
99 class kOmegaSSTSAS
100 :
101  public kOmegaSST<BasicTurbulenceModel>
102 {
103  // Private Member Functions
104 
105  // Disallow default bitwise copy construct and assignment
106  kOmegaSSTSAS(const kOmegaSSTSAS&);
107  void operator=(const kOmegaSSTSAS&);
108 
109 
110 protected:
111 
112  // Protected data
113 
114  // Model constants
121 
122 
123  // Fields
124 
125  //- Run-time selectable delta model
127 
128 
129  // Protected Member Functions
130 
131  //- SAS omega source
132  virtual tmp<fvScalarMatrix> Qsas
133  (
134  const volScalarField::Internal& S2,
137  ) const;
138 
139 
140 public:
142  typedef typename BasicTurbulenceModel::alphaField alphaField;
143  typedef typename BasicTurbulenceModel::rhoField rhoField;
144  typedef typename BasicTurbulenceModel::transportModel transportModel;
145 
146 
147  //- Runtime type information
148  TypeName("kOmegaSSTSAS");
149 
150 
151  // Constructors
152 
153  //- Construct from components
155  (
156  const alphaField& alpha,
157  const rhoField& rho,
158  const volVectorField& U,
159  const surfaceScalarField& alphaRhoPhi,
160  const surfaceScalarField& phi,
161  const transportModel& transport,
162  const word& propertiesName = turbulenceModel::propertiesName,
163  const word& type = typeName
164  );
165 
166 
167  //- Destructor
168  virtual ~kOmegaSSTSAS()
169  {}
170 
171 
172  // Member Functions
173 
174  //- Re-read model coefficients if they have changed
175  virtual bool read();
176 
177  //- Access function to filter width
178  inline const volScalarField& delta() const
179  {
180  return delta_();
181  }
182 };
183 
184 
185 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
186 
187 } // End namespace RASModels
188 } // End namespace Foam
189 
190 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
191 
192 #ifdef NoRepository
193  #include "kOmegaSSTSAS.C"
194 #endif
195 
196 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
197 
198 #endif
199 
200 // ************************************************************************* //
surfaceScalarField & phi
BasicTurbulenceModel::rhoField rhoField
Definition: kOmegaSSTSAS.H:142
U
Definition: pEqn.H:83
autoPtr< Foam::LESdelta > delta_
Run-time selectable delta model.
Definition: kOmegaSSTSAS.H:125
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from string.
Definition: word.H:59
tmp< volScalarField::Internal > beta(const volScalarField::Internal &F1) const
dimensionedScalar sigmaPhi_
Definition: kOmegaSSTSAS.H:118
tmp< volScalarField::Internal > gamma(const volScalarField::Internal &F1) const
fileName::Type type(const fileName &, const bool followLink=true)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:485
virtual bool read()
Re-read model coefficients if they have changed.
Definition: kOmegaSSTSAS.C:180
virtual tmp< fvScalarMatrix > Qsas(const volScalarField::Internal &S2, const volScalarField::Internal &gamma, const volScalarField::Internal &beta) const
SAS omega source.
Definition: kOmegaSSTSAS.C:39
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
A class for managing temporary objects.
Definition: PtrList.H:53
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
virtual ~kOmegaSSTSAS()
Destructor.
Definition: kOmegaSSTSAS.H:167
TypeName("kOmegaSSTSAS")
Runtime type information.
BasicTurbulenceModel::alphaField alphaField
Definition: kOmegaSSTSAS.H:141
const volScalarField & delta() const
Access function to filter width.
Definition: kOmegaSSTSAS.H:177
BasicTurbulenceModel::transportModel transportModel
Definition: kOmegaSSTSAS.H:143
Namespace for OpenFOAM.
Scale-adaptive URAS model based on the k-omega-SST RAS model.
Definition: kOmegaSSTSAS.H:98