mixtureKEpsilon.H
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-2023 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 Class
25  Foam::RASModels::mixtureKEpsilon
26 
27 Description
28  Mixture k-epsilon turbulence model for two-phase gas-liquid systems
29 
30  The basic structure of the model is based on:
31  \verbatim
32  Behzadi, A., Issa, R. I., & Rusche, H. (2004).
33  Modelling of dispersed bubble and droplet flow at high phase fractions.
34  Chemical Engineering Science, 59(4), 759-770.
35  \endverbatim
36 
37  but an effective density for the gas is used in the averaging and an
38  alternative model for bubble-generated turbulence from:
39  \verbatim
40  Lahey Jr, R. T. (2005).
41  The simulation of multidimensional multiphase flows.
42  Nuclear Engineering and Design, 235(10), 1043-1060.
43  \endverbatim
44 
45  The default model coefficients are
46  \verbatim
47  mixtureKEpsilonCoeffs
48  {
49  Cmu 0.09;
50  C1 1.44;
51  C2 1.92;
52  C3 C2;
53  Cp 0.25; // Bubble generated turbulence
54  alphap 0.3; // Gas phase fraction below which
55  // bubble generated turbulence is included
56  sigmak 1.0;
57  sigmaEps 1.3;
58  }
59  \endverbatim
60 
61 SourceFiles
62  mixtureKEpsilon.C
63 
64 \*---------------------------------------------------------------------------*/
65 
66 #ifndef mixtureKEpsilon_H
67 #define mixtureKEpsilon_H
68 
69 #include "RASModel.H"
70 #include "eddyViscosity.H"
71 
72 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
73 
74 namespace Foam
75 {
76 namespace RASModels
77 {
78 
79 /*---------------------------------------------------------------------------*\
80  Class mixtureKEpsilon Declaration
81 \*---------------------------------------------------------------------------*/
82 
83 template<class BasicMomentumTransportModel>
84 class mixtureKEpsilon
85 :
86  public eddyViscosity<RASModel<BasicMomentumTransportModel>>
87 {
88  // Private Data
89 
91  *gasTurbulencePtr_;
92 
93 
94  // Private Member Functions
95 
96  //- Return the turbulence model for the other phase
97  mixtureKEpsilon<BasicMomentumTransportModel>& gasTurbulence() const;
98 
99 
100 protected:
101 
102  // Protected data
103 
104  // Model coefficients
105 
114 
115  // Fields
116 
119 
120  // Mixture fields
121 
126 
127 
128  // Protected Member Functions
129 
130  void initMixtureFields();
131 
132  virtual void correctNut();
133 
134  tmp<volScalarField> Ct2() const;
135 
138  tmp<volScalarField> rhom() const;
139 
141  (
142  const volScalarField& fc,
143  const volScalarField& fd
144  ) const;
145 
147  (
148  const volScalarField& fc,
149  const volScalarField& fd
150  ) const;
151 
153  (
154  const surfaceScalarField& fc,
155  const surfaceScalarField& fd
156  ) const;
157 
159  virtual tmp<fvScalarMatrix> kSource() const;
160  virtual tmp<fvScalarMatrix> epsilonSource() const;
161 
162  //- Return the effective diffusivity for k
163  tmp<volScalarField> DkEff(const volScalarField& nutm) const
164  {
165  return volScalarField::New
166  (
167  "DkEff",
168  nutm/sigmak_
169  );
170  }
171 
172  //- Return the effective diffusivity for epsilon
174  {
175  return volScalarField::New
176  (
177  "DepsilonEff",
178  nutm/sigmaEps_
179  );
180  }
181 
182 public:
183 
184  typedef typename BasicMomentumTransportModel::alphaField alphaField;
185  typedef typename BasicMomentumTransportModel::rhoField rhoField;
186 
187 
188  //- Runtime type information
189  TypeName("mixtureKEpsilon");
190 
191 
192  // Constructors
193 
194  //- Construct from components
196  (
197  const alphaField& alpha,
198  const rhoField& rho,
199  const volVectorField& U,
200  const surfaceScalarField& alphaRhoPhi,
201  const surfaceScalarField& phi,
202  const viscosity& viscosity,
203  const word& type = typeName
204  );
205 
206  //- Disallow default bitwise copy construction
207  mixtureKEpsilon(const mixtureKEpsilon&) = delete;
208 
209 
210  //- Destructor
211  virtual ~mixtureKEpsilon()
212  {}
213 
214 
215  // Member Functions
216 
217  //- Re-read model coefficients if they have changed
218  virtual bool read();
219 
220  //- Return the turbulence kinetic energy
221  virtual tmp<volScalarField> k() const
222  {
223  return k_;
224  }
225 
226  //- Return the turbulence kinetic energy dissipation rate
227  virtual tmp<volScalarField> epsilon() const
228  {
229  return epsilon_;
230  }
231 
232  //- Return the turbulence specific dissipation rate
233  virtual tmp<volScalarField> omega() const
234  {
235  return volScalarField::New
236  (
237  "omega",
238  epsilon_/(Cmu_*k_)
239  );
240  }
241 
242  //- Solve the turbulence equations and correct the turbulence viscosity
243  virtual void correct();
244 
245 
246  // Member Operators
247 
248  //- Disallow default bitwise assignment
249  void operator=(const mixtureKEpsilon&) = delete;
250 };
251 
252 
253 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
254 
255 } // End namespace RASModels
256 } // End namespace Foam
257 
258 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
259 
260 #ifdef NoRepository
261  #include "mixtureKEpsilon.C"
262 #endif
263 
264 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
265 
266 #endif
267 
268 // ************************************************************************* //
Generic GeometricField class.
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &)
Return a temporary field constructed from name,.
Mixture k-epsilon turbulence model for two-phase gas-liquid systems.
BasicMomentumTransportModel::alphaField alphaField
virtual tmp< fvScalarMatrix > epsilonSource() const
tmp< volScalarField > Ct2() const
tmp< volScalarField > rhogEff() const
virtual tmp< volScalarField > k() const
Return the turbulence kinetic energy.
tmp< volScalarField > mixU(const volScalarField &fc, const volScalarField &fd) const
virtual tmp< volScalarField > omega() const
Return the turbulence specific dissipation rate.
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
void operator=(const mixtureKEpsilon &)=delete
Disallow default bitwise assignment.
tmp< volScalarField > rhom() const
tmp< volScalarField > bubbleG() const
autoPtr< volScalarField > rhom_
autoPtr< volScalarField > km_
autoPtr< volScalarField > Ct2_
tmp< volScalarField > DepsilonEff(const volScalarField &nutm) const
Return the effective diffusivity for epsilon.
virtual tmp< fvScalarMatrix > kSource() const
tmp< surfaceScalarField > mixFlux(const surfaceScalarField &fc, const surfaceScalarField &fd) const
tmp< volScalarField > DkEff(const volScalarField &nutm) const
Return the effective diffusivity for k.
tmp< volScalarField > rholEff() const
autoPtr< volScalarField > epsilonm_
tmp< volScalarField > mix(const volScalarField &fc, const volScalarField &fd) const
virtual tmp< volScalarField > epsilon() const
Return the turbulence kinetic energy dissipation rate.
TypeName("mixtureKEpsilon")
Runtime type information.
virtual ~mixtureKEpsilon()
Destructor.
virtual bool read()
Re-read model coefficients if they have changed.
mixtureKEpsilon(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.
BasicMomentumTransportModel::rhoField rhoField
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
Eddy viscosity turbulence model base class.
Definition: eddyViscosity.H:52
BasicMomentumTransportModel::alphaField alphaField
BasicMomentumTransportModel::rhoField rhoField
A class for managing temporary objects.
Definition: tmp.H:55
Abstract base class for all fluid physical properties.
Definition: viscosity.H:50
A class for handling words, derived from string.
Definition: word.H:62
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.
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488