kOmegaSSTSato.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) 2014-2018 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 "kOmegaSSTSato.H"
27 #include "fvOptions.H"
28 #include "twoPhaseSystem.H"
29 
30 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace RASModels
35 {
36 
37 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
38 
39 template<class BasicTurbulenceModel>
40 kOmegaSSTSato<BasicTurbulenceModel>::kOmegaSSTSato
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& propertiesName,
49  const word& type
50 )
51 :
53  (
54  alpha,
55  rho,
56  U,
57  alphaRhoPhi,
58  phi,
59  transport,
60  propertiesName,
61  type
62  ),
63 
64  gasTurbulencePtr_(nullptr),
65 
66  Cmub_
67  (
69  (
70  "Cmub",
71  this->coeffDict_,
72  0.6
73  )
74  )
75 {
76  if (type == typeName)
77  {
78  this->printCoeffs(type);
79  }
80 }
81 
82 
83 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
84 
85 template<class BasicTurbulenceModel>
87 {
89  {
90  Cmub_.readIfPresent(this->coeffDict());
91 
92  return true;
93  }
94  else
95  {
96  return false;
97  }
98 }
99 
100 template<class BasicTurbulenceModel>
102 <
103  typename BasicTurbulenceModel::transportModel
104 >&
106 {
107  if (!gasTurbulencePtr_)
108  {
109  const volVectorField& U = this->U_;
110 
111  const transportModel& liquid = this->transport();
112  const twoPhaseSystem& fluid =
113  refCast<const twoPhaseSystem>(liquid.fluid());
114  const transportModel& gas = fluid.otherPhase(liquid);
115 
116  gasTurbulencePtr_ =
117  &U.db()
119  (
121  (
123  gas.name()
124  )
125  );
126  }
127 
128  return *gasTurbulencePtr_;
129 }
130 
131 
132 template<class BasicTurbulenceModel>
134 (
135  const volScalarField& S2,
136  const volScalarField& F2
137 )
138 {
140  this->gasTurbulence();
141 
143  (
144  pow(this->betaStar_, 0.25)*this->y_*sqrt(this->k_)/this->nu()
145  );
146 
147  this->nut_ =
148  this->a1_*this->k_/max(this->a1_*this->omega_, this->b1_*F2*sqrt(S2))
149  + sqr(1 - exp(-yPlus/16.0))
150  *Cmub_*gasTurbulence.transport().d()*gasTurbulence.alpha()
151  *(mag(this->U_ - gasTurbulence.U()));
152 
153  this->nut_.correctBoundaryConditions();
154  fv::options::New(this->mesh_).correct(this->nut_);
155 
156  BasicTurbulenceModel::correctNut();
157 }
158 
159 
160 template<class BasicTurbulenceModel>
162 {
164 }
165 
166 
167 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
168 
169 } // End namespace RASModels
170 } // End namespace Foam
171 
172 // ************************************************************************* //
Implementation of the k-omega-SST turbulence model for dispersed bubbly flows with Sato (1981) bubble...
surfaceScalarField & phi
multiphaseSystem & fluid
Definition: createFields.H:11
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
const volVectorField & U() const
Access function to velocity field.
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
void correct(GeometricField< Type, fvPatchField, volMesh > &field)
Apply correction to field.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar sqrt(const dimensionedScalar &ds)
BasicTurbulenceModel::transportModel transportModel
#define F2(B, C, D)
Definition: SHA1.C:174
BasicTurbulenceModel::alphaField alphaField
Generic dimensioned Type class.
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
Class which solves the volume fraction equations for two phases.
const transportModel & transport() const
Access function to incompressible transport model.
dimensionedScalar exp(const dimensionedScalar &ds)
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from string.
Definition: word.H:59
static word groupName(Name name, const word &group)
BasicTurbulenceModel::rhoField rhoField
Specialisation for RAS of the generic kOmegaSSTBase base class. For more information, see Description of kOmegaSSTBase.H.
Definition: kOmegaSST.H:57
Templated abstract base class for multiphase compressible turbulence models.
const phaseModel & otherPhase(const phaseModel &phase) const
Return the phase not given as an argument.
fileName::Type type(const fileName &, const bool followLink=true)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:481
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
virtual bool read()
Read model coefficients if they have changed.
Definition: kOmegaSSTSato.C:86
U
Definition: pEqn.H:72
scalar yPlus
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
dimensioned< scalar > mag(const dimensioned< Type > &)
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:361
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
volScalarField & nu
Namespace for OpenFOAM.