LaheyKEpsilon.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::LaheyKEpsilon
26 
27 Group
28  grpRASTurbulence
29 
30 Description
31  Continuous-phase k-epsilon model including bubble-generated turbulence.
32 
33  Reference:
34  \verbatim
35  Lahey Jr, R. T. (2005).
36  The simulation of multidimensional multiphase flows.
37  Nuclear Engineering and Design, 235(10), 1043-1060.
38  \endverbatim
39 
40  The default model coefficients are
41  \verbatim
42  LaheyKEpsilonCoeffs
43  {
44  Cmu 0.09;
45  C1 1.44;
46  C2 1.92;
47  C3 -0.33;
48  sigmak 1.0;
49  sigmaEps 1.3;
50  Cp 0.25;
51  Cmub 0.6;
52  alphaInversion 0.3;
53  }
54  \endverbatim
55 
56 SourceFiles
57  LaheyKEpsilon.C
58 
59 \*---------------------------------------------------------------------------*/
60 
61 #ifndef LaheyKEpsilon_H
62 #define LaheyKEpsilon_H
63 
64 #include "kEpsilon.H"
65 
66 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
67 
68 namespace Foam
69 {
70 namespace RASModels
71 {
72 
73 /*---------------------------------------------------------------------------*\
74  Class LaheyKEpsilon Declaration
75 \*---------------------------------------------------------------------------*/
76 
77 template<class BasicTurbulenceModel>
78 class LaheyKEpsilon
79 :
80  public kEpsilon<BasicTurbulenceModel>
81 {
82  // Private data
83 
85  <
86  typename BasicTurbulenceModel::transportModel
87  > *gasTurbulencePtr_;
88 
89 
90  // Private Member Functions
91 
92  //- Return the turbulence model for the gas phase
94  <
95  typename BasicTurbulenceModel::transportModel
96  >&
97  gasTurbulence() const;
98 
99  // Disallow default bitwise copy construct and assignment
101  void operator=(const LaheyKEpsilon&);
102 
103 
104 protected:
105 
106  // Protected data
107 
108  // Model coefficients
114 
115 
116  // Protected Member Functions
117 
118  virtual void correctNut();
121  virtual tmp<fvScalarMatrix> kSource() const;
122  virtual tmp<fvScalarMatrix> epsilonSource() const;
123 
124 
125 public:
127  typedef typename BasicTurbulenceModel::alphaField alphaField;
128  typedef typename BasicTurbulenceModel::rhoField rhoField;
129  typedef typename BasicTurbulenceModel::transportModel transportModel;
130 
131 
132  //- Runtime type information
133  TypeName("LaheyKEpsilon");
134 
135 
136  // Constructors
137 
138  //- Construct from components
140  (
141  const alphaField& alpha,
142  const rhoField& rho,
143  const volVectorField& U,
144  const surfaceScalarField& alphaRhoPhi,
145  const surfaceScalarField& phi,
146  const transportModel& transport,
147  const word& propertiesName = turbulenceModel::propertiesName,
148  const word& type = typeName
149  );
150 
151 
152  //- Destructor
153  virtual ~LaheyKEpsilon()
154  {}
155 
156 
157  // Member Functions
158 
159  //- Read model coefficients if they have changed
160  virtual bool read();
161 
162  //- Solve the turbulence equations and correct the turbulence viscosity
163  virtual void correct();
164 };
165 
166 
167 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
168 
169 } // End namespace RASModels
170 } // End namespace Foam
171 
172 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
173 
174 #ifdef NoRepository
175  #include "LaheyKEpsilon.C"
176 #endif
177 
178 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
179 
180 #endif
181 
182 // ************************************************************************* //
TypeName("LaheyKEpsilon")
Runtime type information.
surfaceScalarField & phi
U
Definition: pEqn.H:83
BasicTurbulenceModel::rhoField rhoField
dimensionedScalar alphaInversion_
virtual bool read()
Read model coefficients if they have changed.
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Continuous-phase k-epsilon model including bubble-generated turbulence.
Definition: LaheyKEpsilon.H:77
BasicTurbulenceModel::alphaField alphaField
Standard k-epsilon turbulence model for incompressible and compressible flows including rapid distort...
Definition: kEpsilon.H:86
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from string.
Definition: word.H:59
BasicTurbulenceModel::transportModel transportModel
Templated abstract base class for multiphase compressible turbulence models.
tmp< volScalarField > phaseTransferCoeff() const
fileName::Type type(const fileName &, const bool followLink=true)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:485
virtual ~LaheyKEpsilon()
Destructor.
virtual tmp< fvScalarMatrix > kSource() const
tmp< volScalarField > bubbleG() const
virtual tmp< fvScalarMatrix > epsilonSource() const
A class for managing temporary objects.
Definition: PtrList.H:53
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Namespace for OpenFOAM.