WALE.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) 2015-2016 OpenFOAM Foundation
9  Copyright (C) 2019-2021 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 Class
28  Foam::LESModels::WALE
29 
30 Group
31  grpLESTurbulence
32 
33 Description
34  The Wall-adapting local eddy-viscosity (WALE) SGS model.
35 
36  Reference:
37  \verbatim
38  Nicoud, F., & Ducros, F. (1999).
39  Subgrid-scale stress modelling based on the square of the velocity
40  gradient tensor.
41  Flow, Turbulence and Combustion, 62(3), 183-200.
42  \endverbatim
43 
44  The default model coefficients are
45  \verbatim
46  WALECoeffs
47  {
48  Ck 0.094;
49  Ce 1.048;e
50  Cw 0.325;
51  }
52  \endverbatim
53 
54 See also
55  Foam::LESModels::Smagorinsky
56 
57 SourceFiles
58  WALE.C
59 
60 \*---------------------------------------------------------------------------*/
61 
62 #ifndef WALE_H
63 #define WALE_H
64 
65 #include "LESModel.H"
66 #include "LESeddyViscosity.H"
67 
68 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
69 
70 namespace Foam
71 {
72 namespace LESModels
73 {
74 
75 /*---------------------------------------------------------------------------*\
76  Class WALE Declaration
77 \*---------------------------------------------------------------------------*/
78 
79 template<class BasicTurbulenceModel>
80 class WALE
81 :
82  public LESeddyViscosity<BasicTurbulenceModel>
83 {
84  // Private Member Functions
85 
86  //- No copy construct
87  WALE(const WALE&) = delete;
88 
89  //- No copy assignment
90  void operator=(const WALE&) = delete;
91 
92 
93 protected:
94 
95  // Protected data
96 
99 
100 
101  // Protected Member Functions
102 
103  //- Return the deviatoric symmetric part of the square of the given
104  // velocity gradient field
105  tmp<volSymmTensorField> Sd(const volTensorField& gradU) const;
106 
107  //- Return SGS kinetic energy
108  // calculated from the given velocity gradient
109  tmp<volScalarField> k(const volTensorField& gradU) const;
110 
111  //- Update the SGS eddy-viscosity
112  virtual void correctNut();
113 
114 
115 public:
116 
117  typedef typename BasicTurbulenceModel::alphaField alphaField;
118  typedef typename BasicTurbulenceModel::rhoField rhoField;
119  typedef typename BasicTurbulenceModel::transportModel transportModel;
120 
121 
122  //- Runtime type information
123  TypeName("WALE");
124 
125 
126  // Constructors
127 
128  //- Construct from components
129  WALE
130  (
131  const alphaField& alpha,
132  const rhoField& rho,
133  const volVectorField& U,
134  const surfaceScalarField& alphaRhoPhi,
135  const surfaceScalarField& phi,
136  const transportModel& transport,
137  const word& propertiesName = turbulenceModel::propertiesName,
138  const word& type = typeName
139  );
140 
141 
142  //- Destructor
143  virtual ~WALE() = default;
144 
145 
146  // Member Functions
147 
148  //- Read model coefficients if they have changed
149  virtual bool read();
150 
151  //- Return SGS kinetic energy
152  virtual tmp<volScalarField> k() const
153  {
154  return k(fvc::grad(this->U_));
155  }
156 
157  //- Correct Eddy-Viscosity and related properties
158  virtual void correct();
159 };
160 
161 
162 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
163 
164 } // End namespace LESModels
165 } // End namespace Foam
166 
167 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
168 
169 #ifdef NoRepository
170  #include "WALE.C"
171 #endif
172 
173 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
174 
175 #endif
176 
177 // ************************************************************************* //
Foam::LESModels::WALE::transportModel
BasicTurbulenceModel::transportModel transportModel
Definition: WALE.H:118
Foam::fvc::grad
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:54
Foam::LESModels::WALE::Ck_
dimensionedScalar Ck_
Definition: WALE.H:96
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::LESModels::WALE::correctNut
virtual void correctNut()
Update the SGS eddy-viscosity.
Definition: WALE.C:92
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:212
LESeddyViscosity.H
Foam::LESModels::WALE::Cw_
dimensionedScalar Cw_
Definition: WALE.H:97
Foam::turbulenceModel::propertiesName
static const word propertiesName
Default name of the turbulence properties dictionary.
Definition: turbulenceModel.H:100
Foam::LESModels::LESeddyViscosity::transportModel
BasicTurbulenceModel::transportModel transportModel
Definition: LESeddyViscosity.H:75
rho
rho
Definition: readInitialConditions.H:88
Foam::LESModels::WALE::Sd
tmp< volSymmTensorField > Sd(const volTensorField &gradU) const
Return the deviatoric symmetric part of the square of the given.
Definition: WALE.C:43
WALE.C
LESModel.H
Foam::LESModels::WALE::k
virtual tmp< volScalarField > k() const
Return SGS kinetic energy.
Definition: WALE.H:151
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::LESModels::LESeddyViscosity
Eddy viscosity LES SGS model base class.
Definition: LESeddyViscosity.H:58
Foam::LESModels::WALE::TypeName
TypeName("WALE")
Runtime type information.
Foam::dimensioned< scalar >
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::LESModels::WALE::rhoField
BasicTurbulenceModel::rhoField rhoField
Definition: WALE.H:117
Foam::LESModels::WALE::alphaField
BasicTurbulenceModel::alphaField alphaField
Definition: WALE.H:116
Foam::LESModels::LESeddyViscosity::rhoField
BasicTurbulenceModel::rhoField rhoField
Definition: LESeddyViscosity.H:74
U
U
Definition: pEqn.H:72
Foam::LESModels::WALE::read
virtual bool read()
Read model coefficients if they have changed.
Definition: WALE.C:159
Foam::type
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:590
Foam::LESModels::WALE
The Wall-adapting local eddy-viscosity (WALE) SGS model.
Definition: WALE.H:79
Foam::LESModels::WALE::correct
virtual void correct()
Correct Eddy-Viscosity and related properties.
Definition: WALE.C:174
Foam::LESModels::WALE::~WALE
virtual ~WALE()=default
Destructor.
Foam::GeometricField< tensor, fvPatchField, volMesh >
Foam::LESModels::LESeddyViscosity::alphaField
BasicTurbulenceModel::alphaField alphaField
Definition: LESeddyViscosity.H:73