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 | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2013-2016 OpenFOAM Foundation
9 Copyright (C) 2019-2020 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "LaheyKEpsilon.H"
30#include "fvOptions.H"
31#include "twoPhaseSystem.H"
32#include "dragModel.H"
33
34// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35
36namespace Foam
37{
38namespace RASModels
39{
40
41// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
42
43template<class BasicTurbulenceModel>
45(
46 const alphaField& alpha,
47 const rhoField& rho,
48 const volVectorField& U,
49 const surfaceScalarField& alphaRhoPhi,
51 const transportModel& transport,
52 const word& propertiesName,
53 const word& type
54)
55:
56 kEpsilon<BasicTurbulenceModel>
57 (
58 alpha,
59 rho,
60 U,
61 alphaRhoPhi,
62 phi,
63 transport,
64 propertiesName,
65 type
66 ),
67
68 gasTurbulencePtr_(nullptr),
69
70 alphaInversion_
71 (
72 dimensioned<scalar>::getOrAddToDict
73 (
74 "alphaInversion",
75 this->coeffDict_,
76 0.3
77 )
78 ),
79
80 Cp_
81 (
82 dimensioned<scalar>::getOrAddToDict
83 (
84 "Cp",
85 this->coeffDict_,
86 0.25
87 )
88 ),
89
90 C3_
91 (
92 dimensioned<scalar>::getOrAddToDict
93 (
94 "C3",
95 this->coeffDict_,
96 this->C2_.value()
97 )
98 ),
99
100 Cmub_
101 (
102 dimensioned<scalar>::getOrAddToDict
103 (
104 "Cmub",
105 this->coeffDict_,
106 0.6
107 )
108 )
109{
110 if (type == typeName)
111 {
112 this->printCoeffs(type);
113 }
114}
115
116
117// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
118
119template<class BasicTurbulenceModel>
121{
123 {
124 alphaInversion_.readIfPresent(this->coeffDict());
125 Cp_.readIfPresent(this->coeffDict());
126 C3_.readIfPresent(this->coeffDict());
127 Cmub_.readIfPresent(this->coeffDict());
128
129 return true;
130 }
131
132 return false;
133}
134
135
136template<class BasicTurbulenceModel>
138<
139 typename BasicTurbulenceModel::transportModel
140>&
142{
143 if (!gasTurbulencePtr_)
144 {
145 const volVectorField& U = this->U_;
146
147 const transportModel& liquid = this->transport();
148 const twoPhaseSystem& fluid =
149 refCast<const twoPhaseSystem>(liquid.fluid());
150 const transportModel& gas = fluid.otherPhase(liquid);
151
152 gasTurbulencePtr_ =
153 &U.db()
155 (
157 (
159 gas.name()
160 )
161 );
162 }
163
164 return *gasTurbulencePtr_;
165}
166
167
168template<class BasicTurbulenceModel>
170{
172 this->gasTurbulence();
173
174 this->nut_ =
175 this->Cmu_*sqr(this->k_)/this->epsilon_
176 + Cmub_*gasTurbulence.transport().d()*gasTurbulence.alpha()
177 *(mag(this->U_ - gasTurbulence.U()));
178
179 this->nut_.correctBoundaryConditions();
180 fv::options::New(this->mesh_).correct(this->nut_);
181
182 BasicTurbulenceModel::correctNut();
183}
184
185
186template<class BasicTurbulenceModel>
188{
190 this->gasTurbulence();
191
192 const transportModel& liquid = this->transport();
193 const twoPhaseSystem& fluid = refCast<const twoPhaseSystem>(liquid.fluid());
194 const transportModel& gas = fluid.otherPhase(liquid);
195
196 const dragModel& drag = fluid.lookupSubModel<dragModel>(gas, liquid);
197
198 volScalarField magUr(mag(this->U_ - gasTurbulence.U()));
199
200 tmp<volScalarField> bubbleG
201 (
202 Cp_
203 *(
204 pow3(magUr)
205 + pow(drag.CdRe()*liquid.nu()/gas.d(), 4.0/3.0)
206 *pow(magUr, 5.0/3.0)
207 )
208 *gas
209 /gas.d()
210 );
211
212 return bubbleG;
213}
214
215
216template<class BasicTurbulenceModel>
219{
220 const volVectorField& U = this->U_;
221 const alphaField& alpha = this->alpha_;
222 const rhoField& rho = this->rho_;
223
224 const turbulenceModel& gasTurbulence = this->gasTurbulence();
225
226 return
227 (
228 max(alphaInversion_ - alpha, scalar(0))
229 *rho
230 *min(gasTurbulence.epsilon()/gasTurbulence.k(), 1.0/U.time().deltaT())
231 );
232}
233
234
235template<class BasicTurbulenceModel>
237{
238 const alphaField& alpha = this->alpha_;
239 const rhoField& rho = this->rho_;
240
242 this->gasTurbulence();
243
244 const volScalarField phaseTransferCoeff(this->phaseTransferCoeff());
245
246 return
247 alpha*rho*bubbleG()
248 + phaseTransferCoeff*gasTurbulence.k()
249 - fvm::Sp(phaseTransferCoeff, this->k_);
250}
251
252
253template<class BasicTurbulenceModel>
255{
256 const alphaField& alpha = this->alpha_;
257 const rhoField& rho = this->rho_;
258
260 this->gasTurbulence();
261
262 const volScalarField phaseTransferCoeff(this->phaseTransferCoeff());
263
264 return
265 alpha*rho*this->C3_*this->epsilon_*bubbleG()/this->k_
266 + phaseTransferCoeff*gasTurbulence.epsilon()
267 - fvm::Sp(phaseTransferCoeff, this->epsilon_);
268}
269
270
271template<class BasicTurbulenceModel>
273{
275}
276
277
278// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
279
280} // End namespace RASModels
281} // End namespace Foam
282
283// ************************************************************************* //
surfaceScalarField & phi
twoPhaseSystem & fluid
void correctBoundaryConditions()
Correct boundary field.
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
Templated abstract base class for multiphase compressible turbulence models.
Continuous-phase k-epsilon model including bubble-generated turbulence.
Definition: LaheyKEpsilon.H:83
BasicTurbulenceModel::alphaField alphaField
virtual tmp< fvScalarMatrix > epsilonSource() const
BasicTurbulenceModel::rhoField rhoField
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
tmp< volScalarField > bubbleG() const
tmp< volScalarField > phaseTransferCoeff() const
virtual tmp< fvScalarMatrix > kSource() const
BasicTurbulenceModel::transportModel transportModel
virtual bool read()
Read model coefficients if they have changed.
Standard k-epsilon turbulence model for incompressible and compressible flows including rapid distort...
Definition: kEpsilon.H:92
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: kEpsilon.C:222
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
const alphaField & alpha() const
Access function to phase fraction.
const transportModel & transport() const
Access function to incompressible transport model.
Generic dimensioned Type class.
Generic thermophysical properties class for a liquid in which the functions and coefficients for each...
Definition: liquid.H:57
A class for managing temporary objects.
Definition: tmp.H:65
Base-class for all transport models used by the incompressible turbulence models.
Abstract base class for turbulence models (RAS, LES and laminar).
const volVectorField & U() const
Access function to velocity field.
static const word propertiesName
Default name of the turbulence properties dictionary.
virtual tmp< volScalarField > k() const =0
Return the turbulence kinetic energy.
virtual tmp< volScalarField > epsilon() const =0
Return the turbulence kinetic energy dissipation rate.
Class which solves the volume fraction equations for two phases.
A class for handling words, derived from Foam::string.
Definition: word.H:68
U
Definition: pEqn.H:72
zeroField Sp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
Namespace for OpenFOAM.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar pow3(const dimensionedScalar &ds)
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:598
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
volScalarField & alpha
Info<< "Reading strained laminar flame speed field Su\n"<< endl;volScalarField Su(IOobject("Su", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);Info<< "Reading field betav\n"<< endl;volScalarField betav(IOobject("betav", mesh.facesInstance(), mesh, IOobject::MUST_READ, IOobject::NO_WRITE), mesh);Info<< "Reading field Lobs\n"<< endl;volScalarField Lobs(IOobject("Lobs", mesh.facesInstance(), mesh, IOobject::MUST_READ, IOobject::NO_WRITE), mesh);Info<< "Reading field CT\n"<< endl;volSymmTensorField CT(IOobject("CT", mesh.facesInstance(), mesh, IOobject::MUST_READ, IOobject::NO_WRITE), mesh);Info<< "Reading field Nv\n"<< endl;volScalarField Nv(IOobject("Nv", mesh.facesInstance(), mesh, IOobject::MUST_READ, IOobject::NO_WRITE), mesh);Info<< "Reading field nsv\n"<< 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:165