WALE.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::LESModels::WALE
26 
27 Description
28  The Wall-adapting local eddy-viscosity (WALE) SGS model.
29 
30  Reference:
31  \verbatim
32  Nicoud, F., & Ducros, F. (1999).
33  Subgrid-scale stress modelling based on the square of the velocity
34  gradient tensor.
35  Flow, Turbulence and Combustion, 62(3), 183-200.
36  \endverbatim
37 
38  The default model coefficients are
39  \verbatim
40  WALECoeffs
41  {
42  Ck 0.094;
43  Ce 1.048;e
44  Cw 0.325;
45  }
46  \endverbatim
47 
48 See also
49  Foam::LESModels::Smagorinsky
50 
51 SourceFiles
52  WALE.C
53 
54 \*---------------------------------------------------------------------------*/
55 
56 #ifndef WALE_H
57 #define WALE_H
58 
59 #include "LESModel.H"
60 #include "LESeddyViscosity.H"
61 
62 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
63 
64 namespace Foam
65 {
66 namespace LESModels
67 {
68 
69 /*---------------------------------------------------------------------------*\
70  Class WALE Declaration
71 \*---------------------------------------------------------------------------*/
72 
73 template<class BasicTurbulenceModel>
74 class WALE
75 :
76  public LESeddyViscosity<BasicTurbulenceModel>
77 {
78  // Private Member Functions
79 
80  // Disallow default bitwise copy construct and assignment
81  WALE(const WALE&);
82  void operator=(const WALE&);
83 
84 
85 protected:
86 
87  // Protected data
88 
91 
92 
93  // Protected Member Functions
94 
95  //- Return the deviatoric symmetric part of the square of the given
96  // velocity gradient field
97  tmp<volSymmTensorField> Sd(const volTensorField& gradU) const;
98 
99  //- Return SGS kinetic energy
100  // calculated from the given velocity gradient
101  tmp<volScalarField> k(const volTensorField& gradU) const;
102 
103  //- Update the SGS eddy-viscosity
104  virtual void correctNut();
105 
106 
107 public:
109  typedef typename BasicTurbulenceModel::alphaField alphaField;
110  typedef typename BasicTurbulenceModel::rhoField rhoField;
111  typedef typename BasicTurbulenceModel::transportModel transportModel;
112 
113 
114  //- Runtime type information
115  TypeName("WALE");
116 
117 
118  // Constructors
119 
120  //- Construct from components
121  WALE
122  (
123  const alphaField& alpha,
124  const rhoField& rho,
125  const volVectorField& U,
126  const surfaceScalarField& alphaRhoPhi,
127  const surfaceScalarField& phi,
128  const transportModel& transport,
129  const word& propertiesName = turbulenceModel::propertiesName,
130  const word& type = typeName
131  );
132 
133 
134  //- Destructor
135  virtual ~WALE()
136  {}
137 
138 
139  // Member Functions
140 
141  //- Read model coefficients if they have changed
142  virtual bool read();
143 
144  //- Return SGS kinetic energy
145  virtual tmp<volScalarField> k() const
146  {
147  return k(fvc::grad(this->U_));
148  }
149 
150  //- Return sub-grid disipation rate
151  virtual tmp<volScalarField> epsilon() const;
152 
153  //- Correct Eddy-Viscosity and related properties
154  virtual void correct();
155 };
156 
157 
158 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
159 
160 } // End namespace LESModels
161 } // End namespace Foam
162 
163 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
164 
165 #ifdef NoRepository
166  #include "WALE.C"
167 #endif
168 
169 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
170 
171 #endif
172 
173 // ************************************************************************* //
virtual bool read()
Read model coefficients if they have changed.
Definition: WALE.C:156
Eddy viscosity LES SGS model base class.
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
surfaceScalarField & phi
virtual ~WALE()
Destructor.
Definition: WALE.H:134
virtual tmp< volScalarField > epsilon() const
Return sub-grid disipation rate.
Definition: WALE.C:173
BasicTurbulenceModel::rhoField rhoField
Definition: WALE.H:109
virtual void correctNut()
Update the SGS eddy-viscosity.
Definition: WALE.C:89
tmp< volSymmTensorField > Sd(const volTensorField &gradU) const
Return the deviatoric symmetric part of the square of the given.
Definition: WALE.C:40
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from string.
Definition: word.H:59
TypeName("WALE")
Runtime type information.
dimensionedScalar Cw_
Definition: WALE.H:89
dimensionedScalar Ck_
Definition: WALE.H:88
BasicTurbulenceModel::transportModel transportModel
Definition: WALE.H:110
fileName::Type type(const fileName &, const bool followLink=true)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:481
BasicTurbulenceModel::alphaField alphaField
Definition: WALE.H:108
U
Definition: pEqn.H:72
The Wall-adapting local eddy-viscosity (WALE) SGS model.
Definition: WALE.H:73
virtual void correct()
Correct Eddy-Viscosity and related properties.
Definition: WALE.C:196
A class for managing temporary objects.
Definition: PtrList.H:53
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
virtual tmp< volScalarField > k() const
Return SGS kinetic energy.
Definition: WALE.H:144
Namespace for OpenFOAM.