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 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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 
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  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<Switch>("YInfZero", false))
86 {}
87 
88 
89 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
90 
91 template<class YInfType>
93 (
94  const scalar dt,
95  scalarField& availableMass,
96  scalarField& dMass,
97  scalarField& dEnergy,
98  YInfType YInf
99 )
100 {
101  const thermoSingleLayer& film = filmType<thermoSingleLayer>();
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 // ************************************************************************* //
Foam::regionModels::surfaceFilmModels::filmThermoModel::Cp
virtual scalar Cp(const scalar p, const scalar T) const =0
Return specific heat capacity [J/kg/K].
Foam::Switch
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.H:77
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::hs
virtual const volScalarField & hs() const
Return the film sensible enthalpy [J/kg].
Definition: thermoSingleLayer.C:680
Foam::regionModels::surfaceFilmModels::filmThermoModel::hl
virtual scalar hl(const scalar p, const scalar T) const =0
Return latent heat [J/kg].
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::YPrimary
const PtrList< volScalarField > & YPrimary() const
Specie mass fractions [0-1].
Definition: thermoSingleLayerI.H:130
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::Us
virtual const volVectorField & Us() const
Return the film surface velocity [m/s].
Definition: kinematicSingleLayer.C:958
Foam::SLGThermo
Thermo package for (S)olids (L)iquids and (G)ases Takes reference to thermo package,...
Definition: SLGThermo.H:64
Foam::regionModels::surfaceFilmModels::filmThermoModel::W
virtual scalar W() const =0
Return molecular weight [kg/kmol].
Foam::regionModels::surfaceFilmModels::filmThermoModel::name
virtual const word & name() const =0
Return the specie name.
Foam::vtk::Tools::zeroField
vtkSmartPointer< vtkFloatArray > zeroField(const word &name, const label size)
Create named field initialized to zero.
Definition: foamVtkToolsTemplates.C:327
standardPhaseChange.H
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::thermo
const SLGThermo & thermo() const
Return const reference to the SLGThermo object.
Definition: thermoSingleLayerI.H:44
thermo
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::UPrimary
const volVectorField & UPrimary() const
Velocity [m/s].
Definition: kinematicSingleLayerI.H:152
Foam::regionModels::surfaceFilmModels::addToRunTimeSelectionTable
addToRunTimeSelectionTable(surfaceFilmRegionModel, kinematicSingleLayer, mesh)
rho
rho
Definition: readInitialConditions.H:88
Foam::regionModels::surfaceFilmModels::standardPhaseChange::standardPhaseChange
standardPhaseChange(const standardPhaseChange &)=delete
No copy construct.
Foam::regionModels::surfaceFilmModels::filmThermoModel::Tb
virtual scalar Tb(const scalar p) const =0
Return boiling temperature [K].
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
thermoSingleLayer.H
zeroField.H
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::rhoPrimary
const volScalarField & rhoPrimary() const
Density [kg/m3].
Definition: kinematicSingleLayerI.H:164
Foam::regionModels::surfaceFilmModels::thermoSingleLayer
Thermodynamic form of single-cell layer surface film model.
Definition: thermoSingleLayer.H:67
Foam::Field< scalar >
Foam::regionModels::surfaceFilmModels::surfaceFilmRegionModel
Base class for surface film models.
Definition: surfaceFilmRegionModel.H:55
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::filmThermo
const filmThermoModel & filmThermo() const
Film thermo.
Definition: kinematicSingleLayerI.H:176
Foam::regionModels::surfaceFilmModels::standardPhaseChange::Sh
scalar Sh(const scalar Re, const scalar Sc) const
Return Sherwood number as a function of Reynolds and Schmidt numbers.
Definition: standardPhaseChange.C:57
Foam::regionModels::singleLayerRegion::magSf
virtual const volScalarField & magSf() const
Return the face area magnitudes / [m2].
Definition: singleLayerRegion.C:223
delta
scalar delta
Definition: LISASMDCalcMethod2.H:8
Foam::T
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Definition: FieldFieldFunctions.C:58
Foam::regionModels::surfaceFilmModels::standardPhaseChange::correctModel
void correctModel(const scalar dt, scalarField &availableMass, scalarField &dMass, scalarField &dEnergy, YInfType YInf)
Definition: standardPhaseChange.C:93
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::pPrimary
const volScalarField & pPrimary() const
Pressure [Pa].
Definition: kinematicSingleLayerI.H:158
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
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::T
virtual const volScalarField & T() const
Return the film mean temperature [K].
Definition: thermoSingleLayer.C:662
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::muPrimary
const volScalarField & muPrimary() const
Viscosity [Pa.s].
Definition: kinematicSingleLayerI.H:170
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::regionModels::surfaceFilmModels::filmThermoModel::pv
virtual scalar pv(const scalar p, const scalar T) const =0
Return vapour pressure [Pa].
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::regionModels::surfaceFilmModels::filmThermoModel::D
virtual scalar D(const scalar p, const scalar T) const =0
Return diffusivity [m2/s].
Foam::regionModels::surfaceFilmModels::defineTypeNameAndDebug
defineTypeNameAndDebug(kinematicSingleLayer, 0)
Foam::regionModels::surfaceFilmModels::filmThermoModel
Base class for film thermo models.
Definition: filmThermoModel.H:56
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::delta
const volScalarField & delta() const
Return const access to the film thickness [m].
Definition: kinematicSingleLayerI.H:86
Foam::SLGThermo::carrierId
label carrierId(const word &cmptName, bool allowNotFound=false) const
Index of carrier component.
Definition: SLGThermo.C:151
Foam::Re
scalarField Re(const UList< complex > &cf)
Extract real component.
Definition: complexField.C:159
Cp
const volScalarField & Cp
Definition: EEqn.H:7
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::rho
virtual const volScalarField & rho() const
Return the film density [kg/m3].
Definition: kinematicSingleLayer.C:982
Foam::cbrt
dimensionedScalar cbrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:155
Foam::regionModels::surfaceFilmModels::phaseChangeModel
Base class for surface film phase change models.
Definition: phaseChangeModel.H:57