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-------------------------------------------------------------------------------
11License
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
34namespace Foam
35{
36namespace LESModels
37{
38
39// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
40
41template<class BasicTurbulenceModel>
43(
44 const tmp<volTensorField>& gradU
45) const
46{
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
54 (
56 (
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
69template<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
84template<class BasicTurbulenceModel>
86(
87 const alphaField& alpha,
88 const rhoField& rho,
89 const volVectorField& U,
90 const surfaceScalarField& alphaRhoPhi,
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
128template<class BasicTurbulenceModel>
130{
132 {
133 Ck_.readIfPresent(this->coeffDict());
134
135 return true;
136 }
137
138 return false;
139}
140
141
142template<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// ************************************************************************* //
label k
scalar delta
surfaceScalarField & phi
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
Eddy viscosity LES SGS model base class.
The Smagorinsky SGS model.
Definition: Smagorinsky.H:95
BasicTurbulenceModel::alphaField alphaField
Definition: Smagorinsky.H:124
BasicTurbulenceModel::rhoField rhoField
Definition: Smagorinsky.H:125
virtual void correct()
Correct Eddy-Viscosity and related properties.
Definition: Smagorinsky.C:143
virtual void correctNut()
Update the SGS eddy viscosity.
Definition: Smagorinsky.C:70
BasicTurbulenceModel::transportModel transportModel
Definition: Smagorinsky.H:126
virtual tmp< volScalarField > k() const
Return SGS kinetic energy.
Definition: Smagorinsky.H:159
virtual bool read()
Read model coefficients if they have changed.
Definition: Smagorinsky.C:129
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
Generic dimensioned Type class.
A class for managing temporary objects.
Definition: tmp.H:65
A class for handling words, derived from Foam::string.
Definition: word.H:68
U
Definition: pEqn.H:72
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:54
Namespace for OpenFOAM.
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:598
dimensionedScalar sqrt(const dimensionedScalar &ds)
volScalarField & alpha
const dimensionedScalar & D
volScalarField & b
Definition: createFields.H:27