SpalartAllmarasDDES.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) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2022 Upstream CFD GmbH
10  Copyright (C) 2019-2022 OpenCFD Ltd.
11 -------------------------------------------------------------------------------
12 License
13  This file is part of OpenFOAM.
14 
15  OpenFOAM is free software: you can redistribute it and/or modify it
16  under the terms of the GNU General Public License as published by
17  the Free Software Foundation, either version 3 of the License, or
18  (at your option) any later version.
19 
20  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
21  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
22  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
23  for more details.
24 
25  You should have received a copy of the GNU General Public License
26  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
27 
28 Class
29  Foam::LESModels::SpalartAllmarasDDES
30 
31 Group
32  grpDESTurbulence
33 
34 Description
35  SpalartAllmaras DDES turbulence model for incompressible and compressible
36  flows.
37 
38  Reference:
39  \verbatim
40  Spalart, P. R., Deck, S., Shur, M. L., Squires,
41  K. D., Strelets, M. K., & Travin, A. (2006).
42  A new version of detached-eddy simulation,
43  resistant to ambiguous grid densities.
44  Theoretical and computational fluid dynamics, 20(3), 181-195.
45  DOI:10.1007/s00162-006-0015-0
46  \endverbatim
47 
48 SourceFiles
49  SpalartAllmarasDDES.C
50 
51 \*---------------------------------------------------------------------------*/
52 
53 #ifndef Foam_SpalartAllmarasDDES_H
54 #define Foam_SpalartAllmarasDDES_H
55 
56 #include "SpalartAllmarasBase.H"
57 
58 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
59 
60 namespace Foam
61 {
62 namespace LESModels
63 {
64 
65 /*---------------------------------------------------------------------------*\
66  Class SpalartAllmarasDDES Declaration
67 \*---------------------------------------------------------------------------*/
68 
69 template<class BasicTurbulenceModel>
71 :
72  public SpalartAllmarasDES<BasicTurbulenceModel>
73 {
74  // Private Member Functions
75 
76  //- Return the shielding function
77  tmp<volScalarField> fd(const volScalarField& magGradU) const;
78 
79  //- No copy construct
81 
82  //- No copy assignment
83  void operator=(const SpalartAllmarasDDES&) = delete;
84 
85 
86 protected:
87 
88  // Protected Data
89 
90  // Model coefficients
91 
94 
95 
96  // Protected Member Functions
97 
98  //- Return the production term
100  (
101  const volScalarField& chi,
102  const volScalarField& fv1,
103  const volTensorField& gradU,
104  const volScalarField& dTilda
105  ) const;
106 
107  //- Return the length scale
109  (
110  const volScalarField& chi,
111  const volScalarField& fv1,
112  const volTensorField& gradU
113  ) const;
114 
115 
116 public:
117 
118  typedef typename BasicTurbulenceModel::alphaField alphaField;
119  typedef typename BasicTurbulenceModel::rhoField rhoField;
120  typedef typename BasicTurbulenceModel::transportModel transportModel;
121 
122 
123  //- Runtime type information
124  TypeName("SpalartAllmarasDDES");
126 
127  // Constructors
128 
129  //- Construct from components
131  (
132  const alphaField& alpha,
133  const rhoField& rho,
134  const volVectorField& U,
135  const surfaceScalarField& alphaRhoPhi,
136  const surfaceScalarField& phi,
137  const transportModel& transport,
138  const word& propertiesName = turbulenceModel::propertiesName,
139  const word& type = typeName
140  );
141 
142 
143  //- Destructor
144  virtual ~SpalartAllmarasDDES() = default;
145 
146 
147  // Member Functions
148 
149  //- Read from dictionary
150  virtual bool read();
151 
152  //- Return the shielding function
153  virtual tmp<volScalarField> fd() const;
154 };
155 
156 
157 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
158 
159 } // End namespace LESModels
160 } // End namespace Foam
161 
162 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
163 
164 #ifdef NoRepository
165  #include "SpalartAllmarasDDES.C"
166 #endif
167 
168 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
169 
170 #endif
171 
172 // ************************************************************************* //
BasicTurbulenceModel::transportModel transportModel
virtual ~SpalartAllmarasDDES()=default
Destructor.
SpalartAllmarasDES DES turbulence model for incompressible and compressible flows.
BasicTurbulenceModel::alphaField alphaField
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from Foam::string.
Definition: word.H:63
virtual tmp< volScalarField > dTilda(const volScalarField &chi, const volScalarField &fv1, const volTensorField &gradU) const
Return the length scale.
SpalartAllmaras DDES turbulence model for incompressible and compressible flows.
TypeName("SpalartAllmarasDDES")
Runtime type information.
virtual tmp< volScalarField > Stilda(const volScalarField &chi, const volScalarField &fv1, const volTensorField &gradU, const volScalarField &dTilda) const
Return the production term.
BasicTurbulenceModel::rhoField rhoField
U
Definition: pEqn.H:72
virtual bool read()
Read from dictionary.
A class for managing temporary objects.
Definition: HashPtrTable.H:50
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
virtual tmp< volScalarField > fd() const
Return the shielding function.
Namespace for OpenFOAM.