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-2023 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*this->Ck_*this->delta()*(dev(D) && D));
50 
51  return volScalarField::New
52  (
53  this->groupName("k"),
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_ = this->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 viscosity& viscosity,
81  const word& type
82 )
83 :
84  LESeddyViscosity<BasicMomentumTransportModel>
85  (
86  type,
87  alpha,
88  rho,
89  U,
90  alphaRhoPhi,
91  phi,
92  viscosity
93  )
94 {
95  if (type == typeName)
96  {
97  this->printCoeffs(type);
98  }
99 }
100 
101 
102 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
103 
104 template<class BasicMomentumTransportModel>
106 {
108 }
109 
110 
111 template<class BasicMomentumTransportModel>
113 {
115  correctNut();
116 }
117 
118 
119 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
120 
121 } // End namespace LESModels
122 } // End namespace Foam
123 
124 // ************************************************************************* //
static const Foam::dimensionedScalar D("D", Foam::dimTemperature, 257.14)
label k
scalar delta
Generic GeometricField class.
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &)
Return a temporary field constructed from name,.
Eddy viscosity LES SGS model base class.
BasicMomentumTransportModel::alphaField alphaField
virtual bool read()
Read model coefficients if they have changed.
BasicMomentumTransportModel::rhoField rhoField
Smagorinsky(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &viscosity, const word &type=typeName)
Construct from components.
Definition: Smagorinsky.C:74
virtual tmp< volScalarField > k() const
Return SGS kinetic energy.
Definition: Smagorinsky.H:141
virtual void correct()
Correct Eddy-Viscosity and related properties.
Definition: Smagorinsky.C:112
virtual void correctNut()
Update the SGS eddy viscosity.
Definition: Smagorinsky.C:60
virtual bool read()
Read model coefficients if they have changed.
Definition: Smagorinsky.C:105
static autoPtr< dictionary > New(Istream &)
Construct top-level dictionary on freestore from Istream.
Definition: dictionaryIO.C:100
virtual void correct()=0
Solve the turbulence equations and correct the turbulence viscosity.
A class for managing temporary objects.
Definition: tmp.H:55
Abstract base class for all fluid physical properties.
Definition: viscosity.H:50
A class for handling words, derived from string.
Definition: word.H:62
volScalarField & b
Definition: createFields.H:27
U
Definition: pEqn.H:72
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
const dimensionedScalar c
Speed of light in a vacuum.
tmp< VolField< typename outerProduct< vector, Type >::type > > grad(const SurfaceField< Type > &ssf)
Definition: fvcGrad.C:46
Namespace for OpenFOAM.
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar sqrt(const dimensionedScalar &ds)
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488