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