continuousGasKEpsilon.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::continuousGasKEpsilon
26 
27 Description
28  k-epsilon model for the gas-phase in a two-phase system
29  supporting phase-inversion.
30 
31  In the limit that the gas-phase fraction approaches zero a contribution from
32  the other phase is blended into the k and epsilon equations up to the
33  phase-fraction of alphaInversion at which point phase-inversion is
34  considered to have occurred and the model reverts to the pure single-phase
35  form.
36 
37  This model is unpublished and is provided as a stable numerical framework
38  on which a more physical model may be built.
39 
40  The default model coefficients are
41  \verbatim
42  continuousGasKEpsilonCoeffs
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  alphaInversion 0.7;
51  }
52  \endverbatim
53 
54 SourceFiles
55  continuousGasKEpsilon.C
56 
57 \*---------------------------------------------------------------------------*/
58 
59 #ifndef continuousGasKEpsilon_H
60 #define continuousGasKEpsilon_H
61 
62 #include "kEpsilon.H"
63 
64 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
65 
66 namespace Foam
67 {
68 namespace RASModels
69 {
70 
71 /*---------------------------------------------------------------------------*\
72  Class continuousGasKEpsilon Declaration
73 \*---------------------------------------------------------------------------*/
74 
75 template<class BasicTurbulenceModel>
77 :
78  public kEpsilon<BasicTurbulenceModel>
79 {
80  // Private data
81 
82  mutable const turbulenceModel *liquidTurbulencePtr_;
83 
84  volScalarField nutEff_;
85 
86 
87  // Private Member Functions
88 
89  // Disallow default bitwise copy construct and assignment
91  void operator=(const continuousGasKEpsilon&);
92 
93 
94 protected:
95 
96  // Protected data
97 
98  // Model coefficients
99 
101 
102 
103  // Protected Member Functions
104 
105  virtual void correctNut();
107  virtual tmp<fvScalarMatrix> kSource() const;
108  virtual tmp<fvScalarMatrix> epsilonSource() const;
109 
110 
111 public:
113  typedef typename BasicTurbulenceModel::alphaField alphaField;
114  typedef typename BasicTurbulenceModel::rhoField rhoField;
115  typedef typename BasicTurbulenceModel::transportModel transportModel;
116 
117 
118  //- Runtime type information
119  TypeName("continuousGasKEpsilon");
120 
121 
122  // Constructors
123 
124  //- Construct from components
126  (
127  const alphaField& alpha,
128  const rhoField& rho,
129  const volVectorField& U,
130  const surfaceScalarField& alphaRhoPhi,
131  const surfaceScalarField& phi,
132  const transportModel& transport,
133  const word& propertiesName = turbulenceModel::propertiesName,
134  const word& type = typeName
135  );
136 
137 
138  //- Destructor
139  virtual ~continuousGasKEpsilon()
140  {}
141 
142 
143  // Member Functions
144 
145  //- Re-read model coefficients if they have changed
146  virtual bool read();
147 
148  //- Return the turbulence model for the liquid phase
149  const turbulenceModel& liquidTurbulence() const;
150 
151  //- Return the effective viscosity
152  virtual tmp<volScalarField> nuEff() const;
153 
154  //- Return the effective density for the stress
155  virtual tmp<volScalarField> rhoEff() const;
156 
157  //- Return the Reynolds stress tensor
158  virtual tmp<volSymmTensorField> R() const;
159 };
160 
161 
162 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
163 
164 } // End namespace RASModels
165 } // End namespace Foam
166 
167 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
168 
169 #ifdef NoRepository
170  #include "continuousGasKEpsilon.C"
171 #endif
172 
173 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
174 
175 #endif
176 
177 // ************************************************************************* //
BasicTurbulenceModel::transportModel transportModel
surfaceScalarField & phi
virtual tmp< fvScalarMatrix > epsilonSource() const
virtual bool read()
Re-read model coefficients if they have changed.
Abstract base class for turbulence models (RAS, LES and laminar).
k-epsilon model for the gas-phase in a two-phase system supporting phase-inversion.
Standard k-epsilon turbulence model for incompressible and compressible flows including rapid distort...
Definition: kEpsilon.H:83
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from string.
Definition: word.H:59
TypeName("continuousGasKEpsilon")
Runtime type information.
virtual tmp< volScalarField > nuEff() const
Return the effective viscosity.
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
virtual tmp< volSymmTensorField > R() const
Return the Reynolds stress tensor.
BasicTurbulenceModel::alphaField alphaField
U
Definition: pEqn.H:72
tmp< volScalarField > phaseTransferCoeff() const
BasicTurbulenceModel::rhoField rhoField
virtual tmp< volScalarField > rhoEff() const
Return the effective density for the stress.
const turbulenceModel & liquidTurbulence() const
Return the turbulence model for the liquid phase.
A class for managing temporary objects.
Definition: PtrList.H:53
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Namespace for OpenFOAM.