waxSolventEvaporation.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) 2017 OpenFOAM Foundation
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 "waxSolventEvaporation.H"
30 #include "thermoSingleLayer.H"
31 #include "zeroField.H"
32 
33 #include "fvmDdt.H"
34 #include "fvmDiv.H"
35 #include "fvcDiv.H"
36 #include "fvmSup.H"
37 
38 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
39 
40 namespace Foam
41 {
42 namespace regionModels
43 {
44 namespace surfaceFilmModels
45 {
46 
47 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
48 
50 
52 (
53  phaseChangeModel,
55  dictionary
56 );
57 
58 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
59 
61 (
62  const scalar Re,
63  const scalar Sc
64 ) const
65 {
66  if (Re < 5.0E+05)
67  {
68  return 0.664*sqrt(Re)*cbrt(Sc);
69  }
70  else
71  {
72  return 0.037*pow(Re, 0.8)*cbrt(Sc);
73  }
74 }
75 
76 
77 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
78 
79 waxSolventEvaporation::waxSolventEvaporation
80 (
82  const dictionary& dict
83 )
84 :
85  phaseChangeModel(typeName, film, dict),
86  Wwax_
87  (
88  IOobject
89  (
90  typeName + ":Wwax",
91  film.regionMesh().time().constant(),
92  film.regionMesh()
93  ),
94  coeffDict_.get<scalar>("Wwax")
95  ),
96  Wsolvent_
97  (
98  IOobject
99  (
100  typeName + ":Wsolvent",
101  film.regionMesh().time().constant(),
102  film.regionMesh()
103  ),
104  coeffDict_.get<scalar>("Wsolvent")
105  ),
106  Ysolvent0_
107  (
108  IOobject
109  (
110  typeName + ":Ysolvent0",
111  film.regionMesh().time().constant(),
112  film.regionMesh(),
115  )
116  ),
117  Ysolvent_
118  (
119  IOobject
120  (
121  typeName + ":Ysolvent",
122  film.regionMesh().time().timeName(),
123  film.regionMesh(),
126  ),
127  film.regionMesh()
128  ),
129  deltaMin_(coeffDict_.get<scalar>("deltaMin")),
130  L_(coeffDict_.get<scalar>("L")),
131  TbFactor_(coeffDict_.lookupOrDefault<scalar>("TbFactor", 1.1)),
132  YInfZero_(coeffDict_.lookupOrDefault("YInfZero", false)),
133  activityCoeff_
134  (
135  Function1<scalar>::New("activityCoeff", coeffDict_)
136  )
137 {}
138 
139 
140 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
141 
143 {}
144 
145 
146 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
147 
148 template<class YInfType>
150 (
151  const scalar dt,
152  scalarField& availableMass,
153  scalarField& dMass,
154  scalarField& dEnergy,
155  YInfType YInf
156 )
157 {
158  const thermoSingleLayer& film = filmType<thermoSingleLayer>();
159 
160  const volScalarField& delta = film.delta();
161  const volScalarField& deltaRho = film.deltaRho();
162  const surfaceScalarField& phi = film.phi();
163 
164  // Set local thermo properties
165  const SLGThermo& thermo = film.thermo();
166  const filmThermoModel& filmThermo = film.filmThermo();
167  const label vapId = thermo.carrierId(filmThermo.name());
168 
169  // Retrieve fields from film model
170  const scalarField& pInf = film.pPrimary();
171  const scalarField& T = film.T();
172  const scalarField& hs = film.hs();
173  const scalarField& rho = film.rho();
174  const scalarField& rhoInf = film.rhoPrimary();
175  const scalarField& muInf = film.muPrimary();
176  const scalarField& magSf = film.magSf();
177  const vectorField dU(film.UPrimary() - film.Us());
178  const scalarField limMass
179  (
180  max(scalar(0), availableMass - deltaMin_*rho*magSf)
181  );
182 
183  // Molecular weight of vapour [kg/kmol]
184  const scalar Wvap = thermo.carrier().W(vapId);
185 
186  const scalar Wwax = Wwax_.value();
187  const scalar Wsolvent = Wsolvent_.value();
188 
189  volScalarField::Internal evapRateCoeff
190  (
191  IOobject
192  (
193  typeName + ":evapRateCoeff",
194  film.regionMesh().time().timeName(),
195  film.regionMesh(),
198  false
199  ),
200  film.regionMesh(),
202  );
203 
204  volScalarField::Internal evapRateInf
205  (
206  IOobject
207  (
208  typeName + ":evapRateInf",
209  film.regionMesh().time().timeName(),
210  film.regionMesh(),
213  false
214  ),
215  film.regionMesh(),
217  );
218 
219  bool filmPresent = false;
220 
221  forAll(dMass, celli)
222  {
223  if (delta[celli] > deltaMin_)
224  {
225  filmPresent = true;
226 
227  const scalar Ysolvent = Ysolvent_[celli];
228 
229  // Molefraction of solvent in liquid film
230  const scalar Xsolvent
231  (
232  Ysolvent*Wsolvent/((1 - Ysolvent)*Wwax + Ysolvent*Wsolvent)
233  );
234 
235  // Primary region density [kg/m3]
236  const scalar rhoInfc = rhoInf[celli];
237 
238  // Cell pressure [Pa]
239  const scalar pc = pInf[celli];
240 
241  // Calculate the boiling temperature
242  const scalar Tb = filmThermo.Tb(pc);
243 
244  // Local temperature - impose lower limit of 200 K for stability
245  const scalar Tloc = min(TbFactor_*Tb, max(200.0, T[celli]));
246 
247  const scalar pPartialCoeff
248  (
249  filmThermo.pv(pc, Tloc)*activityCoeff_->value(Xsolvent)
250  );
251 
252  scalar XsCoeff = pPartialCoeff/pc;
253 
254  // Vapour phase mole fraction of solvent at interface
255  scalar Xs = XsCoeff*Xsolvent;
256 
257  if (Xs > 1)
258  {
260  << "Solvent vapour pressure > ambient pressure"
261  << endl;
262 
263  XsCoeff /= Xs;
264  Xs = 1;
265  }
266 
267  // Vapour phase mass fraction of solvent at the interface
268  const scalar YsCoeff
269  (
270  XsCoeff/(XsCoeff*Xsolvent*Wsolvent + (1 - Xs)*Wvap)
271  );
272 
273  // Primary region viscosity [Pa.s]
274  const scalar muInfc = muInf[celli];
275 
276  // Reynolds number
277  const scalar Re = rhoInfc*mag(dU[celli])*L_/muInfc;
278 
279  // Vapour diffusivity [m2/s]
280  const scalar Dab = filmThermo.D(pc, Tloc);
281 
282  // Schmidt number
283  const scalar Sc = muInfc/(rhoInfc*(Dab + ROOTVSMALL));
284 
285  // Sherwood number
286  const scalar Sh = this->Sh(Re, Sc);
287 
288  // Mass transfer coefficient [m/s]
289  evapRateCoeff[celli] = rhoInfc*Sh*Dab/(L_ + ROOTVSMALL);
290 
291  // Solvent mass transfer
292  const scalar dm
293  (
294  max
295  (
296  dt*magSf[celli]
297  *evapRateCoeff[celli]*(YsCoeff*Ysolvent - YInf[celli]),
298  0
299  )
300  );
301 
302  if (dm > limMass[celli])
303  {
304  evapRateCoeff[celli] *= limMass[celli]/dm;
305  }
306 
307  evapRateInf[celli] = evapRateCoeff[celli]*YInf[celli];
308  evapRateCoeff[celli] *= YsCoeff;
309 
310  // hVap[celli] = filmThermo.hl(pc, Tloc);
311  }
312  }
313 
314  const dimensionedScalar deltaRho0Bydt
315  (
316  "deltaRho0",
317  deltaRho.dimensions()/dimTime,
318  ROOTVSMALL/dt
319  );
320 
321  volScalarField::Internal impingementRate
322  (
323  max
324  (
325  -film.rhoSp()(),
326  dimensionedScalar(film.rhoSp().dimensions(), Zero)
327  )
328  );
329 
330  if (filmPresent)
331  {
332  // Solve for the solvent mass fraction
333  fvScalarMatrix YsolventEqn
334  (
335  fvm::ddt(deltaRho, Ysolvent_)
336  + fvm::div(phi, Ysolvent_)
337  ==
338  deltaRho0Bydt*Ysolvent_()
339 
340  + evapRateInf
341 
342  // Include the effect of the impinging droplets
343  // added with Ysolvent = Ysolvent0
344  + impingementRate*Ysolvent0_
345 
346  - fvm::Sp
347  (
348  deltaRho0Bydt
349  + evapRateCoeff
350  + film.rhoSp()()
351  + impingementRate,
352  Ysolvent_
353  )
354  );
355 
356  YsolventEqn.relax();
357  YsolventEqn.solve();
358 
359  Ysolvent_.min(1);
360  Ysolvent_.max(0);
361 
362  scalarField dm
363  (
364  dt*magSf*rhoInf*(evapRateCoeff*Ysolvent_ + evapRateInf)
365  );
366 
367  dMass += dm;
368 
369  // Heat is assumed to be removed by heat-transfer to the wall
370  // so the energy remains unchanged by the phase-change.
371  dEnergy += dm*hs;
372 
373  // Latent heat [J/kg]
374  // dEnergy += dm*(hs[celli] + hVap);
375  }
376 }
377 
378 
380 (
381  const scalar dt,
382  scalarField& availableMass,
383  scalarField& dMass,
384  scalarField& dEnergy
385 )
386 {
387  if (YInfZero_)
388  {
389  correctModel(dt, availableMass, dMass, dEnergy, zeroField());
390  }
391  else
392  {
393  const thermoSingleLayer& film = filmType<thermoSingleLayer>();
394  const label vapId = film.thermo().carrierId(film.filmThermo().name());
395  const scalarField& YInf = film.YPrimary()[vapId];
396 
397  correctModel(dt, availableMass, dMass, dEnergy, YInf);
398  }
399 }
400 
401 
402 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
403 
404 } // End namespace surfaceFilmModels
405 } // End namespace regionModels
406 } // End namespace Foam
407 
408 // ************************************************************************* //
Foam::IOobject::NO_WRITE
Definition: IOobject.H:130
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
Foam::IOobject::AUTO_WRITE
Definition: IOobject.H:129
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::hs
virtual const volScalarField & hs() const
Return the film sensible enthalpy [J/kg].
Definition: thermoSingleLayer.C:703
waxSolventEvaporation
Wax solvent mixture evaporation model.
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:987
Foam::SLGThermo
Thermo package for (S)olids (L)iquids and (G)ases Takes reference to thermo package,...
Definition: SLGThermo.H:64
Foam::Zero
static constexpr const zero Zero
Global zero.
Definition: zero.H:128
Foam::dimVelocity
const dimensionSet dimVelocity
Foam::dimDensity
const dimensionSet dimDensity
Foam::regionModels::surfaceFilmModels::waxSolventEvaporation::Sh
scalar Sh(const scalar Re, const scalar Sc) const
Return Sherwood number as a function of Reynolds and Schmidt numbers.
Definition: waxSolventEvaporation.C:61
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:258
Foam::Time::timeName
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:764
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::thermo
const SLGThermo & thermo() const
Return const reference to the SLGThermo object.
Definition: thermoSingleLayerI.H:44
fvcDiv.H
Calculate the divergence of the given field.
thermo
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::rhoSp
volScalarField & rhoSp()
Mass [kg/m2/s].
Definition: kinematicSingleLayerI.H:127
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:337
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::UPrimary
const volVectorField & UPrimary() const
Velocity [m/s].
Definition: kinematicSingleLayerI.H:151
Foam::regionModels::surfaceFilmModels::addToRunTimeSelectionTable
addToRunTimeSelectionTable(surfaceFilmRegionModel, kinematicSingleLayer, mesh)
rho
rho
Definition: readInitialConditions.H:96
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
Foam::Function1
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
Definition: Function1.H:56
fvmDiv.H
Calculate the matrix for the divergence of the given field and flux.
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
thermoSingleLayer.H
zeroField.H
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::rhoPrimary
const volScalarField & rhoPrimary() const
Density [kg/m3].
Definition: kinematicSingleLayerI.H:163
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:54
Foam::regionModels::surfaceFilmModels::thermoSingleLayer
Thermodynamic form of single-cell layer surface film model.
Definition: thermoSingleLayer.H:67
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:62
Foam::fvm::Sp
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
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:175
Foam::fvMatrix::relax
void relax(const scalar alpha)
Relax matrix (for steady-state solution).
Definition: fvMatrix.C:574
Foam::regionModels::singleLayerRegion::magSf
virtual const volScalarField & magSf() const
Return the face area magnitudes / [m2].
Definition: singleLayerRegion.C:222
Foam::regionModels::regionModel::regionMesh
const fvMesh & regionMesh() const
Return the region mesh database.
Definition: regionModelI.H:63
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::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:43
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::pPrimary
const volScalarField & pPrimary() const
Pressure [Pa].
Definition: kinematicSingleLayerI.H:157
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::kinematicSingleLayer::phi
virtual const surfaceScalarField & phi() const
Return the film flux [kg.m/s].
Definition: kinematicSingleLayer.C:1005
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::T
virtual const volScalarField & T() const
Return the film mean temperature [K].
Definition: thermoSingleLayer.C:685
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::muPrimary
const volScalarField & muPrimary() const
Viscosity [Pa.s].
Definition: kinematicSingleLayerI.H:169
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:121
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::dimensioned< scalar >
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].
fvmSup.H
Calculate the matrix for implicit and explicit sources.
Foam::regionModels::surfaceFilmModels::defineTypeNameAndDebug
defineTypeNameAndDebug(kinematicSingleLayer, 0)
Foam::fvm::ddt
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:48
Foam::regionModels::surfaceFilmModels::waxSolventEvaporation::correctModel
void correctModel(const scalar dt, scalarField &availableMass, scalarField &dMass, scalarField &dEnergy, YInfType YInf)
Definition: waxSolventEvaporation.C:150
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:85
Foam::SLGThermo::carrierId
label carrierId(const word &cmptName, bool allowNotFound=false) const
Index of carrier component.
Definition: SLGThermo.C:150
Foam::Re
scalarField Re(const UList< complex > &cf)
Extract real component.
Definition: complexField.C:159
Foam::regionModels::surfaceFilmModels::waxSolventEvaporation::~waxSolventEvaporation
virtual ~waxSolventEvaporation()
Destructor.
Definition: waxSolventEvaporation.C:142
Foam::fvMatrix
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvPatchField.H:76
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:246
waxSolventEvaporation.H
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::rho
virtual const volScalarField & rho() const
Return the film density [kg/m3].
Definition: kinematicSingleLayer.C:1011
Foam::cbrt
dimensionedScalar cbrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:155
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::deltaRho
virtual const volScalarField & deltaRho() const
Return the film thickness*density (helper field) [kg/m3].
Definition: kinematicSingleLayer.C:999
Foam::TimePaths::constant
const word & constant() const
Return constant name.
Definition: TimePathsI.H:88
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::IOobject::NO_READ
Definition: IOobject.H:123
Foam::fvm::div
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:48
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:294
fvmDdt.H
Calulate the matrix for the first temporal derivative.
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54
Foam::IOobject::MUST_READ
Definition: IOobject.H:120
Foam::regionModels::surfaceFilmModels::phaseChangeModel
Base class for surface film phase change models.
Definition: phaseChangeModel.H:57