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-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 "NicenoKEqn.H"
28 #include "twoPhaseSystem.H"
29 #include "dragModel.H"
30 
31 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace LESModels
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  this->Ck_.value()
84  )
85  ),
86 
87  Cmub_
88  (
90  (
91  "Cmub",
92  this->coeffDict_,
93  0.6
94  )
95  )
96 {
97  if (type == typeName)
98  {
99  this->printCoeffs(type);
100  }
101 }
102 
103 
104 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
105 
106 template<class BasicTurbulenceModel>
108 {
110  {
111  alphaInversion_.readIfPresent(this->coeffDict());
112  Cp_.readIfPresent(this->coeffDict());
113  Cmub_.readIfPresent(this->coeffDict());
114 
115  return true;
116  }
117  else
118  {
119  return false;
120  }
121 }
122 
123 
124 template<class BasicTurbulenceModel>
126 <
127  typename BasicTurbulenceModel::transportModel
128 >&
130 {
131  if (!gasTurbulencePtr_)
132  {
133  const volVectorField& U = this->U_;
134 
135  const transportModel& liquid = this->transport();
136  const twoPhaseSystem& fluid =
137  refCast<const twoPhaseSystem>(liquid.fluid());
138  const transportModel& gas = fluid.otherPhase(liquid);
139 
140  gasTurbulencePtr_ =
141  &U.db()
143  (
145  (
147  gas.name()
148  )
149  );
150  }
151 
152  return *gasTurbulencePtr_;
153 }
154 
155 
156 template<class BasicTurbulenceModel>
158 {
160  this->gasTurbulence();
161 
162  this->nut_ =
163  this->Ck_*sqrt(this->k_)*this->delta()
164  + Cmub_*gasTurbulence.transport().d()*gasTurbulence.alpha()
165  *(mag(this->U_ - gasTurbulence.U()));
166 
167  this->nut_.correctBoundaryConditions();
168 }
169 
170 
171 template<class BasicTurbulenceModel>
173 {
175  this->gasTurbulence();
176 
177  const transportModel& liquid = this->transport();
178  const twoPhaseSystem& fluid =
179  refCast<const twoPhaseSystem>(liquid.fluid());
180  const transportModel& gas = fluid.otherPhase(liquid);
181 
182  volScalarField magUr(mag(this->U_ - gasTurbulence.U()));
183 
184  tmp<volScalarField> bubbleG
185  (
186  Cp_*sqr(magUr)*fluid.drag(gas).K()/liquid.rho()
187  );
188 
189  return bubbleG;
190 }
191 
192 
193 template<class BasicTurbulenceModel>
196 {
197  const volVectorField& U = this->U_;
198  const alphaField& alpha = this->alpha_;
199  const rhoField& rho = this->rho_;
200 
201  const turbulenceModel& gasTurbulence = this->gasTurbulence();
202 
203  return
204  (
205  max(alphaInversion_ - alpha, scalar(0))
206  *rho
207  *min
208  (
209  this->Ce_*sqrt(gasTurbulence.k())/this->delta(),
210  1.0/U.time().deltaT()
211  )
212  );
213 }
214 
215 
216 template<class BasicTurbulenceModel>
218 {
219  const alphaField& alpha = this->alpha_;
220  const rhoField& rho = this->rho_;
221 
223  this->gasTurbulence();
224 
225  const volScalarField phaseTransferCoeff(this->phaseTransferCoeff());
226 
227  return
228  alpha*rho*bubbleG()
229  + phaseTransferCoeff*gasTurbulence.k()
230  - fvm::Sp(phaseTransferCoeff, this->k_);
231 }
232 
233 
234 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
235 
236 } // End namespace LESModels
237 } // End namespace Foam
238 
239 // ************************************************************************* //
Foam::IOobject::groupName
static word groupName(Name name, const word &group)
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::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
Foam::LESModels::NicenoKEqn::phaseTransferCoeff
tmp< volScalarField > phaseTransferCoeff() const
Definition: NicenoKEqn.C:195
fluid
multiphaseSystem & fluid
Definition: createFields.H:10
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::LESModels::NicenoKEqn::rhoField
BasicTurbulenceModel::rhoField rhoField
Definition: NicenoKEqn.H:123
Foam::TurbulenceModel< volScalarField, volScalarField, compressibleTurbulenceModel, TransportModel >::alpha
const alphaField & alpha() const
Access function to phase fraction.
Definition: TurbulenceModel.H:148
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
U
U
Definition: pEqn.H:46
Foam::LESModels::kEqn
One equation eddy-viscosity model.
Definition: kEqn.H:74
Foam::fvm::Sp
tmp< fvMatrix< Type > > Sp(const DimensionedField< scalar, volMesh > &, const GeometricField< Type, fvPatchField, volMesh > &)
NicenoKEqn.H
Foam::LESModels::NicenoKEqn::kSource
virtual tmp< fvScalarMatrix > kSource() const
Definition: NicenoKEqn.C:217
Foam::LESModels::NicenoKEqn::alphaField
BasicTurbulenceModel::alphaField alphaField
Definition: NicenoKEqn.H:122
delta
scalar delta
Definition: LISASMDCalcMethod2.H:8
Foam::LESModels::NicenoKEqn::read
virtual bool read()
Read model coefficients if they have changed.
Definition: NicenoKEqn.C:107
Foam::LESModels::NicenoKEqn::bubbleG
tmp< volScalarField > bubbleG() const
Definition: NicenoKEqn.C:172
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
Foam::LESModels::NicenoKEqn::correctNut
virtual void correctNut()
Definition: NicenoKEqn.C:157
Foam::LESModels::NicenoKEqn::NicenoKEqn
NicenoKEqn(const NicenoKEqn &)
rho
rho
Definition: pEqn.H:3
Foam::LESModels::NicenoKEqn::gasTurbulence
const PhaseCompressibleTurbulenceModel< typename BasicTurbulenceModel::transportModel > & gasTurbulence() const
Return the turbulence model for the gas phase.
Definition: NicenoKEqn.C:129
Foam::max
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:49
Foam::LESModels::NicenoKEqn::transportModel
BasicTurbulenceModel::transportModel transportModel
Definition: NicenoKEqn.H:124
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:142
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::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