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 -------------------------------------------------------------------------------
11 License
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 
34 namespace Foam
35 {
36 namespace LESModels
37 {
38 
39 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
40 
41 template<class BasicTurbulenceModel>
43 (
44  const alphaField& alpha,
45  const rhoField& rho,
46  const volVectorField& U,
47  const surfaceScalarField& alphaRhoPhi,
48  const surfaceScalarField& phi,
49  const transportModel& transport,
50  const word& propertiesName,
51  const word& type
52 )
53 :
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  (
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 
87 template<class BasicTurbulenceModel>
89 {
91  {
92  Cmub_.readIfPresent(this->coeffDict());
93 
94  return true;
95  }
96 
97  return false;
98 }
99 
100 
101 template<class BasicTurbulenceModel>
102 const PhaseCompressibleTurbulenceModel
103 <
104  typename BasicTurbulenceModel::transportModel
105 >&
106 SmagorinskyZhang<BasicTurbulenceModel>::gasTurbulence() const
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()
119  .lookupObject<PhaseCompressibleTurbulenceModel<transportModel>>
120  (
122  (
124  gas.name()
125  )
126  );
127  }
128 
129  return *gasTurbulencePtr_;
130 }
131 
132 
133 template<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 // ************************************************************************* //
twoPhaseSystem & fluid
scalar delta
void correct(GeometricField< Type, PatchField, GeoMesh > &field)
Apply correction to field.
BasicTurbulenceModel::transportModel transportModel
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:47
The Smagorinsky SGS model including bubble-generated turbulence.
virtual bool read()
Read model coefficients if they have changed.
Definition: Smagorinsky.C:122
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
const volVectorField & U() const
Access function to velocity field.
dimensionedScalar sqrt(const dimensionedScalar &ds)
label k
Boltzmann constant.
Generic dimensioned Type class.
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:82
The Smagorinsky SGS model.
Definition: Smagorinsky.H:88
const transportModel & transport() const
Access function to incompressible transport model.
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
BasicTurbulenceModel::rhoField rhoField
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from Foam::string.
Definition: word.H:63
Templated abstract base class for multiphase compressible turbulence models.
virtual bool read()
Read model coefficients if they have changed.
U
Definition: pEqn.H:72
static options & New(const fvMesh &mesh)
Construct fvOptions and register to database if not present.
Definition: fvOptions.C:96
virtual void correctNut()
Update the SGS eddy viscosity.
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
BasicTurbulenceModel::alphaField alphaField
Namespace for OpenFOAM.