continuousGasKEqn.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-2022 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 "continuousGasKEqn.H"
27 #include "phaseSystem.H"
28 
29 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace LESModels
34 {
35 
36 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
37 
38 template<class BasicMomentumTransportModel>
40 (
41  const alphaField& alpha,
42  const rhoField& rho,
43  const volVectorField& U,
44  const surfaceScalarField& alphaRhoPhi,
45  const surfaceScalarField& phi,
46  const viscosity& viscosity,
47  const word& type
48 )
49 :
50  kEqn<BasicMomentumTransportModel>
51  (
52  alpha,
53  rho,
54  U,
55  alphaRhoPhi,
56  phi,
57  viscosity,
58  type
59  ),
60 
61  liquidTurbulencePtr_(nullptr),
62 
63  alphaInversion_
64  (
65  dimensioned<scalar>::lookupOrAddToDict
66  (
67  "alphaInversion",
68  this->coeffDict_,
69  0.7
70  )
71  )
72 {
73  if (type == typeName)
74  {
75  this->printCoeffs(type);
76  }
77 }
78 
79 
80 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
81 
82 template<class BasicMomentumTransportModel>
84 {
86  {
87  alphaInversion_.readIfPresent(this->coeffDict());
88 
89  return true;
90  }
91  else
92  {
93  return false;
94  }
95 }
96 
97 
98 template<class BasicMomentumTransportModel>
101 {
102  if (!liquidTurbulencePtr_)
103  {
104  const volVectorField& U = this->U_;
105 
106  const phaseModel& gas = refCast<const phaseModel>(this->properties());
107  const phaseSystem& fluid = gas.fluid();
108  const phaseModel& liquid = fluid.otherPhase(gas);
109 
110  liquidTurbulencePtr_ =
111  &U.db().lookupType<momentumTransportModel>(liquid.name());
112  }
113 
114  return *liquidTurbulencePtr_;
115 }
116 
117 
118 template<class BasicMomentumTransportModel>
121 {
122  const volVectorField& U = this->U_;
123  const alphaField& alpha = this->alpha_;
124  const rhoField& rho = this->rho_;
125 
126  const momentumTransportModel& liquidTurbulence = this->liquidTurbulence();
127 
128  return
129  (
130  max(alphaInversion_ - alpha, scalar(0))
131  *rho
132  *min
133  (
134  this->Ce_*sqrt(liquidTurbulence.k())/this->delta(),
135  1.0/U.time().deltaT()
136  )
137  );
138 }
139 
140 
141 template<class BasicMomentumTransportModel>
144 {
145  const momentumTransportModel& liquidTurbulence = this->liquidTurbulence();
146  const volScalarField phaseTransferCoeff(this->phaseTransferCoeff());
147 
148  return
149  phaseTransferCoeff*liquidTurbulence.k()
150  - fvm::Sp(phaseTransferCoeff, this->k_);
151 }
152 
153 
154 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
155 
156 } // End namespace LESModels
157 } // End namespace Foam
158 
159 // ************************************************************************* //
Generic GeometricField class.
tmp< volScalarField > phaseTransferCoeff() const
virtual tmp< fvScalarMatrix > kSource() const
const momentumTransportModel & liquidTurbulence() const
Return the turbulence model for the liquid phase.
virtual bool read()
Read model coefficients if they have changed.
continuousGasKEqn(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.
One equation eddy-viscosity model.
Definition: kEqn.H:74
BasicMomentumTransportModel::alphaField alphaField
Definition: kEqn.H:97
BasicMomentumTransportModel::rhoField rhoField
Definition: kEqn.H:98
Generic dimensioned Type class.
virtual const word & name() const
Return the name of the liquid.
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 phaseSystem & fluid() const
Return the system to which this phase belongs.
Definition: phaseModel.C:127
Class to represent a system of phases and model interfacial transfers between them.
Definition: phaseSystem.H:73
const phaseModel & otherPhase(const phaseModel &phase) const
Return the phase not given as an argument in a two-phase system.
Definition: phaseSystemI.H:111
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.
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
dimensionedScalar sqrt(const dimensionedScalar &ds)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488