Loading...
Searching...
No Matches
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-2023 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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
30#include "fvOptions.H"
31#include "twoPhaseSystem.H"
32#include "virtualMassModel.H"
33#include "dragModel.H"
34
35// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36
37namespace Foam
39namespace RASModels
40{
41
42// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
43
44template<class BasicTurbulenceModel>
45continuousGasKEpsilon<BasicTurbulenceModel>::continuousGasKEpsilon
46(
47 const alphaField& alpha,
48 const rhoField& rho,
49 const volVectorField& U,
50 const surfaceScalarField& alphaRhoPhi,
52 const transportModel& transport,
53 const word& propertiesName,
54 const word& type
55)
56:
57 kEpsilon<BasicTurbulenceModel>
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 (
74 (
75 IOobject::groupName("nutEff", alphaRhoPhi.group()),
76 this->runTime_.timeName(),
77 this->mesh_,
80 ),
81 this->nut_
82 ),
83
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
103template<class BasicTurbulenceModel>
105{
107 {
108 alphaInversion_.readIfPresent(this->coeffDict());
109
110 return true;
112
113 return false;
114}
115
116
117template<class BasicTurbulenceModel>
119{
121
122 const turbulenceModel& liquidTurbulence = this->liquidTurbulence();
123 const transportModel& gas = this->transport();
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
148template<class BasicTurbulenceModel>
149const turbulenceModel&
151{
152 if (!liquidTurbulencePtr_)
153 {
154 const volVectorField& U = this->U_;
155
156 const transportModel& gas = this->transport();
157 const twoPhaseSystem& fluid =
159 const transportModel& liquid = fluid.otherPhase(gas);
160
161 liquidTurbulencePtr_ =
162 &U.db().lookupObject<turbulenceModel>
163 (
165 (
167 liquid.name()
168 )
169 );
170 }
172 return *liquidTurbulencePtr_;
173}
174
175
176template<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
194 (
196 (
197 IOobject::groupName("nuEff", this->alphaRhoPhi_.group()),
198 blend*this->nut_
199 + (1.0 - blend)*rhoEff()*nutEff_/this->transport().rho()
200 + this->nu()
202 );
203}
204
205
206template<class BasicTurbulenceModel>
209{
210 const transportModel& gas = this->transport();
212 const transportModel& liquid = fluid.otherPhase(gas);
213
214 const virtualMassModel& virtualMass =
215 fluid.lookupSubModel<virtualMassModel>(gas, liquid);
216
218 (
220 (
221 IOobject::groupName("rhoEff", this->alphaRhoPhi_.group()),
222 gas.rho() + (virtualMass.Cvm() + 3.0/20.0)*liquid.rho()
224 );
225}
226
227
228template<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()
247 );
248}
249
250
251template<class BasicTurbulenceModel>
254{
255 const turbulenceModel& liquidTurbulence = this->liquidTurbulence();
256 const volScalarField phaseTransferCoeff(this->phaseTransferCoeff());
257
258 return
260 - fvm::Sp(phaseTransferCoeff, this->k_);
261}
262
263
264template<class BasicTurbulenceModel>
267{
268 const turbulenceModel& liquidTurbulence = this->liquidTurbulence();
269 const volScalarField phaseTransferCoeff(this->phaseTransferCoeff());
270
271 return
274}
275
276
277template<class BasicTurbulenceModel>
280{
281 tmp<volScalarField> tk(this->k());
282
284 (
286 (
288 (
289 IOobject::groupName("R", this->alphaRhoPhi_.group()),
290 this->runTime_.timeName(),
291 this->mesh_,
294 ),
295 ((2.0/3.0)*I)*tk() - (nutEff_)*devTwoSymm(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// ************************************************************************* //
label k
twoPhaseSystem & fluid
@ NO_READ
Nothing to be read.
@ READ_IF_PRESENT
Reading is optional [identical to LAZY_READ].
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
@ AUTO_WRITE
Automatically write from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
static word group(const word &name)
Return group (extension part of name).
Definition IOobject.C:300
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
BasicTurbulenceModel::alphaField alphaField
virtual tmp< fvScalarMatrix > epsilonSource() const
const turbulenceModel & liquidTurbulence() const
Return the turbulence model for the liquid phase.
BasicTurbulenceModel::rhoField rhoField
virtual tmp< volScalarField > nuEff() const
Return the effective viscosity.
tmp< volScalarField > phaseTransferCoeff() const
virtual tmp< volScalarField > rhoEff() const
Return the effective density for the stress.
virtual tmp< volSymmTensorField > R() const
Return the Reynolds stress tensor.
virtual tmp< fvScalarMatrix > kSource() const
BasicTurbulenceModel::transportModel transportModel
virtual bool read()
Re-read model coefficients if they have changed.
volScalarField epsilon_
Definition kEpsilon.H:137
virtual void correctNut()
Definition kEpsilon.C:37
virtual bool read()
Re-read model coefficients if they have changed.
Definition kEpsilon.C:294
Generic dimensioned Type class.
static dimensioned< Type > getOrAddToDict(const word &name, dictionary &dict, const dimensionSet &dims=dimless, const Type &deflt=Type(Zero))
Construct dimensioned from dictionary, with default value.
void correct(GeometricField< Type, PatchField, GeoMesh > &field)
Apply correction to field.
static options & New(const fvMesh &mesh)
Construct fvOptions and register to database if not present otherwise lookup and return.
Definition fvOptions.C:116
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:21
A class for managing temporary objects.
Definition tmp.H:75
Base-class for all transport models used by the incompressible turbulence models.
Abstract base class for turbulence models (RAS, LES and laminar).
static const word propertiesName
Default name of the turbulence properties dictionary.
virtual tmp< volScalarField > k() const =0
Return the turbulence kinetic energy.
virtual tmp< volScalarField > epsilon() const =0
Return the turbulence kinetic energy dissipation rate.
Class which solves the volume fraction equations for two phases.
virtual tmp< volScalarField > Cvm() const =0
Return the virtual mass coefficient.
A class for handling words, derived from Foam::string.
Definition word.H:66
U
Definition pEqn.H:72
word timeName
Definition getTimeIndex.H:3
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition fvcGrad.C:47
zeroField Sp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
Namespace for OpenFOAM.
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
dimensionedScalar exp(const dimensionedScalar &ds)
GeometricField< vector, fvPatchField, volMesh > volVectorField
dimensionedSymmTensor sqr(const dimensionedVector &dv)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
SymmTensor< Cmpt > devTwoSymm(const SymmTensor< Cmpt > &st)
Return the deviatoric part of twice the symmetric part of a SymmTensor.
static const Identity< scalar > I
Definition Identity.H:100
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
GeometricField< symmTensor, fvPatchField, volMesh > volSymmTensorField
volScalarField & nu
volScalarField & alpha