Loading...
Searching...
No Matches
standardPhaseChange.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) 2011-2018 OpenFOAM Foundation
9 Copyright (C) 2020-2022 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "standardPhaseChange.H"
31#include "thermoSingleLayer.H"
32#include "zeroField.H"
33
34// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35
36namespace Foam
37{
38namespace regionModels
39{
41{
42
43// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
44
46
48(
49 phaseChangeModel,
52);
53
54// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
55
57(
58 const scalar Re,
59 const scalar Sc
60) const
61{
62 if (Re < 5.0E+05)
63 {
64 return 0.664*sqrt(Re)*cbrt(Sc);
65 }
66 else
67 {
68 return 0.037*pow(Re, 0.8)*cbrt(Sc);
69 }
70}
71
72
73// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
74
76(
77 surfaceFilmRegionModel& film,
78 const dictionary& dict
79)
80:
81 phaseChangeModel(typeName, film, dict),
82 deltaMin_(coeffDict_.get<scalar>("deltaMin")),
83 L_(coeffDict_.get<scalar>("L")),
84 TbFactor_(coeffDict_.getOrDefault<scalar>("TbFactor", 1.1)),
85 YInfZero_(coeffDict_.getOrDefault("YInfZero", false))
86{}
87
88
89// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
90
91template<class YInfType>
93(
94 const scalar dt,
95 scalarField& availableMass,
96 scalarField& dMass,
97 scalarField& dEnergy,
98 const YInfType& YInf
99)
100{
102
103 // Set local thermo properties
104 const SLGThermo& thermo = film.thermo();
105 const filmThermoModel& filmThermo = film.filmThermo();
106 const label vapId = thermo.carrierId(filmThermo.name());
107
108 // Retrieve fields from film model
109 const scalarField& delta = film.delta();
110 const scalarField& pInf = film.pPrimary();
111 const scalarField& T = film.T();
112 const scalarField& hs = film.hs();
113 const scalarField& rho = film.rho();
114 const scalarField& rhoInf = film.rhoPrimary();
115 const scalarField& muInf = film.muPrimary();
116 const scalarField& magSf = film.magSf();
117 const vectorField dU(film.UPrimary() - film.Us());
118 const scalarField limMass
119 (
120 max(scalar(0), availableMass - deltaMin_*rho*magSf)
121 );
122
123 // Molecular weight of vapour [kg/kmol]
124 const scalar Wvap = thermo.carrier().W(vapId);
125
126 // Molecular weight of liquid [kg/kmol]
127 const scalar Wliq = filmThermo.W();
128
129 forAll(dMass, celli)
130 {
131 scalar dm = 0;
132
133 if (delta[celli] > deltaMin_)
134 {
135 // Cell pressure [Pa]
136 const scalar pc = pInf[celli];
137
138 // Calculate the boiling temperature
139 const scalar Tb = filmThermo.Tb(pc);
140
141 // Local temperature - impose lower limit of 200 K for stability
142 const scalar Tloc = min(TbFactor_*Tb, max(200.0, T[celli]));
143
144 // Saturation pressure [Pa]
145 const scalar pSat = filmThermo.pv(pc, Tloc);
146
147 // Latent heat [J/kg]
148 const scalar hVap = filmThermo.hl(pc, Tloc);
149
150 // Calculate mass transfer
151 if (pSat >= 0.95*pc)
152 {
153 // Boiling
154 const scalar Cp = filmThermo.Cp(pc, Tloc);
155 const scalar Tcorr = max(0.0, T[celli] - Tb);
156 const scalar qCorr = limMass[celli]*Cp*(Tcorr);
157 dm = qCorr/hVap;
158 }
159 else
160 {
161 // Primary region density [kg/m3]
162 const scalar rhoInfc = rhoInf[celli];
163
164 // Primary region viscosity [Pa.s]
165 const scalar muInfc = muInf[celli];
166
167 // Reynolds number
168 const scalar Re = rhoInfc*mag(dU[celli])*L_/muInfc;
169
170 // Vapour mass fraction at interface
171 const scalar Ys = Wliq*pSat/(Wliq*pSat + Wvap*(pc - pSat));
172
173 // Vapour diffusivity [m2/s]
174 const scalar Dab = filmThermo.D(pc, Tloc);
175
176 // Schmidt number
177 const scalar Sc = muInfc/(rhoInfc*(Dab + ROOTVSMALL));
178
179 // Sherwood number
180 const scalar Sh = this->Sh(Re, Sc);
181
182 // Mass transfer coefficient [m/s]
183 const scalar hm = Sh*Dab/(L_ + ROOTVSMALL);
184
185 // Add mass contribution to source
186 dm = dt*magSf[celli]*rhoInfc*hm*(Ys - YInf[celli])/(1.0 - Ys);
187 }
188
189 dMass[celli] += min(limMass[celli], max(dm, 0));
190
191 // Heat is assumed to be removed by heat-transfer to the wall
192 // so the energy remains unchanged by the phase-change.
193 dEnergy[celli] += dm*hs[celli];
194 // dEnergy[celli] += dm*(hs[celli] + hVap);
195 }
196 }
197}
198
199
201(
202 const scalar dt,
203 scalarField& availableMass,
204 scalarField& dMass,
205 scalarField& dEnergy
206)
207{
208 if (YInfZero_)
209 {
210 correctModel(dt, availableMass, dMass, dEnergy, zeroField());
211 }
212 else
213 {
214 const thermoSingleLayer& film = filmType<thermoSingleLayer>();
215 const label vapId = film.thermo().carrierId(film.filmThermo().name());
216 const scalarField& YInf = film.YPrimary()[vapId];
217
218 correctModel(dt, availableMass, dMass, dEnergy, YInf);
219 }
220}
221
222
223// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
224
225} // End namespace surfaceFilmModels
226} // End namespace regionModels
227} // End namespace Foam
228
229// ************************************************************************* //
scalar delta
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
Thermo package for (S)olids (L)iquids and (G)ases Takes reference to thermo package,...
Definition SLGThermo.H:63
const word & name() const
Name function is needed to disambiguate those inherited from regIOobject and dictionary.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
virtual const volScalarField & magSf() const
Return the face area magnitudes / [m2].
const surfaceFilmRegionModel & film() const
Return const access to the film surface film model.
virtual const word & name() const =0
Return the specie name.
virtual scalar Tb(const scalar p) const =0
Return boiling temperature [K].
virtual scalar W() const =0
Return molecular weight [kg/kmol].
virtual scalar pv(const scalar p, const scalar T) const =0
Return vapour pressure [Pa].
virtual scalar hl(const scalar p, const scalar T) const =0
Return latent heat [J/kg].
virtual scalar D(const scalar p, const scalar T) const =0
Return diffusivity [m2/s].
virtual scalar Cp(const scalar p, const scalar T) const =0
Return specific heat capacity [J/kg/K].
Base class for surface film phase change models.
Standard phase change model with modification for boiling.
const scalar TbFactor_
Boiling temperature factor / [].
standardPhaseChange(const standardPhaseChange &)=delete
No copy construct.
const scalar deltaMin_
Minimum film height for model to be active.
void correctModel(const scalar dt, scalarField &availableMass, scalarField &dMass, scalarField &dEnergy, const YInfType &YInf)
scalar Sh(const scalar Re, const scalar Sc) const
Return Sherwood number as a function of Reynolds and Schmidt numbers.
virtual const volScalarField & hs() const =0
Return the film surface temperature [J/kg].
virtual const volScalarField & T() const =0
Return the film mean temperature [K].
virtual const volVectorField & Us() const =0
Return the film surface velocity [m/s].
virtual const volScalarField & rho() const =0
Return the film density [kg/m3].
virtual const volScalarField & delta() const =0
Return the film thickness [m].
Thermodynamic form of single-cell layer surface film model.
const dictionary coeffDict_
Coefficients dictionary.
const dictionary & dict() const
Return const access to the cloud dictionary.
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
A class representing the concept of a field of 0 used to avoid unnecessary manipulations for objects ...
Definition zeroField.H:50
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
const volScalarField & Cp
Definition EEqn.H:7
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
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
Field< vector > vectorField
Specialisation of Field<T> for vector.
scalarField Re(const UList< complex > &cmplx)
Extract real component.
dimensionedScalar cbrt(const dimensionedScalar &ds)
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299