dynamicLagrangian.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) 2019 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::dynamicLagrangian
29 
30 Group
31  grpLESTurbulence
32 
33 Description
34  Dynamic SGS model with Lagrangian averaging
35 
36  Reference:
37  \verbatim
38  Meneveau, C., Lund, T. S., & Cabot, W. H. (1996).
39  A Lagrangian dynamic subgrid-scale model of turbulence.
40  Journal of Fluid Mechanics, 319, 353-385.
41  \endverbatim
42 
43 SourceFiles
44  dynamicLagrangian.C
45 
46 \*---------------------------------------------------------------------------*/
47 
48 #ifndef dynamicLagrangian_H
49 #define dynamicLagrangian_H
50 
51 #include "LESeddyViscosity.H"
52 #include "simpleFilter.H"
53 
54 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
55 
56 namespace Foam
57 {
58 namespace LESModels
59 {
60 
61 /*---------------------------------------------------------------------------*\
62  Class dynamicLagrangian Declaration
63 \*---------------------------------------------------------------------------*/
64 
65 template<class BasicTurbulenceModel>
67 :
68  public LESeddyViscosity<BasicTurbulenceModel>
69 {
70  // Private Member Functions
71 
72  //- No copy construct
73  dynamicLagrangian(const dynamicLagrangian&) = delete;
74 
75  //- No copy assignment
76  void operator=(const dynamicLagrangian&) = delete;
77 
78 
79 protected:
80 
81  // Protected data
82 
85 
87 
91 
94 
95 
96  // Protected Member Functions
97 
98  //- Update sub-grid eddy-viscosity
99  void correctNut(const tmp<volTensorField>& gradU);
100 
101  virtual void correctNut();
102 
103 
104 public:
105 
106  typedef typename BasicTurbulenceModel::alphaField alphaField;
107  typedef typename BasicTurbulenceModel::rhoField rhoField;
108  typedef typename BasicTurbulenceModel::transportModel transportModel;
110  //- Runtime type information
111  TypeName("dynamicLagrangian");
112 
113 
114  // Constructors
115 
116  //- Construct from components
118  (
119  const alphaField& alpha,
120  const rhoField& rho,
121  const volVectorField& U,
122  const surfaceScalarField& alphaRhoPhi,
123  const surfaceScalarField& phi,
124  const transportModel& transport,
125  const word& propertiesName = turbulenceModel::propertiesName,
126  const word& type = typeName
127  );
128 
129 
130  //- Destructor
131  virtual ~dynamicLagrangian() = default;
132 
133 
134  // Member Functions
135 
136  //- Read model coefficients if they have changed
137  virtual bool read();
138 
139  //- Return SGS kinetic energy
140  tmp<volScalarField> k(const tmp<volTensorField>& gradU) const
141  {
142  return
143  pow(2.0*flm_/fmm_, 2.0/3.0)
144  * pow(this->Ce_, -2.0/3.0)
145  * sqr(this->delta())*magSqr(dev(symm(gradU)));
146  }
147 
148  //- Return SGS kinetic energy
149  virtual tmp<volScalarField> k() const
150  {
151  return k(fvc::grad(this->U_));
152  }
153 
154  //- Return the effective diffusivity for k
155  tmp<volScalarField> DkEff() const
156  {
157  return tmp<volScalarField>
158  (
159  new volScalarField("DkEff", this->nut_ + this->nu())
160  );
161  }
163  //- Correct Eddy-Viscosity and related properties
164  virtual void correct();
165 };
166 
167 
168 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
169 
170 } // End namespace LESModels
171 } // End namespace Foam
172 
173 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
174 
175 #ifdef NoRepository
176  #include "dynamicLagrangian.C"
177 #endif
178 
179 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
180 
181 #endif
182 
183 // ************************************************************************* //
Eddy viscosity LES SGS model base class.
Abstract class for LES filters.
Definition: LESfilter.H:53
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:47
tmp< volScalarField > k(const tmp< volTensorField > &gradU) const
Return SGS kinetic energy.
Dynamic SGS model with Lagrangian averaging.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
TypeName("dynamicLagrangian")
Runtime type information.
dimensionedScalar Ce_
Empirical model constant.
Definition: LESModel.H:87
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:81
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
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from Foam::string.
Definition: word.H:63
const volScalarField & delta() const
Access function to filter width.
Definition: LESModel.H:245
virtual bool read()
Read model coefficients if they have changed.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
virtual ~dynamicLagrangian()=default
Destructor.
U
Definition: pEqn.H:72
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
Simple top-hat filter used in dynamic LES models.
Definition: simpleFilter.H:48
BasicTurbulenceModel::rhoField rhoField
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
A class for managing temporary objects.
Definition: HashPtrTable.H:50
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
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
Namespace for OpenFOAM.