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