kOmegaSSTSAS.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2015-2018 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 Description
28  Scale-adaptive URAS model based on the k-omega-SST RAS model.
29 
30  References:
31  \verbatim
32  Egorov, Y., & Menter F.R. (2008).
33  Development and Application of SST-SAS Model in the DESIDER Project.
34  Advances in Hybrid RANS-LES Modelling,
35  Notes on Num. Fluid Mech. And Multidisciplinary Design,
36  Volume 97, 261-270.
37  \endverbatim
38 
39  The model coefficients are
40  \verbatim
41  kOmegaSSTSASCoeffs
42  {
43  // Default SST coefficients
44  alphaK1 0.85;
45  alphaK2 1.0;
46  alphaOmega1 0.5;
47  alphaOmega2 0.856;
48  beta1 0.075;
49  beta2 0.0828;
50  betaStar 0.09;
51  gamma1 5/9;
52  gamma2 0.44;
53  a1 0.31;
54  b1 1.0;
55  c1 10.0;
56  F3 no;
57 
58  // Default SAS coefficients
59  Cs 0.11;
60  kappa 0.41;
61  zeta2 3.51;
62  sigmaPhi 2.0/3.0;
63  C 2;
64 
65  // Delta must be specified for SAS e.g.
66  delta cubeRootVol;
67 
68  cubeRootVolCoeffs
69  {}
70  }
71  \endverbatim
72 
73 SourceFiles
74  kOmegaSSTSAS.C
75 
76 \*---------------------------------------------------------------------------*/
77 
78 #ifndef kOmegaSSTSAS_H
79 #define kOmegaSSTSAS_H
80 
81 #include "kOmegaSST.H"
82 #include "LESdelta.H"
83 
84 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
85 
86 namespace Foam
87 {
88 namespace RASModels
89 {
90 
91 /*---------------------------------------------------------------------------*\
92  Class kOmegaSSTSAS Declaration
93 \*---------------------------------------------------------------------------*/
94 
95 template<class BasicTurbulenceModel>
96 class kOmegaSSTSAS
97 :
98  public kOmegaSST<BasicTurbulenceModel>
99 {
100  // Private Member Functions
101 
102  // Disallow default bitwise copy construct and assignment
103  kOmegaSSTSAS(const kOmegaSSTSAS&);
104  void operator=(const kOmegaSSTSAS&);
105 
106 
107 protected:
108 
109  // Protected data
110 
111  // Model constants
118 
119 
120  // Fields
121 
122  //- Run-time selectable delta model
124 
125 
126  // Protected Member Functions
127 
128  //- SAS omega source
129  virtual tmp<fvScalarMatrix> Qsas
130  (
131  const volScalarField::Internal& S2,
134  ) const;
135 
136 
137 public:
139  typedef typename BasicTurbulenceModel::alphaField alphaField;
140  typedef typename BasicTurbulenceModel::rhoField rhoField;
141  typedef typename BasicTurbulenceModel::transportModel transportModel;
142 
143 
144  //- Runtime type information
145  TypeName("kOmegaSSTSAS");
146 
147 
148  // Constructors
149 
150  //- Construct from components
152  (
153  const alphaField& alpha,
154  const rhoField& rho,
155  const volVectorField& U,
156  const surfaceScalarField& alphaRhoPhi,
157  const surfaceScalarField& phi,
158  const transportModel& transport,
159  const word& propertiesName = turbulenceModel::propertiesName,
160  const word& type = typeName
161  );
162 
163 
164  //- Destructor
165  virtual ~kOmegaSSTSAS()
166  {}
167 
168 
169  // Member Functions
170 
171  //- Re-read model coefficients if they have changed
172  virtual bool read();
173 
174  //- Access function to filter width
175  inline const volScalarField& delta() const
176  {
177  return delta_();
178  }
179 };
180 
181 
182 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
183 
184 } // End namespace RASModels
185 } // End namespace Foam
186 
187 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
188 
189 #ifdef NoRepository
190  #include "kOmegaSSTSAS.C"
191 #endif
192 
193 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
194 
195 #endif
196 
197 // ************************************************************************* //
BasicTurbulenceModel::rhoField rhoField
Definition: kOmegaSSTSAS.H:139
surfaceScalarField & phi
autoPtr< Foam::LESdelta > delta_
Run-time selectable delta model.
Definition: kOmegaSSTSAS.H:122
DimensionedField< scalar, volMesh > Internal
Type of the internal field from which this GeometricField is derived.
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:115
Specialisation for RAS of the generic kOmegaSSTBase base class. For more information, see Description of kOmegaSSTBase.H.
Definition: kOmegaSST.H:57
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:481
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
U
Definition: pEqn.H:72
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:164
TypeName("kOmegaSSTSAS")
Runtime type information.
BasicTurbulenceModel::alphaField alphaField
Definition: kOmegaSSTSAS.H:138
const volScalarField & delta() const
Access function to filter width.
Definition: kOmegaSSTSAS.H:174
BasicTurbulenceModel::transportModel transportModel
Definition: kOmegaSSTSAS.H:140
Namespace for OpenFOAM.
Scale-adaptive URAS model based on the k-omega-SST RAS model.
Definition: kOmegaSSTSAS.H:95