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 | Copyright (C) 2013-2015 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"
28 #include "twoPhaseSystem.H"
29 #include "virtualMassModel.H"
30 
31 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace RASModels
36 {
37 
38 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
39 
40 template<class BasicTurbulenceModel>
42 (
43  const alphaField& alpha,
44  const rhoField& rho,
45  const volVectorField& U,
46  const surfaceScalarField& alphaRhoPhi,
47  const surfaceScalarField& phi,
48  const transportModel& transport,
49  const word& propertiesName,
50  const word& type
51 )
52 :
54  (
55  alpha,
56  rho,
57  U,
58  alphaRhoPhi,
59  phi,
60  transport,
61  propertiesName,
62  type
63  ),
64 
65  liquidTurbulencePtr_(NULL),
66 
67  nutEff_
68  (
69  IOobject
70  (
71  IOobject::groupName("nutEff", U.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 BasicTurbulenceModel>
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 BasicTurbulenceModel>
117 {
119 
120  const turbulenceModel& liquidTurbulence = this->liquidTurbulence();
121  const transportModel& gas = this->transport();
122  const twoPhaseSystem& fluid = refCast<const twoPhaseSystem>(gas.fluid());
123  const transportModel& liquid = fluid.otherPhase(gas);
124 
125  volScalarField thetal(liquidTurbulence.k()/liquidTurbulence.epsilon());
126  volScalarField rhodv(gas.rho() + fluid.virtualMass(gas).Cvm()*liquid.rho());
127  volScalarField thetag((rhodv/(18*liquid.rho()*liquid.nu()))*sqr(gas.d()));
128  volScalarField expThetar
129  (
130  min
131  (
132  exp(min(thetal/thetag, scalar(50))),
133  scalar(1)
134  )
135  );
136  volScalarField omega((1 - expThetar)/(1 + expThetar));
137 
138  nutEff_ = omega*liquidTurbulence.nut();
139 }
140 
141 
142 template<class BasicTurbulenceModel>
143 const turbulenceModel&
145 {
146  if (!liquidTurbulencePtr_)
147  {
148  const volVectorField& U = this->U_;
149 
150  const transportModel& gas = this->transport();
151  const twoPhaseSystem& fluid =
152  refCast<const twoPhaseSystem>(gas.fluid());
153  const transportModel& liquid = fluid.otherPhase(gas);
154 
155  liquidTurbulencePtr_ =
156  &U.db().lookupObject<turbulenceModel>
157  (
159  (
161  liquid.name()
162  )
163  );
164  }
165 
166  return *liquidTurbulencePtr_;
167 }
168 
169 
170 template<class BasicTurbulenceModel>
173 {
174  volScalarField blend
175  (
176  max
177  (
178  min
179  (
180  (this->alpha_ - scalar(0.5))/(alphaInversion_ - 0.5),
181  scalar(1)
182  ),
183  scalar(0)
184  )
185  );
186 
187  return tmp<volScalarField>
188  (
189  new volScalarField
190  (
191  IOobject::groupName("nuEff", this->U_.group()),
192  blend*this->nut_
193  + (1.0 - blend)*rhoEff()*nutEff_/this->transport().rho()
194  + this->nu()
195  )
196  );
197 }
198 
199 
200 template<class BasicTurbulenceModel>
203 {
204  const transportModel& gas = this->transport();
205  const twoPhaseSystem& fluid = refCast<const twoPhaseSystem>(gas.fluid());
206  const transportModel& liquid = fluid.otherPhase(gas);
207 
208  return tmp<volScalarField>
209  (
210  new volScalarField
211  (
212  IOobject::groupName("rhoEff", this->U_.group()),
213  gas.rho() + (fluid.virtualMass(gas).Cvm() + 3.0/20.0)*liquid.rho()
214  )
215  );
216 }
217 
218 
219 template<class BasicTurbulenceModel>
222 {
223  const volVectorField& U = this->U_;
224  const alphaField& alpha = this->alpha_;
225  const rhoField& rho = this->rho_;
226 
227  const turbulenceModel& liquidTurbulence = this->liquidTurbulence();
228 
229  return
230  (
231  max(alphaInversion_ - alpha, scalar(0))
232  *rho
233  *min
234  (
235  liquidTurbulence.epsilon()/liquidTurbulence.k(),
236  1.0/U.time().deltaT()
237  )
238  );
239 }
240 
241 
242 template<class BasicTurbulenceModel>
245 {
246  const turbulenceModel& liquidTurbulence = this->liquidTurbulence();
247  const volScalarField phaseTransferCoeff(this->phaseTransferCoeff());
248 
249  return
250  phaseTransferCoeff*liquidTurbulence.k()
251  - fvm::Sp(phaseTransferCoeff, this->k_);
252 }
253 
254 
255 template<class BasicTurbulenceModel>
258 {
259  const turbulenceModel& liquidTurbulence = this->liquidTurbulence();
260  const volScalarField phaseTransferCoeff(this->phaseTransferCoeff());
261 
262  return
263  phaseTransferCoeff*liquidTurbulence.epsilon()
264  - fvm::Sp(phaseTransferCoeff, this->epsilon_);
265 }
266 
267 
268 template<class BasicTurbulenceModel>
271 {
272  tmp<volScalarField> tk(this->k());
273 
275  (
277  (
278  IOobject
279  (
280  IOobject::groupName("R", this->U_.group()),
281  this->runTime_.timeName(),
282  this->mesh_,
285  ),
286  ((2.0/3.0)*I)*tk() - (nutEff_)*dev(twoSymm(fvc::grad(this->U_))),
287  tk().boundaryField().types()
288  )
289  );
290 }
291 
292 
293 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
294 
295 } // End namespace RASModels
296 } // End namespace Foam
297 
298 // ************************************************************************* //
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Foam::IOobject::groupName
static word groupName(Name name, const word &group)
Foam::fvc::grad
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::twoPhaseSystem
Class which solves the volume fraction equations for two phases.
Definition: twoPhaseSystem.H:51
Foam::IOobject::AUTO_WRITE
@ AUTO_WRITE
Definition: IOobject.H:117
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::turbulenceModel::nut
virtual tmp< volScalarField > nut() const =0
Return the turbulence viscosity.
fluid
multiphaseSystem & fluid
Definition: createFields.H:10
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:216
Foam::RASModels::continuousGasKEpsilon::rhoEff
virtual tmp< volScalarField > rhoEff() const
Return the effective density for the stress.
Definition: continuousGasKEpsilon.C:202
Foam::turbulenceModel::propertiesName
static const word propertiesName
Default name of the turbulence properties dictionary.
Definition: turbulenceModel.H:97
Foam::RASModels::continuousGasKEpsilon::transportModel
BasicTurbulenceModel::transportModel transportModel
Definition: continuousGasKEpsilon.H:117
Foam::exp
dimensionedScalar exp(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:252
Foam::IOobject::NO_WRITE
@ NO_WRITE
Definition: IOobject.H:118
U
U
Definition: pEqn.H:46
Foam::fvm::Sp
tmp< fvMatrix< Type > > Sp(const DimensionedField< scalar, volMesh > &, const GeometricField< Type, fvPatchField, volMesh > &)
Foam::RASModels::continuousGasKEpsilon::continuousGasKEpsilon
continuousGasKEpsilon(const continuousGasKEpsilon &)
nu
volScalarField & nu
Definition: readMechanicalProperties.H:179
Foam::RASModels::continuousGasKEpsilon::rhoField
BasicTurbulenceModel::rhoField rhoField
Definition: continuousGasKEpsilon.H:116
Foam::IOobject::NO_READ
@ NO_READ
Definition: IOobject.H:111
Foam::RASModels::continuousGasKEpsilon::alphaField
BasicTurbulenceModel::alphaField alphaField
Definition: continuousGasKEpsilon.H:115
Foam::RASModels::continuousGasKEpsilon::liquidTurbulence
const turbulenceModel & liquidTurbulence() const
Return the turbulence model for the liquid phase.
Definition: continuousGasKEpsilon.C:144
Foam::RASModels::continuousGasKEpsilon::read
virtual bool read()
Re-read model coefficients if they have changed.
Definition: continuousGasKEpsilon.C:100
Foam::RASModels::continuousGasKEpsilon::correctNut
virtual void correctNut()
Definition: continuousGasKEpsilon.C:116
Foam::I
static const sphericalTensor I(1)
Foam::turbulenceModel
Abstract base class for turbulence models (RAS, LES and laminar).
Definition: turbulenceModel.H:60
Foam::RASModels::continuousGasKEpsilon::nuEff
virtual tmp< volScalarField > nuEff() const
Return the effective viscosity.
Definition: continuousGasKEpsilon.C:172
Foam::RASModels::kEpsilon::correctNut
virtual void correctNut()
Definition: kEpsilon.C:40
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
boundaryField
cellIbMask *cellIbMaskExt *faceIbMask *cellIbMask boundaryField().evaluateCoupled()
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:41
Foam::turbulenceModel::k
virtual tmp< volScalarField > k() const =0
Return the turbulence kinetic energy.
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::RASModels::continuousGasKEpsilon::phaseTransferCoeff
tmp< volScalarField > phaseTransferCoeff() const
Definition: continuousGasKEpsilon.C:221
continuousGasKEpsilon.H
rho
rho
Definition: pEqn.H:3
Foam::max
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:49
Foam::RASModels::continuousGasKEpsilon::epsilonSource
virtual tmp< fvScalarMatrix > epsilonSource() const
Definition: continuousGasKEpsilon.C:257
Foam::RASModels::continuousGasKEpsilon::R
virtual tmp< volSymmTensorField > R() const
Return the Reynolds stress tensor.
Definition: continuousGasKEpsilon.C:270
k
label k
Boltzmann constant.
Definition: LISASMDCalcMethod2.H:41
Foam::type
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:588
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::RASModels::kEpsilon
Standard k-epsilon turbulence model for incompressible and compressible flows including rapid distort...
Definition: kEpsilon.H:86
Foam::IOobject::READ_IF_PRESENT
@ READ_IF_PRESENT
Definition: IOobject.H:110
Foam::RASModels::continuousGasKEpsilon::kSource
virtual tmp< fvScalarMatrix > kSource() const
Definition: continuousGasKEpsilon.C:244
Foam::turbulenceModel::epsilon
virtual tmp< volScalarField > epsilon() const =0
Return the turbulence kinetic energy dissipation rate.
Foam::min
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::twoSymm
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:93
Foam::dev
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:104