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-2020 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 "fvOptions.H"
28 
29 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace LESModels
34 {
35 
36 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
37 
38 template<class BasicMomentumTransportModel>
40 (
41  const alphaField& alpha,
42  const rhoField& rho,
43  const volVectorField& U,
44  const surfaceScalarField& alphaRhoPhi,
45  const surfaceScalarField& phi,
46  const transportModel& transport,
47  const word& type
48 )
49 :
51  (
52  alpha,
53  rho,
54  U,
55  alphaRhoPhi,
56  phi,
57  transport,
58  type
59  ),
60 
61  gasTurbulencePtr_(nullptr),
62 
63  Cmub_
64  (
66  (
67  "Cmub",
68  this->coeffDict_,
69  0.6
70  )
71  )
72 {
73  if (type == typeName)
74  {
75  this->printCoeffs(type);
76  }
77 }
78 
79 
80 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
81 
82 template<class BasicMomentumTransportModel>
84 {
86  {
87  Cmub_.readIfPresent(this->coeffDict());
88 
89  return true;
90  }
91  else
92  {
93  return false;
94  }
95 }
96 
97 
98 template<class BasicMomentumTransportModel>
100 <
101  typename BasicMomentumTransportModel::transportModel
102 >&
104 {
105  if (!gasTurbulencePtr_)
106  {
107  const volVectorField& U = this->U_;
108 
109  const transportModel& liquid = this->transport();
110  const phaseSystem& fluid = liquid.fluid();
111  const transportModel& gas = fluid.otherPhase(liquid);
112 
113  gasTurbulencePtr_ =
114  &U.db().lookupObject
115  <
117  >
118  (
120  (
121  momentumTransportModel::typeName,
122  gas.name()
123  )
124  );
125  }
126 
127  return *gasTurbulencePtr_;
128 }
129 
130 
131 template<class BasicMomentumTransportModel>
133 {
135  gasTurbulence = this->gasTurbulence();
136 
137  volScalarField k(this->k(fvc::grad(this->U_)));
138 
139  this->nut_ =
140  this->Ck_*sqrt(k)*this->delta()
141  + Cmub_*gasTurbulence.transport().d()*gasTurbulence.alpha()
142  *(mag(this->U_ - gasTurbulence.U()));
143 
144  this->nut_.correctBoundaryConditions();
145  fv::options::New(this->mesh_).correct(this->nut_);
146 }
147 
148 
149 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
150 
151 } // End namespace LESModels
152 } // End namespace Foam
153 
154 // ************************************************************************* //
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
const transportModel & transport() const
Access function to incompressible transport model.
The Smagorinsky SGS model including bubble-generated turbulence.
void correct(GeometricField< Type, fvPatchField, volMesh > &field)
Apply correction to field.
dimensionedScalar sqrt(const dimensionedScalar &ds)
label k
Boltzmann constant.
Generic dimensioned Type class.
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
The Smagorinsky SGS model.
Definition: Smagorinsky.H:86
phi
Definition: pEqn.H:104
Class to represent a system of phases and model interfacial transfers between them.
Definition: phaseSystem.H:67
A class for handling words, derived from string.
Definition: word.H:59
static word groupName(Name name, const word &group)
const phaseModel & otherPhase(const phaseModel &phase) const
Return the phase not given as an argument in a two-phase system.
Definition: phaseSystemI.H:112
phaseSystem & fluid
Definition: createFields.H:11
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
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 > &)
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
static options & New(const fvMesh &mesh)
Construct fvOptions and register to datbase if not present.
Definition: fvOptions.C:101
const objectRegistry & db() const
Return the local objectRegistry.
Definition: IOobject.C:322
virtual void correctNut()
Update the SGS eddy viscosity.
Namespace for OpenFOAM.