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-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 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  sigmak 1.0;
54  sigmaEps 1.3;
55  }
56  \endverbatim
57 
58 SourceFiles
59  mixtureKEpsilon.C
60 
61 \*---------------------------------------------------------------------------*/
62 
63 #ifndef mixtureKEpsilon_H
64 #define mixtureKEpsilon_H
65 
66 #include "RASModel.H"
67 #include "eddyViscosity.H"
68 
69 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
70 
71 namespace Foam
72 {
73 namespace RASModels
74 {
75 
76 /*---------------------------------------------------------------------------*\
77  Class mixtureKEpsilon Declaration
78 \*---------------------------------------------------------------------------*/
79 
80 template<class BasicTurbulenceModel>
81 class mixtureKEpsilon
82 :
83  public eddyViscosity<RASModel<BasicTurbulenceModel>>
84 {
85  // Private data
86 
87  mutable mixtureKEpsilon<BasicTurbulenceModel> *liquidTurbulencePtr_;
88 
89 
90  // Private Member Functions
91 
92  // Disallow default bitwise copy construct and assignment
94  void operator=(const mixtureKEpsilon&);
95 
96  //- Return the turbulence model for the other phase
97  mixtureKEpsilon<BasicTurbulenceModel>& liquidTurbulence() const;
98 
99 
100 protected:
101 
102  // Protected data
103 
104  // Model coefficients
113 
114  // Fields
118 
119  // Mixture fields
125 
126 
127  // Protected Member Functions
128 
130 
131  void correctInletOutlet
132  (
133  volScalarField& vsf,
134  const volScalarField& refVsf
135  ) const;
136 
137  void initMixtureFields();
138 
139  virtual void correctNut();
140 
141  tmp<volScalarField> Ct2() const;
142 
145  tmp<volScalarField> rhom() const;
146 
148  (
149  const volScalarField& fc,
150  const volScalarField& fd
151  ) const;
152 
154  (
155  const volScalarField& fc,
156  const volScalarField& fd
157  ) const;
158 
160  (
161  const surfaceScalarField& fc,
162  const surfaceScalarField& fd
163  ) const;
164 
166  virtual tmp<fvScalarMatrix> kSource() const;
167  virtual tmp<fvScalarMatrix> epsilonSource() const;
168 
169  //- Return the effective diffusivity for k
170  tmp<volScalarField> DkEff(const volScalarField& nutm) const
171  {
172  return tmp<volScalarField>
173  (
174  new volScalarField
175  (
176  "DkEff",
177  nutm/sigmak_
178  )
179  );
180  }
181 
182  //- Return the effective diffusivity for epsilon
184  {
185  return tmp<volScalarField>
186  (
187  new volScalarField
188  (
189  "DepsilonEff",
190  nutm/sigmaEps_
191  )
192  );
193  }
194 
195 public:
197  typedef typename BasicTurbulenceModel::alphaField alphaField;
198  typedef typename BasicTurbulenceModel::rhoField rhoField;
199  typedef typename BasicTurbulenceModel::transportModel transportModel;
200 
201 
202  //- Runtime type information
203  TypeName("mixtureKEpsilon");
204 
205 
206  // Constructors
207 
208  //- Construct from components
210  (
211  const alphaField& alpha,
212  const rhoField& rho,
213  const volVectorField& U,
214  const surfaceScalarField& alphaRhoPhi,
215  const surfaceScalarField& phi,
216  const transportModel& transport,
217  const word& propertiesName = turbulenceModel::propertiesName,
218  const word& type = typeName
219  );
220 
221 
222  //- Destructor
223  virtual ~mixtureKEpsilon()
224  {}
225 
226 
227  // Member Functions
228 
229  //- Re-read model coefficients if they have changed
230  virtual bool read();
231 
232  //- Return the turbulence kinetic energy
233  virtual tmp<volScalarField> k() const
234  {
235  return k_;
236  }
237 
238  //- Return the turbulence kinetic energy dissipation rate
239  virtual tmp<volScalarField> epsilon() const
240  {
241  return epsilon_;
242  }
243 
244  //- Solve the turbulence equations and correct the turbulence viscosity
245  virtual void correct();
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
tmp< volScalarField > bubbleG() const
tmp< volScalarField > rhogEff() const
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
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
fileName::Type type(const fileName &, const bool followLink=true)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:481
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.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
A class for managing temporary objects.
Definition: PtrList.H:53
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
void correctInletOutlet(volScalarField &vsf, const volScalarField &refVsf) const
tmp< volScalarField > rholEff() const
Namespace for OpenFOAM.
virtual ~mixtureKEpsilon()
Destructor.