KocamustafaogullariIshiiNucleationSite.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) 2019-2023 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 
29 #include "phaseSystem.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace wallBoilingModels
36 {
37 namespace nucleationSiteModels
38 {
41  (
45  );
46 }
47 }
48 }
49 
50 
51 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
52 
53 template<class ScalarFieldType>
55 Foam::wallBoilingModels::nucleationSiteModels::
56 KocamustafaogullariIshiiNucleationSite::calculate
57 (
58  const ScalarFieldType& Tsatw,
59  const ScalarFieldType& L,
60  const ScalarFieldType& dDep,
61  const ScalarFieldType& Tw,
62  const ScalarFieldType& rhoLiquid,
63  const ScalarFieldType& rhoVapour,
64  const ScalarFieldType& sigma
65 ) const
66 {
67  const ScalarFieldType rhoM((rhoLiquid - rhoVapour)/rhoVapour);
68 
69  const dimensionedScalar zeroT_(dimTemperature, 0);
70  auto zeroT = coefficient<ScalarFieldType>::value(zeroT_);
71 
73 
74  // eq. (32)
75  const ScalarFieldType f
76  (
77  2.157e-7*pow(rhoM, -3.2)*pow(1 + 0.0049*rhoM, 4.13)
78  );
79 
80  // eq. (17)
81  const ScalarFieldType rRc
82  (
83  dDep*max(Tw - Tsatw, zeroT)*rhoVapour*L/(4*sigma*Tsatw)
84  );
85 
86  return Cn/sqr(dDep)*pow(rRc, 4.4)*f;
87 }
88 
89 
90 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
91 
94 (
95  const dictionary& dict
96 )
97 :
99  Cn_(dimensionedScalar::lookupOrDefault("Cn", dict, dimless, 1))
100 {}
101 
102 
105 (
107 )
108 :
110  Cn_(model.Cn_)
111 {}
112 
113 
114 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
115 
119 {}
120 
121 
122 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
123 
127 (
128  const phaseModel& liquid,
129  const phaseModel& vapour,
130  const label patchi,
131  const scalarField& Tl,
132  const scalarField& Tsatw,
133  const scalarField& L,
134  const scalarField& dDep,
135  const scalarField& fDep
136 ) const
137 {
138  const scalarField& Tw =
139  liquid.thermo().T().boundaryField()[patchi];
140 
141  return
142  calculate
143  (
144  Tsatw,
145  L,
146  dDep,
147  Tw,
148  static_cast<const scalarField&>
149  (
150  liquid.rho().boundaryField()[patchi]
151  ),
152  static_cast<const scalarField&>
153  (
154  vapour.rho().boundaryField()[patchi]
155  ),
156  liquid.fluid().sigma(phaseInterfaceKey(liquid, vapour), patchi)()
157  );
158 }
159 
160 
164 (
165  const phaseModel& liquid,
166  const phaseModel& vapour,
167  const phaseModel& solid,
168  const volScalarField& Tf,
169  const volScalarField& Tsatw,
170  const volScalarField& L,
171  const volScalarField& dDep,
172  const volScalarField& fDep
173 ) const
174 {
175  return calculate
176  (
177  Tsatw,
178  L,
179  dDep,
180  Tf,
181  liquid.rho(),
182  vapour.rho(),
183  liquid.fluid().sigma(phaseInterfaceKey(liquid, vapour))()
184  );
185 }
186 
187 
190 {
192  writeKeyword(os, "Cn") << Cn_ << token::END_STATEMENT << nl;
193 }
194 
195 
196 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
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 keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
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
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
@ END_STATEMENT
Definition: token.H:105
Base class for nucleation site density models.
virtual void write(Ostream &os) const
Write to stream.
virtual tmp< scalarField > nucleationSiteDensity(const phaseModel &liquid, const phaseModel &vapour, const label patchi, const scalarField &Tl, const scalarField &Tsatw, const scalarField &L, const scalarField &dDep, const scalarField &fDep) const
Calculate and return the nucleation-site density.
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, PatchField, GeoMesh > &distance)
Calculate distance data from patches.
defineTypeNameAndDebug(KocamustafaogullariIshiiNucleationSite, 0)
addToRunTimeSelectionTable(nucleationSiteModel, KocamustafaogullariIshiiNucleationSite, dictionary)
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
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimless
const dimensionSet dimTemperature
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
Ostream & writeKeyword(Foam::Ostream &os, const keyType &kw)
Write the keyword to the Ostream with the current level of indentation.
Definition: keyType.C:155
static const char nl
Definition: Ostream.H:260
labelList f(nPoints)
dictionary dict