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-2025 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 
28 #include "phaseSystem.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace wallBoilingModels
36 {
37 namespace departureDiameterModels
38 {
41  (
45  );
46 }
47 }
48 }
49 
50 
51 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
52 
53 template<class ScalarFieldType>
55 Foam::wallBoilingModels::departureDiameterModels::
56 KocamustafaogullariIshiiDepartureDiameter::calculate
57 (
58  const fvMesh& mesh,
59  const ScalarFieldType& Tl,
60  const ScalarFieldType& Tsatw,
61  const ScalarFieldType& L,
62  const ScalarFieldType& rhoLiquid,
63  const ScalarFieldType& rhoVapour,
64  const ScalarFieldType& sigma
65 ) const
66 {
68  (
70  );
71 
72  auto phi = coefficient<ScalarFieldType>::value(phi_);
73 
74  return
75  0.0012*pow((rhoLiquid - rhoVapour)/rhoVapour, 0.9)*0.0208*phi
76  *sqrt(sigma/(mag(g)*(rhoLiquid - rhoVapour)));
77 }
78 
79 
80 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
81 
85 (
86  const dictionary& dict
87 )
88 :
90  phi_("phi", dimless, dict)
91 {}
92 
93 
97 (
99 )
100 :
102  phi_(model.phi_)
103 {}
104 
105 
106 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
107 
111 {}
112 
113 
114 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
115 
119 (
120  const phaseModel& liquid,
121  const phaseModel& vapour,
122  const label patchi,
123  const scalarField& Tl,
124  const scalarField& Tsatw,
125  const scalarField& L
126 ) const
127 {
128 
129  return
130  calculate
131  (
132  liquid.mesh(),
133  Tl,
134  Tsatw,
135  L,
136  static_cast<const scalarField&>
137  (
138  liquid.rho().boundaryField()[patchi]
139  ),
140  static_cast<const scalarField&>
141  (
142  vapour.rho().boundaryField()[patchi]
143  ),
144  liquid.fluid().sigma(phaseInterfaceKey(liquid, vapour), patchi)()
145  );
146 }
147 
148 
152 (
153  const phaseModel& liquid,
154  const phaseModel& vapour,
155  const phaseModel& solid,
156  const volScalarField::Internal& Tf,
157  const volScalarField::Internal& Tsatw,
158  const volScalarField::Internal& L
159 ) const
160 {
161  return calculate
162  (
163  liquid.mesh(),
164  liquid.thermo().T()(),
165  Tsatw,
166  L,
167  liquid.rho()(),
168  vapour.rho()(),
169  liquid.fluid().sigma(phaseInterfaceKey(liquid, vapour))()()
170  );
171 }
172 
173 
177 (
178  const phaseModel& liquid,
179  const phaseModel& vapour,
180  const phaseModel& solid,
181  const volScalarField& Tf,
182  const volScalarField& Tsatw,
183  const volScalarField& L
184 ) const
185 {
186  return calculate
187  (
188  liquid.mesh(),
189  liquid.thermo().T(),
190  Tsatw,
191  L,
192  liquid.rho(),
193  vapour.rho(),
194  liquid.fluid().sigma(phaseInterfaceKey(liquid, vapour))()
195  );
196 }
197 
198 
201 {
203  writeEntry(os, "phi", phi_);
204 }
205 
206 
207 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Generic GeometricField class.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:96
Generic thermophysical properties class for a liquid in which the functions and coefficients for each...
Definition: liquid.H:53
scalar rho(scalar p, scalar T) const
Liquid density [kg/m^3].
Definition: liquidI.H:26
scalar sigma(scalar p, scalar T) const
Surface tension [N/m].
Definition: liquidI.H:104
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type and name.
Word-pair based class used for keying interface models in hash tables.
virtual const volScalarField & rho() const =0
Return the density field.
A class for managing temporary objects.
Definition: tmp.H:55
Base class for bubble departure diameter models.
virtual void write(Ostream &os) const
Write to stream.
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 for a wall patch.
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
label patchi
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m^2/K^4].
label calculate(const fvMesh &mesh, const labelHashSet &patchIDs, const scalar minFaceFraction, GeometricField< scalar, GeoMesh > &distance)
Calculate distance data from patches.
addToRunTimeSelectionTable(departureDiameterModel, KocamustafaogullariIshiiDepartureDiameter, dictionary)
defineTypeNameAndDebug(KocamustafaogullariIshiiDepartureDiameter, 0)
Namespace for OpenFOAM.
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
const dimensionSet dimless
void mag(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
void pow(LagrangianPatchField< typename powProduct< Type, r >::type > &f, const LagrangianPatchField< Type > &f1)
void sqrt(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
dictionary dict