LiquidEvaporationBoil.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-2016 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 "LiquidEvaporationBoil.H"
29 #include "specie.H"
30 #include "mathematicalConstants.H"
31 
32 using namespace Foam::constant::mathematical;
33 
34 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
35 
36 template<class CloudType>
38 (
39  const label celli
40 ) const
41 {
42  scalarField Xc(this->owner().thermo().carrier().Y().size());
43 
44  forAll(Xc, i)
45  {
46  Xc[i] =
47  this->owner().thermo().carrier().Y()[i][celli]
48  /this->owner().thermo().carrier().W(i);
49  }
50 
51  return Xc/sum(Xc);
52 }
53 
54 
55 template<class CloudType>
57 (
58  const scalar Re,
59  const scalar Sc
60 ) const
61 {
62  return 2.0 + 0.6*Foam::sqrt(Re)*cbrt(Sc);
63 }
64 
65 
66 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
67 
68 template<class CloudType>
70 (
71  const dictionary& dict,
72  CloudType& owner
73 )
74 :
75  PhaseChangeModel<CloudType>(dict, owner, typeName),
76  liquids_(owner.thermo().liquids()),
77  activeLiquids_(this->coeffDict().lookup("activeLiquids")),
78  liqToCarrierMap_(activeLiquids_.size(), -1),
79  liqToLiqMap_(activeLiquids_.size(), -1)
80 {
81  if (activeLiquids_.size() == 0)
82  {
84  << "Evaporation model selected, but no active liquids defined"
85  << nl << endl;
86  }
87  else
88  {
89  Info<< "Participating liquid species:" << endl;
90 
91  // Determine mapping between liquid and carrier phase species
92  forAll(activeLiquids_, i)
93  {
94  Info<< " " << activeLiquids_[i] << endl;
95  liqToCarrierMap_[i] =
96  owner.composition().carrierId(activeLiquids_[i]);
97  }
98 
99  // Determine mapping between model active liquids and global liquids
100  const label idLiquid = owner.composition().idLiquid();
101  forAll(activeLiquids_, i)
102  {
103  liqToLiqMap_[i] =
104  owner.composition().localId(idLiquid, activeLiquids_[i]);
105  }
106  }
107 }
108 
109 
110 template<class CloudType>
112 (
114 )
115 :
117  liquids_(pcm.owner().thermo().liquids()),
118  activeLiquids_(pcm.activeLiquids_),
119  liqToCarrierMap_(pcm.liqToCarrierMap_),
120  liqToLiqMap_(pcm.liqToLiqMap_)
121 {}
122 
123 
124 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
125 
126 template<class CloudType>
128 {}
129 
130 
131 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
132 
133 template<class CloudType>
135 (
136  const scalar dt,
137  const label celli,
138  const scalar Re,
139  const scalar Pr,
140  const scalar d,
141  const scalar nu,
142  const scalar T,
143  const scalar Ts,
144  const scalar pc,
145  const scalar Tc,
146  const scalarField& X,
147  scalarField& dMassPC
148 ) const
149 {
150  // immediately evaporate mass that has reached critical condition
151  if ((liquids_.Tc(X) - T) < SMALL)
152  {
153  if (debug)
154  {
156  << "Parcel reached critical conditions: "
157  << "evaporating all available mass" << endl;
158  }
159 
160  forAll(activeLiquids_, i)
161  {
162  const label lid = liqToLiqMap_[i];
163  dMassPC[lid] = GREAT;
164  }
165 
166  return;
167  }
168 
169  // droplet surface pressure assumed to surface vapour pressure
170  scalar ps = liquids_.pv(pc, Ts, X);
171 
172  // vapour density at droplet surface [kg/m3]
173  scalar rhos = ps*liquids_.W(X)/(RR*Ts);
174 
175  // construct carrier phase species volume fractions for cell, celli
176  const scalarField XcMix(calcXc(celli));
177 
178  // carrier thermo properties
179  scalar Hsc = 0.0;
180  scalar Hc = 0.0;
181  scalar Cpc = 0.0;
182  scalar kappac = 0.0;
183  forAll(this->owner().thermo().carrier().Y(), i)
184  {
185  scalar Yc = this->owner().thermo().carrier().Y()[i][celli];
186  Hc += Yc*this->owner().thermo().carrier().Ha(i, pc, Tc);
187  Hsc += Yc*this->owner().thermo().carrier().Ha(i, ps, Ts);
188  Cpc += Yc*this->owner().thermo().carrier().Cp(i, ps, Ts);
189  kappac += Yc*this->owner().thermo().carrier().kappa(i, ps, Ts);
190  }
191 
192  // calculate mass transfer of each specie in liquid
193  forAll(activeLiquids_, i)
194  {
195  const label gid = liqToCarrierMap_[i];
196  const label lid = liqToLiqMap_[i];
197 
198  // boiling temperature at cell pressure for liquid species lid [K]
199  const scalar TBoil = liquids_.properties()[lid].pvInvert(pc);
200 
201  // limit droplet temperature to boiling/critical temperature
202  const scalar Td = min(T, 0.999*TBoil);
203 
204  // saturation pressure for liquid species lid [Pa]
205  const scalar pSat = liquids_.properties()[lid].pv(pc, Td);
206 
207  // carrier phase concentration
208  const scalar Xc = XcMix[gid];
209 
210 
211  if (Xc*pc > pSat)
212  {
213  // saturated vapour - no phase change
214  }
215  else
216  {
217  // vapour diffusivity [m2/s]
218  const scalar Dab = liquids_.properties()[lid].D(ps, Ts);
219 
220  // Schmidt number
221  const scalar Sc = nu/(Dab + ROOTVSMALL);
222 
223  // Sherwood number
224  const scalar Sh = this->Sh(Re, Sc);
225 
226 
227  if (pSat > 0.999*pc)
228  {
229  // boiling
230 
231  const scalar deltaT = max(T - TBoil, 0.5);
232 
233  // vapour heat of formation
234  const scalar hv = liquids_.properties()[lid].hl(pc, Td);
235 
236  // empirical heat transfer coefficient W/m2/K
237  scalar alphaS = 0.0;
238  if (deltaT < 5.0)
239  {
240  alphaS = 760.0*pow(deltaT, 0.26);
241  }
242  else if (deltaT < 25.0)
243  {
244  alphaS = 27.0*pow(deltaT, 2.33);
245  }
246  else
247  {
248  alphaS = 13800.0*pow(deltaT, 0.39);
249  }
250 
251  // flash-boil vaporisation rate
252  const scalar Gf = alphaS*deltaT*pi*sqr(d)/hv;
253 
254  // model constants
255  // NOTE: using Sherwood number instead of Nusselt number
256  const scalar A = (Hc - Hsc)/hv;
257  const scalar B = pi*kappac/Cpc*d*Sh;
258 
259  scalar G = 0.0;
260  if (A > 0.0)
261  {
262  // heat transfer from the surroundings contributes
263  // to the vaporisation process
264  scalar Gr = 1e-5;
265 
266  for (label i=0; i<50; i++)
267  {
268  scalar GrDash = Gr;
269 
270  G = B/(1.0 + Gr)*log(1.0 + A*(1.0 + Gr));
271  Gr = Gf/G;
272 
273  if (mag(Gr - GrDash)/GrDash < 1e-3)
274  {
275  break;
276  }
277  }
278  }
279 
280  dMassPC[lid] += (G + Gf)*dt;
281  }
282  else
283  {
284  // evaporation
285 
286  // surface molar fraction - Raoult's Law
287  const scalar Xs = X[lid]*pSat/pc;
288 
289  // molar ratio
290  const scalar Xr = (Xs - Xc)/max(SMALL, 1.0 - Xs);
291 
292  if (Xr > 0)
293  {
294  // mass transfer [kg]
295  dMassPC[lid] += pi*d*Sh*Dab*rhos*log(1.0 + Xr)*dt;
296  }
297  }
298  }
299  }
300 }
301 
302 
303 template<class CloudType>
305 (
306  const label idc,
307  const label idl,
308  const scalar p,
309  const scalar T
310 ) const
311 {
312  scalar dh = 0;
313 
314  scalar TDash = T;
315  if (liquids_.properties()[idl].pv(p, T) >= 0.999*p)
316  {
317  TDash = liquids_.properties()[idl].pvInvert(p);
318  }
319 
320  typedef PhaseChangeModel<CloudType> parent;
321  switch (parent::enthalpyTransfer_)
322  {
323  case (parent::etLatentHeat):
324  {
325  dh = liquids_.properties()[idl].hl(p, TDash);
326  break;
327  }
328  case (parent::etEnthalpyDifference):
329  {
330  scalar hc = this->owner().composition().carrier().Ha(idc, p, TDash);
331  scalar hp = liquids_.properties()[idl].h(p, TDash);
332 
333  dh = hc - hp;
334  break;
335  }
336  default:
337  {
339  << "Unknown enthalpyTransfer type" << abort(FatalError);
340  }
341  }
342 
343  return dh;
344 }
345 
346 
347 template<class CloudType>
349 (
350  const scalarField& X
351 ) const
352 {
353  return liquids_.Tpt(X);
354 }
355 
356 
357 template<class CloudType>
359 (
360  const scalar p,
361  const scalarField& X
362 ) const
363 {
364  return liquids_.pvInvert(p, X);
365 }
366 
367 
368 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::constant::thermodynamic::RR
const scalar RR
Universal gas constant: default in [J/(kmol K)].
Definition: thermodynamicConstants.C:46
Foam::PhaseChangeModel
Templated phase change model class.
Definition: ReactingCloud.H:61
mathematicalConstants.H
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::LiquidEvaporationBoil::calculate
virtual void calculate(const scalar dt, const label celli, const scalar Re, const scalar Pr, const scalar d, const scalar nu, const scalar T, const scalar Ts, const scalar pc, const scalar Tc, const scalarField &X, scalarField &dMassPC) const
Update model.
Definition: LiquidEvaporationBoil.C:135
Foam::constant::universal::G
const dimensionedScalar G
Newtonian constant of gravitation.
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:59
Foam::LiquidEvaporationBoil::liqToLiqMap_
List< label > liqToLiqMap_
Mapping between local and global liquid species.
Definition: LiquidEvaporationBoil.H:80
thermo
psiReactionThermo & thermo
Definition: createFields.H:28
B
static const Foam::dimensionedScalar B("", Foam::dimless, 18.678)
A
static const Foam::dimensionedScalar A("", Foam::dimPressure, 611.21)
Foam::constant::mathematical::e
constexpr scalar e(M_E)
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:350
Foam::LiquidEvaporationBoil::dh
virtual scalar dh(const label idc, const label idl, const scalar p, const scalar T) const
Return the enthalpy per unit mass.
Definition: LiquidEvaporationBoil.C:305
specie.H
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::LiquidEvaporationBoil::Tvap
virtual scalar Tvap(const scalarField &X) const
Return vapourisation temperature.
Definition: LiquidEvaporationBoil.C:349
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::LiquidEvaporationBoil::~LiquidEvaporationBoil
virtual ~LiquidEvaporationBoil()
Destructor.
Definition: LiquidEvaporationBoil.C:127
nu
volScalarField & nu
Definition: readMechanicalProperties.H:176
Foam::LiquidEvaporationBoil::TMax
virtual scalar TMax(const scalar p, const scalarField &X) const
Return maximum/limiting temperature.
Definition: LiquidEvaporationBoil.C:359
Foam::LiquidEvaporationBoil::liqToCarrierMap_
List< label > liqToCarrierMap_
Mapping between liquid and carrier species.
Definition: LiquidEvaporationBoil.H:77
Foam::LiquidEvaporationBoil::activeLiquids_
List< word > activeLiquids_
List of active liquid names.
Definition: LiquidEvaporationBoil.H:74
Foam::Field< scalar >
LiquidEvaporationBoil.H
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
Pr
dimensionedScalar Pr("Pr", dimless, laminarTransport)
Foam::T
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Definition: FieldFieldFunctions.C:58
Foam::LiquidEvaporationBoil::LiquidEvaporationBoil
LiquidEvaporationBoil(const dictionary &dict, CloudType &cloud)
Construct from dictionary.
Definition: LiquidEvaporationBoil.C:70
Foam::LiquidEvaporationBoil
Liquid evaporation model.
Definition: LiquidEvaporationBoil.H:62
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::DSMCCloud
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:71
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
Foam::LiquidEvaporationBoil::calcXc
tmp< scalarField > calcXc(const label celli) const
Calculate the carrier phase component volume fractions at celli.
Definition: LiquidEvaporationBoil.C:38
Foam::log
dimensionedScalar log(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:262
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:137
Y
PtrList< volScalarField > & Y
Definition: createFieldRefs.H:7
Foam::Xr
spatialTransform Xr(const vector &a, const scalar omega)
Rotational spatial transformation tensor about axis a by omega radians.
Definition: spatialTransformI.H:275
Foam::PhaseChangeModel::Sh
scalar Sh() const
Sherwood number.
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:372
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::constant::mathematical
Mathematical constants.
Foam::constant::mathematical::pi
constexpr scalar pi(M_PI)
Foam::nl
constexpr char nl
Definition: Ostream.H:385
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::Re
scalarField Re(const UList< complex > &cf)
Extract real component.
Definition: complexField.C:159
Foam::sum
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
Definition: DimensionedFieldFunctions.C:327
Foam::cbrt
dimensionedScalar cbrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:155
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:298