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-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 "LaheyKEpsilon.H"
28 #include "twoPhaseSystem.H"
29 #include "dragModel.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  gasTurbulencePtr_(NULL),
66 
67  alphaInversion_
68  (
70  (
71  "alphaInversion",
72  this->coeffDict_,
73  0.3
74  )
75  ),
76 
77  Cp_
78  (
80  (
81  "Cp",
82  this->coeffDict_,
83  0.25
84  )
85  ),
86 
87  C3_
88  (
90  (
91  "C3",
92  this->coeffDict_,
93  this->C2_.value()
94  )
95  ),
96 
97  Cmub_
98  (
100  (
101  "Cmub",
102  this->coeffDict_,
103  0.6
104  )
105  )
106 {
107  if (type == typeName)
108  {
109  this->printCoeffs(type);
110  }
111 }
112 
113 
114 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
115 
116 template<class BasicTurbulenceModel>
118 {
120  {
121  alphaInversion_.readIfPresent(this->coeffDict());
122  Cp_.readIfPresent(this->coeffDict());
123  C3_.readIfPresent(this->coeffDict());
124  Cmub_.readIfPresent(this->coeffDict());
125 
126  return true;
127  }
128  else
129  {
130  return false;
131  }
132 }
133 
134 
135 template<class BasicTurbulenceModel>
137 <
138  typename BasicTurbulenceModel::transportModel
139 >&
141 {
142  if (!gasTurbulencePtr_)
143  {
144  const volVectorField& U = this->U_;
145 
146  const transportModel& liquid = this->transport();
147  const twoPhaseSystem& fluid =
148  refCast<const twoPhaseSystem>(liquid.fluid());
149  const transportModel& gas = fluid.otherPhase(liquid);
150 
151  gasTurbulencePtr_ =
152  &U.db()
154  (
156  (
158  gas.name()
159  )
160  );
161  }
162 
163  return *gasTurbulencePtr_;
164 }
165 
166 
167 template<class BasicTurbulenceModel>
169 {
171  this->gasTurbulence();
172 
173  this->nut_ =
174  this->Cmu_*sqr(this->k_)/this->epsilon_
175  + Cmub_*gasTurbulence.transport().d()*gasTurbulence.alpha()
176  *(mag(this->U_ - gasTurbulence.U()));
177 
178  this->nut_.correctBoundaryConditions();
179 }
180 
181 
182 template<class BasicTurbulenceModel>
184 {
186  this->gasTurbulence();
187 
188  const transportModel& liquid = this->transport();
189  const twoPhaseSystem& fluid = refCast<const twoPhaseSystem>(liquid.fluid());
190  const transportModel& gas = fluid.otherPhase(liquid);
191 
192  volScalarField magUr(mag(this->U_ - gasTurbulence.U()));
193 
194  tmp<volScalarField> bubbleG
195  (
196  Cp_
197  *(
198  pow3(magUr)
199  + pow(fluid.drag(gas).CdRe()*liquid.nu()/gas.d(), 4.0/3.0)
200  *pow(magUr, 5.0/3.0)
201  )
202  *gas
203  /gas.d()
204  );
205 
206  return bubbleG;
207 }
208 
209 
210 template<class BasicTurbulenceModel>
213 {
214  const volVectorField& U = this->U_;
215  const alphaField& alpha = this->alpha_;
216  const rhoField& rho = this->rho_;
217 
218  const turbulenceModel& gasTurbulence = this->gasTurbulence();
219 
220  return
221  (
222  max(alphaInversion_ - alpha, scalar(0))
223  *rho
224  *min(gasTurbulence.epsilon()/gasTurbulence.k(), 1.0/U.time().deltaT())
225  );
226 }
227 
228 
229 template<class BasicTurbulenceModel>
231 {
232  const alphaField& alpha = this->alpha_;
233  const rhoField& rho = this->rho_;
234 
236  this->gasTurbulence();
237 
238  const volScalarField phaseTransferCoeff(this->phaseTransferCoeff());
239 
240  return
241  alpha*rho*bubbleG()
242  + phaseTransferCoeff*gasTurbulence.k()
243  - fvm::Sp(phaseTransferCoeff, this->k_);
244 }
245 
246 
247 template<class BasicTurbulenceModel>
249 {
250  const alphaField& alpha = this->alpha_;
251  const rhoField& rho = this->rho_;
252 
254  this->gasTurbulence();
255 
256  const volScalarField phaseTransferCoeff(this->phaseTransferCoeff());
257 
258  return
259  alpha*rho*this->C3_*this->epsilon_*bubbleG()/this->k_
260  + phaseTransferCoeff*gasTurbulence.epsilon()
261  - fvm::Sp(phaseTransferCoeff, this->epsilon_);
262 }
263 
264 
265 template<class BasicTurbulenceModel>
267 {
269 }
270 
271 
272 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
273 
274 } // End namespace RASModels
275 } // End namespace Foam
276 
277 // ************************************************************************* //
Foam::RASModels::kEpsilon::correct
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: kEpsilon.C:223
Foam::IOobject::groupName
static word groupName(Name name, const word &group)
LaheyKEpsilon.H
Foam::RASModels::LaheyKEpsilon::transportModel
BasicTurbulenceModel::transportModel transportModel
Definition: LaheyKEpsilon.H:128
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::RASModels::LaheyKEpsilon::phaseTransferCoeff
tmp< volScalarField > phaseTransferCoeff() const
Definition: LaheyKEpsilon.C:212
Foam::PhaseCompressibleTurbulenceModel
Templated abstract base class for multiphase compressible turbulence models.
Definition: PhaseCompressibleTurbulenceModel.H:51
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
fluid
multiphaseSystem & fluid
Definition: createFields.H:10
Foam::RASModels::LaheyKEpsilon::gasTurbulence
const PhaseCompressibleTurbulenceModel< typename BasicTurbulenceModel::transportModel > & gasTurbulence() const
Return the turbulence model for the gas phase.
Definition: LaheyKEpsilon.C:140
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::turbulenceModel::propertiesName
static const word propertiesName
Default name of the turbulence properties dictionary.
Definition: turbulenceModel.H:97
Foam::RASModels::LaheyKEpsilon::bubbleG
tmp< volScalarField > bubbleG() const
Definition: LaheyKEpsilon.C:183
Foam::TurbulenceModel< volScalarField, volScalarField, compressibleTurbulenceModel, TransportModel >::alpha
const alphaField & alpha() const
Access function to phase fraction.
Definition: TurbulenceModel.H:148
Foam::RASModels::LaheyKEpsilon::kSource
virtual tmp< fvScalarMatrix > kSource() const
Definition: LaheyKEpsilon.C:230
Foam::RASModels::LaheyKEpsilon::correctNut
virtual void correctNut()
Definition: LaheyKEpsilon.C:168
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
Foam::RASModels::LaheyKEpsilon::epsilonSource
virtual tmp< fvScalarMatrix > epsilonSource() const
Definition: LaheyKEpsilon.C:248
U
U
Definition: pEqn.H:46
Foam::fvm::Sp
tmp< fvMatrix< Type > > Sp(const DimensionedField< scalar, volMesh > &, const GeometricField< Type, fvPatchField, volMesh > &)
Foam::pow3
dimensionedScalar pow3(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:87
Foam::RASModels::LaheyKEpsilon::rhoField
BasicTurbulenceModel::rhoField rhoField
Definition: LaheyKEpsilon.H:127
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:73
Foam::turbulenceModel
Abstract base class for turbulence models (RAS, LES and laminar).
Definition: turbulenceModel.H:60
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
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::TurbulenceModel< volScalarField, volScalarField, compressibleTurbulenceModel, TransportModel >::transport
const transportModel & transport() const
Access function to incompressible transport model.
Definition: TurbulenceModel.H:154
rho
rho
Definition: pEqn.H:3
Foam::RASModels::LaheyKEpsilon::correct
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: LaheyKEpsilon.C:266
Foam::max
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::RASModels::LaheyKEpsilon::read
virtual bool read()
Read model coefficients if they have changed.
Definition: LaheyKEpsilon.C:117
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:49
Foam::RASModels::LaheyKEpsilon::LaheyKEpsilon
LaheyKEpsilon(const LaheyKEpsilon &)
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::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::turbulenceModel::U
const volVectorField & U() const
Access function to velocity field.
Definition: turbulenceModel.H:142
Foam::RASModels::LaheyKEpsilon::alphaField
BasicTurbulenceModel::alphaField alphaField
Definition: LaheyKEpsilon.H:126