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-2020 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 BasicMomentumTransportModel>
96 class kOmegaSSTSAS
97 :
98  public kOmegaSST<BasicMomentumTransportModel>
99 {
100 protected:
101 
102  // Protected data
103 
104  // Model constants
111 
112 
113  // Fields
114 
115  //- Run-time selectable delta model
117 
118 
119  // Protected Member Functions
120 
121  //- SAS omega source
122  virtual tmp<fvScalarMatrix> Qsas
123  (
124  const volScalarField::Internal& S2,
127  ) const;
128 
129 
130 public:
132  typedef typename BasicMomentumTransportModel::alphaField alphaField;
133  typedef typename BasicMomentumTransportModel::rhoField rhoField;
134  typedef typename BasicMomentumTransportModel::transportModel transportModel;
135 
136 
137  //- Runtime type information
138  TypeName("kOmegaSSTSAS");
139 
140 
141  // Constructors
142 
143  //- Construct from components
145  (
146  const alphaField& alpha,
147  const rhoField& rho,
148  const volVectorField& U,
149  const surfaceScalarField& alphaRhoPhi,
150  const surfaceScalarField& phi,
151  const transportModel& transport,
152  const word& type = typeName
153  );
154 
155  //- Disallow default bitwise copy construction
156  kOmegaSSTSAS(const kOmegaSSTSAS&) = delete;
157 
158 
159  //- Destructor
160  virtual ~kOmegaSSTSAS()
161  {}
162 
163 
164  // Member Functions
165 
166  //- Re-read model coefficients if they have changed
167  virtual bool read();
168 
169  //- Access function to filter width
170  inline const volScalarField& delta() const
171  {
172  return delta_();
173  }
174 
175 
176  // Member Operators
177 
178  //- Disallow default bitwise assignment
179  void operator=(const kOmegaSSTSAS&) = delete;
180 };
181 
182 
183 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
184 
185 } // End namespace RASModels
186 } // End namespace Foam
187 
188 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
189 
190 #ifdef NoRepository
191  #include "kOmegaSSTSAS.C"
192 #endif
193 
194 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
195 
196 #endif
197 
198 // ************************************************************************* //
autoPtr< Foam::LESdelta > delta_
Run-time selectable delta model.
Definition: kOmegaSSTSAS.H:115
DimensionedField< scalar, volMesh > Internal
Type of the internal field from which this GeometricField is derived.
phi
Definition: pEqn.H:104
BasicMomentumTransportModel::rhoField rhoField
Definition: kOmegaSSTSAS.H:132
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:108
Specialisation for RAS of the generic kOmegaSSTBase base class. For more information, see Description of kOmegaSSTBase.H.
Definition: kOmegaSST.H:57
kOmegaSSTSAS(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &type=typeName)
Construct from components.
Definition: kOmegaSSTSAS.C:97
tmp< volScalarField::Internal > gamma(const volScalarField::Internal &F1) const
BasicMomentumTransportModel::alphaField alphaField
Definition: kOmegaSSTSAS.H:131
void operator=(const kOmegaSSTSAS &)=delete
Disallow default bitwise assignment.
virtual bool read()
Re-read model coefficients if they have changed.
Definition: kOmegaSSTSAS.C:178
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
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
A class for managing temporary objects.
Definition: PtrList.H:53
virtual ~kOmegaSSTSAS()
Destructor.
Definition: kOmegaSSTSAS.H:159
TypeName("kOmegaSSTSAS")
Runtime type information.
const volScalarField & delta() const
Access function to filter width.
Definition: kOmegaSSTSAS.H:169
Namespace for OpenFOAM.
BasicMomentumTransportModel::transportModel transportModel
Definition: kOmegaSSTSAS.H:133
Scale-adaptive URAS model based on the k-omega-SST RAS model.
Definition: kOmegaSSTSAS.H:95