KocamustafaogullariIshii.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) 2016-2019 OpenFOAM Foundation
9  Copyright (C) 2018-2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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 
33 #include "ThermalDiffusivity.H"
35 #include "phaseSystem.H"
36 
37 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38 
39 namespace Foam
40 {
41 namespace wallBoilingModels
42 {
43 namespace departureDiameterModels
44 {
47  (
51  );
52 }
53 }
54 }
55 
56 
57 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
58 
59 Foam::wallBoilingModels::departureDiameterModels::
60 KocamustafaogullariIshii::KocamustafaogullariIshii
61 (
62  const dictionary& dict
63 )
64 :
66  phi_(dict.get<scalar>("phi"))
67 {}
68 
69 
70 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
71 
75 (
76  const phaseModel& liquid,
77  const phaseModel& vapor,
78  const label patchi,
79  const scalarField& Tl,
80  const scalarField& Tsatw,
81  const scalarField& L
82 ) const
83 {
84  // Gravitational acceleration
85  const auto& g =
86  liquid.mesh().time().lookupObject<uniformDimensionedVectorField>("g");
87 
88  const scalarField rhoLiquid(liquid.thermo().rho(patchi));
89  const scalarField rhoVapor(vapor.thermo().rho(patchi));
90 
91  const scalarField rhoM((rhoLiquid - rhoVapor)/rhoVapor);
92 
93  const tmp<volScalarField>& tsigma =
94  liquid.fluid().sigma(phasePairKey(liquid.name(), vapor.name()));
95 
96  const volScalarField& sigma = tsigma();
97  const fvPatchScalarField& sigmaw = sigma.boundaryField()[patchi];
98 
99  return
100  0.0012*pow(rhoM, 0.9)*0.0208*phi_
101  *sqrt(sigmaw/(mag(g.value())*(rhoLiquid - rhoVapor)));
102 }
103 
104 
107 {
109  os.writeEntry("phi", phi_);
110 }
111 
112 
113 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
Generic thermophysical properties class for a liquid in which the functions and coefficients for each...
Definition: liquid.H:50
dictionary dict
addToRunTimeSelectionTable(departureDiameterModel, KocamustafaogullariIshii, dictionary)
Base class for bubble departure diameter models for boiling flows.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
const vector L(dict.get< vector >("L"))
Various UniformDimensionedField types.
virtual tmp< scalarField > dDeparture(const phaseModel &liquid, const phaseModel &vapor, const label patchi, const scalarField &Tl, const scalarField &Tsatw, const scalarField &L) const
Calculate and return the departure diameter field.
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensionedScalar sigma("sigma", dimMass/sqr(dimTime), transportProperties)
Macros for easy insertion into run-time selection tables.
scalar rho(scalar p, scalar T) const
Liquid density [kg/m^3].
Definition: liquidI.H:21
defineTypeNameAndDebug(KocamustafaogullariIshii, 0)
An ordered or unorder pair of phase names. Typically specified as follows.
Definition: phasePairKey.H:61
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
const uniformDimensionedVectorField & g
OBJstream os(runTime.globalPath()/outputName)
virtual tmp< volScalarField > rho() const
Density [kg/m^3].
Definition: rhoThermo.C:217
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
A correlation for bubble departure diameter modelling based on Kocamustafaogullari-Ishii (1983) for b...
const word & name() const
Definition: phaseModel.H:166
A class for managing temporary objects.
Definition: HashPtrTable.H:50
scalar sigma(scalar p, scalar T) const
Surface tension [N/m].
Definition: liquidI.H:87
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phaseModel.H:53
virtual void write(Ostream &os) const
Write.
virtual const rhoThermo & thermo() const =0
Return the thermophysical model.
Namespace for OpenFOAM.