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