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-2021 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 transportModel& transport,
48  const word& type
49 )
50 :
52  (
53  alpha,
54  rho,
55  U,
56  alphaRhoPhi,
57  phi,
58  transport,
59  type
60  ),
61 
62  gasTurbulencePtr_(nullptr),
63 
64  Cmub_
65  (
67  (
68  "Cmub",
69  this->coeffDict_,
70  0.6
71  )
72  )
73 {
74  if (type == typeName)
75  {
76  this->printCoeffs(type);
77  }
78 }
79 
80 
81 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
82 
83 template<class BasicMomentumTransportModel>
85 {
87  {
88  Cmub_.readIfPresent(this->coeffDict());
89 
90  return true;
91  }
92  else
93  {
94  return false;
95  }
96 }
97 
98 
99 template<class BasicMomentumTransportModel>
101 <
102  typename BasicMomentumTransportModel::transportModel
103 >&
105 {
106  if (!gasTurbulencePtr_)
107  {
108  const volVectorField& U = this->U_;
109 
110  const phaseModel& liquid = refCast<const phaseModel>(this->transport());
111  const phaseSystem& fluid = liquid.fluid();
112  const phaseModel& gas = fluid.otherPhase(liquid);
113 
114  gasTurbulencePtr_ =
115  &U.db().lookupObject
116  <
118  >
119  (
121  (
122  momentumTransportModel::typeName,
123  gas.name()
124  )
125  );
126  }
127 
128  return *gasTurbulencePtr_;
129 }
130 
131 
132 template<class BasicMomentumTransportModel>
134 {
136  gasTurbulence = this->gasTurbulence();
137 
138  const phaseModel& liquid = refCast<const phaseModel>(this->transport());
139  const phaseSystem& fluid = liquid.fluid();
140  const phaseModel& gas = fluid.otherPhase(liquid);
141 
142  volScalarField k(this->k(fvc::grad(this->U_)));
143 
144  this->nut_ =
145  this->Ck_*sqrt(k)*this->delta()
146  + Cmub_*gas.d()*gasTurbulence.alpha()
147  *(mag(this->U_ - gasTurbulence.U()));
148 
149  this->nut_.correctBoundaryConditions();
150  fvConstraints::New(this->mesh_).constrain(this->nut_);
151 }
152 
153 
154 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
155 
156 } // End namespace LESModels
157 } // End namespace Foam
158 
159 // ************************************************************************* //
Generic thermophysical properties class for a liquid in which the functions and coefficients for each...
Definition: liquid.H:50
Templated abstract base class for multiphase compressible turbulence models.
scalar delta
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
BasicMomentumTransportModel::rhoField rhoField
The Smagorinsky SGS model including bubble-generated turbulence.
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
dimensionedScalar sqrt(const dimensionedScalar &ds)
const word & name() const
Definition: phaseModel.H:110
label k
Boltzmann constant.
Generic dimensioned Type class.
tmp< volScalarField > d() const
Return the Sauter-mean diameter.
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
The Smagorinsky SGS model.
Definition: Smagorinsky.H:86
Class to represent a system of phases and model interfacial transfers between them.
Definition: phaseSystem.H:70
A class for handling words, derived from string.
Definition: word.H:59
static word groupName(Name name, const word &group)
phaseSystem & fluid
Definition: createFields.H:11
phi
Definition: correctPhi.H:3
static autoPtr< dictionary > New(Istream &)
Construct top-level dictionary on freestore from Istream.
Definition: dictionaryIO.C:96
SmagorinskyZhang(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &type=typeName)
Construct from components.
virtual bool read()
Read model coefficients if they have changed.
const volVectorField & U() const
Access function to velocity field.
BasicMomentumTransportModel::alphaField alphaField
U
Definition: pEqn.H:72
const phaseSystem & fluid() const
Return the system to which this phase belongs.
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
BasicMomentumTransportModel::transportModel transportModel
dimensioned< scalar > mag(const dimensioned< Type > &)
const objectRegistry & db() const
Return the local objectRegistry.
Definition: IOobject.C:312
virtual void correctNut()
Update the SGS eddy viscosity.
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phaseModel.H:53
Namespace for OpenFOAM.