standardPhaseChange.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2021 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "standardPhaseChange.H"
27 #include "thermoSingleLayer.H"
28 #include "liquidThermo.H"
29 #include "basicSpecieMixture.H"
30 #include "zeroField.H"
31 
33 
34 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38 namespace regionModels
39 {
40 namespace surfaceFilmModels
41 {
42 
43 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
44 
45 defineTypeNameAndDebug(standardPhaseChange, 0);
46 
48 (
49  phaseChangeModel,
50  standardPhaseChange,
51  dictionary
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 (
78  const dictionary& dict
79 )
80 :
81  speciePhaseChange(typeName, film, dict),
82  deltaMin_(coeffDict_.lookup<scalar>("deltaMin")),
83  L_(coeffDict_.lookup<scalar>("L")),
84  TbFactor_(coeffDict_.lookupOrDefault<scalar>("TbFactor", 1.1)),
85  YInfZero_(coeffDict_.lookupOrDefault<Switch>("YInfZero", false))
86 {}
87 
88 
89 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
90 
92 {}
93 
94 
95 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
96 
97 template<class YInfType>
99 (
100  const scalar dt,
101  scalarField& availableMass,
102  scalarField& dMass,
103  scalarField& dEnergy,
104  YInfType YInf
105 )
106 {
107  const thermoSingleLayer& film = filmType<thermoSingleLayer>();
108 
109  // Set local liquidThermo properties
110  const liquidProperties& liquidThermo =
111  refCast<const heRhoThermopureMixtureliquidProperties>(film.thermo())
112  .cellThermoMixture(0).properties();
113 
114  const basicSpecieMixture& primarySpecieThermo =
115  refCast<const basicSpecieMixture>(film.primaryThermo());
116 
117  // Retrieve fields from film model
118  const scalarField& delta = film.delta();
119  const scalarField& pInf = film.pPrimary();
120  const scalarField& T = film.thermo().T();
121  const scalarField& rho = film.rho();
122  const scalarField& rhoInf = film.rhoPrimary();
123  const scalarField& muInf = film.muPrimary();
124  const scalarField& magSf = film.magSf();
125  const vectorField dU(film.UPrimary() - film.Us());
126  const scalarField limMass
127  (
128  max(scalar(0), availableMass - deltaMin_*rho*magSf)
129  );
130 
131  // Molecular weight of vapour [kg/kmol]
132  const scalar Wvap = this->Wvap();
133 
134  // Molecular weight of liquid [kg/kmol]
135  const scalar Wliq = liquidThermo.W();
136 
137  forAll(dMass, celli)
138  {
139  scalar dm = 0;
140 
141  if (delta[celli] > deltaMin_)
142  {
143  // Cell pressure [Pa]
144  const scalar pc = pInf[celli];
145 
146  // Calculate the boiling temperature
147  const scalar Tb = liquidThermo.pvInvert(pc);
148 
149  // Local surface temperature at which evaporation takes place
150  // impose lower limit of 200 K for stability
151  const scalar Tloc = min(TbFactor_*Tb, max(200.0, T[celli]));
152 
153  // Saturation pressure [Pa]
154  const scalar pSat = liquidThermo.pv(pc, Tloc);
155 
156  // Latent heat [J/kg]
157  const scalar hVap = liquidThermo.hl(pc, Tloc);
158 
159  // Calculate mass transfer
160  if (pSat >= 0.95*pc)
161  {
162  // Boiling
163  const scalar Cp = liquidThermo.Cp(pc, Tloc);
164  const scalar Tcorr = max(0.0, T[celli] - Tb);
165  const scalar qCorr = limMass[celli]*Cp*(Tcorr);
166  dm = qCorr/hVap;
167  }
168  else
169  {
170  // Primary region density [kg/m^3]
171  const scalar rhoInfc = rhoInf[celli];
172 
173  // Primary region viscosity [Pa.s]
174  const scalar muInfc = muInf[celli];
175 
176  // Reynolds number
177  const scalar Re = rhoInfc*mag(dU[celli])*L_/muInfc;
178 
179  // Vapour mass fraction at interface
180  const scalar Ys = Wliq*pSat/(Wliq*pSat + Wvap*(pc - pSat));
181 
182  // Vapour diffusivity [m^2/s]
183  const scalar Dab = liquidThermo.D(pc, Tloc);
184 
185  // Schmidt number
186  const scalar Sc = muInfc/(rhoInfc*(Dab + rootVSmall));
187 
188  // Sherwood number
189  const scalar Sh = this->Sh(Re, Sc);
190 
191  // Mass transfer coefficient [m/s]
192  const scalar hm = Sh*Dab/(L_ + rootVSmall);
193 
194  // Add mass contribution to source
195  dm = dt*magSf[celli]*rhoInfc*hm*(Ys - YInf[celli])/(1.0 - Ys);
196  }
197 
198  dm = min(limMass[celli], max(dm, 0));
199 
200  dMass[celli] += dm;
201 
202  // Assume that the vapour transferred to the primary region is
203  // already at temperature Tloc so that all heat required for
204  // the phase-change is provided by the film
205  dEnergy[celli] += dm*primarySpecieThermo.Hs(vapId(), pc, Tloc);
206  }
207  }
208 }
209 
210 
212 (
213  const scalar dt,
214  scalarField& availableMass,
215  scalarField& dMass,
216  scalarField& dEnergy
217 )
218 {
219  if (YInfZero_)
220  {
221  correctModel(dt, availableMass, dMass, dEnergy, zeroField());
222  }
223  else
224  {
225  const thermoSingleLayer& film = filmType<thermoSingleLayer>();
226  const scalarField& YInf = film.YPrimary()[vapId()];
227 
228  correctModel(dt, availableMass, dMass, dEnergy, YInf);
229  }
230 }
231 
232 
233 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
234 
235 } // End namespace surfaceFilmModels
236 } // End namespace regionModels
237 } // End namespace Foam
238 
239 // ************************************************************************* //
scalar delta
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
dimensionedScalar sqrt(const dimensionedScalar &ds)
A simple wrapper around bool so that it can be read as a word: true/false, on/off, yes/no, y/n, t/f, or none/any.
Definition: Switch.H:60
const volScalarField & rhoPrimary() const
Density [kg/m^3].
Specialisation of basicMixture for a mixture consisting of a number for molecular species...
tmp< volVectorField::Internal > Us() const
Return the film surface velocity [m/s].
Macros for easy insertion into run-time selection tables.
addToRunTimeSelectionTable(ejectionModel, BrunDrippingEjection, dictionary)
const volScalarField & muPrimary() const
Viscosity [Pa.s].
const volVectorField & UPrimary() const
Velocity [m/s].
scalar Sh(const scalar Re, const scalar Sc) const
Return Sherwood number as a function of Reynolds and Schmidt numbers.
A class representing the concept of a field of 0 used to avoid unnecessary manipulations for objects ...
Definition: zeroField.H:50
standardPhaseChange(surfaceFilmRegionModel &film, const dictionary &dict)
Construct from surface film model.
const scalar TbFactor_
Boiling temperature factor / [].
dimensionedScalar cbrt(const dimensionedScalar &ds)
The thermophysical properties of a liquid.
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
const PtrList< volScalarField > & YPrimary() const
Specie mass fractions [0-1].
const volScalarField & pPrimary() const
Pressure [Pa].
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
const volScalarField & delta() const
Return const access to the film thickness [m].
void correctModel(const scalar dt, scalarField &availableMass, scalarField &dMass, scalarField &dEnergy, YInfType YInf)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
const volScalarField::Internal & magSf() const
Return the face area magnitudes [m^2].
virtual const volScalarField & T() const =0
Temperature [K].
virtual scalar Hs(const label speciei, const scalar p, const scalar T) const =0
Sensible enthalpy [J/kg].
scalar Cp(const scalar p, const scalar T) const
Definition: EtoHthermo.H:2
volScalarField rhoInf(IOobject("rho", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), mesh, rhoInfValue)
dimensioned< scalar > mag(const dimensioned< Type > &)
const volScalarField & rho() const
Return the film density [kg/m^3].
const scalar deltaMin_
Minimum film height for model to be active.
Specie phase change model abstract base class.
const fluidThermo & primaryThermo() const
Return const reference to the primary region thermo object.
defineTypeNameAndDebug(kinematicSingleLayer, 0)
Thermodynamic form of single-cell layer surface film model.
Namespace for OpenFOAM.