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 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  typeName + ":Wwax",
92  film.regionMesh().time().constant(),
93  film.regionMesh().thisDb()
94  ),
95  coeffDict_.get<scalar>("Wwax")
96  ),
97  Wsolvent_
98  (
99  IOobject
100  (
101  typeName + ":Wsolvent",
102  film.regionMesh().time().constant(),
103  film.regionMesh().thisDb()
104  ),
105  coeffDict_.get<scalar>("Wsolvent")
106  ),
107  Ysolvent0_
108  (
109  IOobject
110  (
111  typeName + ":Ysolvent0",
112  film.regionMesh().time().constant(),
113  film.regionMesh().thisDb(),
114  IOobject::MUST_READ,
115  IOobject::NO_WRITE
116  )
117  ),
118  Ysolvent_
119  (
120  IOobject
121  (
122  typeName + ":Ysolvent",
123  film.regionMesh().time().timeName(),
124  film.regionMesh().thisDb(),
125  IOobject::MUST_READ,
126  IOobject::AUTO_WRITE
127  ),
128  film.regionMesh()
129  ),
130  deltaMin_(coeffDict_.get<scalar>("deltaMin")),
131  L_(coeffDict_.get<scalar>("L")),
132  TbFactor_(coeffDict_.getOrDefault<scalar>("TbFactor", 1.1)),
133  YInfZero_(coeffDict_.getOrDefault("YInfZero", false)),
134  activityCoeff_
135  (
136  Function1<scalar>::New("activityCoeff", coeffDict_, &film.regionMesh())
137  )
138 {}
139 
140 
141 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
142 
144 {}
145 
146 
147 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
148 
149 template<class YInfType>
151 (
152  const scalar dt,
153  scalarField& availableMass,
154  scalarField& dMass,
155  scalarField& dEnergy,
156  const YInfType& YInf
157 )
158 {
159  const thermoSingleLayer& film = filmType<thermoSingleLayer>();
160 
161  const volScalarField& delta = film.delta();
162  const volScalarField& deltaRho = film.deltaRho();
163  const surfaceScalarField& phi = film.phi();
164 
165  // Set local thermo properties
166  const SLGThermo& thermo = film.thermo();
167  const filmThermoModel& filmThermo = film.filmThermo();
168  const label vapId = thermo.carrierId(filmThermo.name());
169 
170  // Retrieve fields from film model
171  const scalarField& pInf = film.pPrimary();
172  const scalarField& T = film.T();
173  const scalarField& hs = film.hs();
174  const scalarField& rho = film.rho();
175  const scalarField& rhoInf = film.rhoPrimary();
176  const scalarField& muInf = film.muPrimary();
177  const scalarField& magSf = film.magSf();
178  const vectorField dU(film.UPrimary() - film.Us());
179  const scalarField limMass
180  (
181  max(scalar(0), availableMass - deltaMin_*rho*magSf)
182  );
183 
184  // Molecular weight of vapour [kg/kmol]
185  const scalar Wvap = thermo.carrier().W(vapId);
186 
187  const scalar Wwax = Wwax_.value();
188  const scalar Wsolvent = Wsolvent_.value();
189 
190  volScalarField::Internal evapRateCoeff
191  (
192  IOobject
193  (
194  typeName + ":evapRateCoeff",
196  film.regionMesh().thisDb(),
200  ),
201  film.regionMesh(),
203  );
204 
205  volScalarField::Internal evapRateInf
206  (
207  IOobject
208  (
209  typeName + ":evapRateInf",
211  film.regionMesh().thisDb(),
215  ),
216  film.regionMesh(),
218  );
219 
220  bool filmPresent = false;
221 
222  forAll(dMass, celli)
223  {
224  if (delta[celli] > deltaMin_)
225  {
226  filmPresent = true;
227 
228  const scalar Ysolvent = Ysolvent_[celli];
229 
230  // Molefraction of solvent in liquid film
231  const scalar Xsolvent
232  (
233  Ysolvent*Wsolvent/((1 - Ysolvent)*Wwax + Ysolvent*Wsolvent)
234  );
235 
236  // Primary region density [kg/m3]
237  const scalar rhoInfc = rhoInf[celli];
238 
239  // Cell pressure [Pa]
240  const scalar pc = pInf[celli];
241 
242  // Calculate the boiling temperature
243  const scalar Tb = filmThermo.Tb(pc);
244 
245  // Local temperature - impose lower limit of 200 K for stability
246  const scalar Tloc = min(TbFactor_*Tb, max(200.0, T[celli]));
247 
248  const scalar pPartialCoeff
249  (
250  filmThermo.pv(pc, Tloc)*activityCoeff_->value(Xsolvent)
251  );
252 
253  scalar XsCoeff = pPartialCoeff/pc;
254 
255  // Vapour phase mole fraction of solvent at interface
256  scalar Xs = XsCoeff*Xsolvent;
257 
258  if (Xs > 1)
259  {
261  << "Solvent vapour pressure > ambient pressure"
262  << endl;
263 
264  XsCoeff /= Xs;
265  Xs = 1;
266  }
267 
268  // Vapour phase mass fraction of solvent at the interface
269  const scalar YsCoeff
270  (
271  XsCoeff/(XsCoeff*Xsolvent*Wsolvent + (1 - Xs)*Wvap)
272  );
273 
274  // Primary region viscosity [Pa.s]
275  const scalar muInfc = muInf[celli];
276 
277  // Reynolds number
278  const scalar Re = rhoInfc*mag(dU[celli])*L_/muInfc;
279 
280  // Vapour diffusivity [m2/s]
281  const scalar Dab = filmThermo.D(pc, Tloc);
282 
283  // Schmidt number
284  const scalar Sc = muInfc/(rhoInfc*(Dab + ROOTVSMALL));
285 
286  // Sherwood number
287  const scalar Sh = this->Sh(Re, Sc);
288 
289  // Mass transfer coefficient [m/s]
290  evapRateCoeff[celli] = rhoInfc*Sh*Dab/(L_ + ROOTVSMALL);
291 
292  // Solvent mass transfer
293  const scalar dm
294  (
295  max
296  (
297  dt*magSf[celli]
298  *evapRateCoeff[celli]*(YsCoeff*Ysolvent - YInf[celli]),
299  0
300  )
301  );
302 
303  if (dm > limMass[celli])
304  {
305  evapRateCoeff[celli] *= limMass[celli]/dm;
306  }
307 
308  evapRateInf[celli] = evapRateCoeff[celli]*YInf[celli];
309  evapRateCoeff[celli] *= YsCoeff;
310 
311  // hVap[celli] = filmThermo.hl(pc, Tloc);
312  }
313  }
314 
315  const dimensionedScalar deltaRho0Bydt
316  (
317  "deltaRho0",
318  deltaRho.dimensions()/dimTime,
319  ROOTVSMALL/dt
320  );
321 
322  volScalarField::Internal impingementRate
323  (
324  max
325  (
326  -film.rhoSp()(),
327  dimensionedScalar(film.rhoSp().dimensions(), Zero)
328  )
329  );
330 
331  if (filmPresent)
332  {
333  // Solve for the solvent mass fraction
334  fvScalarMatrix YsolventEqn
335  (
336  fvm::ddt(deltaRho, Ysolvent_)
338  ==
339  deltaRho0Bydt*Ysolvent_()
340 
341  + evapRateInf
342 
343  // Include the effect of the impinging droplets
344  // added with Ysolvent = Ysolvent0
345  + impingementRate*Ysolvent0_
346 
347  - fvm::Sp
348  (
349  deltaRho0Bydt
350  + evapRateCoeff
351  + film.rhoSp()()
352  + impingementRate,
353  Ysolvent_
354  )
355  );
356 
357  YsolventEqn.relax();
358  YsolventEqn.solve();
359 
360  Ysolvent_.clamp_range(zero_one{});
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 // ************************************************************************* //
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.
Ignore writing from objectRegistry::writeObject()
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:360
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.
virtual const objectRegistry & thisDb() const
Return the object registry - resolve conflict polyMesh/lduMesh.
Definition: fvMesh.H:376
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
word timeName
Definition: getTimeIndex.H:3
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.
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
Definition: Time.C:714
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.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: areaFieldsFwd.H:42
Nothing to be read.
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:172
defineTypeNameAndDebug(kinematicSingleLayer, 0)
void correctModel(const scalar dt, scalarField &availableMass, scalarField &dMass, scalarField &dEnergy, const YInfType &YInf)
Do not request registration (bool: false)
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].