continuousGasKEqn.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 "continuousGasKEqn.H"
28 #include "twoPhaseSystem.H"
29 
30 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace LESModels
35 {
36 
37 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
38 
39 template<class BasicTurbulenceModel>
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  liquidTurbulencePtr_(NULL),
65 
66  alphaInversion_
67  (
69  (
70  "alphaInversion",
71  this->coeffDict_,
72  0.7
73  )
74  )
75 {
76  if (type == typeName)
77  {
78  this->printCoeffs(type);
79  }
80 }
81 
82 
83 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
84 
85 template<class BasicTurbulenceModel>
87 {
89  {
90  alphaInversion_.readIfPresent(this->coeffDict());
91 
92  return true;
93  }
94  else
95  {
96  return false;
97  }
98 }
99 
100 
101 template<class BasicTurbulenceModel>
102 const turbulenceModel&
104 {
105  if (!liquidTurbulencePtr_)
106  {
107  const volVectorField& U = this->U_;
108 
109  const transportModel& gas = this->transport();
110  const twoPhaseSystem& fluid =
111  refCast<const twoPhaseSystem>(gas.fluid());
112  const transportModel& liquid = fluid.otherPhase(gas);
113 
114  liquidTurbulencePtr_ =
115  &U.db().lookupObject<turbulenceModel>
116  (
118  (
120  liquid.name()
121  )
122  );
123  }
124 
125  return *liquidTurbulencePtr_;
126 }
127 
128 
129 template<class BasicTurbulenceModel>
132 {
133  const volVectorField& U = this->U_;
134  const alphaField& alpha = this->alpha_;
135  const rhoField& rho = this->rho_;
136 
137  const turbulenceModel& liquidTurbulence = this->liquidTurbulence();
138 
139  return
140  (
141  max(alphaInversion_ - alpha, scalar(0))
142  *rho
143  *min
144  (
145  this->Ce_*sqrt(liquidTurbulence.k())/this->delta(),
146  1.0/U.time().deltaT()
147  )
148  );
149 }
150 
151 
152 template<class BasicTurbulenceModel>
155 {
156  const turbulenceModel& liquidTurbulence = this->liquidTurbulence();
157  const volScalarField phaseTransferCoeff(this->phaseTransferCoeff());
158 
159  return
160  phaseTransferCoeff*liquidTurbulence.k()
161  - fvm::Sp(phaseTransferCoeff, this->k_);
162 }
163 
164 
165 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
166 
167 } // End namespace LESModels
168 } // End namespace Foam
169 
170 // ************************************************************************* //
Foam::LESModels::continuousGasKEqn::alphaField
BasicTurbulenceModel::alphaField alphaField
Definition: continuousGasKEqn.H:109
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::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
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::LESModels::continuousGasKEqn::kSource
virtual tmp< fvScalarMatrix > kSource() const
Definition: continuousGasKEqn.C:154
Foam::turbulenceModel::propertiesName
static const word propertiesName
Default name of the turbulence properties dictionary.
Definition: turbulenceModel.H:97
Foam::LESModels::continuousGasKEqn::read
virtual bool read()
Read model coefficients if they have changed.
Definition: continuousGasKEqn.C:86
U
U
Definition: pEqn.H:46
Foam::LESModels::kEqn
One equation eddy-viscosity model.
Definition: kEqn.H:74
continuousGasKEqn.H
Foam::fvm::Sp
tmp< fvMatrix< Type > > Sp(const DimensionedField< scalar, volMesh > &, const GeometricField< Type, fvPatchField, volMesh > &)
Foam::LESModels::continuousGasKEqn::phaseTransferCoeff
tmp< volScalarField > phaseTransferCoeff() const
Definition: continuousGasKEqn.C:131
Foam::LESModels::continuousGasKEqn::liquidTurbulence
const turbulenceModel & liquidTurbulence() const
Return the turbulence model for the liquid phase.
Definition: continuousGasKEqn.C:103
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::LESModels::continuousGasKEqn::continuousGasKEqn
continuousGasKEqn(const continuousGasKEqn &)
rho
rho
Definition: pEqn.H:3
Foam::max
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::LESModels::continuousGasKEqn::rhoField
BasicTurbulenceModel::rhoField rhoField
Definition: continuousGasKEqn.H:110
Foam::LESModels::continuousGasKEqn::transportModel
BasicTurbulenceModel::transportModel transportModel
Definition: continuousGasKEqn.H:111
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 > &)