Smagorinsky.C
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-2018 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 volScalarField::New
51  (
52  IOobject::groupName("k", this->alphaRhoPhi_.group()),
53  sqr((-b + sqrt(sqr(b) + 4*a*c))/(2*a))
54  );
55 }
56 
57 
58 template<class BasicTurbulenceModel>
60 {
61  volScalarField k(this->k(fvc::grad(this->U_)));
62 
63  this->nut_ = Ck_*this->delta()*sqrt(k);
64  this->nut_.correctBoundaryConditions();
65  fv::options::New(this->mesh_).correct(this->nut_);
66 
67  BasicTurbulenceModel::correctNut();
68 }
69 
70 
71 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
72 
73 template<class BasicTurbulenceModel>
75 (
76  const alphaField& alpha,
77  const rhoField& rho,
78  const volVectorField& U,
79  const surfaceScalarField& alphaRhoPhi,
80  const surfaceScalarField& phi,
81  const transportModel& transport,
82  const word& propertiesName,
83  const word& type
84 )
85 :
87  (
88  type,
89  alpha,
90  rho,
91  U,
92  alphaRhoPhi,
93  phi,
94  transport,
95  propertiesName
96  ),
97 
98  Ck_
99  (
101  (
102  "Ck",
103  this->coeffDict_,
104  0.094
105  )
106  )
107 {
108  if (type == typeName)
109  {
110  this->printCoeffs(type);
111  }
112 }
113 
114 
115 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
116 
117 template<class BasicTurbulenceModel>
119 {
121  {
122  Ck_.readIfPresent(this->coeffDict());
123 
124  return true;
125  }
126  else
127  {
128  return false;
129  }
130 }
131 
132 
133 template<class BasicTurbulenceModel>
135 {
136  volScalarField k(this->k(fvc::grad(this->U_)));
137 
138  return volScalarField::New
139  (
140  IOobject::groupName("epsilon", this->alphaRhoPhi_.group()),
141  this->Ce_*k*sqrt(k)/this->delta()
142  );
143 }
144 
145 
146 template<class BasicTurbulenceModel>
148 {
150  correctNut();
151 }
152 
153 
154 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
155 
156 } // End namespace LESModels
157 } // End namespace Foam
158 
159 // ************************************************************************* //
Eddy viscosity LES SGS model base class.
scalar delta
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
virtual tmp< volScalarField > epsilon() const
Return sub-grid disipation rate.
Definition: Smagorinsky.C:134
surfaceScalarField & phi
virtual bool read()
Read model coefficients if they have changed.
Definition: Smagorinsky.C:118
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:59
label k
Boltzmann constant.
Generic dimensioned Type class.
BasicTurbulenceModel::alphaField alphaField
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Mesh &, const dimensionSet &, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
Return a temporary field constructed from name, mesh, dimensionSet.
Smagorinsky(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName=turbulenceModel::propertiesName, const word &type=typeName)
Construct from components.
Definition: Smagorinsky.C:75
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 void correct()
Correct Eddy-Viscosity and related properties.
Definition: Smagorinsky.C:147
virtual tmp< volScalarField > k() const
Return SGS kinetic energy.
Definition: Smagorinsky.H:148
U
Definition: pEqn.H:72
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
const dimensionedScalar c
Speed of light in a vacuum.
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
BasicTurbulenceModel::rhoField rhoField
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
static options & New(const fvMesh &mesh)
Construct fvOptions and register to datbase if not present.
Definition: fvOptions.C:101
Namespace for OpenFOAM.