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-2020 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 BasicMomentumTransportModel>
77 :
78  public kEpsilon<BasicMomentumTransportModel>
79 {
80  // Private Data
81 
82  mutable const momentumTransportModel *liquidTurbulencePtr_;
83 
84  volScalarField nutEff_;
85 
86 
87 protected:
88 
89  // Protected data
90 
91  // Model coefficients
92 
94 
95 
96  // Protected Member Functions
97 
98  virtual void correctNut();
100  virtual tmp<fvScalarMatrix> kSource() const;
101  virtual tmp<fvScalarMatrix> epsilonSource() const;
102 
103 
104 public:
106  typedef typename BasicMomentumTransportModel::alphaField alphaField;
107  typedef typename BasicMomentumTransportModel::rhoField rhoField;
108  typedef typename BasicMomentumTransportModel::transportModel transportModel;
109 
110 
111  //- Runtime type information
112  TypeName("continuousGasKEpsilon");
113 
114 
115  // Constructors
116 
117  //- Construct from components
119  (
120  const alphaField& alpha,
121  const rhoField& rho,
122  const volVectorField& U,
123  const surfaceScalarField& alphaRhoPhi,
124  const surfaceScalarField& phi,
125  const transportModel& transport,
126  const word& type = typeName
127  );
128 
129  //- Disallow default bitwise copy construction
131 
132 
133  //- Destructor
134  virtual ~continuousGasKEpsilon()
135  {}
136 
137 
138  // Member Functions
139 
140  //- Re-read model coefficients if they have changed
141  virtual bool read();
142 
143  //- Return the turbulence model for the liquid phase
145 
146  //- Return the effective viscosity
147  virtual tmp<volScalarField> nuEff() const;
148 
149  //- Return the effective density for the stress
150  virtual tmp<volScalarField> rhoEff() const;
151 
152  //- Return the stress tensor [m^2/s^2]
153  virtual tmp<volSymmTensorField> sigma() const;
154 
155 
156  // Member Operators
157 
158  //- Disallow default bitwise assignment
159  void operator=(const continuousGasKEpsilon&) = delete;
160 };
161 
162 
163 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
164 
165 } // End namespace RASModels
166 } // End namespace Foam
167 
168 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
169 
170 #ifdef NoRepository
171  #include "continuousGasKEpsilon.C"
172 #endif
173 
174 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
175 
176 #endif
177 
178 // ************************************************************************* //
virtual tmp< fvScalarMatrix > epsilonSource() const
virtual bool read()
Re-read model coefficients if they have changed.
BasicMomentumTransportModel::transportModel transportModel
phi
Definition: pEqn.H:104
BasicMomentumTransportModel::rhoField rhoField
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
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
Abstract base class for turbulence models (RAS, LES and laminar).
continuousGasKEpsilon(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &type=typeName)
Construct from components.
const momentumTransportModel & liquidTurbulence() const
Return the turbulence model for the liquid phase.
U
Definition: pEqn.H:72
tmp< volScalarField > phaseTransferCoeff() const
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
virtual tmp< volSymmTensorField > sigma() const
Return the stress tensor [m^2/s^2].
void operator=(const continuousGasKEpsilon &)=delete
Disallow default bitwise assignment.
virtual tmp< volScalarField > rhoEff() const
Return the effective density for the stress.
BasicMomentumTransportModel::alphaField alphaField
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
A class for managing temporary objects.
Definition: PtrList.H:53
Namespace for OpenFOAM.