Smagorinsky.C
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 \*---------------------------------------------------------------------------*/
25 
26 #include "Smagorinsky.H"
27 #include "fvOptions.H"
28 
29 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace LESModels
34 {
35 
36 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
37 
38 template<class BasicTurbulenceModel>
39 tmp<volScalarField> Smagorinsky<BasicTurbulenceModel>::k
40 (
41  const tmp<volTensorField>& gradU
42 ) const
43 {
44  volSymmTensorField D(symm(gradU));
45 
46  volScalarField a(this->Ce_/this->delta());
47  volScalarField b((2.0/3.0)*tr(D));
48  volScalarField c(2*Ck_*this->delta()*(dev(D) && D));
49 
50  return tmp<volScalarField>
51  (
52  new volScalarField
53  (
54  IOobject
55  (
56  IOobject::groupName("k", this->U_.group()),
57  this->runTime_.timeName(),
58  this->mesh_
59  ),
60  sqr((-b + sqrt(sqr(b) + 4*a*c))/(2*a))
61  )
62  );
63 }
64 
65 
66 template<class BasicTurbulenceModel>
68 {
69  volScalarField k(this->k(fvc::grad(this->U_)));
70 
71  this->nut_ = Ck_*this->delta()*sqrt(k);
72  this->nut_.correctBoundaryConditions();
73  fv::options::New(this->mesh_).correct(this->nut_);
74 
75  BasicTurbulenceModel::correctNut();
76 }
77 
78 
79 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
80 
81 template<class BasicTurbulenceModel>
83 (
84  const alphaField& alpha,
85  const rhoField& rho,
86  const volVectorField& U,
87  const surfaceScalarField& alphaRhoPhi,
88  const surfaceScalarField& phi,
89  const transportModel& transport,
90  const word& propertiesName,
91  const word& type
92 )
93 :
95  (
96  type,
97  alpha,
98  rho,
99  U,
100  alphaRhoPhi,
101  phi,
102  transport,
103  propertiesName
104  ),
105 
106  Ck_
107  (
109  (
110  "Ck",
111  this->coeffDict_,
112  0.094
113  )
114  )
115 {
116  if (type == typeName)
117  {
118  this->printCoeffs(type);
119  }
120 }
121 
122 
123 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
124 
125 template<class BasicTurbulenceModel>
127 {
129  {
130  Ck_.readIfPresent(this->coeffDict());
131 
132  return true;
133  }
134  else
135  {
136  return false;
137  }
138 }
139 
140 
141 template<class BasicTurbulenceModel>
143 {
144  volScalarField k(this->k(fvc::grad(this->U_)));
145 
146  return tmp<volScalarField>
147  (
148  new volScalarField
149  (
150  IOobject
151  (
152  IOobject::groupName("epsilon", this->U_.group()),
153  this->runTime_.timeName(),
154  this->mesh_
155  ),
156  this->Ce_*k*sqrt(k)/this->delta()
157  )
158  );
159 }
160 
161 
162 template<class BasicTurbulenceModel>
164 {
166  correctNut();
167 }
168 
169 
170 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
171 
172 } // End namespace LESModels
173 } // End namespace Foam
174 
175 // ************************************************************************* //
Eddy viscosity LES SGS model base class.
scalar delta
surfaceScalarField & phi
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
U
Definition: pEqn.H:83
virtual bool read()
Read model coefficients if they have changed.
Definition: Smagorinsky.C:126
BasicTurbulenceModel::transportModel transportModel
void correct(GeometricField< Type, fvPatchField, volMesh > &field)
Apply correction to field.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar sqrt(const dimensionedScalar &ds)
virtual void correctNut()
Update the SGS eddy viscosity.
Definition: Smagorinsky.C:67
label k
Boltzmann constant.
Generic dimensioned Type class.
BasicTurbulenceModel::alphaField alphaField
The Smagorinsky SGS model.
Definition: Smagorinsky.H:89
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
virtual void correct()=0
Solve the turbulence equations and correct the turbulence viscosity.
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
A class for handling words, derived from string.
Definition: word.H:59
static word groupName(Name name, const word &group)
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
virtual tmp< volScalarField > epsilon() const
Return sub-grid disipation rate.
Definition: Smagorinsky.C:142
virtual void correct()
Correct Eddy-Viscosity and related properties.
Definition: Smagorinsky.C:163
virtual tmp< volScalarField > k() const
Return SGS kinetic energy.
Definition: Smagorinsky.H:155
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:461
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
const dimensionedScalar c
Speed of light in a vacuum.
BasicTurbulenceModel::rhoField rhoField
A class for managing temporary objects.
Definition: PtrList.H:54
static options & New(const fvMesh &mesh)
Construct fvOptions and register to datbase if not present.
Definition: fvOptions.C:103
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Namespace for OpenFOAM.