mixtureKEpsilon.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2013-2016 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 Group
28  grpRASTurbulence
29 
30 Description
31  Mixture k-epsilon turbulence model for two-phase gas-liquid systems
32 
33  The basic structure of the model is based on:
34  \verbatim
35  Behzadi, A., Issa, R. I., & Rusche, H. (2004).
36  Modelling of dispersed bubble and droplet flow at high phase fractions.
37  Chemical Engineering Science, 59(4), 759-770.
38  \endverbatim
39 
40  but an effective density for the gas is used in the averaging and an
41  alternative model for bubble-generated turbulence from:
42  \verbatim
43  Lahey Jr, R. T. (2005).
44  The simulation of multidimensional multiphase flows.
45  Nuclear Engineering and Design, 235(10), 1043-1060.
46  \endverbatim
47 
48  The default model coefficients are
49  \verbatim
50  mixtureKEpsilonCoeffs
51  {
52  Cmu 0.09;
53  C1 1.44;
54  C2 1.92;
55  C3 C2;
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 BasicTurbulenceModel>
84 class mixtureKEpsilon
85 :
86  public eddyViscosity<RASModel<BasicTurbulenceModel>>
87 {
88  // Private data
89 
90  mutable mixtureKEpsilon<BasicTurbulenceModel> *liquidTurbulencePtr_;
91 
92 
93  // Private Member Functions
94 
95  // Disallow default bitwise copy construct and assignment
97  void operator=(const mixtureKEpsilon&);
98 
99  //- Return the turbulence model for the other phase
100  mixtureKEpsilon<BasicTurbulenceModel>& liquidTurbulence() const;
101 
102 
103 protected:
104 
105  // Protected data
106 
107  // Model coefficients
116 
117  // Fields
121 
122  // Mixture fields
128 
129 
130  // Protected Member Functions
131 
133 
134  void correctInletOutlet
135  (
136  volScalarField& vsf,
137  const volScalarField& refVsf
138  ) const;
139 
140  void initMixtureFields();
141 
142  virtual void correctNut();
143 
144  tmp<volScalarField> Ct2() const;
145 
148  tmp<volScalarField> rhom() const;
149 
151  (
152  const volScalarField& fc,
153  const volScalarField& fd
154  ) const;
155 
157  (
158  const volScalarField& fc,
159  const volScalarField& fd
160  ) const;
161 
163  (
164  const surfaceScalarField& fc,
165  const surfaceScalarField& fd
166  ) const;
167 
169  virtual tmp<fvScalarMatrix> kSource() const;
170  virtual tmp<fvScalarMatrix> epsilonSource() const;
171 
172  //- Return the effective diffusivity for k
173  tmp<volScalarField> DkEff(const volScalarField& nutm) const
174  {
175  return tmp<volScalarField>
176  (
177  new volScalarField
178  (
179  "DkEff",
180  nutm/sigmak_
181  )
182  );
183  }
184 
185  //- Return the effective diffusivity for epsilon
187  {
188  return tmp<volScalarField>
189  (
190  new volScalarField
191  (
192  "DepsilonEff",
193  nutm/sigmaEps_
194  )
195  );
196  }
197 
198 public:
200  typedef typename BasicTurbulenceModel::alphaField alphaField;
201  typedef typename BasicTurbulenceModel::rhoField rhoField;
202  typedef typename BasicTurbulenceModel::transportModel transportModel;
203 
204 
205  //- Runtime type information
206  TypeName("mixtureKEpsilon");
207 
208 
209  // Constructors
210 
211  //- Construct from components
213  (
214  const alphaField& alpha,
215  const rhoField& rho,
216  const volVectorField& U,
217  const surfaceScalarField& alphaRhoPhi,
218  const surfaceScalarField& phi,
219  const transportModel& transport,
220  const word& propertiesName = turbulenceModel::propertiesName,
221  const word& type = typeName
222  );
223 
224 
225  //- Destructor
226  virtual ~mixtureKEpsilon()
227  {}
228 
229 
230  // Member Functions
231 
232  //- Re-read model coefficients if they have changed
233  virtual bool read();
234 
235  //- Return the turbulence kinetic energy
236  virtual tmp<volScalarField> k() const
237  {
238  return k_;
239  }
240 
241  //- Return the turbulence kinetic energy dissipation rate
242  virtual tmp<volScalarField> epsilon() const
243  {
244  return epsilon_;
245  }
246 
247  //- Solve the turbulence equations and correct the turbulence viscosity
248  virtual void correct();
249 };
250 
251 
252 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
253 
254 } // End namespace RASModels
255 } // End namespace Foam
256 
257 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
258 
259 #ifdef NoRepository
260  #include "mixtureKEpsilon.C"
261 #endif
262 
263 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
264 
265 #endif
266 
267 // ************************************************************************* //
virtual tmp< volScalarField > epsilon() const
Return the turbulence kinetic energy dissipation rate.
surfaceScalarField & phi
BasicTurbulenceModel::alphaField alphaField
tmp< surfaceScalarField > mixFlux(const surfaceScalarField &fc, const surfaceScalarField &fd) const
tmp< volScalarField > mixU(const volScalarField &fc, const volScalarField &fd) const
U
Definition: pEqn.H:83
BasicTurbulenceModel::rhoField rhoField
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:52
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:485
tmp< volScalarField > Ct2() const
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.