continuousGasKEpsilon.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2013-2017 OpenFOAM Foundation
9  Copyright (C) 2019-2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "continuousGasKEpsilon.H"
30 #include "fvOptions.H"
31 #include "twoPhaseSystem.H"
32 #include "virtualMassModel.H"
33 #include "dragModel.H"
34 
35 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39 namespace RASModels
40 {
41 
42 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
43 
44 template<class BasicTurbulenceModel>
45 continuousGasKEpsilon<BasicTurbulenceModel>::continuousGasKEpsilon
46 (
47  const alphaField& alpha,
48  const rhoField& rho,
49  const volVectorField& U,
50  const surfaceScalarField& alphaRhoPhi,
51  const surfaceScalarField& phi,
52  const transportModel& transport,
53  const word& propertiesName,
54  const word& type
55 )
56 :
58  (
59  alpha,
60  rho,
61  U,
62  alphaRhoPhi,
63  phi,
64  transport,
65  propertiesName,
66  type
67  ),
68 
69  liquidTurbulencePtr_(nullptr),
70 
71  nutEff_
72  (
73  IOobject
74  (
75  IOobject::groupName("nutEff", alphaRhoPhi.group()),
76  this->runTime_.timeName(),
77  this->mesh_,
80  ),
81  this->nut_
82  ),
83 
84  alphaInversion_
85  (
87  (
88  "alphaInversion",
89  this->coeffDict_,
90  0.7
91  )
92  )
93 {
94  if (type == typeName)
95  {
96  this->printCoeffs(type);
97  }
98 }
99 
100 
101 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
102 
103 template<class BasicTurbulenceModel>
105 {
107  {
108  alphaInversion_.readIfPresent(this->coeffDict());
109 
110  return true;
111  }
112 
113  return false;
114 }
115 
116 
117 template<class BasicTurbulenceModel>
119 {
121 
122  const turbulenceModel& liquidTurbulence = this->liquidTurbulence();
123  const transportModel& gas = this->transport();
124  const twoPhaseSystem& fluid = refCast<const twoPhaseSystem>(gas.fluid());
125  const transportModel& liquid = fluid.otherPhase(gas);
126 
127  const virtualMassModel& virtualMass =
128  fluid.lookupSubModel<virtualMassModel>(gas, liquid);
129 
130  volScalarField thetal(liquidTurbulence.k()/liquidTurbulence.epsilon());
131  volScalarField rhodv(gas.rho() + virtualMass.Cvm()*liquid.rho());
132  volScalarField thetag((rhodv/(18*liquid.rho()*liquid.nu()))*sqr(gas.d()));
133  volScalarField expThetar
134  (
135  min
136  (
137  exp(min(thetal/thetag, scalar(50))),
138  scalar(1)
139  )
140  );
141  volScalarField omega((1 - expThetar)/(1 + expThetar));
142 
143  nutEff_ = omega*liquidTurbulence.nut();
144  fv::options::New(this->mesh_).correct(nutEff_);
145 }
146 
147 
148 template<class BasicTurbulenceModel>
149 const turbulenceModel&
151 {
152  if (!liquidTurbulencePtr_)
153  {
154  const volVectorField& U = this->U_;
155 
156  const transportModel& gas = this->transport();
157  const twoPhaseSystem& fluid =
158  refCast<const twoPhaseSystem>(gas.fluid());
159  const transportModel& liquid = fluid.otherPhase(gas);
160 
161  liquidTurbulencePtr_ =
162  &U.db().lookupObject<turbulenceModel>
163  (
165  (
167  liquid.name()
168  )
169  );
170  }
171 
172  return *liquidTurbulencePtr_;
173 }
174 
175 
176 template<class BasicTurbulenceModel>
179 {
180  volScalarField blend
181  (
182  max
183  (
184  min
185  (
186  (this->alpha_ - scalar(0.5))/(alphaInversion_ - 0.5),
187  scalar(1)
188  ),
189  scalar(0)
190  )
191  );
192 
193  return tmp<volScalarField>
194  (
195  new volScalarField
196  (
197  IOobject::groupName("nuEff", this->alphaRhoPhi_.group()),
198  blend*this->nut_
199  + (1.0 - blend)*rhoEff()*nutEff_/this->transport().rho()
200  + this->nu()
201  )
202  );
203 }
204 
205 
206 template<class BasicTurbulenceModel>
209 {
210  const transportModel& gas = this->transport();
211  const twoPhaseSystem& fluid = refCast<const twoPhaseSystem>(gas.fluid());
212  const transportModel& liquid = fluid.otherPhase(gas);
213 
214  const virtualMassModel& virtualMass =
215  fluid.lookupSubModel<virtualMassModel>(gas, liquid);
216 
217  return tmp<volScalarField>
218  (
219  new volScalarField
220  (
221  IOobject::groupName("rhoEff", this->alphaRhoPhi_.group()),
222  gas.rho() + (virtualMass.Cvm() + 3.0/20.0)*liquid.rho()
223  )
224  );
225 }
226 
227 
228 template<class BasicTurbulenceModel>
231 {
232  const volVectorField& U = this->U_;
233  const alphaField& alpha = this->alpha_;
234  const rhoField& rho = this->rho_;
235 
236  const turbulenceModel& liquidTurbulence = this->liquidTurbulence();
237 
238  return
239  (
240  max(alphaInversion_ - alpha, scalar(0))
241  *rho
242  *min
243  (
244  liquidTurbulence.epsilon()/liquidTurbulence.k(),
245  1.0/U.time().deltaT()
246  )
247  );
248 }
249 
250 
251 template<class BasicTurbulenceModel>
254 {
255  const turbulenceModel& liquidTurbulence = this->liquidTurbulence();
256  const volScalarField phaseTransferCoeff(this->phaseTransferCoeff());
257 
258  return
259  phaseTransferCoeff*liquidTurbulence.k()
260  - fvm::Sp(phaseTransferCoeff, this->k_);
261 }
262 
263 
264 template<class BasicTurbulenceModel>
267 {
268  const turbulenceModel& liquidTurbulence = this->liquidTurbulence();
269  const volScalarField phaseTransferCoeff(this->phaseTransferCoeff());
270 
271  return
272  phaseTransferCoeff*liquidTurbulence.epsilon()
273  - fvm::Sp(phaseTransferCoeff, this->epsilon_);
274 }
275 
276 
277 template<class BasicTurbulenceModel>
280 {
281  tmp<volScalarField> tk(this->k());
282 
284  (
286  (
287  IOobject
288  (
289  IOobject::groupName("R", this->alphaRhoPhi_.group()),
290  this->runTime_.timeName(),
291  this->mesh_,
294  ),
295  ((2.0/3.0)*I)*tk() - (nutEff_)*dev(twoSymm(fvc::grad(this->U_))),
296  tk().boundaryField().types()
297  )
298  );
299 }
300 
301 
302 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
303 
304 } // End namespace RASModels
305 } // End namespace Foam
306 
307 // ************************************************************************* //
Foam::IOobject::NO_WRITE
Definition: IOobject.H:195
Foam::virtualMassModel
Definition: virtualMassModel.H:55
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
Foam::IOobject::AUTO_WRITE
Definition: IOobject.H:194
Foam::fvc::grad
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:54
fvOptions.H
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::twoPhaseSystem
Class which solves the volume fraction equations for two phases.
Definition: twoPhaseSystem.H:53
Foam::fv::optionList::correct
void correct(GeometricField< Type, fvPatchField, volMesh > &field)
Apply correction to field.
Definition: fvOptionListTemplates.C:355
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::turbulenceModel::nut
virtual tmp< volScalarField > nut() const =0
Return the turbulence viscosity.
Foam::fv::options::New
static options & New(const fvMesh &mesh)
Construct fvOptions and register to database if not present.
Definition: fvOptions.C:103
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:212
Foam::RASModels::continuousGasKEpsilon::rhoEff
virtual tmp< volScalarField > rhoEff() const
Return the effective density for the stress.
Definition: continuousGasKEpsilon.C:208
Foam::turbulenceModel::propertiesName
static const word propertiesName
Default name of the turbulence properties dictionary.
Definition: turbulenceModel.H:100
Foam::liquid
Generic thermophysical properties class for a liquid in which the functions and coefficients for each...
Definition: liquid.H:54
Foam::RASModels::continuousGasKEpsilon::transportModel
BasicTurbulenceModel::transportModel transportModel
Definition: continuousGasKEpsilon.H:122
Foam::exp
dimensionedScalar exp(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:261
fluid
twoPhaseSystem & fluid
Definition: setRegionFluidFields.H:3
rho
rho
Definition: readInitialConditions.H:88
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
nu
volScalarField & nu
Definition: readMechanicalProperties.H:176
Foam::RASModels::continuousGasKEpsilon::rhoField
BasicTurbulenceModel::rhoField rhoField
Definition: continuousGasKEpsilon.H:121
Foam::fvm::Sp
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
Foam::virtualMassModel::Cvm
virtual tmp< volScalarField > Cvm() const =0
Return the virtual mass coefficient.
Foam::RASModels::continuousGasKEpsilon::alphaField
BasicTurbulenceModel::alphaField alphaField
Definition: continuousGasKEpsilon.H:120
Foam::IOobject::READ_IF_PRESENT
Definition: IOobject.H:187
Foam::RASModels::continuousGasKEpsilon::liquidTurbulence
const turbulenceModel & liquidTurbulence() const
Return the turbulence model for the liquid phase.
Definition: continuousGasKEpsilon.C:150
Foam::RASModels::continuousGasKEpsilon::read
virtual bool read()
Re-read model coefficients if they have changed.
Definition: continuousGasKEpsilon.C:104
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::RASModels::continuousGasKEpsilon::correctNut
virtual void correctNut()
Definition: continuousGasKEpsilon.C:118
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
Foam::turbulenceModel
Abstract base class for turbulence models (RAS, LES and laminar).
Definition: turbulenceModel.H:63
Foam::RASModels::continuousGasKEpsilon::nuEff
virtual tmp< volScalarField > nuEff() const
Return the effective viscosity.
Definition: continuousGasKEpsilon.C:178
Foam::RASModels::kEpsilon::correctNut
virtual void correctNut()
Definition: kEpsilon.C:43
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:42
Foam::turbulenceModel::k
virtual tmp< volScalarField > k() const =0
Return the turbulence kinetic energy.
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::RASModels::continuousGasKEpsilon::phaseTransferCoeff
tmp< volScalarField > phaseTransferCoeff() const
Definition: continuousGasKEpsilon.C:230
continuousGasKEpsilon.H
Foam::liquid::rho
scalar rho(scalar p, scalar T) const
Liquid density [kg/m^3].
Definition: liquidI.H:28
U
U
Definition: pEqn.H:72
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::RASModels::continuousGasKEpsilon::epsilonSource
virtual tmp< fvScalarMatrix > epsilonSource() const
Definition: continuousGasKEpsilon.C:266
Foam::RASModels::continuousGasKEpsilon::R
virtual tmp< volSymmTensorField > R() const
Return the Reynolds stress tensor.
Definition: continuousGasKEpsilon.C:279
Foam::type
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:590
k
label k
Boltzmann constant.
Definition: LISASMDCalcMethod2.H:41
Foam::IOobject::groupName
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
Foam::GeometricField< vector, fvPatchField, volMesh >
Foam::RASModels::kEpsilon
Standard k-epsilon turbulence model for incompressible and compressible flows including rapid distort...
Definition: kEpsilon.H:89
Foam::IOobject::NO_READ
Definition: IOobject.H:188
Foam::RASModels::continuousGasKEpsilon::kSource
virtual tmp< fvScalarMatrix > kSource() const
Definition: continuousGasKEpsilon.C:253
Foam::turbulenceModel::epsilon
virtual tmp< volScalarField > epsilon() const =0
Return the turbulence kinetic energy dissipation rate.
Foam::twoSymm
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:95
Foam::I
static const Identity< scalar > I
Definition: Identity.H:95
Foam::dev
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:106