SmagorinskyZhang.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) 2013-2016 OpenFOAM Foundation
9 Copyright (C) 2019-2020 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 "SmagorinskyZhang.H"
30#include "fvOptions.H"
31
32// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33
34namespace Foam
35{
36namespace LESModels
37{
38
39// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
40
41template<class BasicTurbulenceModel>
43(
44 const alphaField& alpha,
45 const rhoField& rho,
46 const volVectorField& U,
47 const surfaceScalarField& alphaRhoPhi,
49 const transportModel& transport,
50 const word& propertiesName,
51 const word& type
52)
53:
54 Smagorinsky<BasicTurbulenceModel>
55 (
56 alpha,
57 rho,
58 U,
59 alphaRhoPhi,
60 phi,
61 transport,
62 propertiesName,
63 type
64 ),
65
66 gasTurbulencePtr_(nullptr),
67
68 Cmub_
69 (
70 dimensioned<scalar>::getOrAddToDict
71 (
72 "Cmub",
73 this->coeffDict_,
74 0.6
75 )
76 )
77{
78 if (type == typeName)
79 {
80 this->printCoeffs(type);
81 }
82}
83
84
85// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
86
87template<class BasicTurbulenceModel>
89{
91 {
92 Cmub_.readIfPresent(this->coeffDict());
93
94 return true;
95 }
96
97 return false;
98}
99
100
101template<class BasicTurbulenceModel>
103<
104 typename BasicTurbulenceModel::transportModel
105>&
107{
108 if (!gasTurbulencePtr_)
109 {
110 const volVectorField& U = this->U_;
111
112 const transportModel& liquid = this->transport();
113 const twoPhaseSystem& fluid =
114 refCast<const twoPhaseSystem>(liquid.fluid());
115 const transportModel& gas = fluid.otherPhase(liquid);
116
117 gasTurbulencePtr_ =
118 &U.db()
120 (
122 (
124 gas.name()
125 )
126 );
127 }
128
129 return *gasTurbulencePtr_;
130}
131
132
133template<class BasicTurbulenceModel>
135{
137 this->gasTurbulence();
138
139 volScalarField k(this->k(fvc::grad(this->U_)));
140
141 this->nut_ =
142 this->Ck_*sqrt(k)*this->delta()
143 + Cmub_*gasTurbulence.transport().d()*gasTurbulence.alpha()
144 *(mag(this->U_ - gasTurbulence.U()));
145
146 this->nut_.correctBoundaryConditions();
147 fv::options::New(this->mesh_).correct(this->nut_);
148
149 BasicTurbulenceModel::correctNut();
150}
151
152
153// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
154
155} // End namespace LESModels
156} // End namespace Foam
157
158// ************************************************************************* //
label k
scalar delta
surfaceScalarField & phi
twoPhaseSystem & fluid
void correctBoundaryConditions()
Correct boundary field.
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
The Smagorinsky SGS model including bubble-generated turbulence.
BasicTurbulenceModel::alphaField alphaField
BasicTurbulenceModel::rhoField rhoField
virtual void correctNut()
Update the SGS eddy viscosity.
BasicTurbulenceModel::transportModel transportModel
virtual bool read()
Read model coefficients if they have changed.
The Smagorinsky SGS model.
Definition: Smagorinsky.H:95
Templated abstract base class for multiphase compressible turbulence models.
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
const alphaField & alpha() const
Access function to phase fraction.
const transportModel & transport() const
Access function to incompressible transport model.
Generic dimensioned Type class.
Generic thermophysical properties class for a liquid in which the functions and coefficients for each...
Definition: liquid.H:57
Base-class for all transport models used by the incompressible turbulence models.
const volVectorField & U() const
Access function to velocity field.
static const word propertiesName
Default name of the turbulence properties dictionary.
Class which solves the volume fraction equations for two phases.
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.
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)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
volScalarField & alpha