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"
27 #include "fvOptions.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>
41 continuousGasKEpsilon<BasicTurbulenceModel>::continuousGasKEpsilon
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  fv::options::New(this->mesh_).correct(nutEff_);
140 }
141 
142 
143 template<class BasicTurbulenceModel>
144 const turbulenceModel&
146 {
147  if (!liquidTurbulencePtr_)
148  {
149  const volVectorField& U = this->U_;
150 
151  const transportModel& gas = this->transport();
152  const twoPhaseSystem& fluid =
153  refCast<const twoPhaseSystem>(gas.fluid());
154  const transportModel& liquid = fluid.otherPhase(gas);
155 
156  liquidTurbulencePtr_ =
158  (
160  (
162  liquid.name()
163  )
164  );
165  }
166 
167  return *liquidTurbulencePtr_;
168 }
169 
170 
171 template<class BasicTurbulenceModel>
174 {
175  volScalarField blend
176  (
177  max
178  (
179  min
180  (
181  (this->alpha_ - scalar(0.5))/(alphaInversion_ - 0.5),
182  scalar(1)
183  ),
184  scalar(0)
185  )
186  );
187 
188  return tmp<volScalarField>
189  (
190  new volScalarField
191  (
192  IOobject::groupName("nuEff", this->U_.group()),
193  blend*this->nut_
194  + (1.0 - blend)*rhoEff()*nutEff_/this->transport().rho()
195  + this->nu()
196  )
197  );
198 }
199 
200 
201 template<class BasicTurbulenceModel>
204 {
205  const transportModel& gas = this->transport();
206  const twoPhaseSystem& fluid = refCast<const twoPhaseSystem>(gas.fluid());
207  const transportModel& liquid = fluid.otherPhase(gas);
208 
209  return tmp<volScalarField>
210  (
211  new volScalarField
212  (
213  IOobject::groupName("rhoEff", this->U_.group()),
214  gas.rho() + (fluid.virtualMass(gas).Cvm() + 3.0/20.0)*liquid.rho()
215  )
216  );
217 }
218 
219 
220 template<class BasicTurbulenceModel>
223 {
224  const volVectorField& U = this->U_;
225  const alphaField& alpha = this->alpha_;
226  const rhoField& rho = this->rho_;
227 
228  const turbulenceModel& liquidTurbulence = this->liquidTurbulence();
229 
230  return
231  (
232  max(alphaInversion_ - alpha, scalar(0))
233  *rho
234  *min
235  (
236  liquidTurbulence.epsilon()/liquidTurbulence.k(),
237  1.0/U.time().deltaT()
238  )
239  );
240 }
241 
242 
243 template<class BasicTurbulenceModel>
246 {
247  const turbulenceModel& liquidTurbulence = this->liquidTurbulence();
248  const volScalarField phaseTransferCoeff(this->phaseTransferCoeff());
249 
250  return
251  phaseTransferCoeff*liquidTurbulence.k()
252  - fvm::Sp(phaseTransferCoeff, this->k_);
253 }
254 
255 
256 template<class BasicTurbulenceModel>
259 {
260  const turbulenceModel& liquidTurbulence = this->liquidTurbulence();
261  const volScalarField phaseTransferCoeff(this->phaseTransferCoeff());
262 
263  return
264  phaseTransferCoeff*liquidTurbulence.epsilon()
265  - fvm::Sp(phaseTransferCoeff, this->epsilon_);
266 }
267 
268 
269 template<class BasicTurbulenceModel>
272 {
273  tmp<volScalarField> tk(this->k());
274 
276  (
278  (
279  IOobject
280  (
281  IOobject::groupName("R", this->U_.group()),
282  this->runTime_.timeName(),
283  this->mesh_,
286  ),
287  ((2.0/3.0)*I)*tk() - (nutEff_)*dev(twoSymm(fvc::grad(this->U_))),
288  tk().boundaryField().types()
289  )
290  );
291 }
292 
293 
294 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
295 
296 } // End namespace RASModels
297 } // End namespace Foam
298 
299 // ************************************************************************* //
BasicTurbulenceModel::transportModel transportModel
surfaceScalarField & phi
GeometricField< symmTensor, fvPatchField, volMesh > volSymmTensorField
Definition: volFieldsFwd.H:58
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
U
Definition: pEqn.H:83
const virtualMassModel & virtualMass(const phaseModel &phase) const
Return the virtual mass model for the given phase.
multiphaseSystem & fluid
Definition: createFields.H:10
const objectRegistry & db() const
Return the local objectRegistry.
Definition: IOobject.C:221
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
tmp< fvMatrix< Type > > Sp(const DimensionedField< scalar, volMesh > &, const GeometricField< Type, fvPatchField, volMesh > &)
void correct(GeometricField< Type, fvPatchField, volMesh > &field)
Apply correction to field.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
virtual bool read()
Re-read model coefficients if they have changed.
label k
Boltzmann constant.
Generic dimensioned Type class.
tmp< volScalarField > phaseTransferCoeff() const
virtual tmp< fvScalarMatrix > epsilonSource() const
Abstract base class for turbulence models (RAS, LES and laminar).
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
Class which solves the volume fraction equations for two phases.
virtual tmp< volScalarField > k() const =0
Return the turbulence kinetic energy.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
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:86
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from string.
Definition: word.H:59
const Time & time() const
Return time.
Definition: IOobject.C:227
static word groupName(Name name, const word &group)
virtual tmp< volScalarField > rhoEff() const
Return the effective density for the stress.
virtual tmp< volScalarField > nuEff() const
Return the effective viscosity.
const phaseModel & otherPhase(const phaseModel &phase) const
Constant access the phase not given as an argument.
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
virtual tmp< volScalarField > Cvm() const =0
Return the virtual mass coefficient.
BasicTurbulenceModel::alphaField alphaField
dimensionedScalar deltaT() const
Return time step.
Definition: TimeStateI.H:53
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:461
virtual tmp< volScalarField > epsilon() const =0
Return the turbulence kinetic energy dissipation rate.
const turbulenceModel & liquidTurbulence() const
Return the turbulence model for the liquid phase.
virtual tmp< fvScalarMatrix > kSource() const
virtual void correctNut()
Definition: kEpsilon.C:40
BasicTurbulenceModel::rhoField rhoField
virtual tmp< volScalarField > nut() const =0
Return the turbulence viscosity.
A class for managing temporary objects.
Definition: PtrList.H:54
static options & New(const fvMesh &mesh)
Construct fvOptions and register to datbase if not present.
Definition: fvOptions.C:103
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
volScalarField & nu
Namespace for OpenFOAM.
virtual tmp< volSymmTensorField > R() const
Return the Reynolds stress tensor.