Bromley.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) 2018-2021 OpenCFD Ltd
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
28#include "Bromley.H"
31#include "constants.H"
32
33using namespace Foam::constant;
34
35// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36
37namespace Foam
38{
39namespace wallBoilingModels
40{
41namespace filmBoilingModels
42{
45 (
47 Bromley,
49 );
50}
51}
52}
53
54
55// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
56
58(
59 const dictionary& dict
60)
61:
63 Cn_(dict.getOrDefault<scalar>("Cn", 0.62)),
64 emissivity_(dict.getOrDefault<scalar>("emissivity", 1)),
65 L_(dict.get<scalar>("L"))
66{}
67
68
69// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
70
73(
74 const phaseModel& liquid,
75 const phaseModel& vapor,
76 const label patchi,
77 const scalarField& Tl,
78 const scalarField& Tsatw,
79 const scalarField& L
80) const
81{
82 const fvPatchScalarField& Tw =
83 liquid.thermo().T().boundaryField()[patchi];
84 const auto& g =
85 liquid.mesh().time().lookupObject<uniformDimensionedVectorField>("g");
86
87 const labelUList& cells = liquid.mesh().boundary()[patchi].faceCells();
88
89 const scalarField& pw = liquid.thermo().p().boundaryField()[patchi];
90
91 tmp<scalarField> trhoVapor = vapor.thermo().rhoEoS(pw, Tsatw, cells);
92 const scalarField& rhoVapor = trhoVapor.ref();
93
94 tmp<scalarField> trhoLiq = liquid.thermo().rhoEoS(pw, Tsatw, cells);
95 const scalarField& rhoLiq = trhoLiq.ref();
96
97
98 const scalarField kappaVapor(vapor.kappa(patchi));
99
100 tmp<scalarField> tCp = vapor.thermo().Cp(pw, Tsatw, cells);
101 const scalarField& CpVapor = tCp();
102
103 const scalarField muVapor(vapor.mu(patchi));
104
105 const scalarField htcRad
106 (
107 emissivity_*physicoChemical::sigma.value()*(pow4(Tw) - pow4(Tsatw))
108 / max((Tw - Tsatw), scalar(1e-4))
109 );
110
111 return
112 Cn_*pow025
113 (
114 pow3(kappaVapor)
115 *rhoVapor*(rhoLiq - rhoVapor)*mag(g.value())
116 *(L + scalar(0.4)*CpVapor*max((Tw-Tsatw), scalar(0)))
117 /(L_*muVapor*max((Tw-Tsatw), scalar(1e-4)))
118 ) + scalar(0.75)*htcRad;
119}
120
121
123(
124 Ostream& os
125) const
126{
128 os.writeEntry("Cn", Cn_);
129 os.writeEntry("L", L_);
130 os.writeEntry("emissivity", emissivity_);
131}
132
133
134// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
const uniformDimensionedVectorField & g
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
virtual tmp< scalarField > rhoEoS(const scalarField &p, const scalarField &T, const labelList &cells) const =0
Density from pressure and temperature from EoS.
virtual tmp< volScalarField > Cp() const =0
Heat capacity at constant pressure [J/kg/K].
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
virtual bool write()
Write the output fields.
Generic thermophysical properties class for a liquid in which the functions and coefficients for each...
Definition: liquid.H:57
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phaseModel.H:61
virtual const rhoThermo & thermo() const =0
Return the thermophysical model.
virtual tmp< volScalarField > mu() const =0
Return the laminar dynamic viscosity.
const dimensionedScalar & kappa() const
Definition: phaseModel.H:161
A class for managing temporary objects.
Definition: tmp.H:65
T & ref() const
Definition: tmpI.H:227
A correlation for boiling film modelling based on Bromley (1950) for boiling flows.
Definition: Bromley.H:119
virtual tmp< scalarField > htcFilmBoil(const phaseModel &liquid, const phaseModel &vapor, const label patchi, const scalarField &Tl, const scalarField &Tsatw, const scalarField &L) const
Calculate and return the film boiling correlation.
Definition: Bromley.C:73
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
const tmp< volScalarField > & tCp
Definition: EEqn.H:4
OBJstream os(runTime.globalPath()/outputName)
const cellShapeList & cells
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m2/K4].
Different types of constants.
Namespace for OpenFOAM.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
dimensionedScalar pow3(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
dimensionedScalar pow4(const dimensionedScalar &ds)
dimensionedScalar pow025(const dimensionedScalar &ds)
dictionary dict
volScalarField & e
Definition: createFields.H:11
const vector L(dict.get< vector >("L"))