NicenoKEqn.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 "NicenoKEqn.H"
27 #include "fvModels.H"
28 #include "fvConstraints.H"
29 #include "phaseSystem.H"
30 #include "dragModel.H"
31 
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace LESModels
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  kEqn<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(), this->Ck_.value()),
68  Cmub_("Cmub", this->coeffDict(), 0.6)
69 {}
70 
71 
72 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
73 
74 template<class BasicMomentumTransportModel>
76 {
78  {
79  alphaInversion_.readIfPresent(this->coeffDict());
80  Cp_.readIfPresent(this->coeffDict());
81  Cmub_.readIfPresent(this->coeffDict());
82 
83  return true;
84  }
85  else
86  {
87  return false;
88  }
89 }
90 
91 
92 template<class BasicMomentumTransportModel>
95 {
96  if (!gasTurbulencePtr_)
97  {
98  const volVectorField& U = this->U_;
99 
100  const phaseModel& liquid =
101  refCast<const phaseModel>(this->properties());
102  const phaseSystem& fluid = liquid.fluid();
103  const phaseModel& gas = fluid.otherPhase(liquid);
104 
105  gasTurbulencePtr_ =
106  &U.db().lookupObject
107  <
109  >
110  (
112  (
113  momentumTransportModel::typeName,
114  gas.name()
115  )
116  );
117  }
118 
119  return *gasTurbulencePtr_;
120 }
121 
122 
123 template<class BasicMomentumTransportModel>
125 {
126  const phaseCompressibleMomentumTransportModel& gasTurbulence =
127  this->gasTurbulence();
128 
129  const phaseModel& liquid = refCast<const phaseModel>(this->properties());
130  const phaseSystem& fluid = liquid.fluid();
131  const phaseModel& gas = fluid.otherPhase(liquid);
132 
133  this->nut_ =
134  this->Ck_*sqrt(this->k_)*this->delta()
135  + Cmub_*gas.d()*gasTurbulence.alpha()
136  *(mag(this->U_ - gasTurbulence.U()));
137 
138  this->nut_.correctBoundaryConditions();
139  fvConstraints::New(this->mesh_).constrain(this->nut_);
140 }
141 
142 
143 template<class BasicMomentumTransportModel>
145 {
146  const phaseCompressibleMomentumTransportModel& gasTurbulence =
147  this->gasTurbulence();
148 
149  const phaseModel& liquid = refCast<const phaseModel>(this->properties());
150  const phaseSystem& fluid = liquid.fluid();
151  const phaseModel& gas = fluid.otherPhase(liquid);
152 
156 
157  volScalarField magUr(mag(this->U_ - gasTurbulence.U()));
158 
159  tmp<volScalarField> bubbleG
160  (
161  Cp_*sqr(magUr)*drag.K()/liquid.rho()
162  );
163 
164  return bubbleG;
165 }
166 
167 
168 template<class BasicMomentumTransportModel>
171 {
172  const volVectorField& U = this->U_;
173  const alphaField& alpha = this->alpha_;
174  const rhoField& rho = this->rho_;
175 
176  const momentumTransportModel& gasTurbulence = this->gasTurbulence();
177 
178  return
179  (
180  max(alphaInversion_ - alpha, scalar(0))
181  *rho
182  *min
183  (
184  this->Ce_*sqrt(gasTurbulence.k())/this->delta(),
185  1.0/U.time().deltaT()
186  )
187  );
188 }
189 
190 
191 template<class BasicMomentumTransportModel>
193 {
194  const alphaField& alpha = this->alpha_;
195  const rhoField& rho = this->rho_;
196 
197  const phaseCompressibleMomentumTransportModel& gasTurbulence =
198  this->gasTurbulence();
199 
200  const volScalarField phaseTransferCoeff(this->phaseTransferCoeff());
201 
202  return
203  alpha*rho*bubbleG()
204  + phaseTransferCoeff*gasTurbulence.k()
205  - fvm::Sp(phaseTransferCoeff, this->k_);
206 }
207 
208 
209 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
210 
211 } // End namespace LESModels
212 } // End namespace Foam
213 
214 // ************************************************************************* //
scalar delta
Generic GeometricField class.
static word groupName(Name name, const word &group)
One-equation SGS model for the continuous phase in a two-phase system including bubble-generated turb...
Definition: NicenoKEqn.H:75
tmp< volScalarField > bubbleG() const
Definition: NicenoKEqn.C:144
NicenoKEqn(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: NicenoKEqn.C:43
tmp< volScalarField > phaseTransferCoeff() const
Definition: NicenoKEqn.C:170
virtual void correctNut()
Definition: NicenoKEqn.C:124
virtual tmp< fvScalarMatrix > kSource() const
Definition: NicenoKEqn.C:192
virtual bool read()
Read model coefficients if they have changed.
Definition: NicenoKEqn.C:75
One equation eddy-viscosity model.
Definition: kEqn.H:74
BasicMomentumTransportModel::alphaField alphaField
Definition: kEqn.H:97
BasicMomentumTransportModel::rhoField rhoField
Definition: kEqn.H:98
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.
Generic thermophysical properties class for a liquid in which the functions and coefficients for each...
Definition: liquid.H:53
scalar rho(scalar p, scalar T) const
Liquid density [kg/m^3].
Definition: liquidI.H:26
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.
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)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
void sqr(LagrangianPatchField< typename outerProduct< Type, Type >::type > &f, const LagrangianPatchField< Type > &f1)
void sqrt(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488