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  Copyright (C) 2020-2023 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 "waxSolventEvaporation.H"
31 #include "thermoSingleLayer.H"
32 #include "zeroField.H"
33 
34 #include "fvmDdt.H"
35 #include "fvmDiv.H"
36 #include "fvcDiv.H"
37 #include "fvmSup.H"
38 
39 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
40 
41 namespace Foam
42 {
43 namespace regionModels
44 {
45 namespace surfaceFilmModels
46 {
47 
48 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
49 
50 defineTypeNameAndDebug(waxSolventEvaporation, 0);
51 
53 (
54  phaseChangeModel,
57 );
58 
59 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
60 
62 (
63  const scalar Re,
64  const scalar Sc
65 ) const
66 {
67  if (Re < 5.0E+05)
68  {
69  return 0.664*sqrt(Re)*cbrt(Sc);
70  }
71  else
72  {
73  return 0.037*pow(Re, 0.8)*cbrt(Sc);
74  }
75 }
76 
77 
78 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
79 
81 (
82  surfaceFilmRegionModel& film,
83  const dictionary& dict
84 )
85 :
86  phaseChangeModel(typeName, film, dict),
87  Wwax_
88  (
89  IOobject
90  (
91  IOobject::scopedName(typeName, "Wwax"),
92  film.regionMesh().time().constant(),
93  film.regionMesh().thisDb(),
94  IOobject::NO_READ,
95  IOobject::NO_WRITE,
96  IOobject::REGISTER
97  ),
98  coeffDict_.get<scalar>("Wwax")
99  ),
100  Wsolvent_
101  (
102  IOobject
103  (
104  IOobject::scopedName(typeName, "Wsolvent"),
105  film.regionMesh().time().constant(),
106  film.regionMesh().thisDb(),
107  IOobject::NO_READ,
108  IOobject::NO_WRITE,
109  IOobject::REGISTER
110  ),
111  coeffDict_.get<scalar>("Wsolvent")
112  ),
113  Ysolvent0_
114  (
115  IOobject
116  (
117  IOobject::scopedName(typeName, "Ysolvent0"),
118  film.regionMesh().time().constant(),
119  film.regionMesh().thisDb(),
120  IOobject::MUST_READ,
121  IOobject::NO_WRITE,
122  IOobject::REGISTER
123  )
124  ),
125  Ysolvent_
126  (
127  IOobject
128  (
129  IOobject::scopedName(typeName, "Ysolvent"),
130  film.regionMesh().time().timeName(),
131  film.regionMesh().thisDb(),
132  IOobject::MUST_READ,
133  IOobject::AUTO_WRITE,
134  IOobject::REGISTER
135  ),
136  film.regionMesh()
137  ),
138  deltaMin_(coeffDict_.get<scalar>("deltaMin")),
139  L_(coeffDict_.get<scalar>("L")),
140  TbFactor_(coeffDict_.getOrDefault<scalar>("TbFactor", 1.1)),
141  YInfZero_(coeffDict_.getOrDefault("YInfZero", false)),
142  activityCoeff_
143  (
144  Function1<scalar>::New("activityCoeff", coeffDict_, &film.regionMesh())
145  )
146 {}
147 
148 
149 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
150 
152 {}
153 
154 
155 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
156 
157 template<class YInfType>
159 (
160  const scalar dt,
161  scalarField& availableMass,
162  scalarField& dMass,
163  scalarField& dEnergy,
164  const YInfType& YInf
165 )
166 {
167  const thermoSingleLayer& film = filmType<thermoSingleLayer>();
168 
169  const volScalarField& delta = film.delta();
170  const volScalarField& deltaRho = film.deltaRho();
171  const surfaceScalarField& phi = film.phi();
172 
173  // Set local thermo properties
174  const SLGThermo& thermo = film.thermo();
175  const filmThermoModel& filmThermo = film.filmThermo();
176  const label vapId = thermo.carrierId(filmThermo.name());
177 
178  // Retrieve fields from film model
179  const scalarField& pInf = film.pPrimary();
180  const scalarField& T = film.T();
181  const scalarField& hs = film.hs();
182  const scalarField& rho = film.rho();
183  const scalarField& rhoInf = film.rhoPrimary();
184  const scalarField& muInf = film.muPrimary();
185  const scalarField& magSf = film.magSf();
186  const vectorField dU(film.UPrimary() - film.Us());
187  const scalarField limMass
188  (
189  max(scalar(0), availableMass - deltaMin_*rho*magSf)
190  );
191 
192  // Molecular weight of vapour [kg/kmol]
193  const scalar Wvap = thermo.carrier().W(vapId);
194 
195  const scalar Wwax = Wwax_.value();
196  const scalar Wsolvent = Wsolvent_.value();
197 
198  auto tevapRateCoeff = volScalarField::Internal::New
199  (
200  IOobject::scopedName(typeName, "evapRateCoeff"),
201  film.regionMesh(),
203  );
204  auto& evapRateCoeff = tevapRateCoeff.ref();
205 
206  auto tevapRateInf = volScalarField::Internal::New
207  (
208  IOobject::scopedName(typeName, "evapRateInf"),
209  film.regionMesh(),
211  );
212  auto& evapRateInf = tevapRateInf.ref();
213 
214  bool filmPresent = false;
215 
216  forAll(dMass, celli)
217  {
218  if (delta[celli] > deltaMin_)
219  {
220  filmPresent = true;
221 
222  const scalar Ysolvent = Ysolvent_[celli];
223 
224  // Molefraction of solvent in liquid film
225  const scalar Xsolvent
226  (
227  Ysolvent*Wsolvent/((1 - Ysolvent)*Wwax + Ysolvent*Wsolvent)
228  );
229 
230  // Primary region density [kg/m3]
231  const scalar rhoInfc = rhoInf[celli];
232 
233  // Cell pressure [Pa]
234  const scalar pc = pInf[celli];
235 
236  // Calculate the boiling temperature
237  const scalar Tb = filmThermo.Tb(pc);
238 
239  // Local temperature - impose lower limit of 200 K for stability
240  const scalar Tloc = min(TbFactor_*Tb, max(200.0, T[celli]));
241 
242  const scalar pPartialCoeff
243  (
244  filmThermo.pv(pc, Tloc)*activityCoeff_->value(Xsolvent)
245  );
246 
247  scalar XsCoeff = pPartialCoeff/pc;
248 
249  // Vapour phase mole fraction of solvent at interface
250  scalar Xs = XsCoeff*Xsolvent;
251 
252  if (Xs > 1)
253  {
255  << "Solvent vapour pressure > ambient pressure"
256  << endl;
257 
258  XsCoeff /= Xs;
259  Xs = 1;
260  }
261 
262  // Vapour phase mass fraction of solvent at the interface
263  const scalar YsCoeff
264  (
265  XsCoeff/(XsCoeff*Xsolvent*Wsolvent + (1 - Xs)*Wvap)
266  );
267 
268  // Primary region viscosity [Pa.s]
269  const scalar muInfc = muInf[celli];
270 
271  // Reynolds number
272  const scalar Re = rhoInfc*mag(dU[celli])*L_/muInfc;
273 
274  // Vapour diffusivity [m2/s]
275  const scalar Dab = filmThermo.D(pc, Tloc);
276 
277  // Schmidt number
278  const scalar Sc = muInfc/(rhoInfc*(Dab + ROOTVSMALL));
279 
280  // Sherwood number
281  const scalar Sh = this->Sh(Re, Sc);
282 
283  // Mass transfer coefficient [m/s]
284  evapRateCoeff[celli] = rhoInfc*Sh*Dab/(L_ + ROOTVSMALL);
285 
286  // Solvent mass transfer
287  const scalar dm
288  (
289  max
290  (
291  dt*magSf[celli]
292  *evapRateCoeff[celli]*(YsCoeff*Ysolvent - YInf[celli]),
293  0
294  )
295  );
296 
297  if (dm > limMass[celli])
298  {
299  evapRateCoeff[celli] *= limMass[celli]/dm;
300  }
301 
302  evapRateInf[celli] = evapRateCoeff[celli]*YInf[celli];
303  evapRateCoeff[celli] *= YsCoeff;
304 
305  // hVap[celli] = filmThermo.hl(pc, Tloc);
306  }
307  }
308 
309  const dimensionedScalar deltaRho0Bydt
310  (
311  "deltaRho0",
312  deltaRho.dimensions()/dimTime,
313  ROOTVSMALL/dt
314  );
315 
316  volScalarField::Internal impingementRate
317  (
318  max
319  (
320  -film.rhoSp()(),
321  dimensionedScalar(film.rhoSp().dimensions(), Zero)
322  )
323  );
324 
325  if (filmPresent)
326  {
327  // Solve for the solvent mass fraction
328  fvScalarMatrix YsolventEqn
329  (
330  fvm::ddt(deltaRho, Ysolvent_)
332  ==
333  deltaRho0Bydt*Ysolvent_()
334 
335  + evapRateInf
336 
337  // Include the effect of the impinging droplets
338  // added with Ysolvent = Ysolvent0
339  + impingementRate*Ysolvent0_
340 
341  - fvm::Sp
342  (
343  deltaRho0Bydt
344  + evapRateCoeff
345  + film.rhoSp()()
346  + impingementRate,
347  Ysolvent_
348  )
349  );
350 
351  YsolventEqn.relax();
352  YsolventEqn.solve();
353 
354  Ysolvent_.clamp_range(zero_one{});
355 
356  scalarField dm
357  (
358  dt*magSf*rhoInf*(evapRateCoeff*Ysolvent_ + evapRateInf)
359  );
360 
361  dMass += dm;
362 
363  // Heat is assumed to be removed by heat-transfer to the wall
364  // so the energy remains unchanged by the phase-change.
365  dEnergy += dm*hs;
366 
367  // Latent heat [J/kg]
368  // dEnergy += dm*(hs[celli] + hVap);
369  }
370 }
371 
372 
374 (
375  const scalar dt,
376  scalarField& availableMass,
377  scalarField& dMass,
378  scalarField& dEnergy
379 )
380 {
381  if (YInfZero_)
382  {
383  correctModel(dt, availableMass, dMass, dEnergy, zeroField());
384  }
385  else
386  {
387  const thermoSingleLayer& film = filmType<thermoSingleLayer>();
388  const label vapId = film.thermo().carrierId(film.filmThermo().name());
389  const scalarField& YInf = film.YPrimary()[vapId];
390 
391  correctModel(dt, availableMass, dMass, dEnergy, YInf);
392  }
393 }
394 
395 
396 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
397 
398 } // End namespace surfaceFilmModels
399 } // End namespace regionModels
400 } // End namespace Foam
401 
402 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
scalar delta
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
void clamp_range(const dimensioned< MinMax< Type >> &range)
Clamp field values (in-place) to the specified range.
const Type & value() const noexcept
Return const reference to value.
dictionary dict
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:37
virtual const word & name() const =0
Return the specie name.
virtual const volVectorField & Us() const =0
Return the film surface velocity [m/s].
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
virtual const volScalarField & hs() const =0
Return the film surface temperature [J/kg].
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
autoPtr< Function1< scalar > > activityCoeff_
Activity coefficient as a function of solvent mole fraction.
uniformDimensionedScalarField Wwax_
Molecular weight of wax [kg/kmol].
dimensionedScalar sqrt(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
vtkSmartPointer< vtkFloatArray > zeroField(const word &name, const label size)
Create named field initialized to zero.
const scalar deltaMin_
Minimum film height for model to be active.
virtual const volScalarField & rho() const =0
Return the film density [kg/m3].
addToRunTimeSelectionTable(surfaceFilmRegionModel, kinematicSingleLayer, mesh)
Macros for easy insertion into run-time selection tables.
waxSolventEvaporation(const waxSolventEvaporation &)=delete
No copy construct.
const surfaceFilmRegionModel & film() const
Return const access to the film surface film model.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
word timeName
Definition: getTimeIndex.H:3
static word scopedName(const std::string &scope, const word &name)
Create scope:name or scope_name string.
Definition: IOobjectI.H:40
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensionedScalar cbrt(const dimensionedScalar &ds)
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:41
Calculate the matrix for the first temporal derivative.
virtual const volScalarField & delta() const =0
Return the film thickness [m].
const fvMesh & regionMesh() const
Return the region mesh database.
Definition: regionModelI.H:26
Thermo package for (S)olids (L)iquids and (G)ases Takes reference to thermo package, and provides:
Definition: SLGThermo.H:60
const word & name() const
Name function is needed to disambiguate those inherited from regIOobject and dictionary.
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
zeroField Sp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
virtual scalar D(const scalar p, const scalar T) const =0
Return diffusivity [m2/s].
virtual scalar Tb(const scalar p) const =0
Return boiling temperature [K].
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
Calculate the divergence of the given field.
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:41
uniformDimensionedScalarField Wsolvent_
Molecular weight of liquid [kg/kmol].
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
const dimensionSet dimDensity
scalarField Re(const UList< complex > &cmplx)
Extract real component.
Definition: complexField.C:207
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Calculate the matrix for the divergence of the given field and flux.
virtual const volScalarField & magSf() const
Return the face area magnitudes / [m2].
#define WarningInFunction
Report a warning using Foam::Warning.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
static tmp< DimensionedField< Type, GeoMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const Field< Type > &iField)
Return tmp field (NO_READ, NO_WRITE) from name, mesh, dimensions, copy of internal field...
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
scalar Sh(const scalar Re, const scalar Sc) const
Return Sherwood number as a function of Reynolds and Schmidt numbers.
uniformDimensionedScalarField Ysolvent0_
Initial solvent mass-fraction.
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:180
defineTypeNameAndDebug(kinematicSingleLayer, 0)
void correctModel(const scalar dt, scalarField &availableMass, scalarField &dMass, scalarField &dEnergy, const YInfType &YInf)
Calculate the finiteVolume matrix for implicit and explicit sources.
virtual scalar pv(const scalar p, const scalar T) const =0
Return vapour pressure [Pa].
Thermodynamic form of single-cell layer surface film model.
Namespace for OpenFOAM.
const dimensionSet & dimensions() const noexcept
Return dimensions.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127
const dimensionSet dimVelocity
virtual const volScalarField & T() const =0
Return the film mean temperature [K].