Tatsumoto.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 -------------------------------------------------------------------------------
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 "Tatsumoto.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(Tatsumoto, 0);
44  (
45  CHFSubCoolModel,
46  Tatsumoto,
47  dictionary
48  );
49 }
50 }
51 }
52 
53 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
54 
55 Foam::wallBoilingModels::CHFModels::Tatsumoto::Tatsumoto
56 (
57  const dictionary& dict
58 )
59 :
60  CHFSubCoolModel(),
61  K_(dict.getOrDefault<scalar>("K", 0.23))
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 labelUList& cells = liquid.mesh().boundary()[patchi].faceCells();
79  const scalarField& pw = liquid.thermo().p().boundaryField()[patchi];
80 
81  tmp<scalarField> trhoVapor = vapor.thermo().rhoEoS(pw, Tsatw, cells);
82  const scalarField& rhoVapor = trhoVapor.ref();
83 
84  tmp<scalarField> trhoLiq = liquid.thermo().rhoEoS(pw, Tsatw, cells);
85  const scalarField& rhoLiq = trhoLiq.ref();
86 
87  tmp<scalarField> tCp = liquid.thermo().Cp(pw, Tsatw, cells);
88  const scalarField& Cp = tCp();
89 
90  return
91  1 + K_*pow(rhoVapor/rhoLiq, 0.8)*Cp*max(Tsatw - Tl, scalar(0))/L;
92 }
93 
94 
96 (
97  Ostream& os
98 ) const
99 {
101  os.writeEntry("K", K_);
102 }
103 
104 
105 // ************************************************************************* //
virtual tmp< scalarField > rhoEoS(const scalarField &p, const scalarField &T, const labelList &cells) const =0
Density from pressure and temperature from EoS.
Generic thermophysical properties class for a liquid in which the functions and coefficients for each...
Definition: liquid.H:50
dictionary dict
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
const vector L(dict.get< vector >("L"))
const tmp< volScalarField > & tCp
Definition: EEqn.H:4
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition: tmpI.H:235
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
Various UniformDimensionedField types.
virtual void write(Ostream &os) const
Write.
Macros for easy insertion into run-time selection tables.
const cellShapeList & cells
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: Tatsumoto.C:62
scalar Cp(scalar p, scalar T) const
Liquid heat capacity [J/(kg K)].
Definition: liquidI.H:39
addToRunTimeSelectionTable(CHFModel, Zuber, dictionary)
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
const volScalarField & Cp
Definition: EEqn.H:7
OBJstream os(runTime.globalPath()/outputName)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
A class for managing temporary objects.
Definition: HashPtrTable.H:50
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phaseModel.H:53
virtual const rhoThermo & thermo() const =0
Return the thermophysical model.
virtual void write(Ostream &os) const
Write.
Definition: Tatsumoto.C:89
Namespace for OpenFOAM.