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-2021 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 BasicMomentumTransportModel>
61 :
62  public LESeddyViscosity<BasicMomentumTransportModel>
63 {
64 protected:
65 
66  // Protected data
67 
70 
72 
76 
79 
80 
81  // Protected Member Functions
82 
83  //- Update sub-grid eddy-viscosity
84  void correctNut(const tmp<volTensorField>& gradU);
85 
86  virtual void correctNut();
87 
88 
89 public:
90 
91  typedef typename BasicMomentumTransportModel::alphaField alphaField;
92  typedef typename BasicMomentumTransportModel::rhoField rhoField;
93 
94  //- Runtime type information
95  TypeName("dynamicLagrangian");
96 
97 
98  // Constructors
99 
100  //- Construct from components
102  (
103  const alphaField& alpha,
104  const rhoField& rho,
105  const volVectorField& U,
106  const surfaceScalarField& alphaRhoPhi,
107  const surfaceScalarField& phi,
108  const viscosity& viscosity,
109  const word& type = typeName
110  );
111 
112  //- Disallow default bitwise copy construction
113  dynamicLagrangian(const dynamicLagrangian&) = delete;
114 
115 
116  //- Destructor
117  virtual ~dynamicLagrangian()
118  {}
119 
120 
121  // Member Functions
122 
123  //- Read model coefficients if they have changed
124  virtual bool read();
125 
126  //- Return SGS kinetic energy
127  tmp<volScalarField> k(const tmp<volTensorField>& gradU) const
128  {
129  return
130  pow(2.0*flm_/fmm_, 2.0/3.0)
131  * pow(this->Ce_, -2.0/3.0)
132  * sqr(this->delta())*magSqr(dev(symm(gradU)));
133  }
134 
135  //- Return SGS kinetic energy
136  virtual tmp<volScalarField> k() const
137  {
138  return k(fvc::grad(this->U_));
139  }
140 
141  //- Return the effective diffusivity for k
142  tmp<volScalarField> DkEff() const
143  {
144  return volScalarField::New
145  (
146  "DkEff",
147  this->nut_ + this->nu()
148  );
149  }
150 
151  //- Correct Eddy-Viscosity and related properties
152  virtual void correct();
153 
154 
155  // Member Operators
156 
157  //- Disallow default bitwise assignment
158  void operator=(const dynamicLagrangian&) = delete;
159 };
160 
161 
162 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
163 
164 } // End namespace LESModels
165 } // End namespace Foam
166 
167 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
168 
169 #ifdef NoRepository
170  #include "dynamicLagrangian.C"
171 #endif
172 
173 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
174 
175 #endif
176 
177 // ************************************************************************* //
scalar delta
Generic GeometricField class.
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &)
Return a temporary field constructed from name,.
Eddy viscosity LES SGS model base class.
Dynamic SGS model with Lagrangian averaging.
BasicMomentumTransportModel::alphaField alphaField
tmp< volScalarField > DkEff() const
Return the effective diffusivity for k.
virtual tmp< volScalarField > k() const
Return SGS kinetic energy.
TypeName("dynamicLagrangian")
Runtime type information.
virtual void correct()
Correct Eddy-Viscosity and related properties.
void operator=(const dynamicLagrangian &)=delete
Disallow default bitwise assignment.
dynamicLagrangian(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &viscosity, const word &type=typeName)
Construct from components.
virtual bool read()
Read model coefficients if they have changed.
BasicMomentumTransportModel::rhoField rhoField
Abstract class for LES filters.
Definition: LESfilter.H:55
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
Simple top-hat filter used in dynamic LES models.
Definition: simpleFilter.H:54
A class for managing temporary objects.
Definition: tmp.H:55
Abstract base class for all fluid physical properties.
Definition: viscosity.H:50
A class for handling words, derived from string.
Definition: word.H:62
U
Definition: pEqn.H:72
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
tmp< VolField< typename outerProduct< vector, Type >::type > > grad(const SurfaceField< Type > &ssf)
Definition: fvcGrad.C:46
Namespace for OpenFOAM.
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensioned< scalar > magSqr(const dimensioned< Type > &)
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488