HuaXu.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 "HuaXu.H"
31 #include "phasePairKey.H"
32 #include "phaseSystem.H"
33 
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38 namespace wallBoilingModels
39 {
40 namespace CHFModels
41 {
42  defineTypeNameAndDebug(HuaXu, 0);
44  (
46  HuaXu,
47  dictionary
48  );
49 }
50 }
51 }
52 
53 
54 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
55 
57 (
58  const dictionary& dict
59 )
60 :
62  Kburn_(dict.getOrDefault<scalar>("Kburn", 1.5))
63 {}
64 
65 
66 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
67 
69 {}
70 
71 
72 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
73 
76 (
77  const phaseModel& liquid,
78  const phaseModel& vapor,
79  const label patchi,
80  const scalarField& Tl,
81  const scalarField& Tsatw,
82  const scalarField& L
83 ) const
84 {
86  liquid.mesh().time().lookupObject<uniformDimensionedVectorField>("g");
87 
88  const scalarField alphaLiq(liquid.alpha(patchi));
89 
90  const scalarField rhoVapor(vapor.thermo().rho(patchi));
91  const scalarField rhoLiq(liquid.thermo().rho(patchi));
92  tmp<volScalarField> tCp = liquid.thermo().Cp();
93  const volScalarField& Cp = tCp();
94  const fvPatchScalarField& Cpw = Cp.boundaryField()[patchi];
95 
96  const phasePairKey pair(liquid.name(), vapor.name());
97  const scalarField sigma
98  (
99  liquid.fluid().sigma(pair)().boundaryField()[patchi]
100  );
101 
102  const scalarField Pe
103  (
104  pow(sigma, 0.75)
105  /
106  (
107  alphaLiq
108  * pow(mag(g.value())*(rhoLiq-rhoVapor), 0.25)
109  * sqrt(rhoVapor)
110  )
111  );
112 
113  const scalarField Ja
114  (
115  rhoLiq*Cpw*max(Tsatw - Tl, scalar(0))/(rhoVapor*L)
116  );
117 
118  return
119  Kburn_*(1 + 0.345*Ja/pow(Pe, 0.25));
120 }
121 
122 
124 (
125  Ostream& os
126 ) const
127 {
129  os.writeEntry("Kburn", Kburn_);
130 }
131 
132 
133 // ************************************************************************* //
Foam::fvPatchField< scalar >
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::wallBoilingModels::CHFSubCoolModel
Definition: CHFSubCoolModel.H:56
Foam::wallBoilingModels::CHFModels::defineTypeNameAndDebug
defineTypeNameAndDebug(Zuber, 0)
Foam::wallBoilingModels::CHFModels::addToRunTimeSelectionTable
addToRunTimeSelectionTable(CHFModel, Zuber, dictionary)
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::phaseSystem::sigma
tmp< volScalarField > sigma(const phasePairKey &key) const
Return the surface tension coefficient for a pair.
Definition: phaseSystem.C:319
Foam::basicThermo::Cp
virtual tmp< volScalarField > Cp() const =0
Heat capacity at constant pressure [J/kg/K].
tCp
const tmp< volScalarField > & tCp
Definition: EEqn.H:4
Foam::wallBoilingModels::CHFSubCoolModel::write
virtual void write(Ostream &os) const
Definition: CHFSubCoolModel.C:74
Foam::phaseModel::thermo
virtual const rhoThermo & thermo() const =0
Access const to phase thermo.
CHFSubCoolModel
Base class for nucleation site density models.
Foam::wallBoilingModels::CHFModels::HuaXu::write
virtual void write(Ostream &os) const
Definition: HuaXu.C:124
wallBoilingModels
Critical heat flux (CHF) correlation.
Foam::wallBoilingModels::CHFModels::HuaXu::~HuaXu
virtual ~HuaXu()
Destructor.
Definition: HuaXu.C:68
Foam::UniformDimensionedField< vector >
Foam::rhoThermo::rho
virtual tmp< volScalarField > rho() const
Density [kg/m^3].
Definition: rhoThermo.C:224
Foam::Field< scalar >
Foam::phasePairKey
Definition: phasePairKey.H:57
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::wallBoilingModels::CHFModels::HuaXu::HuaXu
HuaXu(const dictionary &dict)
Construct from a dictionary.
Definition: HuaXu.C:57
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.
g
const uniformDimensionedVectorField & g
Definition: createFluidFields.H:24
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::wallBoilingModels::CHFModels::HuaXu::CHFSubCool
virtual tmp< scalarField > CHFSubCool(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: HuaXu.C:76
Foam::phaseModel::name
const word & name() const
Definition: phaseModel.H:140
uniformDimensionedFields.H
Foam::phaseModel::alpha
const volScalarField & alpha() const
Thermal diffusivity for enthalpy of mixture [kg/m/s].
Definition: phaseModel.C:221
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Cp
const volScalarField & Cp
Definition: EEqn.H:7
Foam::phaseModel::fluid
const phaseSystem & fluid() const
Return the system to which this phase belongs.
Definition: phaseModel.C:100
Foam::Ostream::writeEntry
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:232
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
sigma
dimensionedScalar sigma("sigma", dimMass/sqr(dimTime), transportProperties)
Foam::GeometricField< scalar, fvPatchField, volMesh >
HuaXu.H