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-2020 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  typedef typename BasicMomentumTransportModel::transportModel transportModel;
94 
95  //- Runtime type information
96  TypeName("dynamicLagrangian");
97 
98 
99  // Constructors
100 
101  //- Construct from components
103  (
104  const alphaField& alpha,
105  const rhoField& rho,
106  const volVectorField& U,
107  const surfaceScalarField& alphaRhoPhi,
108  const surfaceScalarField& phi,
109  const transportModel& transport,
110  const word& type = typeName
111  );
112 
113  //- Disallow default bitwise copy construction
114  dynamicLagrangian(const dynamicLagrangian&) = delete;
115 
116 
117  //- Destructor
118  virtual ~dynamicLagrangian()
119  {}
120 
121 
122  // Member Functions
123 
124  //- Read model coefficients if they have changed
125  virtual bool read();
126 
127  //- Return SGS kinetic energy
128  tmp<volScalarField> k(const tmp<volTensorField>& gradU) const
129  {
130  return
131  pow(2.0*flm_/fmm_, 2.0/3.0)
132  * pow(this->Ce_, -2.0/3.0)
133  * sqr(this->delta())*magSqr(dev(symm(gradU)));
134  }
135 
136  //- Return SGS kinetic energy
137  virtual tmp<volScalarField> k() const
138  {
139  return k(fvc::grad(this->U_));
140  }
141 
142  //- Return the effective diffusivity for k
143  tmp<volScalarField> DkEff() const
144  {
145  return volScalarField::New
146  (
147  "DkEff",
148  this->nut_ + this->nu()
149  );
150  }
151 
152  //- Correct Eddy-Viscosity and related properties
153  virtual void correct();
154 
155 
156  // Member Operators
157 
158  //- Disallow default bitwise assignment
159  void operator=(const dynamicLagrangian&) = delete;
160 };
161 
162 
163 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
164 
165 } // End namespace LESModels
166 } // End namespace Foam
167 
168 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
169 
170 #ifdef NoRepository
171  #include "dynamicLagrangian.C"
172 #endif
173 
174 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
175 
176 #endif
177 
178 // ************************************************************************* //
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
Dynamic SGS model with Lagrangian averaging.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
TypeName("dynamicLagrangian")
Runtime type information.
dynamicLagrangian(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &type=typeName)
Construct from components.
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Internal &, const PtrList< fvPatchField< scalar >> &)
Return a temporary field constructed from name,.
tmp< volScalarField > DkEff() const
Return the effective diffusivity for k.
phi
Definition: pEqn.H:104
virtual tmp< volScalarField > k() const
Return SGS kinetic energy.
BasicMomentumTransportModel::transportModel transportModel
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
A class for handling words, derived from string.
Definition: word.H:59
const volScalarField & delta() const
Access function to filter width.
Definition: LESModel.H:190
virtual ~dynamicLagrangian()
Destructor.
dimensioned< scalar > magSqr(const dimensioned< Type > &)
virtual bool read()
Read model coefficients if they have changed.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
U
Definition: pEqn.H:72
BasicMomentumTransportModel::rhoField rhoField
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
Simple top-hat filter used in dynamic LES models.
Definition: simpleFilter.H:50
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
A class for managing temporary objects.
Definition: PtrList.H:53
virtual void correct()
Correct Eddy-Viscosity and related properties.
BasicMomentumTransportModel::alphaField alphaField
void operator=(const dynamicLagrangian &)=delete
Disallow default bitwise assignment.
Namespace for OpenFOAM.