KocamustafaogullariIshiiDepartureDiameter.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2016-2020 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 
31 #include "phaseSystem.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 namespace wallBoilingModels
38 {
39 namespace departureDiameterModels
40 {
41  defineTypeNameAndDebug(KocamustafaogullariIshiiDepartureDiameter, 0);
43  (
44  departureDiameterModel,
45  KocamustafaogullariIshiiDepartureDiameter,
46  dictionary
47  );
48 }
49 }
50 }
51 
52 
53 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
54 
58 (
59  const dictionary& dict
60 )
61 :
62  departureDiameterModel(),
63  phi_(dict.lookup<scalar>("phi"))
64 {}
65 
66 
70 (
71  const KocamustafaogullariIshiiDepartureDiameter& model
72 )
73 :
75  phi_(model.phi_)
76 {}
77 
78 
79 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
80 
84 {}
85 
86 
87 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
88 
92 (
93  const phaseModel& liquid,
94  const phaseModel& vapor,
95  const label patchi,
96  const scalarField& Tl,
97  const scalarField& Tsatw,
98  const scalarField& L
99 ) const
100 {
101  // Gravitational acceleration
103  liquid.mesh().lookupObject<uniformDimensionedVectorField>("g");
104 
105  const scalarField rhoLiquid(liquid.thermo().rho(patchi));
106  const scalarField rhoVapor(vapor.thermo().rho(patchi));
107 
108  const scalarField rhoM((rhoLiquid - rhoVapor)/rhoVapor);
109 
110  const scalarField sigmaw
111  (
112  liquid.fluid().sigma(phasePairKey(liquid.name(), vapor.name()), patchi)
113  );
114 
115  return
116  0.0012*pow(rhoM, 0.9)*0.0208*phi_
117  *sqrt(sigmaw/(mag(g.value())*(rhoLiquid - rhoVapor)));
118 }
119 
120 
123 {
125  writeEntry(os, "phi", phi_);
126 }
127 
128 
129 // ************************************************************************* //
dictionary dict
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
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.
UniformDimensionedField< vector > uniformDimensionedVectorField
dimensionedScalar sqrt(const dimensionedScalar &ds)
Macros for easy insertion into run-time selection tables.
virtual void write(Ostream &os) const
stressControl lookup("compactNormalStress") >> compactNormalStress
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
defineTypeNameAndDebug(combustionModel, 0)
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensioned< scalar > mag(const dimensioned< Type > &)
A class for managing temporary objects.
Definition: PtrList.H:53
KocamustafaogullariIshiiDepartureDiameter(const dictionary &dict)
Construct from a dictionary.
Namespace for OpenFOAM.