continuousGasKEqn.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::LESModels::continuousGasKEqn
26 
27 Description
28  One-equation SGS 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-equation up to the phase-fraction of
33  alphaInversion at which point phase-inversion is considered to have occurred
34  and the model reverts to the pure single-phase form.
35 
36  This model is unpublished and is provided as a stable numerical framework
37  on which a more physical model may be built.
38 
39  The default model coefficients are
40  \verbatim
41  continuousKEqnCoeffs
42  {
43  Ck 0.094;
44  Ce 1.048;
45  alphaInversion 0.7;
46  }
47  \endverbatim
48 
49 SourceFiles
50  continuousGasKEqn.C
51 
52 \*---------------------------------------------------------------------------*/
53 
54 #ifndef continuousGasKEqn_H
55 #define continuousGasKEqn_H
56 
57 #include "kEqn.H"
58 
59 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
60 
61 namespace Foam
62 {
63 namespace LESModels
64 {
65 
66 /*---------------------------------------------------------------------------*\
67  Class continuousGasKEqn Declaration
68 \*---------------------------------------------------------------------------*/
69 
70 template<class BasicMomentumTransportModel>
72 :
73  public kEqn<BasicMomentumTransportModel>
74 {
75  // Private Data
76 
77  mutable const momentumTransportModel *liquidTurbulencePtr_;
78 
79 
80 protected:
81 
82  // Protected data
83 
84  // Model coefficients
85 
87 
88 
89  // Protected Member Functions
90 
91  //- Return the turbulence model for the liquid phase
93 
95  virtual tmp<fvScalarMatrix> kSource() const;
96 
97 
98 public:
99 
100  typedef typename BasicMomentumTransportModel::alphaField alphaField;
101  typedef typename BasicMomentumTransportModel::rhoField rhoField;
102  typedef typename BasicMomentumTransportModel::transportModel transportModel;
103 
104 
105  //- Runtime type information
106  TypeName("continuousGasKEqn");
107 
108 
109  // Constructors
110 
111  //- Construct from components
113  (
114  const alphaField& alpha,
115  const rhoField& rho,
116  const volVectorField& U,
117  const surfaceScalarField& alphaRhoPhi,
118  const surfaceScalarField& phi,
119  const transportModel& transport,
120  const word& type = typeName
121  );
122 
123  //- Disallow default bitwise copy construction
124  continuousGasKEqn(const continuousGasKEqn&) = delete;
125 
126 
127  //- Destructor
128  virtual ~continuousGasKEqn()
129  {}
130 
131 
132  // Member Functions
133 
134  //- Read model coefficients if they have changed
135  virtual bool read();
136 
137 
138  // Member Operators
139 
140  //- Disallow default bitwise assignment
141  void operator=(const continuousGasKEqn&) = delete;
142 };
143 
144 
145 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
146 
147 } // End namespace LESModels
148 } // End namespace Foam
149 
150 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
151 
152 #ifdef NoRepository
153  #include "continuousGasKEqn.C"
154 #endif
155 
156 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
157 
158 #endif
159 
160 // ************************************************************************* //
One equation eddy-viscosity model.
Definition: kEqn.H:71
BasicMomentumTransportModel::alphaField alphaField
phi
Definition: pEqn.H:104
continuousGasKEqn(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.
A class for handling words, derived from string.
Definition: word.H:59
virtual tmp< fvScalarMatrix > kSource() const
tmp< volScalarField > phaseTransferCoeff() const
virtual bool read()
Read model coefficients if they have changed.
One-equation SGS model for the gas-phase in a two-phase system supporting phase-inversion.
Abstract base class for turbulence models (RAS, LES and laminar).
TypeName("continuousGasKEqn")
Runtime type information.
U
Definition: pEqn.H:72
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
virtual ~continuousGasKEqn()
Destructor.
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
BasicMomentumTransportModel::transportModel transportModel
A class for managing temporary objects.
Definition: PtrList.H:53
void operator=(const continuousGasKEqn &)=delete
Disallow default bitwise assignment.
Namespace for OpenFOAM.
BasicMomentumTransportModel::rhoField rhoField