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-2021 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 "fvModels.H"
28 #include "fvConstraints.H"
29 
30 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace LESModels
35 {
36 
37 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
38 
39 template<class BasicMomentumTransportModel>
41 (
42  const tmp<volTensorField>& gradU
43 ) const
44 {
45  volSymmTensorField D(symm(gradU));
46 
47  volScalarField a(this->Ce_/this->delta());
48  volScalarField b((2.0/3.0)*tr(D));
49  volScalarField c(2*Ck_*this->delta()*(dev(D) && D));
50 
51  return volScalarField::New
52  (
53  IOobject::groupName("k", this->alphaRhoPhi_.group()),
54  sqr((-b + sqrt(sqr(b) + 4*a*c))/(2*a))
55  );
56 }
57 
58 
59 template<class BasicMomentumTransportModel>
61 {
62  volScalarField k(this->k(fvc::grad(this->U_)));
63 
64  this->nut_ = Ck_*this->delta()*sqrt(k);
65  this->nut_.correctBoundaryConditions();
66  fvConstraints::New(this->mesh_).constrain(this->nut_);
67 }
68 
69 
70 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
71 
72 template<class BasicMomentumTransportModel>
74 (
75  const alphaField& alpha,
76  const rhoField& rho,
77  const volVectorField& U,
78  const surfaceScalarField& alphaRhoPhi,
79  const surfaceScalarField& phi,
80  const transportModel& transport,
81  const word& type
82 )
83 :
85  (
86  type,
87  alpha,
88  rho,
89  U,
90  alphaRhoPhi,
91  phi,
92  transport
93  ),
94 
95  Ck_
96  (
98  (
99  "Ck",
100  this->coeffDict_,
101  0.094
102  )
103  )
104 {
105  if (type == typeName)
106  {
107  this->printCoeffs(type);
108  }
109 }
110 
111 
112 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
113 
114 template<class BasicMomentumTransportModel>
116 {
118  {
119  Ck_.readIfPresent(this->coeffDict());
120 
121  return true;
122  }
123  else
124  {
125  return false;
126  }
127 }
128 
129 
130 template<class BasicMomentumTransportModel>
132 {
133  volScalarField k(this->k(fvc::grad(this->U_)));
134 
135  return volScalarField::New
136  (
137  IOobject::groupName("epsilon", this->alphaRhoPhi_.group()),
138  this->Ce_*k*sqrt(k)/this->delta()
139  );
140 }
141 
142 
143 template<class BasicMomentumTransportModel>
145 {
147  correctNut();
148 }
149 
150 
151 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
152 
153 } // End namespace LESModels
154 } // End namespace Foam
155 
156 // ************************************************************************* //
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:131
BasicMomentumTransportModel::rhoField rhoField
virtual bool read()
Read model coefficients if they have changed.
Definition: Smagorinsky.C:115
dimensionedSymmTensor sqr(const dimensionedVector &dv)
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Internal &, const PtrList< fvPatchField< scalar >> &)
Return a temporary field constructed from name,.
dimensionedScalar sqrt(const dimensionedScalar &ds)
const dimensionedScalar b
Wien displacement law constant: default SI units: [m K].
Definition: createFields.H:27
virtual void correctNut()
Update the SGS eddy viscosity.
Definition: Smagorinsky.C:60
label k
Boltzmann constant.
Generic dimensioned Type class.
const dimensionedScalar c
Speed of light in a vacuum.
virtual void correct()=0
Solve the turbulence equations and correct the turbulence viscosity.
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
Smagorinsky(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &type=typeName)
Construct from components.
Definition: Smagorinsky.C:74
A class for handling words, derived from string.
Definition: word.H:59
static word groupName(Name name, const word &group)
phi
Definition: correctPhi.H:3
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
virtual void correct()
Correct Eddy-Viscosity and related properties.
Definition: Smagorinsky.C:144
static autoPtr< dictionary > New(Istream &)
Construct top-level dictionary on freestore from Istream.
Definition: dictionaryIO.C:96
virtual tmp< volScalarField > k() const
Return SGS kinetic energy.
Definition: Smagorinsky.H:147
BasicMomentumTransportModel::alphaField alphaField
U
Definition: pEqn.H:72
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
BasicMomentumTransportModel::transportModel transportModel
A class for managing temporary objects.
Definition: PtrList.H:53
Namespace for OpenFOAM.