LaheyKEpsilon.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-2016 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 "LaheyKEpsilon.H"
27 #include "fvOptions.H"
28 #include "twoPhaseSystem.H"
29 
30 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace RASModels
35 {
36 
37 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
38 
39 template<class BasicTurbulenceModel>
40 LaheyKEpsilon<BasicTurbulenceModel>::LaheyKEpsilon
41 (
42  const alphaField& alpha,
43  const rhoField& rho,
44  const volVectorField& U,
45  const surfaceScalarField& alphaRhoPhi,
46  const surfaceScalarField& phi,
47  const transportModel& transport,
48  const word& propertiesName,
49  const word& type
50 )
51 :
53  (
54  alpha,
55  rho,
56  U,
57  alphaRhoPhi,
58  phi,
59  transport,
60  propertiesName,
61  type
62  ),
63 
64  gasTurbulencePtr_(NULL),
65 
66  alphaInversion_
67  (
69  (
70  "alphaInversion",
71  this->coeffDict_,
72  0.3
73  )
74  ),
75 
76  Cp_
77  (
79  (
80  "Cp",
81  this->coeffDict_,
82  0.25
83  )
84  ),
85 
86  C3_
87  (
89  (
90  "C3",
91  this->coeffDict_,
92  this->C2_.value()
93  )
94  ),
95 
96  Cmub_
97  (
99  (
100  "Cmub",
101  this->coeffDict_,
102  0.6
103  )
104  )
105 {
106  if (type == typeName)
107  {
108  this->printCoeffs(type);
109  }
110 }
111 
112 
113 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
114 
115 template<class BasicTurbulenceModel>
117 {
119  {
120  alphaInversion_.readIfPresent(this->coeffDict());
121  Cp_.readIfPresent(this->coeffDict());
122  C3_.readIfPresent(this->coeffDict());
123  Cmub_.readIfPresent(this->coeffDict());
124 
125  return true;
126  }
127  else
128  {
129  return false;
130  }
131 }
132 
133 
134 template<class BasicTurbulenceModel>
136 <
137  typename BasicTurbulenceModel::transportModel
138 >&
140 {
141  if (!gasTurbulencePtr_)
142  {
143  const volVectorField& U = this->U_;
144 
145  const transportModel& liquid = this->transport();
146  const twoPhaseSystem& fluid =
147  refCast<const twoPhaseSystem>(liquid.fluid());
148  const transportModel& gas = fluid.otherPhase(liquid);
149 
150  gasTurbulencePtr_ =
151  &U.db()
153  (
155  (
157  gas.name()
158  )
159  );
160  }
161 
162  return *gasTurbulencePtr_;
163 }
164 
165 
166 template<class BasicTurbulenceModel>
168 {
170  this->gasTurbulence();
171 
172  this->nut_ =
173  this->Cmu_*sqr(this->k_)/this->epsilon_
174  + Cmub_*gasTurbulence.transport().d()*gasTurbulence.alpha()
175  *(mag(this->U_ - gasTurbulence.U()));
176 
177  this->nut_.correctBoundaryConditions();
178  fv::options::New(this->mesh_).correct(this->nut_);
179 
180  BasicTurbulenceModel::correctNut();
181 }
182 
183 
184 template<class BasicTurbulenceModel>
186 {
188  this->gasTurbulence();
189 
190  const transportModel& liquid = this->transport();
191  const twoPhaseSystem& fluid = refCast<const twoPhaseSystem>(liquid.fluid());
192  const transportModel& gas = fluid.otherPhase(liquid);
193 
194  volScalarField magUr(mag(this->U_ - gasTurbulence.U()));
195 
196  tmp<volScalarField> bubbleG
197  (
198  Cp_
199  *(
200  pow3(magUr)
201  + pow(fluid.drag(gas).CdRe()*liquid.nu()/gas.d(), 4.0/3.0)
202  *pow(magUr, 5.0/3.0)
203  )
204  *gas
205  /gas.d()
206  );
207 
208  return bubbleG;
209 }
210 
211 
212 template<class BasicTurbulenceModel>
215 {
216  const volVectorField& U = this->U_;
217  const alphaField& alpha = this->alpha_;
218  const rhoField& rho = this->rho_;
219 
220  const turbulenceModel& gasTurbulence = this->gasTurbulence();
221 
222  return
223  (
224  max(alphaInversion_ - alpha, scalar(0))
225  *rho
226  *min(gasTurbulence.epsilon()/gasTurbulence.k(), 1.0/U.time().deltaT())
227  );
228 }
229 
230 
231 template<class BasicTurbulenceModel>
233 {
234  const alphaField& alpha = this->alpha_;
235  const rhoField& rho = this->rho_;
236 
238  this->gasTurbulence();
239 
240  const volScalarField phaseTransferCoeff(this->phaseTransferCoeff());
241 
242  return
243  alpha*rho*bubbleG()
244  + phaseTransferCoeff*gasTurbulence.k()
245  - fvm::Sp(phaseTransferCoeff, this->k_);
246 }
247 
248 
249 template<class BasicTurbulenceModel>
251 {
252  const alphaField& alpha = this->alpha_;
253  const rhoField& rho = this->rho_;
254 
256  this->gasTurbulence();
257 
258  const volScalarField phaseTransferCoeff(this->phaseTransferCoeff());
259 
260  return
261  alpha*rho*this->C3_*this->epsilon_*bubbleG()/this->k_
262  + phaseTransferCoeff*gasTurbulence.epsilon()
263  - fvm::Sp(phaseTransferCoeff, this->epsilon_);
264 }
265 
266 
267 template<class BasicTurbulenceModel>
269 {
271 }
272 
273 
274 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
275 
276 } // End namespace RASModels
277 } // End namespace Foam
278 
279 // ************************************************************************* //
surfaceScalarField & phi
const transportModel & transport() const
Access function to incompressible transport model.
U
Definition: pEqn.H:83
BasicTurbulenceModel::rhoField rhoField
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 > &)
virtual bool read()
Read model coefficients if they have changed.
void correct(GeometricField< Type, fvPatchField, volMesh > &field)
Apply correction to field.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
const volVectorField & U() const
Access function to velocity field.
Generic dimensioned Type class.
tmp< volScalarField > bubbleG() const
Continuous-phase k-epsilon model including bubble-generated turbulence.
Definition: LaheyKEpsilon.H:77
Abstract base class for turbulence models (RAS, LES and laminar).
Class which solves the volume fraction equations for two phases.
virtual tmp< volScalarField > k() const =0
Return the turbulence kinetic energy.
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
BasicTurbulenceModel::alphaField alphaField
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< fvScalarMatrix > epsilonSource() const
BasicTurbulenceModel::transportModel transportModel
tmp< volScalarField > phaseTransferCoeff() const
const phaseModel & otherPhase(const phaseModel &phase) const
Constant access the phase not given as an argument.
Templated abstract base class for multiphase compressible turbulence models.
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedScalar pow3(const dimensionedScalar &ds)
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.
dimensioned< scalar > mag(const dimensioned< Type > &)
const dragModel & drag(const phaseModel &phase) const
Return the drag model for the given phase.
A class for managing temporary objects.
Definition: PtrList.H:54
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: kEpsilon.C:221
virtual tmp< volScalarField > CdRe() const =0
Drag coefficient.
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: [].
virtual tmp< fvScalarMatrix > kSource() const
Namespace for OpenFOAM.