LaheyKEpsilon.C
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-2024 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 \*---------------------------------------------------------------------------*/
25 
26 #include "LaheyKEpsilon.H"
27 #include "fvModels.H"
28 #include "fvConstraints.H"
29 #include "phaseSystem.H"
30 #include "dispersedDragModel.H"
31 
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace RASModels
37 {
38 
39 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
40 
41 template<class BasicMomentumTransportModel>
43 (
44  const alphaField& alpha,
45  const rhoField& rho,
46  const volVectorField& U,
47  const surfaceScalarField& alphaRhoPhi,
48  const surfaceScalarField& phi,
49  const viscosity& viscosity,
50  const word& type
51 )
52 :
53  kEpsilon<BasicMomentumTransportModel>
54  (
55  alpha,
56  rho,
57  U,
58  alphaRhoPhi,
59  phi,
60  viscosity,
61  type
62  ),
63 
64  gasTurbulencePtr_(nullptr),
65 
66  alphaInversion_("alphaInversion", this->coeffDict(), 0.3),
67  Cp_("Cp", this->coeffDict(), 0.25),
68  C4_("C4", this->coeffDict(), this->C2_.value()),
69  Cmub_("Cmub", this->coeffDict(), 0.6)
70 {}
71 
72 
73 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
74 
75 template<class BasicMomentumTransportModel>
77 {
79  {
80  alphaInversion_.readIfPresent(this->coeffDict());
81  Cp_.readIfPresent(this->coeffDict());
82  C4_.readIfPresent(this->coeffDict());
83  Cmub_.readIfPresent(this->coeffDict());
84 
85  return true;
86  }
87  else
88  {
89  return false;
90  }
91 }
92 
93 
94 template<class BasicMomentumTransportModel>
97 {
98  if (!gasTurbulencePtr_)
99  {
100  const volVectorField& U = this->U_;
101 
102  const phaseModel& liquid =
103  refCast<const phaseModel>(this->properties());
104  const phaseSystem& fluid = liquid.fluid();
105  const phaseModel& gas = fluid.otherPhase(liquid);
106 
107  gasTurbulencePtr_ =
108  &U.db().lookupObject
109  <
111  >
112  (
114  (
115  momentumTransportModel::typeName,
116  gas.name()
117  )
118  );
119  }
120 
121  return *gasTurbulencePtr_;
122 }
123 
124 
125 template<class BasicMomentumTransportModel>
127 {
128  const phaseCompressibleMomentumTransportModel& gasTurbulence =
129  this->gasTurbulence();
130 
131  const phaseModel& liquid = refCast<const phaseModel>(this->properties());
132  const phaseSystem& fluid = liquid.fluid();
133  const phaseModel& gas = fluid.otherPhase(liquid);
134 
135  this->nut_ =
136  this->Cmu_*sqr(this->k_)/this->epsilon_
137  + Cmub_*gas.d()*gasTurbulence.alpha()
138  *(mag(this->U_ - gasTurbulence.U()));
139 
140  this->nut_.correctBoundaryConditions();
141  fvConstraints::New(this->mesh_).constrain(this->nut_);
142 }
143 
144 
145 template<class BasicMomentumTransportModel>
147 {
148  const phaseCompressibleMomentumTransportModel& gasTurbulence =
149  this->gasTurbulence();
150 
151  const phaseModel& liquid = refCast<const phaseModel>(this->properties());
152  const phaseSystem& fluid = liquid.fluid();
153  const phaseModel& gas = fluid.otherPhase(liquid);
154 
158 
159  volScalarField magUr(mag(this->U_ - gasTurbulence.U()));
160 
161  tmp<volScalarField> bubbleG
162  (
163  Cp_
164  *(
165  pow3(magUr)
166  + pow(drag.CdRe()*liquid.fluidThermo().nu()/gas.d(), 4.0/3.0)
167  *pow(magUr, 5.0/3.0)
168  )
169  *gas
170  /gas.d()
171  );
172 
173  return bubbleG;
174 }
175 
176 
177 template<class BasicMomentumTransportModel>
180 {
181  const volVectorField& U = this->U_;
182  const alphaField& alpha = this->alpha_;
183  const rhoField& rho = this->rho_;
184 
185  const momentumTransportModel& gasTurbulence = this->gasTurbulence();
186 
187  return
188  (
189  max(alphaInversion_ - alpha, scalar(0))
190  *rho
191  *min(gasTurbulence.epsilon()/gasTurbulence.k(), 1.0/U.time().deltaT())
192  );
193 }
194 
195 
196 template<class BasicMomentumTransportModel>
198 {
199  const alphaField& alpha = this->alpha_;
200  const rhoField& rho = this->rho_;
201 
202  const phaseCompressibleMomentumTransportModel& gasTurbulence =
203  this->gasTurbulence();
204 
205  const volScalarField phaseTransferCoeff(this->phaseTransferCoeff());
206 
207  return
208  alpha*rho*bubbleG()
209  + phaseTransferCoeff*gasTurbulence.k()
210  - fvm::Sp(phaseTransferCoeff, this->k_);
211 }
212 
213 
214 template<class BasicMomentumTransportModel>
217 {
218  const alphaField& alpha = this->alpha_;
219  const rhoField& rho = this->rho_;
220 
221  const phaseCompressibleMomentumTransportModel& gasTurbulence =
222  this->gasTurbulence();
223 
224  const volScalarField phaseTransferCoeff(this->phaseTransferCoeff());
225 
226  return
227  alpha*rho*this->C4_*this->epsilon_*bubbleG()/this->k_
228  + phaseTransferCoeff*gasTurbulence.epsilon()
229  - fvm::Sp(phaseTransferCoeff, this->epsilon_);
230 }
231 
232 
233 template<class BasicMomentumTransportModel>
235 {
237 }
238 
239 
240 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
241 
242 } // End namespace RASModels
243 } // End namespace Foam
244 
245 // ************************************************************************* //
Generic GeometricField class.
static word groupName(Name name, const word &group)
Continuous-phase k-epsilon model including bubble-generated turbulence.
Definition: LaheyKEpsilon.H:78
virtual tmp< fvScalarMatrix > epsilonSource() const
Source term for the epsilon equation.
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
LaheyKEpsilon(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &viscosity, const word &type=typeName)
Construct from components.
Definition: LaheyKEpsilon.C:43
tmp< volScalarField > bubbleG() const
tmp< volScalarField > phaseTransferCoeff() const
virtual void correctNut()
Correct the eddy-viscosity nut.
virtual tmp< fvScalarMatrix > kSource() const
Source term for the k equation.
virtual bool read()
Read model coefficients if they have changed.
Definition: LaheyKEpsilon.C:76
Standard k-epsilon turbulence model for incompressible and compressible flows including rapid distort...
Definition: kEpsilon.H:86
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: kEpsilon.C:174
static autoPtr< dictionary > New(Istream &)
Construct top-level dictionary on freestore from Istream.
Definition: dictionaryIO.C:103
Class to represent a interface between phases where one phase is considered dispersed within the othe...
Base class for Lagrangian drag models.
BasicMomentumTransportModel::alphaField alphaField
BasicMomentumTransportModel::rhoField rhoField
Generic thermophysical properties class for a liquid in which the functions and coefficients for each...
Definition: liquid.H:53
Abstract base class for turbulence models (RAS, LES and laminar).
virtual tmp< volScalarField > k() const =0
Return the turbulence kinetic energy.
const volVectorField & U() const
Access function to velocity field.
virtual tmp< volScalarField > epsilon() const =0
Return the turbulence kinetic energy dissipation rate.
Templated abstract base class for multiphase compressible turbulence models.
const alphaField & alpha() const
Access function to phase fraction.
tmp< volScalarField > d() const
Return the Sauter-mean diameter.
Definition: phaseModel.C:181
const word & name() const
Return the name of this phase.
Definition: phaseModel.C:145
Class to represent a system of phases.
Definition: phaseSystem.H:74
const phaseModel & otherPhase(const phaseModel &phase) const
Return the phase not given as an argument in a two-phase system.
Definition: phaseSystemI.H:174
const ModelType & lookupInterfacialModel(const phaseInterface &interface) const
Return a sub model for an interface.
A class for managing temporary objects.
Definition: tmp.H:55
Abstract base class for all fluid physical properties.
Definition: viscosity.H:50
A class for handling words, derived from string.
Definition: word.H:62
U
Definition: pEqn.H:72
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const VolField< Type > &)
Namespace for OpenFOAM.
void mag(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
void pow(LagrangianPatchField< typename powProduct< Type, r >::type > &f, const LagrangianPatchField< Type > &f1)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
void pow3(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
void sqr(LagrangianPatchField< typename outerProduct< Type, Type >::type > &f, const LagrangianPatchField< Type > &f1)
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488