LaheyKEpsilon.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-2021 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 "fvModels.H"
28 #include "fvConstraints.H"
29 #include "phaseSystem.H"
30 #include "dragModel.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  gasTurbulencePtr_(nullptr),
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 BasicMomentumTransportModel>
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 BasicMomentumTransportModel>
136 <
137  typename BasicMomentumTransportModel::transportModel
138 >&
140 {
141  if (!gasTurbulencePtr_)
142  {
143  const volVectorField& U = this->U_;
144 
145  const phaseModel& liquid = refCast<const phaseModel>(this->transport());
146  const phaseSystem& fluid = liquid.fluid();
147  const phaseModel& gas = fluid.otherPhase(liquid);
148 
149  gasTurbulencePtr_ =
150  &U.db().lookupObject
151  <
153  >
154  (
156  (
157  momentumTransportModel::typeName,
158  gas.name()
159  )
160  );
161  }
162 
163  return *gasTurbulencePtr_;
164 }
165 
166 
167 template<class BasicMomentumTransportModel>
169 {
171  gasTurbulence = this->gasTurbulence();
172 
173  const phaseModel& liquid = refCast<const phaseModel>(this->transport());
174  const phaseSystem& fluid = liquid.fluid();
175  const phaseModel& gas = fluid.otherPhase(liquid);
176 
177  this->nut_ =
178  this->Cmu_*sqr(this->k_)/this->epsilon_
179  + Cmub_*gas.d()*gasTurbulence.alpha()
180  *(mag(this->U_ - gasTurbulence.U()));
181 
182  this->nut_.correctBoundaryConditions();
183  fvConstraints::New(this->mesh_).constrain(this->nut_);
184 }
185 
186 
187 template<class BasicMomentumTransportModel>
189 {
191  gasTurbulence = this->gasTurbulence();
192 
193  const phaseModel& liquid = refCast<const phaseModel>(this->transport());
194  const phaseSystem& fluid = liquid.fluid();
195  const phaseModel& gas = fluid.otherPhase(liquid);
196 
197  const dragModel& drag = fluid.lookupSubModel<dragModel>(gas, liquid);
198 
199  volScalarField magUr(mag(this->U_ - gasTurbulence.U()));
200 
201  tmp<volScalarField> bubbleG
202  (
203  Cp_
204  *(
205  pow3(magUr)
206  + pow(drag.CdRe()*liquid.thermo().nu()/gas.d(), 4.0/3.0)
207  *pow(magUr, 5.0/3.0)
208  )
209  *gas
210  /gas.d()
211  );
212 
213  return bubbleG;
214 }
215 
216 
217 template<class BasicMomentumTransportModel>
220 {
221  const volVectorField& U = this->U_;
222  const alphaField& alpha = this->alpha_;
223  const rhoField& rho = this->rho_;
224 
225  const momentumTransportModel& gasTurbulence = this->gasTurbulence();
226 
227  return
228  (
229  max(alphaInversion_ - alpha, scalar(0))
230  *rho
231  *min(gasTurbulence.epsilon()/gasTurbulence.k(), 1.0/U.time().deltaT())
232  );
233 }
234 
235 
236 template<class BasicMomentumTransportModel>
238 {
239  const alphaField& alpha = this->alpha_;
240  const rhoField& rho = this->rho_;
241 
243  gasTurbulence = this->gasTurbulence();
244 
245  const volScalarField phaseTransferCoeff(this->phaseTransferCoeff());
246 
247  return
248  alpha*rho*bubbleG()
249  + phaseTransferCoeff*gasTurbulence.k()
250  - fvm::Sp(phaseTransferCoeff, this->k_);
251 }
252 
253 
254 template<class BasicMomentumTransportModel>
257 {
258  const alphaField& alpha = this->alpha_;
259  const rhoField& rho = this->rho_;
260 
262  gasTurbulence = this->gasTurbulence();
263 
264  const volScalarField phaseTransferCoeff(this->phaseTransferCoeff());
265 
266  return
267  alpha*rho*this->C3_*this->epsilon_*bubbleG()/this->k_
268  + phaseTransferCoeff*gasTurbulence.epsilon()
269  - fvm::Sp(phaseTransferCoeff, this->epsilon_);
270 }
271 
272 
273 template<class BasicMomentumTransportModel>
275 {
277 }
278 
279 
280 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
281 
282 } // End namespace RASModels
283 } // End namespace Foam
284 
285 // ************************************************************************* //
Generic thermophysical properties class for a liquid in which the functions and coefficients for each...
Definition: liquid.H:50
Templated abstract base class for multiphase compressible turbulence models.
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
virtual bool read()
Read model coefficients if they have changed.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
const word & name() const
Definition: phaseModel.H:110
Generic dimensioned Type class.
tmp< volScalarField > d() const
Return the Sauter-mean diameter.
LaheyKEpsilon(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.
Definition: LaheyKEpsilon.C:43
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.
Continuous-phase k-epsilon model including bubble-generated turbulence.
Definition: LaheyKEpsilon.H:74
Info<< "Reading strained laminar flame speed field Su\"<< endl;volScalarField Su(IOobject("Su", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);Info<< "Reading field betav\"<< endl;volScalarField betav(IOobject("betav", mesh.facesInstance(), mesh, IOobject::MUST_READ, IOobject::NO_WRITE), mesh);Info<< "Reading field Lobs\"<< endl;volScalarField Lobs(IOobject("Lobs", mesh.facesInstance(), mesh, IOobject::MUST_READ, IOobject::NO_WRITE), mesh);Info<< "Reading field CT\"<< endl;volSymmTensorField CT(IOobject("CT", mesh.facesInstance(), mesh, IOobject::MUST_READ, IOobject::NO_WRITE), mesh);Info<< "Reading field Nv\"<< endl;volScalarField Nv(IOobject("Nv", mesh.facesInstance(), mesh, IOobject::MUST_READ, IOobject::NO_WRITE), mesh);Info<< "Reading field nsv\"<< endl;volSymmTensorField nsv(IOobject("nsv", mesh.facesInstance(), mesh, IOobject::MUST_READ, IOobject::NO_WRITE), mesh);IOdictionary PDRProperties(IOobject("PDRProperties", runTime.constant(), mesh, IOobject::MUST_READ_IF_MODIFIED, IOobject::NO_WRITE));autoPtr< PDRDragModel > drag
Definition: createFields.H:192
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:70
A class for handling words, derived from string.
Definition: word.H:59
const rhoThermo & thermo() const
Return const-access to phase rhoThermo.
Definition: phaseModel.H:121
static word groupName(Name name, const word &group)
BasicMomentumTransportModel::alphaField alphaField
Definition: kEpsilon.H:115
phaseSystem & fluid
Definition: createFields.H:11
phi
Definition: correctPhi.H:3
static autoPtr< dictionary > New(Istream &)
Construct top-level dictionary on freestore from Istream.
Definition: dictionaryIO.C:96
tmp< volScalarField > phaseTransferCoeff() const
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
virtual tmp< volScalarField > k() const =0
Return the turbulence kinetic energy.
Abstract base class for turbulence models (RAS, LES and laminar).
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
virtual tmp< fvScalarMatrix > kSource() const
dimensionedScalar pow3(const dimensionedScalar &ds)
const volVectorField & U() const
Access function to velocity field.
U
Definition: pEqn.H:72
const Time & time() const
Return time.
Definition: IOobject.C:318
const phaseSystem & fluid() const
Return the system to which this phase belongs.
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
BasicMomentumTransportModel::transportModel transportModel
Definition: kEpsilon.H:117
virtual tmp< volScalarField > nu() const
Kinematic viscosity of mixture [m^2/s].
Definition: fluidThermo.C:84
tmp< volScalarField > bubbleG() const
dimensioned< scalar > mag(const dimensioned< Type > &)
virtual tmp< fvScalarMatrix > epsilonSource() const
A class for managing temporary objects.
Definition: PtrList.H:53
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: kEpsilon.C:218
virtual tmp< volScalarField > CdRe() const =0
Drag coefficient.
const objectRegistry & db() const
Return the local objectRegistry.
Definition: IOobject.C:312
virtual tmp< volScalarField > epsilon() const =0
Return the turbulence kinetic energy dissipation rate.
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phaseModel.H:53
dimensionedScalar deltaT() const
Return time step.
Definition: TimeStateI.H:53
Namespace for OpenFOAM.