SmagorinskyZhang.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) 2013-2024 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 "SmagorinskyZhang.H"
27 #include "fvModels.H"
28 #include "fvConstraints.H"
29 
30 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace LESModels
35 {
36 
37 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
38 
39 template<class BasicMomentumTransportModel>
41 (
42  const alphaField& alpha,
43  const rhoField& rho,
44  const volVectorField& U,
45  const surfaceScalarField& alphaRhoPhi,
46  const surfaceScalarField& phi,
47  const viscosity& viscosity,
48  const word& type
49 )
50 :
51  Smagorinsky<BasicMomentumTransportModel>
52  (
53  alpha,
54  rho,
55  U,
56  alphaRhoPhi,
57  phi,
58  viscosity,
59  type
60  ),
61 
62  gasTurbulencePtr_(nullptr),
63 
64  Cmub_("Cmub", this->coeffDict(), 0.6)
65 {}
66 
67 
68 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
69 
70 template<class BasicMomentumTransportModel>
72 {
74  {
75  Cmub_.readIfPresent(this->coeffDict());
76 
77  return true;
78  }
79  else
80  {
81  return false;
82  }
83 }
84 
85 
86 template<class BasicMomentumTransportModel>
89 {
90  if (!gasTurbulencePtr_)
91  {
92  const volVectorField& U = this->U_;
93 
94  const phaseModel& liquid =
95  refCast<const phaseModel>(this->properties());
96  const phaseSystem& fluid = liquid.fluid();
97  const phaseModel& gas = fluid.otherPhase(liquid);
98 
99  gasTurbulencePtr_ =
100  &U.db().lookupObject
101  <
103  >
104  (
106  (
107  momentumTransportModel::typeName,
108  gas.name()
109  )
110  );
111  }
112 
113  return *gasTurbulencePtr_;
114 }
115 
116 
117 template<class BasicMomentumTransportModel>
119 {
120  const phaseCompressibleMomentumTransportModel& gasTurbulence =
121  this->gasTurbulence();
122 
123  const phaseModel& liquid = refCast<const phaseModel>(this->properties());
124  const phaseSystem& fluid = liquid.fluid();
125  const phaseModel& gas = fluid.otherPhase(liquid);
126 
127  volScalarField k(this->k(fvc::grad(this->U_)));
128 
129  this->nut_ =
130  this->Ck_*sqrt(k)*this->delta()
131  + Cmub_*gas.d()*gasTurbulence.alpha()
132  *(mag(this->U_ - gasTurbulence.U()));
133 
134  this->nut_.correctBoundaryConditions();
135  fvConstraints::New(this->mesh_).constrain(this->nut_);
136 }
137 
138 
139 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
140 
141 } // End namespace LESModels
142 } // End namespace Foam
143 
144 // ************************************************************************* //
label k
scalar delta
Generic GeometricField class.
static word groupName(Name name, const word &group)
BasicMomentumTransportModel::alphaField alphaField
BasicMomentumTransportModel::rhoField rhoField
The Smagorinsky SGS model including bubble-generated turbulence.
SmagorinskyZhang(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.
virtual void correctNut()
Update the SGS eddy viscosity.
virtual bool read()
Read model coefficients if they have changed.
The Smagorinsky SGS model.
Definition: Smagorinsky.H:89
static autoPtr< dictionary > New(Istream &)
Construct top-level dictionary on freestore from Istream.
Definition: dictionaryIO.C:103
Generic thermophysical properties class for a liquid in which the functions and coefficients for each...
Definition: liquid.H:53
const volVectorField & U() const
Access function to velocity field.
Templated abstract base class for multiphase compressible turbulence models.
const alphaField & alpha() const
Access function to phase fraction.
tmp< volScalarField > d() const
Return the Sauter-mean diameter.
Definition: phaseModel.C:181
const word & name() const
Return the name of this phase.
Definition: phaseModel.C:145
Class to represent a system of phases.
Definition: phaseSystem.H:74
const phaseModel & otherPhase(const phaseModel &phase) const
Return the phase not given as an argument in a two-phase system.
Definition: phaseSystemI.H:174
Abstract base class for all fluid physical properties.
Definition: viscosity.H:50
A class for handling words, derived from string.
Definition: word.H:62
U
Definition: pEqn.H:72
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
tmp< VolField< typename outerProduct< vector, Type >::type > > grad(const SurfaceField< Type > &ssf)
Definition: fvcGrad.C:46
Namespace for OpenFOAM.
void mag(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
void sqrt(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488