NicenoKEqn.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::NicenoKEqn
26 
27 Description
28  One-equation SGS model for the continuous phase in a two-phase system
29  including bubble-generated turbulence.
30 
31  Reference:
32  \verbatim
33  Niceno, B., Dhotre, M. T., & Deen, N. G. (2008).
34  One-equation sub-grid scale (SGS) modelling for
35  Euler–Euler large eddy simulation (EELES) of dispersed bubbly flow.
36  Chemical Engineering Science, 63(15), 3923-3931.
37  \endverbatim
38 
39  The default model coefficients are:
40  \verbatim
41  NicenoKEqnCoeffs
42  {
43  Ck 0.094;
44  Ce 1.048;
45  alphaInversion 0.3;
46  Cp Ck;
47  Cmub 0.6;
48  }
49  \endverbatim
50 
51 SourceFiles
52  NicenoKEqn.C
53 
54 \*---------------------------------------------------------------------------*/
55 
56 #ifndef NicenoKEqn_H
57 #define NicenoKEqn_H
58 
59 #include "kEqn.H"
60 
61 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
62 
63 namespace Foam
64 {
65 namespace LESModels
66 {
67 
68 /*---------------------------------------------------------------------------*\
69  Class NicenoKEqn Declaration
70 \*---------------------------------------------------------------------------*/
71 
72 template<class BasicMomentumTransportModel>
73 class NicenoKEqn
74 :
75  public kEqn<BasicMomentumTransportModel>
76 {
77  // Private Data
78 
80  <
81  typename BasicMomentumTransportModel::transportModel
82  > *gasTurbulencePtr_;
83 
84 
85  // Private Member Functions
86 
87  //- Return the turbulence model for the gas phase
89  <
90  typename BasicMomentumTransportModel::transportModel
91  >&
92  gasTurbulence() const;
93 
94 
95 protected:
96 
97  // Protected data
98 
99  // Model coefficients
104 
105 
106  // Protected Member Functions
107 
108  virtual void correctNut();
111  virtual tmp<fvScalarMatrix> kSource() const;
112 
113 
114 public:
116  typedef typename BasicMomentumTransportModel::alphaField alphaField;
117  typedef typename BasicMomentumTransportModel::rhoField rhoField;
118  typedef typename BasicMomentumTransportModel::transportModel transportModel;
119 
120 
121  //- Runtime type information
122  TypeName("NicenoKEqn");
123 
124 
125  // Constructors
126 
127  //- Construct from components
128  NicenoKEqn
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& type = typeName
137  );
138 
139  //- Disallow default bitwise copy construction
140  NicenoKEqn(const NicenoKEqn&) = delete;
141 
142 
143  //- Destructor
144  virtual ~NicenoKEqn()
145  {}
146 
147 
148  // Member Functions
149 
150  //- Read model coefficients if they have changed
151  virtual bool read();
152 
153 
154  // Member Operators
155 
156  //- Disallow default bitwise assignment
157  void operator=(const NicenoKEqn&) = delete;
158 };
159 
160 
161 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
162 
163 } // End namespace LESModels
164 } // End namespace Foam
165 
166 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
167 
168 #ifdef NoRepository
169  #include "NicenoKEqn.C"
170 #endif
171 
172 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
173 
174 #endif
175 
176 // ************************************************************************* //
Templated abstract base class for multiphase compressible turbulence models.
BasicMomentumTransportModel::rhoField rhoField
Definition: NicenoKEqn.H:116
dimensionedScalar Cmub_
Definition: NicenoKEqn.H:102
dimensionedScalar Cp_
Definition: NicenoKEqn.H:101
One equation eddy-viscosity model.
Definition: kEqn.H:71
phi
Definition: pEqn.H:104
void operator=(const NicenoKEqn &)=delete
Disallow default bitwise assignment.
virtual ~NicenoKEqn()
Destructor.
Definition: NicenoKEqn.H:143
virtual bool read()
Read model coefficients if they have changed.
Definition: NicenoKEqn.C:105
A class for handling words, derived from string.
Definition: word.H:59
NicenoKEqn(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.
Definition: NicenoKEqn.C:42
tmp< volScalarField > bubbleG() const
Definition: NicenoKEqn.C:172
tmp< volScalarField > phaseTransferCoeff() const
Definition: NicenoKEqn.C:196
BasicMomentumTransportModel::alphaField alphaField
Definition: NicenoKEqn.H:115
virtual tmp< fvScalarMatrix > kSource() const
Definition: NicenoKEqn.C:218
U
Definition: pEqn.H:72
One-equation SGS model for the continuous phase in a two-phase system including bubble-generated turb...
Definition: NicenoKEqn.H:72
TypeName("NicenoKEqn")
Runtime type information.
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
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
virtual void correctNut()
Definition: NicenoKEqn.C:156
BasicMomentumTransportModel::transportModel transportModel
Definition: NicenoKEqn.H:117
Namespace for OpenFOAM.
dimensionedScalar alphaInversion_
Definition: NicenoKEqn.H:100