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-2019 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  sigmak 1.0;
55  sigmaEps 1.3;
56  }
57  \endverbatim
58 
59 SourceFiles
60  mixtureKEpsilon.C
61 
62 \*---------------------------------------------------------------------------*/
63 
64 #ifndef mixtureKEpsilon_H
65 #define mixtureKEpsilon_H
66 
67 #include "RASModel.H"
68 #include "eddyViscosity.H"
69 
70 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
71 
72 namespace Foam
73 {
74 namespace RASModels
75 {
76 
77 /*---------------------------------------------------------------------------*\
78  Class mixtureKEpsilon Declaration
79 \*---------------------------------------------------------------------------*/
80 
81 template<class BasicTurbulenceModel>
82 class mixtureKEpsilon
83 :
84  public eddyViscosity<RASModel<BasicTurbulenceModel>>
85 {
86  // Private Data
87 
88  mutable mixtureKEpsilon<BasicTurbulenceModel> *liquidTurbulencePtr_;
89 
90 
91  // Private Member Functions
92 
93  //- Return the turbulence model for the other phase
94  mixtureKEpsilon<BasicTurbulenceModel>& liquidTurbulence() const;
95 
96 
97 protected:
98 
99  // Protected data
100 
101  // Model coefficients
110 
111  // Fields
115 
116  // Mixture fields
122 
123 
124  // Protected Member Functions
125 
127 
128  void correctInletOutlet
129  (
130  volScalarField& vsf,
131  const volScalarField& refVsf
132  ) const;
133 
134  void initMixtureFields();
135 
136  virtual void correctNut();
137 
138  tmp<volScalarField> Ct2() const;
139 
142  tmp<volScalarField> rhom() const;
143 
145  (
146  const volScalarField& fc,
147  const volScalarField& fd
148  ) const;
149 
151  (
152  const volScalarField& fc,
153  const volScalarField& fd
154  ) const;
155 
157  (
158  const surfaceScalarField& fc,
159  const surfaceScalarField& fd
160  ) const;
161 
163  virtual tmp<fvScalarMatrix> kSource() const;
164  virtual tmp<fvScalarMatrix> epsilonSource() const;
165 
166  //- Return the effective diffusivity for k
167  tmp<volScalarField> DkEff(const volScalarField& nutm) const
168  {
169  return volScalarField::New
170  (
171  "DkEff",
172  nutm/sigmak_
173  );
174  }
175 
176  //- Return the effective diffusivity for epsilon
178  {
179  return volScalarField::New
180  (
181  "DepsilonEff",
182  nutm/sigmaEps_
183  );
184  }
185 
186 public:
188  typedef typename BasicTurbulenceModel::alphaField alphaField;
189  typedef typename BasicTurbulenceModel::rhoField rhoField;
190  typedef typename BasicTurbulenceModel::transportModel transportModel;
191 
192 
193  //- Runtime type information
194  TypeName("mixtureKEpsilon");
195 
196 
197  // Constructors
198 
199  //- Construct from components
201  (
202  const alphaField& alpha,
203  const rhoField& rho,
204  const volVectorField& U,
205  const surfaceScalarField& alphaRhoPhi,
206  const surfaceScalarField& phi,
207  const transportModel& transport,
208  const word& propertiesName = turbulenceModel::propertiesName,
209  const word& type = typeName
210  );
211 
212  //- Disallow default bitwise copy construction
213  mixtureKEpsilon(const mixtureKEpsilon&) = delete;
214 
215 
216  //- Destructor
217  virtual ~mixtureKEpsilon()
218  {}
219 
220 
221  // Member Functions
222 
223  //- Re-read model coefficients if they have changed
224  virtual bool read();
225 
226  //- Return the turbulence kinetic energy
227  virtual tmp<volScalarField> k() const
228  {
229  return k_;
230  }
231 
232  //- Return the turbulence kinetic energy dissipation rate
233  virtual tmp<volScalarField> epsilon() const
234  {
235  return epsilon_;
236  }
237 
238  //- Solve the turbulence equations and correct the turbulence viscosity
239  virtual void correct();
240 
241 
242  // Member Operators
243 
244  //- Disallow default bitwise assignment
245  void operator=(const mixtureKEpsilon&) = delete;
246 };
247 
248 
249 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
250 
251 } // End namespace RASModels
252 } // End namespace Foam
253 
254 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
255 
256 #ifdef NoRepository
257  #include "mixtureKEpsilon.C"
258 #endif
259 
260 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
261 
262 #endif
263 
264 // ************************************************************************* //
virtual tmp< volScalarField > epsilon() const
Return the turbulence kinetic energy dissipation rate.
BasicTurbulenceModel::alphaField alphaField
tmp< surfaceScalarField > mixFlux(const surfaceScalarField &fc, const surfaceScalarField &fd) const
tmp< volScalarField > mixU(const volScalarField &fc, const volScalarField &fd) const
BasicTurbulenceModel::rhoField rhoField
surfaceScalarField & phi
autoPtr< volScalarField > epsilonm_
TypeName("mixtureKEpsilon")
Runtime type information.
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
tmp< volScalarField > DepsilonEff(const volScalarField &nutm) const
Return the effective diffusivity for epsilon.
BasicTurbulenceModel::transportModel transportModel
autoPtr< volScalarField > Ct2_
tmp< volScalarField > mix(const volScalarField &fc, const volScalarField &fd) const
autoPtr< volScalarField > rhom_
virtual tmp< volScalarField > k() const
Return the turbulence kinetic energy.
wordList epsilonBoundaryTypes(const volScalarField &epsilon) const
autoPtr< volScalarField > km_
Eddy viscosity turbulence model base class.
Definition: eddyViscosity.H:49
tmp< volScalarField > rhom() const
mixtureKEpsilon(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName=turbulenceModel::propertiesName, const word &type=typeName)
Construct from components.
tmp< volScalarField > bubbleG() const
tmp< volScalarField > rhogEff() const
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Mesh &, const dimensionSet &, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
Return a temporary field constructed from name, mesh, dimensionSet.
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from string.
Definition: word.H:59
virtual bool read()
Re-read model coefficients if they have changed.
virtual tmp< fvScalarMatrix > epsilonSource() const
virtual tmp< fvScalarMatrix > kSource() const
void operator=(const mixtureKEpsilon &)=delete
Disallow default bitwise assignment.
tmp< volScalarField > Ct2() const
U
Definition: pEqn.H:72
Mixture k-epsilon turbulence model for two-phase gas-liquid systems.
tmp< volScalarField > DkEff(const volScalarField &nutm) const
Return the effective diffusivity for k.
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
A class for managing temporary objects.
Definition: PtrList.H:53
void correctInletOutlet(volScalarField &vsf, const volScalarField &refVsf) const
tmp< volScalarField > rholEff() const
Namespace for OpenFOAM.
virtual ~mixtureKEpsilon()
Destructor.