NicenoKEqn.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 "NicenoKEqn.H"
27 #include "fvOptions.H"
28 #include "twoPhaseSystem.H"
29 
30 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace LESModels
35 {
36 
37 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
38 
39 template<class BasicTurbulenceModel>
40 NicenoKEqn<BasicTurbulenceModel>::NicenoKEqn
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  this->Ck_.value()
83  )
84  ),
85 
86  Cmub_
87  (
89  (
90  "Cmub",
91  this->coeffDict_,
92  0.6
93  )
94  )
95 {
96  if (type == typeName)
97  {
98  this->printCoeffs(type);
99  }
100 }
101 
102 
103 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
104 
105 template<class BasicTurbulenceModel>
107 {
109  {
110  alphaInversion_.readIfPresent(this->coeffDict());
111  Cp_.readIfPresent(this->coeffDict());
112  Cmub_.readIfPresent(this->coeffDict());
113 
114  return true;
115  }
116  else
117  {
118  return false;
119  }
120 }
121 
122 
123 template<class BasicTurbulenceModel>
125 <
126  typename BasicTurbulenceModel::transportModel
127 >&
129 {
130  if (!gasTurbulencePtr_)
131  {
132  const volVectorField& U = this->U_;
133 
134  const transportModel& liquid = this->transport();
135  const twoPhaseSystem& fluid =
136  refCast<const twoPhaseSystem>(liquid.fluid());
137  const transportModel& gas = fluid.otherPhase(liquid);
138 
139  gasTurbulencePtr_ =
140  &U.db()
142  (
144  (
146  gas.name()
147  )
148  );
149  }
150 
151  return *gasTurbulencePtr_;
152 }
153 
154 
155 template<class BasicTurbulenceModel>
157 {
159  this->gasTurbulence();
160 
161  this->nut_ =
162  this->Ck_*sqrt(this->k_)*this->delta()
163  + Cmub_*gasTurbulence.transport().d()*gasTurbulence.alpha()
164  *(mag(this->U_ - gasTurbulence.U()));
165 
166  this->nut_.correctBoundaryConditions();
167  fv::options::New(this->mesh_).correct(this->nut_);
168 
169  BasicTurbulenceModel::correctNut();
170 }
171 
172 
173 template<class BasicTurbulenceModel>
175 {
177  this->gasTurbulence();
178 
179  const transportModel& liquid = this->transport();
180  const twoPhaseSystem& fluid =
181  refCast<const twoPhaseSystem>(liquid.fluid());
182  const transportModel& gas = fluid.otherPhase(liquid);
183 
184  volScalarField magUr(mag(this->U_ - gasTurbulence.U()));
185 
186  tmp<volScalarField> bubbleG
187  (
188  Cp_*sqr(magUr)*fluid.drag(gas).K()/liquid.rho()
189  );
190 
191  return bubbleG;
192 }
193 
194 
195 template<class BasicTurbulenceModel>
198 {
199  const volVectorField& U = this->U_;
200  const alphaField& alpha = this->alpha_;
201  const rhoField& rho = this->rho_;
202 
203  const turbulenceModel& gasTurbulence = this->gasTurbulence();
204 
205  return
206  (
207  max(alphaInversion_ - alpha, scalar(0))
208  *rho
209  *min
210  (
211  this->Ce_*sqrt(gasTurbulence.k())/this->delta(),
212  1.0/U.time().deltaT()
213  )
214  );
215 }
216 
217 
218 template<class BasicTurbulenceModel>
220 {
221  const alphaField& alpha = this->alpha_;
222  const rhoField& rho = this->rho_;
223 
225  this->gasTurbulence();
226 
227  const volScalarField phaseTransferCoeff(this->phaseTransferCoeff());
228 
229  return
230  alpha*rho*bubbleG()
231  + phaseTransferCoeff*gasTurbulence.k()
232  - fvm::Sp(phaseTransferCoeff, this->k_);
233 }
234 
235 
236 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
237 
238 } // End namespace LESModels
239 } // End namespace Foam
240 
241 // ************************************************************************* //
scalar delta
surfaceScalarField & phi
tmp< volScalarField > bubbleG() const
Definition: NicenoKEqn.C:174
const transportModel & transport() const
Access function to incompressible transport model.
U
Definition: pEqn.H:83
BasicTurbulenceModel::alphaField alphaField
Definition: NicenoKEqn.H:122
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)
const volVectorField & U() const
Access function to velocity field.
dimensionedScalar sqrt(const dimensionedScalar &ds)
Generic dimensioned Type class.
BasicTurbulenceModel::rhoField rhoField
Definition: NicenoKEqn.H:123
Abstract base class for turbulence models (RAS, LES and laminar).
One equation eddy-viscosity model.
Definition: kEqn.H:74
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.
static const word propertiesName
Default name of the turbulence properties dictionary.
virtual bool read()
Read model coefficients if they have changed.
Definition: NicenoKEqn.C:106
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 > kSource() const
Definition: NicenoKEqn.C:219
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 deltaT() const
Return time step.
Definition: TimeStateI.H:53
One-equation SGS model for the continuous phase in a two-phase system including bubble-generated turb...
Definition: NicenoKEqn.H:75
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:461
tmp< volScalarField > phaseTransferCoeff() const
Definition: NicenoKEqn.C:197
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
static options & New(const fvMesh &mesh)
Construct fvOptions and register to datbase if not present.
Definition: fvOptions.C:103
BasicTurbulenceModel::transportModel transportModel
Definition: NicenoKEqn.H:124
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
virtual tmp< volScalarField > K(const volScalarField &Ur) const =0
The dragfunction K used in the momentum eq.
virtual void correctNut()
Definition: NicenoKEqn.C:156
Namespace for OpenFOAM.