continuousGasKEpsilon.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-2021 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 "continuousGasKEpsilon.H"
27 #include "fvModels.H"
28 #include "fvConstraints.H"
29 #include "phaseSystem.H"
30 #include "dragModel.H"
31 #include "virtualMassModel.H"
32 
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 namespace RASModels
38 {
39 
40 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
41 
42 template<class BasicMomentumTransportModel>
44 (
45  const alphaField& alpha,
46  const rhoField& rho,
47  const volVectorField& U,
48  const surfaceScalarField& alphaRhoPhi,
49  const surfaceScalarField& phi,
50  const transportModel& transport,
51  const word& type
52 )
53 :
55  (
56  alpha,
57  rho,
58  U,
59  alphaRhoPhi,
60  phi,
61  transport,
62  type
63  ),
64 
65  liquidTurbulencePtr_(nullptr),
66 
67  nutEff_
68  (
69  IOobject
70  (
71  IOobject::groupName("nutEff", alphaRhoPhi.group()),
72  this->runTime_.timeName(),
73  this->mesh_,
76  ),
77  this->nut_
78  ),
79 
80  alphaInversion_
81  (
83  (
84  "alphaInversion",
85  this->coeffDict_,
86  0.7
87  )
88  )
89 {
90  if (type == typeName)
91  {
92  this->printCoeffs(type);
93  }
94 }
95 
96 
97 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
98 
99 template<class BasicMomentumTransportModel>
101 {
103  {
104  alphaInversion_.readIfPresent(this->coeffDict());
105 
106  return true;
107  }
108  else
109  {
110  return false;
111  }
112 }
113 
114 
115 template<class BasicMomentumTransportModel>
117 {
119 
120  const momentumTransportModel& liquidTurbulence = this->liquidTurbulence();
121 
122  const phaseModel& gas = refCast<const phaseModel>(this->transport());
123  const phaseSystem& fluid = gas.fluid();
124  const phaseModel& liquid = fluid.otherPhase(gas);
125 
126  const virtualMassModel& virtualMass =
127  fluid.lookupSubModel<virtualMassModel>(gas, liquid);
128 
129  volScalarField thetal(liquidTurbulence.k()/liquidTurbulence.epsilon());
130  volScalarField rhodv(gas.rho() + virtualMass.Cvm()*liquid.rho());
131  volScalarField thetag
132  (
133  (rhodv/(18*liquid.rho()*liquid.thermo().nu()))*sqr(gas.d())
134  );
135  volScalarField expThetar
136  (
137  min
138  (
139  exp(min(thetal/thetag, scalar(50))),
140  scalar(1)
141  )
142  );
143  volScalarField omega((1 - expThetar)/(1 + expThetar));
144 
145  nutEff_ = omega*liquidTurbulence.nut();
146  fvConstraints::New(this->mesh_).constrain(nutEff_);
147 }
148 
149 
150 template<class BasicMomentumTransportModel>
153 {
154  if (!liquidTurbulencePtr_)
155  {
156  const volVectorField& U = this->U_;
157 
158  const phaseModel& gas = refCast<const phaseModel>(this->transport());
159  const phaseSystem& fluid = gas.fluid();
160  const phaseModel& liquid = fluid.otherPhase(gas);
161 
162  liquidTurbulencePtr_ =
164  (
166  (
167  momentumTransportModel::typeName,
168  liquid.name()
169  )
170  );
171  }
172 
173  return *liquidTurbulencePtr_;
174 }
175 
176 
177 template<class BasicMomentumTransportModel>
180 {
181  volScalarField blend
182  (
183  max
184  (
185  min
186  (
187  (this->alpha_ - scalar(0.5))/(alphaInversion_ - 0.5),
188  scalar(1)
189  ),
190  scalar(0)
191  )
192  );
193 
194  return volScalarField::New
195  (
196  IOobject::groupName("nuEff", this->alphaRhoPhi_.group()),
197  blend*this->nut_
198  + (1.0 - blend)*rhoEff()*nutEff_/this->rho_
199  + this->nu()
200  );
201 }
202 
203 
204 template<class BasicMomentumTransportModel>
207 {
208  const phaseModel& gas = refCast<const phaseModel>(this->transport());
209  const phaseSystem& fluid = gas.fluid();
210  const phaseModel& liquid = fluid.otherPhase(gas);
211 
212  const virtualMassModel& virtualMass =
213  fluid.lookupSubModel<virtualMassModel>(gas, liquid);
214 
215  return volScalarField::New
216  (
217  IOobject::groupName("rhoEff", this->alphaRhoPhi_.group()),
218  gas.rho() + (virtualMass.Cvm() + 3.0/20.0)*liquid.rho()
219  );
220 }
221 
222 
223 template<class BasicMomentumTransportModel>
226 {
227  const volVectorField& U = this->U_;
228  const alphaField& alpha = this->alpha_;
229  const rhoField& rho = this->rho_;
230 
231  const momentumTransportModel& liquidTurbulence = this->liquidTurbulence();
232 
233  return
234  (
235  max(alphaInversion_ - alpha, scalar(0))
236  *rho
237  *min
238  (
239  liquidTurbulence.epsilon()/liquidTurbulence.k(),
240  1.0/U.time().deltaT()
241  )
242  );
243 }
244 
245 
246 template<class BasicMomentumTransportModel>
249 {
250  const momentumTransportModel& liquidTurbulence = this->liquidTurbulence();
251  const volScalarField phaseTransferCoeff(this->phaseTransferCoeff());
252 
253  return
254  phaseTransferCoeff*liquidTurbulence.k()
255  - fvm::Sp(phaseTransferCoeff, this->k_);
256 }
257 
258 
259 template<class BasicMomentumTransportModel>
262 {
263  const momentumTransportModel& liquidTurbulence = this->liquidTurbulence();
264  const volScalarField phaseTransferCoeff(this->phaseTransferCoeff());
265 
266  return
267  phaseTransferCoeff*liquidTurbulence.epsilon()
268  - fvm::Sp(phaseTransferCoeff, this->epsilon_);
269 }
270 
271 
272 template<class BasicMomentumTransportModel>
275 {
276  tmp<volScalarField> tk(this->k());
277 
279  (
280  IOobject::groupName("R", this->alphaRhoPhi_.group()),
281  ((2.0/3.0)*I)*tk() - (nutEff_)*dev(twoSymm(fvc::grad(this->U_))),
282  tk().boundaryField().types()
283  );
284 }
285 
286 
287 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
288 
289 } // End namespace RASModels
290 } // End namespace Foam
291 
292 // ************************************************************************* //
Generic thermophysical properties class for a liquid in which the functions and coefficients for each...
Definition: liquid.H:50
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
virtual tmp< fvScalarMatrix > epsilonSource() const
dimensionedSymmTensor sqr(const dimensionedVector &dv)
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Internal &, const PtrList< fvPatchField< scalar >> &)
Return a temporary field constructed from name,.
virtual tmp< volScalarField > rho() const =0
Return the density field.
const word & name() const
Definition: phaseModel.H:110
virtual bool read()
Re-read model coefficients if they have changed.
label k
Boltzmann constant.
Generic dimensioned Type class.
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
dimensionedScalar exp(const dimensionedScalar &ds)
static const Identity< scalar > I
Definition: Identity.H:93
Standard k-epsilon turbulence model for incompressible and compressible flows including rapid distort...
Definition: kEpsilon.H:83
Class to represent a system of phases and model interfacial transfers between them.
Definition: phaseSystem.H:70
A class for handling words, derived from string.
Definition: word.H:59
const rhoThermo & thermo() const
Return const-access to phase rhoThermo.
Definition: phaseModel.H:121
static word groupName(Name name, const word &group)
BasicMomentumTransportModel::alphaField alphaField
Definition: kEpsilon.H:115
phaseSystem & fluid
Definition: createFields.H:11
phi
Definition: correctPhi.H:3
virtual tmp< volScalarField > nuEff() const
Return the effective viscosity.
virtual tmp< fvScalarMatrix > kSource() const
virtual tmp< volScalarField > nut() const =0
Return the turbulence viscosity.
static autoPtr< dictionary > New(Istream &)
Construct top-level dictionary on freestore from Istream.
Definition: dictionaryIO.C:96
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
virtual tmp< volScalarField > Cvm() const =0
Return the virtual mass coefficient.
virtual tmp< volScalarField > k() const =0
Return the turbulence kinetic energy.
Abstract base class for turbulence models (RAS, LES and laminar).
continuousGasKEpsilon(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.
const momentumTransportModel & liquidTurbulence() const
Return the turbulence model for the liquid phase.
U
Definition: pEqn.H:72
tmp< volScalarField > phaseTransferCoeff() const
const Time & time() const
Return time.
Definition: IOobject.C:318
const phaseSystem & fluid() const
Return the system to which this phase belongs.
BasicMomentumTransportModel::rhoField rhoField
Definition: kEpsilon.H:116
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
virtual tmp< volSymmTensorField > sigma() const
Return the stress tensor [m^2/s^2].
BasicMomentumTransportModel::transportModel transportModel
Definition: kEpsilon.H:117
virtual tmp< volScalarField > nu() const
Kinematic viscosity of mixture [m^2/s].
Definition: fluidThermo.C:84
virtual void correctNut()
Definition: kEpsilon.C:41
virtual tmp< volScalarField > rhoEff() const
Return the effective density for the stress.
A class for managing temporary objects.
Definition: PtrList.H:53
const objectRegistry & db() const
Return the local objectRegistry.
Definition: IOobject.C:312
virtual tmp< volScalarField > epsilon() const =0
Return the turbulence kinetic energy dissipation rate.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phaseModel.H:53
dimensionedScalar deltaT() const
Return time step.
Definition: TimeStateI.H:53
Namespace for OpenFOAM.