LaheyKEpsilon.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-2022 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 Description
28  Continuous-phase k-epsilon model including bubble-generated turbulence.
29 
30  Reference:
31  \verbatim
32  Lahey Jr, R. T. (2005).
33  The simulation of multidimensional multiphase flows.
34  Nuclear Engineering and Design, 235(10), 1043-1060.
35  \endverbatim
36 
37  The default model coefficients are
38  \verbatim
39  LaheyKEpsilonCoeffs
40  {
41  Cmu 0.09;
42  C1 1.44;
43  C2 1.92;
44  C3 0;
45  C4 1.92;
46  sigmak 1.0;
47  sigmaEps 1.3;
48  Cp 0.25;
49  Cmub 0.6;
50  alphaInversion 0.3;
51  }
52  \endverbatim
53 
54 SourceFiles
55  LaheyKEpsilon.C
56 
57 \*---------------------------------------------------------------------------*/
58 
59 #ifndef LaheyKEpsilon_H
60 #define LaheyKEpsilon_H
61 
62 #include "kEpsilon.H"
63 
64 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
65 
66 namespace Foam
67 {
68 namespace RASModels
69 {
70 
71 /*---------------------------------------------------------------------------*\
72  Class LaheyKEpsilon Declaration
73 \*---------------------------------------------------------------------------*/
74 
75 template<class BasicMomentumTransportModel>
76 class LaheyKEpsilon
77 :
78  public kEpsilon<BasicMomentumTransportModel>
79 {
80  // Private Data
81 
83  *gasTurbulencePtr_;
84 
85 
86  // Private Member Functions
87 
88  //- Return the turbulence model for the gas phase
89  const phaseCompressibleMomentumTransportModel& gasTurbulence() const;
90 
91 
92 protected:
93 
94  // Protected data
95 
96  // Model coefficients
97 
102 
103 
104  // Protected Member Functions
105 
106  virtual void correctNut();
109  virtual tmp<fvScalarMatrix> kSource() const;
110  virtual tmp<fvScalarMatrix> epsilonSource() const;
111 
112 
113 public:
114 
115  typedef typename BasicMomentumTransportModel::alphaField alphaField;
116  typedef typename BasicMomentumTransportModel::rhoField rhoField;
117 
118 
119  //- Runtime type information
120  TypeName("LaheyKEpsilon");
121 
122 
123  // Constructors
124 
125  //- Construct from components
127  (
128  const alphaField& alpha,
129  const rhoField& rho,
130  const volVectorField& U,
131  const surfaceScalarField& alphaRhoPhi,
132  const surfaceScalarField& phi,
133  const viscosity& viscosity,
134  const word& type = typeName
135  );
136 
137  //- Disallow default bitwise copy construction
138  LaheyKEpsilon(const LaheyKEpsilon&) = delete;
139 
140 
141  //- Destructor
142  virtual ~LaheyKEpsilon()
143  {}
144 
145 
146  // Member Functions
147 
148  //- Read model coefficients if they have changed
149  virtual bool read();
150 
151  //- Solve the turbulence equations and correct the turbulence viscosity
152  virtual void correct();
153 
154 
155  // Member Operators
156 
157  //- Disallow default bitwise assignment
158  void operator=(const LaheyKEpsilon&) = delete;
159 };
160 
161 
162 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
163 
164 } // End namespace RASModels
165 } // End namespace Foam
166 
167 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
168 
169 #ifdef NoRepository
170  #include "LaheyKEpsilon.C"
171 #endif
172 
173 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
174 
175 #endif
176 
177 // ************************************************************************* //
Generic GeometricField class.
Continuous-phase k-epsilon model including bubble-generated turbulence.
Definition: LaheyKEpsilon.H:78
BasicMomentumTransportModel::alphaField alphaField
virtual tmp< fvScalarMatrix > epsilonSource() const
void operator=(const LaheyKEpsilon &)=delete
Disallow default bitwise assignment.
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
LaheyKEpsilon(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.
Definition: LaheyKEpsilon.C:43
dimensionedScalar alphaInversion_
Definition: LaheyKEpsilon.H:97
tmp< volScalarField > bubbleG() const
tmp< volScalarField > phaseTransferCoeff() const
virtual tmp< fvScalarMatrix > kSource() const
virtual ~LaheyKEpsilon()
Destructor.
TypeName("LaheyKEpsilon")
Runtime type information.
virtual bool read()
Read model coefficients if they have changed.
BasicMomentumTransportModel::rhoField rhoField
Standard k-epsilon turbulence model for incompressible and compressible flows including rapid distort...
Definition: kEpsilon.H:86
BasicMomentumTransportModel::alphaField alphaField
BasicMomentumTransportModel::rhoField rhoField
Templated abstract base class for multiphase compressible turbulence models.
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