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