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-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 "kOmegaSSTSato.H"
27 #include "fvModels.H"
28 #include "fvConstraints.H"
29 #include "phaseSystem.H"
30 
31 
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace RASModels
37 {
38 
39 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
40 
41 template<class BasicMomentumTransportModel>
43 (
44  const alphaField& alpha,
45  const rhoField& rho,
46  const volVectorField& U,
47  const surfaceScalarField& alphaRhoPhi,
48  const surfaceScalarField& phi,
49  const viscosity& viscosity,
50  const word& type
51 )
52 :
53  kOmegaSST<BasicMomentumTransportModel>
54  (
55  alpha,
56  rho,
57  U,
58  alphaRhoPhi,
59  phi,
60  viscosity,
61  type
62  ),
63 
64  phase_(refCast<const phaseModel>(viscosity)),
65 
66  hasDispersedPhaseNames_(this->coeffDict().found("dispersedPhases")),
67 
68  dispersedPhaseNames_
69  (
70  this->coeffDict().lookupOrDefault("dispersedPhases", hashedWordList())
71  ),
72 
73  Cmub_("Cmub", this->coeffDict(), 0.6)
74 {}
75 
76 
77 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
78 
79 template<class BasicMomentumTransportModel>
81 {
83  {
84  Cmub_.readIfPresent(this->coeffDict());
85 
86  return true;
87  }
88  else
89  {
90  return false;
91  }
92 }
93 
94 
95 template<class BasicMomentumTransportModel>
98 {
99  UPtrList<const phaseModel> dispersedPhases;
100 
101  const phaseSystem& fluid = phase_.fluid();
102 
103  if (hasDispersedPhaseNames_)
104  {
105  dispersedPhases.resize(dispersedPhaseNames_.size());
106 
107  forAll(dispersedPhaseNames_, dispersedPhasei)
108  {
109  dispersedPhases.set
110  (
111  dispersedPhasei,
112  &fluid.phases()[dispersedPhaseNames_[dispersedPhasei]]
113  );
114  }
115  }
116  else
117  {
118  dispersedPhases.resize(fluid.movingPhases().size() - 1);
119 
120  label dispersedPhasei = 0;
121 
122  forAll(fluid.movingPhases(), movingPhasei)
123  {
124  const phaseModel& otherPhase = fluid.movingPhases()[movingPhasei];
125 
126  if (&otherPhase != &phase_)
127  {
128  dispersedPhases.set
129  (
130  dispersedPhasei ++,
131  &otherPhase
132  );
133  }
134  }
135  }
136 
137  return dispersedPhases;
138 }
139 
140 
141 template<class BasicMomentumTransportModel>
143 (
144  const volScalarField& S2,
145  const volScalarField& F2
146 )
147 {
149  (
150  pow(this->betaStar_, 0.25)*this->y()*sqrt(this->k_)/this->nu()
151  );
152 
153  this->nut_ =
154  this->a1_*this->k_/max(this->a1_*this->omega_, this->b1_*F2*sqrt(S2));
155 
156  UPtrList<const phaseModel> dispersedPhases(getDispersedPhases());
157  forAllIter(UPtrList<const phaseModel>, dispersedPhases, phaseIter)
158  {
159  this->nut_ +=
160  sqr(1 - exp(-yPlus/16.0))
161  *Cmub_*phaseIter().d()*phaseIter()
162  *(mag(this->U_ - phaseIter().U()));
163  }
164 
165  this->nut_.correctBoundaryConditions();
166  fvConstraints::New(this->mesh_).constrain(this->nut_);
167 }
168 
169 
170 template<class BasicMomentumTransportModel>
172 {
174 }
175 
176 
177 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
178 
179 } // End namespace RASModels
180 } // End namespace Foam
181 
182 // ************************************************************************* //
scalar y
#define F2(B, C, D)
Definition: SHA1.C:174
bool found
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:458
Generic GeometricField class.
Implementation of the k-omega-SST turbulence model for dispersed bubbly flows with Sato (1981) bubble...
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
kOmegaSSTSato(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.
Definition: kOmegaSSTSato.C:43
virtual bool read()
Read model coefficients if they have changed.
Definition: kOmegaSSTSato.C:80
Specialisation for RAS of the generic kOmegaSSTBase base class. For more information,...
Definition: kOmegaSST.H:64
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: UPtrList.H:66
bool set(const label) const
Is element set.
Definition: UPtrListI.H:87
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
void resize(const label)
Reset size of UPtrList. This can only be used to set the size.
Definition: UPtrListI.H:71
static autoPtr< dictionary > New(Istream &)
Construct top-level dictionary on freestore from Istream.
Definition: dictionaryIO.C:103
A wordList with hashed indices for faster lookup by name.
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
BasicMomentumTransportModel::alphaField alphaField
BasicMomentumTransportModel::rhoField rhoField
Class to represent a system of phases.
Definition: phaseSystem.H:74
const phaseModelPartialList & movingPhases() const
Return the models for phases that are moving.
Definition: phaseSystemI.H:118
const phaseModelList & phases() const
Return the phase models.
Definition: phaseSystemI.H:104
Abstract base class for all fluid physical properties.
Definition: viscosity.H:50
A class for handling words, derived from string.
Definition: word.H:62
const scalar yPlus
U
Definition: pEqn.H:72
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
Namespace for OpenFOAM.
dimensionedScalar exp(const dimensionedScalar &ds)
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
To & refCast(From &r)
Reference type cast template function.
Definition: typeInfo.H:134
void mag(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
void pow(LagrangianPatchField< typename powProduct< Type, r >::type > &f, const LagrangianPatchField< Type > &f1)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
void sqr(LagrangianPatchField< typename outerProduct< Type, Type >::type > &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