kOmegaSSTSato.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) 2014-2017 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 "kOmegaSSTSato.H"
30 #include "fvOptions.H"
31 #include "twoPhaseSystem.H"
32 
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 namespace RASModels
38 {
39 
40 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
41 
42 template<class BasicTurbulenceModel>
43 kOmegaSSTSato<BasicTurbulenceModel>::kOmegaSSTSato
44 (
45  const alphaField& alpha,
46  const rhoField& rho,
47  const volVectorField& U,
48  const surfaceScalarField& alphaRhoPhi,
49  const surfaceScalarField& phi,
50  const transportModel& transport,
51  const word& propertiesName,
52  const word& type
53 )
54 :
56  (
57  alpha,
58  rho,
59  U,
60  alphaRhoPhi,
61  phi,
62  transport,
63  propertiesName,
64  type
65  ),
66 
67  gasTurbulencePtr_(nullptr),
68 
69  Cmub_
70  (
72  (
73  "Cmub",
74  this->coeffDict_,
75  0.6
76  )
77  )
78 {
79  if (type == typeName)
80  {
81  this->printCoeffs(type);
82  }
83 }
84 
85 
86 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
87 
88 template<class BasicTurbulenceModel>
90 {
92  {
93  Cmub_.readIfPresent(this->coeffDict());
94 
95  return true;
96  }
97 
98  return false;
99 }
100 
101 
102 template<class BasicTurbulenceModel>
104 <
105  typename BasicTurbulenceModel::transportModel
106 >&
108 {
109  if (!gasTurbulencePtr_)
110  {
111  const volVectorField& U = this->U_;
112 
113  const transportModel& liquid = this->transport();
114  const twoPhaseSystem& fluid =
115  refCast<const twoPhaseSystem>(liquid.fluid());
116  const transportModel& gas = fluid.otherPhase(liquid);
117 
118  gasTurbulencePtr_ =
119  &U.db()
121  (
123  (
125  gas.name()
126  )
127  );
128  }
129 
130  return *gasTurbulencePtr_;
131 }
132 
133 
134 template<class BasicTurbulenceModel>
136 (
137  const volScalarField& S2
138 )
139 {
141  this->gasTurbulence();
142 
144  (
145  pow(this->betaStar_, 0.25)*this->y_*sqrt(this->k_)/this->nu()
146  );
147 
148  this->nut_ =
149  this->a1_*this->k_
150  /max
151  (
152  this->a1_*this->omega_,
153  this->b1_*this->F23()*sqrt(S2)
154  )
155  + sqr(1 - exp(-yPlus/16.0))
156  *Cmub_*gasTurbulence.transport().d()*gasTurbulence.alpha()
157  *(mag(this->U_ - gasTurbulence.U()));
158 
159  this->nut_.correctBoundaryConditions();
160  fv::options::New(this->mesh_).correct(this->nut_);
161 
162  BasicTurbulenceModel::correctNut();
163 }
164 
165 
166 template<class BasicTurbulenceModel>
168 {
170 }
171 
172 
173 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
174 
175 } // End namespace RASModels
176 } // End namespace Foam
177 
178 // ************************************************************************* //
Foam::RASModels::kOmegaSSTSato::rhoField
BasicTurbulenceModel::rhoField rhoField
Definition: kOmegaSSTSato.H:169
fvOptions.H
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::twoPhaseSystem
Class which solves the volume fraction equations for two phases.
Definition: twoPhaseSystem.H:53
Foam::fv::optionList::correct
void correct(GeometricField< Type, fvPatchField, volMesh > &field)
Apply correction to field.
Definition: fvOptionListTemplates.C:355
Foam::RASModels::kOmegaSST
Implementation of the k-omega-SST turbulence model for incompressible and compressible flows.
Definition: kOmegaSST.H:129
Foam::PhaseCompressibleTurbulenceModel
Templated abstract base class for multiphase compressible turbulence models.
Definition: phaseCompressibleTurbulenceModelFwd.H:44
Foam::fv::options::New
static options & New(const fvMesh &mesh)
Construct fvOptions and register to database if not present.
Definition: fvOptions.C:103
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:212
Foam::RASModels::kOmegaSSTSato::read
virtual bool read()
Read model coefficients if they have changed.
Definition: kOmegaSSTSato.C:89
Foam::turbulenceModel::propertiesName
static const word propertiesName
Default name of the turbulence properties dictionary.
Definition: turbulenceModel.H:100
Foam::liquid
Generic thermophysical properties class for a liquid in which the functions and coefficients for each...
Definition: liquid.H:54
Foam::exp
dimensionedScalar exp(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:261
fluid
twoPhaseSystem & fluid
Definition: setRegionFluidFields.H:3
rho
rho
Definition: readInitialConditions.H:88
Foam::RASModels::kOmegaSSTSato::correct
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: kOmegaSSTSato.C:167
nu
volScalarField & nu
Definition: readMechanicalProperties.H:176
kOmegaSSTSato.H
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:75
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
Foam::RASModels::kOmegaSSTSato::transportModel
BasicTurbulenceModel::transportModel transportModel
Definition: kOmegaSSTSato.H:170
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:42
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::RASModels::kOmegaSSTSato
Implementation of the k-omega-SST turbulence model for dispersed bubbly flows with Sato (1981) bubble...
Definition: kOmegaSSTSato.H:123
Foam::kOmegaSSTBase< eddyViscosity< RASModel< BasicTurbulenceModel > > >::correct
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: kOmegaSSTBase.C:484
Foam::transportModel
Base-class for all transport models used by the incompressible turbulence models.
Definition: transportModel.H:53
U
U
Definition: pEqn.H:72
Foam::RASModels::kOmegaSSTSato::alphaField
BasicTurbulenceModel::alphaField alphaField
Definition: kOmegaSSTSato.H:168
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::type
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:590
yPlus
scalar yPlus
Definition: evaluateNearWall.H:16
Foam::RASModels::kOmegaSST::correctNut
virtual void correctNut()
Definition: kOmegaSST.C:54
Foam::IOobject::groupName
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
Foam::GeometricField< vector, fvPatchField, volMesh >