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-2020 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 "fvOptions.H"
28 #include "phaseSystem.H"
29 #include "dragModel.H"
30 #include "virtualMassModel.H"
31 
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace RASModels
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 transportModel& transport,
50  const word& type
51 )
52 :
54  (
55  alpha,
56  rho,
57  U,
58  alphaRhoPhi,
59  phi,
60  transport,
61  type
62  ),
63 
64  liquidTurbulencePtr_(nullptr),
65 
66  nutEff_
67  (
68  IOobject
69  (
70  IOobject::groupName("nutEff", alphaRhoPhi.group()),
71  this->runTime_.timeName(),
72  this->mesh_,
75  ),
76  this->nut_
77  ),
78 
79  alphaInversion_
80  (
82  (
83  "alphaInversion",
84  this->coeffDict_,
85  0.7
86  )
87  )
88 {
89  if (type == typeName)
90  {
91  this->printCoeffs(type);
92  }
93 }
94 
95 
96 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
97 
98 template<class BasicMomentumTransportModel>
100 {
102  {
103  alphaInversion_.readIfPresent(this->coeffDict());
104 
105  return true;
106  }
107  else
108  {
109  return false;
110  }
111 }
112 
113 
114 template<class BasicMomentumTransportModel>
116 {
118 
119  const momentumTransportModel& liquidTurbulence = this->liquidTurbulence();
120  const transportModel& gas = this->transport();
121  const phaseSystem& fluid = gas.fluid();
122  const transportModel& liquid = fluid.otherPhase(gas);
123 
124  const virtualMassModel& virtualMass =
125  fluid.lookupSubModel<virtualMassModel>(gas, liquid);
126 
127  volScalarField thetal(liquidTurbulence.k()/liquidTurbulence.epsilon());
128  volScalarField rhodv(gas.rho() + virtualMass.Cvm()*liquid.rho());
129  volScalarField thetag
130  (
131  (rhodv/(18*liquid.rho()*liquid.thermo().nu()))*sqr(gas.d())
132  );
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 BasicMomentumTransportModel>
151 {
152  if (!liquidTurbulencePtr_)
153  {
154  const volVectorField& U = this->U_;
155 
156  const transportModel& gas = this->transport();
157  const phaseSystem& fluid = gas.fluid();
158  const transportModel& liquid = fluid.otherPhase(gas);
159 
160  liquidTurbulencePtr_ =
162  (
164  (
165  momentumTransportModel::typeName,
166  liquid.name()
167  )
168  );
169  }
170 
171  return *liquidTurbulencePtr_;
172 }
173 
174 
175 template<class BasicMomentumTransportModel>
178 {
179  volScalarField blend
180  (
181  max
182  (
183  min
184  (
185  (this->alpha_ - scalar(0.5))/(alphaInversion_ - 0.5),
186  scalar(1)
187  ),
188  scalar(0)
189  )
190  );
191 
192  return volScalarField::New
193  (
194  IOobject::groupName("nuEff", this->alphaRhoPhi_.group()),
195  blend*this->nut_
196  + (1.0 - blend)*rhoEff()*nutEff_/this->transport().rho()
197  + this->nu()
198  );
199 }
200 
201 
202 template<class BasicMomentumTransportModel>
205 {
206  const transportModel& gas = this->transport();
207  const phaseSystem& fluid = gas.fluid();
208  const transportModel& liquid = fluid.otherPhase(gas);
209 
210  const virtualMassModel& virtualMass =
211  fluid.lookupSubModel<virtualMassModel>(gas, liquid);
212 
213  return volScalarField::New
214  (
215  IOobject::groupName("rhoEff", this->alphaRhoPhi_.group()),
216  gas.rho() + (virtualMass.Cvm() + 3.0/20.0)*liquid.rho()
217  );
218 }
219 
220 
221 template<class BasicMomentumTransportModel>
224 {
225  const volVectorField& U = this->U_;
226  const alphaField& alpha = this->alpha_;
227  const rhoField& rho = this->rho_;
228 
229  const momentumTransportModel& liquidTurbulence = this->liquidTurbulence();
230 
231  return
232  (
233  max(alphaInversion_ - alpha, scalar(0))
234  *rho
235  *min
236  (
237  liquidTurbulence.epsilon()/liquidTurbulence.k(),
238  1.0/U.time().deltaT()
239  )
240  );
241 }
242 
243 
244 template<class BasicMomentumTransportModel>
247 {
248  const momentumTransportModel& liquidTurbulence = this->liquidTurbulence();
249  const volScalarField phaseTransferCoeff(this->phaseTransferCoeff());
250 
251  return
252  phaseTransferCoeff*liquidTurbulence.k()
253  - fvm::Sp(phaseTransferCoeff, this->k_);
254 }
255 
256 
257 template<class BasicMomentumTransportModel>
260 {
261  const momentumTransportModel& liquidTurbulence = this->liquidTurbulence();
262  const volScalarField phaseTransferCoeff(this->phaseTransferCoeff());
263 
264  return
265  phaseTransferCoeff*liquidTurbulence.epsilon()
266  - fvm::Sp(phaseTransferCoeff, this->epsilon_);
267 }
268 
269 
270 template<class BasicMomentumTransportModel>
273 {
274  tmp<volScalarField> tk(this->k());
275 
277  (
278  IOobject::groupName("R", this->alphaRhoPhi_.group()),
279  ((2.0/3.0)*I)*tk() - (nutEff_)*dev(twoSymm(fvc::grad(this->U_))),
280  tk().boundaryField().types()
281  );
282 }
283 
284 
285 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
286 
287 } // End namespace RASModels
288 } // End namespace Foam
289 
290 // ************************************************************************* //
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
void correct(GeometricField< Type, fvPatchField, volMesh > &field)
Apply correction to field.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const modelType & lookupSubModel(const phasePair &key) const
Return a sub model between a phase pair.
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Internal &, const PtrList< fvPatchField< scalar >> &)
Return a temporary field constructed from name,.
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)
phi
Definition: pEqn.H:104
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:67
A class for handling words, derived from string.
Definition: word.H:59
static word groupName(Name name, const word &group)
const phaseModel & otherPhase(const phaseModel &phase) const
Return the phase not given as an argument in a two-phase system.
Definition: phaseSystemI.H:112
BasicMomentumTransportModel::alphaField alphaField
Definition: kEpsilon.H:115
phaseSystem & fluid
Definition: createFields.H:11
virtual tmp< volScalarField > nuEff() const
Return the effective viscosity.
virtual tmp< fvScalarMatrix > kSource() const
virtual tmp< volScalarField > nut() const =0
Return the turbulence viscosity.
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:328
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 void correctNut()
Definition: kEpsilon.C:40
virtual tmp< volScalarField > rhoEff() const
Return the effective density for the stress.
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
A class for managing temporary objects.
Definition: PtrList.H:53
static options & New(const fvMesh &mesh)
Construct fvOptions and register to datbase if not present.
Definition: fvOptions.C:101
const objectRegistry & db() const
Return the local objectRegistry.
Definition: IOobject.C:322
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
dimensionedScalar deltaT() const
Return time step.
Definition: TimeStateI.H:53
Namespace for OpenFOAM.