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 | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2013-2016 OpenFOAM Foundation
9 Copyright (C) 2019 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27Class
28 Foam::RASModels::continuousGasKEpsilon
29
30Group
31 grpRASTurbulence
32
33Description
34 k-epsilon model for the gas-phase in a two-phase system
35 supporting phase-inversion.
36
37 In the limit that the gas-phase fraction approaches zero a contribution from
38 the other phase is blended into the k and epsilon equations up to the
39 phase-fraction of alphaInversion at which point phase-inversion is
40 considered to have occurred and the model reverts to the pure single-phase
41 form.
42
43 This model is unpublished and is provided as a stable numerical framework
44 on which a more physical model may be built.
45
46 The default model coefficients are
47 \verbatim
48 continuousGasKEpsilonCoeffs
49 {
50 Cmu 0.09;
51 C1 1.44;
52 C2 1.92;
53 C3 -0.33;
54 sigmak 1.0;
55 sigmaEps 1.3;
56 alphaInversion 0.7;
57 }
58 \endverbatim
59
60SourceFiles
61 continuousGasKEpsilon.C
62
63\*---------------------------------------------------------------------------*/
64
65#ifndef continuousGasKEpsilon_H
66#define continuousGasKEpsilon_H
67
68#include "kEpsilon.H"
69
70// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
71
72namespace Foam
73{
74namespace RASModels
75{
76
77/*---------------------------------------------------------------------------*\
78 Class continuousGasKEpsilon Declaration
79\*---------------------------------------------------------------------------*/
80
81template<class BasicTurbulenceModel>
83:
84 public kEpsilon<BasicTurbulenceModel>
85{
86 // Private data
87
88 mutable const turbulenceModel *liquidTurbulencePtr_;
89
90 volScalarField nutEff_;
91
92
93 // Private Member Functions
94
95 //- No copy construct
97
98 //- No copy assignment
99 void operator=(const continuousGasKEpsilon&) = delete;
100
101
102protected:
103
104 // Protected data
105
106 // Model coefficients
109
110
111 // Protected Member Functions
112
113 virtual void correctNut();
115 virtual tmp<fvScalarMatrix> kSource() const;
116 virtual tmp<fvScalarMatrix> epsilonSource() const;
117
118
119public:
121 typedef typename BasicTurbulenceModel::alphaField alphaField;
122 typedef typename BasicTurbulenceModel::rhoField rhoField;
123 typedef typename BasicTurbulenceModel::transportModel transportModel;
124
125
126 //- Runtime type information
127 TypeName("continuousGasKEpsilon");
128
129
130 // Constructors
131
132 //- Construct from components
134 (
135 const alphaField& alpha,
136 const rhoField& rho,
137 const volVectorField& U,
138 const surfaceScalarField& alphaRhoPhi,
139 const surfaceScalarField& phi,
140 const transportModel& transport,
141 const word& propertiesName = turbulenceModel::propertiesName,
142 const word& type = typeName
143 );
144
145
146 //- Destructor
147 virtual ~continuousGasKEpsilon() = default;
148
149
150 // Member Functions
151
152 //- Re-read model coefficients if they have changed
153 virtual bool read();
154
155 //- Return the turbulence model for the liquid phase
156 const turbulenceModel& liquidTurbulence() const;
157
158 //- Return the effective viscosity
159 virtual tmp<volScalarField> nuEff() const;
160
161 //- Return the effective density for the stress
162 virtual tmp<volScalarField> rhoEff() const;
163
164 //- Return the Reynolds stress tensor
165 virtual tmp<volSymmTensorField> R() const;
166};
167
168
169// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
170
171} // End namespace RASModels
172} // End namespace Foam
173
174// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
175
176#ifdef NoRepository
177 #include "continuousGasKEpsilon.C"
178#endif
179
180// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
181
182#endif
183
184// ************************************************************************* //
surfaceScalarField & phi
k-epsilon model for the gas-phase in a two-phase system supporting phase-inversion.
BasicTurbulenceModel::alphaField alphaField
virtual tmp< fvScalarMatrix > epsilonSource() const
const turbulenceModel & liquidTurbulence() const
Return the turbulence model for the liquid phase.
BasicTurbulenceModel::rhoField rhoField
virtual tmp< volScalarField > nuEff() const
Return the effective viscosity.
TypeName("continuousGasKEpsilon")
Runtime type information.
tmp< volScalarField > phaseTransferCoeff() const
virtual tmp< volScalarField > rhoEff() const
Return the effective density for the stress.
virtual tmp< volSymmTensorField > R() const
Return the Reynolds stress tensor.
virtual ~continuousGasKEpsilon()=default
Destructor.
virtual tmp< fvScalarMatrix > kSource() const
BasicTurbulenceModel::transportModel transportModel
virtual bool read()
Re-read model coefficients if they have changed.
Standard k-epsilon turbulence model for incompressible and compressible flows including rapid distort...
Definition: kEpsilon.H:92
A class for managing temporary objects.
Definition: tmp.H:65
Abstract base class for turbulence models (RAS, LES and laminar).
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from Foam::string.
Definition: word.H:68
U
Definition: pEqn.H:72
Namespace for OpenFOAM.
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:598
volScalarField & alpha
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition: typeInfo.H:73