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-2020 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 BasicMomentumTransportModel>
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 BasicMomentumTransportModel>
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 
68 
69 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
70 
71 template<class BasicMomentumTransportModel>
73 (
74  const alphaField& alpha,
75  const rhoField& rho,
76  const volVectorField& U,
77  const surfaceScalarField& alphaRhoPhi,
78  const surfaceScalarField& phi,
79  const transportModel& transport,
80  const word& type
81 )
82 :
84  (
85  type,
86  alpha,
87  rho,
88  U,
89  alphaRhoPhi,
90  phi,
91  transport
92  ),
93 
94  Ck_
95  (
97  (
98  "Ck",
99  this->coeffDict_,
100  0.094
101  )
102  )
103 {
104  if (type == typeName)
105  {
106  this->printCoeffs(type);
107  }
108 }
109 
110 
111 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
112 
113 template<class BasicMomentumTransportModel>
115 {
117  {
118  Ck_.readIfPresent(this->coeffDict());
119 
120  return true;
121  }
122  else
123  {
124  return false;
125  }
126 }
127 
128 
129 template<class BasicMomentumTransportModel>
131 {
132  volScalarField k(this->k(fvc::grad(this->U_)));
133 
134  return volScalarField::New
135  (
136  IOobject::groupName("epsilon", this->alphaRhoPhi_.group()),
137  this->Ce_*k*sqrt(k)/this->delta()
138  );
139 }
140 
141 
142 template<class BasicMomentumTransportModel>
144 {
146  correctNut();
147 }
148 
149 
150 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
151 
152 } // End namespace LESModels
153 } // End namespace Foam
154 
155 // ************************************************************************* //
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:130
BasicMomentumTransportModel::rhoField rhoField
virtual bool read()
Read model coefficients if they have changed.
Definition: Smagorinsky.C:114
void correct(GeometricField< Type, fvPatchField, volMesh > &field)
Apply correction to field.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
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)
virtual void correctNut()
Update the SGS eddy viscosity.
Definition: Smagorinsky.C:59
label k
Boltzmann constant.
Generic dimensioned Type class.
const dimensionedScalar & c
Speed of light in a vacuum.
phi
Definition: pEqn.H:104
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:73
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:143
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
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.