SpalartAllmarasIDDES.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) 2011-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::LESModels::SpalartAllmarasIDDES
26 
27 Description
28  SpalartAllmaras IDDES turbulence model for incompressible and compressible
29  flows
30 
31  Reference:
32  \verbatim
33  Shur, M. L., Spalart, P. R., Strelets, M. K., & Travin, A. K. (2008).
34  A hybrid RANS-LES approach with delayed-DES and wall-modelled LES
35  capabilities.
36  International Journal of Heat and Fluid Flow, 29(6), 1638-1649.
37  \endverbatim
38 
39 SourceFiles
40  SpalartAllmarasIDDES.C
41 
42 \*---------------------------------------------------------------------------*/
43 
44 #ifndef SpalartAllmarasIDDES_H
45 #define SpalartAllmarasIDDES_H
46 
47 #include "SpalartAllmarasDES.H"
48 #include "IDDESDelta.H"
49 
50 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
51 
52 namespace Foam
53 {
54 namespace LESModels
55 {
56 
57 /*---------------------------------------------------------------------------*\
58  Class SpalartAllmarasIDDES Declaration
59 \*---------------------------------------------------------------------------*/
60 
61 template<class BasicMomentumTransportModel>
63 :
64  public SpalartAllmarasDES<BasicMomentumTransportModel>
65 {
66  // Private Data
67 
68  // Model constants
69 
70  dimensionedScalar fwStar_;
73 
74  // Fields
75 
76  const IDDESDelta& IDDESDelta_;
77 
78 
79  // Private Member Functions
80 
81  tmp<volScalarField::Internal> IDDESalpha() const;
82 
84  (
85  const volScalarField::Internal& magGradU
86  ) const;
87 
89  (
90  const volScalarField::Internal& magGradU
91  ) const;
92 
94  (
95  const volScalarField::Internal& nur,
96  const volScalarField::Internal& magGradU
97  ) const;
98 
99  //- Delay function
101  (
102  const volScalarField::Internal& magGradU
103  ) const;
104 
105 
106 protected:
107 
108  // Protected Member Functions
109 
110  using IOobject::modelName;
111 
112  //- Length scale
114  (
117  const volTensorField::Internal& gradU
118  ) const;
119 
120 
121 public:
123  typedef typename BasicMomentumTransportModel::alphaField alphaField;
124  typedef typename BasicMomentumTransportModel::rhoField rhoField;
125  typedef typename BasicMomentumTransportModel::transportModel transportModel;
126 
127 
128  //- Runtime type information
129  TypeName("SpalartAllmarasIDDES");
130 
131 
132  // Constructors
133 
134  //- Construct from components
136  (
137  const alphaField& alpha,
138  const rhoField& rho,
139  const volVectorField& U,
140  const surfaceScalarField& alphaRhoPhi,
141  const surfaceScalarField& phi,
142  const transportModel& transport,
143  const word& type = typeName
144  );
145 
146  //- Disallow default bitwise copy construction
148 
149 
150  //- Destructor
151  virtual ~SpalartAllmarasIDDES()
152  {}
153 
154 
155  // Member Functions
156 
157  //- Read model coefficients if they have changed
158  virtual bool read();
159 
160 
161  // Member Operators
162 
163  //- Disallow default bitwise assignment
164  void operator=(const SpalartAllmarasIDDES&) = delete;
165 };
166 
167 
168 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
169 
170 } // End namespace LESModels
171 } // End namespace Foam
172 
173 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
174 
175 #ifdef NoRepository
176  #include "SpalartAllmarasIDDES.C"
177 #endif
178 
179 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
180 
181 #endif
182 
183 // ************************************************************************* //
SpalartAllmarasDES DES turbulence model for incompressible and compressible flows.
BasicMomentumTransportModel::alphaField alphaField
DimensionedField< scalar, volMesh > Internal
Type of the internal field from which this GeometricField is derived.
virtual bool read()
Read model coefficients if they have changed.
SpalartAllmaras IDDES turbulence model for incompressible and compressible flows. ...
phi
Definition: pEqn.H:104
BasicMomentumTransportModel::rhoField rhoField
A class for handling words, derived from string.
Definition: word.H:59
tmp< volScalarField > fv1(const volScalarField &chi) const
IDDESDelta used by the IDDES (improved low Re Spalart-Allmaras DES model) The min and max delta are c...
Definition: IDDESDelta.H:52
U
Definition: pEqn.H:72
TypeName("SpalartAllmarasIDDES")
Runtime type information.
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
SpalartAllmarasIDDES(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.
static word modelName(Name name, const word &model)
Return the name of the object within the given model.
virtual tmp< volScalarField::Internal > dTilda(const volScalarField::Internal &chi, const volScalarField::Internal &fv1, const volTensorField::Internal &gradU) const
Length scale.
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
BasicMomentumTransportModel::transportModel transportModel
tmp< volScalarField > chi() const
void operator=(const SpalartAllmarasIDDES &)=delete
Disallow default bitwise assignment.
Namespace for OpenFOAM.