Loading...
Searching...
No Matches
BreenWestwater.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) 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 "BreenWestwater.H"
31#include "phasePairKey.H"
32#include "phaseSystem.H"
34// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36namespace Foam
37{
38namespace wallBoilingModels
39{
40namespace filmBoilingModels
41{
44 (
48 );
49}
50}
51}
52
53// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
54
55Foam::wallBoilingModels::filmBoilingModels::BreenWestwater::BreenWestwater
56(
57 const dictionary& dict
58)
59:
61 Cn_(dict.getOrDefault<scalar>("Cn", 0.37))
62{}
63
64
65// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
66
69(
70 const phaseModel& liquid,
71 const phaseModel& vapor,
72 const label patchi,
73 const scalarField& Tl,
74 const scalarField& Tsatw,
75 const scalarField& L
76) const
77{
78 const fvPatchScalarField& Tw =
79 liquid.thermo().T().boundaryField()[patchi];
80
81 const auto& g =
82 liquid.mesh().time().lookupObject<uniformDimensionedVectorField>("g");
83
84 const labelUList& cells = liquid.mesh().boundary()[patchi].faceCells();
85 const scalarField& pw = liquid.thermo().p().boundaryField()[patchi];
86
87 tmp<scalarField> trhoVapor = vapor.thermo().rhoEoS(pw, Tsatw, cells);
88 const scalarField& rhoVapor = trhoVapor.ref();
89
90 tmp<scalarField> trhoLiq = liquid.thermo().rhoEoS(pw, Tsatw, cells);
91 const scalarField& rhoLiq = trhoLiq.ref();
92
93
94 const scalarField kappaLiquid(liquid.kappa(patchi));
95
96 tmp<scalarField> tCp = vapor.thermo().Cp(pw, Tsatw, cells);
97 const scalarField& CpVapor = tCp();
98
99 const scalarField nuLiquid(liquid.nu(patchi));
100
101 const scalarField Leff
102 (
103 L*sqr(1 + 0.34*CpVapor*max((Tw-Tsatw), scalar(0))/L)
104 );
105
106 const scalarField rhoDiff(rhoLiq - rhoVapor);
107
108 const phasePairKey pair(liquid.name(), vapor.name());
109 const scalarField sigma
110 (
111 liquid.fluid().sigma(pair)().boundaryField()[patchi]
112 );
113
114 return
115 Cn_
116 /pow(sigma/mag(g.value())/rhoDiff, 1/8)
117 /pow
118 (
119 nuLiquid*max((Tw-Tsatw), scalar(1e-3))
120 /(pow3(kappaLiquid)*rhoVapor*Leff*mag(g.value())*rhoDiff)
121 , 0.25
122 );
123}
124
125
127(
128 Ostream& os
129) const
130{
132 os.writeEntry("Cn", Cn_);
133}
134
135
136// ************************************************************************* //
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:59
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:133
Base class for film boiling models.
Generic thermophysical properties class for a liquid in which the functions and coefficients for each...
Definition liquid.H:53
scalar sigma(scalar p, scalar T) const
Surface tension [N/m].
Definition liquidI.H:87
scalar kappa(scalar p, scalar T) const
Liquid thermal conductivity [W/(m K)].
Definition liquidI.H:75
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition phaseModel.H:56
const word & name() const
Definition phaseModel.H:166
virtual const rhoThermo & thermo() const =0
Return the thermophysical model.
An ordered or unorder pair of phase names. Typically specified as follows.
A class for managing temporary objects.
Definition tmp.H:75
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition tmpI.H:235
filmBoilingModel()=default
Default construct.
virtual void write(Ostream &os) const
Write.
Boiling film correlation. A correlation for boiling film modelling based on Breen & Westwater (1965) ...
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.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
const tmp< volScalarField > & tCp
Definition EEqn.H:4
OBJstream os(runTime.globalPath()/outputName)
const cellShapeList & cells
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:40
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar pow3(const dimensionedScalar &ds)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
UniformDimensionedField< vector > uniformDimensionedVectorField
UList< label > labelUList
A UList of labels.
Definition UList.H:75
fvPatchField< scalar > fvPatchScalarField
dictionary dict
volScalarField & e
dimensionedScalar sigma("sigma", dimMass/sqr(dimTime), transportProperties)
Various UniformDimensionedField types.
const vector L(dict.get< vector >("L"))