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