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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2017 OpenFOAM Foundation
9  Copyright (C) 2019-2021 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "Smagorinsky.H"
30 #include "fvOptions.H"
31 
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace LESModels
37 {
38 
39 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
40 
41 template<class BasicTurbulenceModel>
43 (
44  const tmp<volTensorField>& gradU
45 ) const
46 {
47  volSymmTensorField D(symm(gradU));
48 
49  volScalarField a(this->Ce_/this->delta());
50  volScalarField b((2.0/3.0)*tr(D));
51  volScalarField c(2*Ck_*this->delta()*(dev(D) && D));
52 
53  return tmp<volScalarField>
54  (
55  new volScalarField
56  (
57  IOobject
58  (
59  IOobject::groupName("k", this->alphaRhoPhi_.group()),
60  this->runTime_.timeName(),
61  this->mesh_
62  ),
63  sqr((-b + sqrt(sqr(b) + 4*a*c))/(2*a))
64  )
65  );
66 }
67 
68 
69 template<class BasicTurbulenceModel>
71 {
72  volScalarField k(this->k(fvc::grad(this->U_)));
73 
74  this->nut_ = Ck_*this->delta()*sqrt(k);
75  this->nut_.correctBoundaryConditions();
76  fv::options::New(this->mesh_).correct(this->nut_);
77 
78  BasicTurbulenceModel::correctNut();
79 }
80 
81 
82 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
83 
84 template<class BasicTurbulenceModel>
86 (
87  const alphaField& alpha,
88  const rhoField& rho,
89  const volVectorField& U,
90  const surfaceScalarField& alphaRhoPhi,
91  const surfaceScalarField& phi,
92  const transportModel& transport,
93  const word& propertiesName,
94  const word& type
95 )
96 :
97  LESeddyViscosity<BasicTurbulenceModel>
98  (
99  type,
100  alpha,
101  rho,
102  U,
103  alphaRhoPhi,
104  phi,
105  transport,
106  propertiesName
107  ),
108 
109  Ck_
110  (
111  dimensioned<scalar>::getOrAddToDict
112  (
113  "Ck",
114  this->coeffDict_,
115  0.094
116  )
117  )
118 {
119  if (type == typeName)
120  {
121  this->printCoeffs(type);
122  }
123 }
124 
125 
126 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
127 
128 template<class BasicTurbulenceModel>
130 {
132  {
133  Ck_.readIfPresent(this->coeffDict());
134 
135  return true;
136  }
137 
138  return false;
139 }
140 
141 
142 template<class BasicTurbulenceModel>
144 {
145  if (!this->turbulence_)
146  {
147  return;
148  }
149 
151  correctNut();
152 }
153 
154 
155 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
156 
157 } // End namespace LESModels
158 } // End namespace Foam
159 
160 // ************************************************************************* //
Eddy viscosity LES SGS model base class.
scalar delta
void correct(GeometricField< Type, PatchField, GeoMesh > &field)
Apply correction to field.
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:47
virtual bool read()
Read model coefficients if they have changed.
Definition: Smagorinsky.C:122
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar sqrt(const dimensionedScalar &ds)
virtual void printCoeffs(const word &type)
Print model coefficients.
Definition: LESModel.C:27
virtual void correctNut()
Update the SGS eddy viscosity.
Definition: Smagorinsky.C:63
label k
Boltzmann constant.
Generic dimensioned Type class.
The Smagorinsky SGS model.
Definition: Smagorinsky.H:88
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
virtual void correct()=0
Solve the turbulence equations and correct the turbulence viscosity.
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
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 Foam::string.
Definition: word.H:63
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
virtual void correct()
Correct Eddy-Viscosity and related properties.
Definition: Smagorinsky.C:136
virtual bool read()
Read model coefficients if they have changed.
virtual tmp< volScalarField > k() const
Return SGS kinetic energy.
Definition: Smagorinsky.H:174
U
Definition: pEqn.H:72
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
Base-class for all transport models used by the incompressible turbulence models. ...
const dimensionedScalar c
Speed of light in a vacuum.
const dimensionedScalar & D
A class for managing temporary objects.
Definition: HashPtrTable.H:50
static options & New(const fvMesh &mesh)
Construct fvOptions and register to database if not present.
Definition: fvOptions.C:96
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
Namespace for OpenFOAM.