dynamicLagrangian.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-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::dynamicLagrangian
26 
27 Description
28  Dynamic SGS model with Lagrangian averaging
29 
30  Reference:
31  \verbatim
32  Meneveau, C., Lund, T. S., & Cabot, W. H. (1996).
33  A Lagrangian dynamic subgrid-scale model of turbulence.
34  Journal of Fluid Mechanics, 319, 353-385.
35  \endverbatim
36 
37 SourceFiles
38  dynamicLagrangian.C
39 
40 \*---------------------------------------------------------------------------*/
41 
42 #ifndef dynamicLagrangian_H
43 #define dynamicLagrangian_H
44 
45 #include "LESeddyViscosity.H"
46 #include "simpleFilter.H"
47 
48 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
49 
50 namespace Foam
51 {
52 namespace LESModels
53 {
54 
55 /*---------------------------------------------------------------------------*\
56  Class dynamicLagrangian Declaration
57 \*---------------------------------------------------------------------------*/
58 
59 template<class BasicTurbulenceModel>
61 :
62  public LESeddyViscosity<BasicTurbulenceModel>
63 {
64  // Private Member Functions
65 
66  // Disallow default bitwise copy construct and assignment
68  void operator=(const dynamicLagrangian&);
69 
70 
71 protected:
72 
73  // Protected data
74 
77 
79 
83 
86 
87 
88  // Protected Member Functions
89 
90  //- Update sub-grid eddy-viscosity
91  void correctNut(const tmp<volTensorField>& gradU);
92 
93  virtual void correctNut();
94 
95 
96 public:
97 
98  typedef typename BasicTurbulenceModel::alphaField alphaField;
99  typedef typename BasicTurbulenceModel::rhoField rhoField;
100  typedef typename BasicTurbulenceModel::transportModel transportModel;
101 
102  //- Runtime type information
103  TypeName("dynamicLagrangian");
104 
105 
106  // Constructors
107 
108  //- Construct from components
110  (
111  const alphaField& alpha,
112  const rhoField& rho,
113  const volVectorField& U,
114  const surfaceScalarField& alphaRhoPhi,
115  const surfaceScalarField& phi,
116  const transportModel& transport,
117  const word& propertiesName = turbulenceModel::propertiesName,
118  const word& type = typeName
119  );
120 
121 
122  //- Destructor
123  virtual ~dynamicLagrangian()
124  {}
125 
126 
127  // Member Functions
128 
129  //- Read model coefficients if they have changed
130  virtual bool read();
131 
132  //- Return SGS kinetic energy
133  tmp<volScalarField> k(const tmp<volTensorField>& gradU) const
134  {
135  return
136  pow(2.0*flm_/fmm_, 2.0/3.0)
137  * pow(this->Ce_, -2.0/3.0)
138  * sqr(this->delta())*magSqr(dev(symm(gradU)));
139  }
140 
141  //- Return SGS kinetic energy
142  virtual tmp<volScalarField> k() const
143  {
144  return k(fvc::grad(this->U_));
145  }
146 
147  //- Return the effective diffusivity for k
148  tmp<volScalarField> DkEff() const
149  {
150  return tmp<volScalarField>
151  (
152  new volScalarField("DkEff", this->nut_ + this->nu())
153  );
154  }
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 "dynamicLagrangian.C"
170 #endif
171 
172 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
173 
174 #endif
175 
176 // ************************************************************************* //
Eddy viscosity LES SGS model base class.
Abstract class for LES filters.
Definition: LESfilter.H:54
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
surfaceScalarField & phi
Dynamic SGS model with Lagrangian averaging.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
TypeName("dynamicLagrangian")
Runtime type information.
tmp< volScalarField > DkEff() const
Return the effective diffusivity for k.
virtual tmp< volScalarField > k() const
Return SGS kinetic energy.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from string.
Definition: word.H:59
const volScalarField & delta() const
Access function to filter width.
Definition: LESModel.H:201
virtual ~dynamicLagrangian()
Destructor.
dimensioned< scalar > magSqr(const dimensioned< Type > &)
virtual bool read()
Read model coefficients if they have changed.
fileName::Type type(const fileName &, const bool followLink=true)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:481
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
U
Definition: pEqn.H:72
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
Simple top-hat filter used in dynamic LES models.
Definition: simpleFilter.H:50
BasicTurbulenceModel::rhoField rhoField
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
A class for managing temporary objects.
Definition: PtrList.H:53
virtual void correct()
Correct Eddy-Viscosity and related properties.
BasicTurbulenceModel::alphaField alphaField
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
volScalarField & nu
BasicTurbulenceModel::transportModel transportModel
Namespace for OpenFOAM.