KocamustafaogullariIshii.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-2019 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 
30 #include "ThermalDiffusivity.H"
32 #include "phaseSystem.H"
33 
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38 namespace wallBoilingModels
39 {
40 namespace departureDiameterModels
41 {
42  defineTypeNameAndDebug(KocamustafaogullariIshii, 0);
44  (
45  departureDiameterModel,
46  KocamustafaogullariIshii,
47  dictionary
48  );
49 }
50 }
51 }
52 
53 
54 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
55 
58 (
59  const dictionary& dict
60 )
61 :
62  departureDiameterModel(),
63  phi_(readScalar(dict.lookup("phi")))
64 {}
65 
66 
67 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
68 
71 {}
72 
73 
74 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
75 
79 (
80  const phaseModel& liquid,
81  const phaseModel& vapor,
82  const label patchi,
83  const scalarField& Tl,
84  const scalarField& Tsatw,
85  const scalarField& L
86 ) const
87 {
88  // Gravitational acceleration
90  liquid.mesh().lookupObject<uniformDimensionedVectorField>("g");
91 
92  const scalarField rhoLiquid(liquid.thermo().rho(patchi));
93  const scalarField rhoVapor(vapor.thermo().rho(patchi));
94 
95  const scalarField rhoM((rhoLiquid - rhoVapor)/rhoVapor);
96 
97  const tmp<volScalarField>& tsigma
98  (
99  liquid.fluid().sigma(phasePairKey(liquid.name(), vapor.name()))
100  );
101  const volScalarField& sigma = tsigma();
102  const fvPatchScalarField& sigmaw = sigma.boundaryField()[patchi];
103 
104  return
105  0.0012*pow(rhoM, 0.9)*0.0208*phi_
106  *sqrt(sigmaw/(mag(g.value())*(rhoLiquid - rhoVapor)));
107 }
108 
109 
111 KocamustafaogullariIshii::write(Ostream& os) const
112 {
114  writeEntry(os, "phi", phi_);
115 }
116 
117 
118 // ************************************************************************* //
KocamustafaogullariIshii(const dictionary &dict)
Construct from a dictionary.
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.
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
UniformDimensionedField< vector > uniformDimensionedVectorField
dimensionedScalar sqrt(const dimensionedScalar &ds)
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
virtual void write(Ostream &os) const
stressControl lookup("compactNormalStress") >> compactNormalStress
fvPatchField< scalar > fvPatchScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if successful.
Definition: doubleScalar.H:68
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)
label patchi
dimensioned< scalar > mag(const dimensioned< Type > &)
A class for managing temporary objects.
Definition: PtrList.H:53
Namespace for OpenFOAM.