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-2020 OpenCFD Ltd
9 -------------------------------------------------------------------------------
10 License
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 
33 using namespace Foam::constant;
34 
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39 namespace wallBoilingModels
40 {
41 namespace 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 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
70 
72 {}
73 
74 
75 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
76 
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 
89  const fvPatchScalarField& Tw =
90  liquid.thermo().T().boundaryField()[patchi];
92  liquid.mesh().time().lookupObject<uniformDimensionedVectorField>("g");
93 
94  const fvPatchScalarField& rhoVaporw
95  (
96  vapor.thermo().rho()().boundaryField()[patchi]
97  );
98 
99  const scalarField rhoLiq(liquid.thermo().rho(patchi));
100  const scalarField kappaVapor(vapor.kappa(patchi));
101 
102  tmp<volScalarField> tCp = vapor.thermo().Cp();
103  const volScalarField& Cp = tCp();
104  const scalarField& CpVapor = Cp.boundaryField()[patchi];
105 
106  const scalarField muVapor(vapor.mu(patchi));
107  //const scalarField dbVapor(vapor.d()().boundaryField()[patchi]);
108 
109  const scalarField htcRad
110  (
111  emissivity_*physicoChemical::sigma.value()*(pow4(Tw) - pow4(Tsatw))
112  / max((Tw - Tsatw), scalar(1e-4))
113  );
114 
115  return
116  Cn_*pow
117  (
118  pow3(kappaVapor)
119  *rhoVaporw*(rhoLiq - rhoVaporw)*mag(g.value())
120  *(L + 0.4*CpVapor*max((Tw-Tsatw), scalar(0)))
121  /(L_*muVapor*max((Tw-Tsatw), scalar(1e-4))),
122  0.25
123  ) + 0.75*htcRad;
124 }
125 
126 
128 (
129  Ostream& os
130 ) const
131 {
133  os.writeEntry("Cn", Cn_);
134  os.writeEntry("L", L_);
135  os.writeEntry("emissivity", emissivity_);
136 }
137 
138 
139 // ************************************************************************* //
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
Foam::fvPatchField< scalar >
Foam::phaseModel::mu
virtual tmp< volScalarField > mu() const
Return the mixture dymanic viscosity.
Definition: phaseModel.C:297
Foam::phaseModel
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phaseModel.H:54
L
const vector L(dict.get< vector >("L"))
Foam::phaseModel::kappa
const dimensionedScalar & kappa() const
Definition: phaseModel.H:157
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::basicThermo::Cp
virtual tmp< volScalarField > Cp() const =0
Heat capacity at constant pressure [J/kg/K].
Foam::constant
Different types of constants.
Definition: atomicConstants.C:38
Foam::wallBoilingModels::filmBoilingModel
Definition: filmBoilingModel.H:56
tCp
const tmp< volScalarField > & tCp
Definition: EEqn.H:4
Foam::phaseModel::thermo
virtual const rhoThermo & thermo() const =0
Access const to phase thermo.
Foam::wallBoilingModels::filmBoilingModels::Bromley::Bromley
Bromley(const dictionary &dict)
Construct from a dictionary.
Definition: Bromley.C:58
wallBoilingModels
Critical heat flux (CHF) correlation.
Foam::pow4
dimensionedScalar pow4(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:100
Foam::UniformDimensionedField< vector >
Bromley.H
Foam::rhoThermo::rho
virtual tmp< volScalarField > rho() const
Density [kg/m^3].
Definition: rhoThermo.C:224
Foam::Field< scalar >
Foam::pow3
dimensionedScalar pow3(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:89
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:75
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::wallBoilingModels::filmBoilingModels::Bromley::htcFilmBoil
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 nucleation-site density.
Definition: Bromley.C:79
g
const uniformDimensionedVectorField & g
Definition: createFluidFields.H:24
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::wallBoilingModels::filmBoilingModels::Bromley
Definition: Bromley.H:60
constants.H
Foam::wallBoilingModels::filmBoilingModels::Bromley::~Bromley
virtual ~Bromley()
Destructor.
Definition: Bromley.C:71
uniformDimensionedFields.H
Foam::basicThermo::T
virtual const volScalarField & T() const
Temperature [K].
Definition: basicThermo.C:524
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Cp
const volScalarField & Cp
Definition: EEqn.H:7
Foam::constant::physicoChemical::sigma
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m2/K4].
Foam::Ostream::writeEntry
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:232
Foam::vtk::write
void write(vtk::formatter &fmt, const Type &val, const label n=1)
Component-wise write of a value (N times)
Definition: foamVtkOutputTemplates.C:35
Foam::wallBoilingModels::filmBoilingModels::Bromley::write
virtual void write(Ostream &os) const
Definition: Bromley.C:128
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::GeometricField::boundaryField
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Definition: GeometricFieldI.H:62