LiquidEvaporation.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2018 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "LiquidEvaporation.H"
27 #include "specie.H"
28 #include "mathematicalConstants.H"
29 
30 using namespace Foam::constant::mathematical;
31 
32 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
33 
34 template<class CloudType>
36 (
37  const label celli
38 ) const
39 {
40  scalarField Xc(this->owner().thermo().carrier().Y().size());
41 
42  forAll(Xc, i)
43  {
44  Xc[i] =
45  this->owner().thermo().carrier().Y()[i][celli]
46  /this->owner().thermo().carrier().W(i);
47  }
48 
49  return Xc/sum(Xc);
50 }
51 
52 
53 template<class CloudType>
55 (
56  const scalar Re,
57  const scalar Sc
58 ) const
59 {
60  return 2.0 + 0.6*Foam::sqrt(Re)*cbrt(Sc);
61 }
62 
63 
64 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
65 
66 template<class CloudType>
68 (
69  const dictionary& dict,
70  CloudType& owner
71 )
72 :
73  PhaseChangeModel<CloudType>(dict, owner, typeName),
74  liquids_(owner.thermo().liquids()),
75  activeLiquids_(this->coeffDict().lookup("activeLiquids")),
76  liqToCarrierMap_(activeLiquids_.size(), -1),
77  liqToLiqMap_(activeLiquids_.size(), -1)
78 {
79  if (activeLiquids_.size() == 0)
80  {
82  << "Evaporation model selected, but no active liquids defined"
83  << nl << endl;
84  }
85  else
86  {
87  Info<< "Participating liquid species:" << endl;
88 
89  // Determine mapping between liquid and carrier phase species
90  forAll(activeLiquids_, i)
91  {
92  Info<< " " << activeLiquids_[i] << endl;
93  liqToCarrierMap_[i] =
94  owner.composition().carrierId(activeLiquids_[i]);
95  }
96 
97  // Determine mapping between model active liquids and global liquids
98  const label idLiquid = owner.composition().idLiquid();
99  forAll(activeLiquids_, i)
100  {
101  liqToLiqMap_[i] =
102  owner.composition().localId(idLiquid, activeLiquids_[i]);
103  }
104  }
105 }
106 
107 
108 template<class CloudType>
110 (
112 )
113 :
115  liquids_(pcm.owner().thermo().liquids()),
116  activeLiquids_(pcm.activeLiquids_),
117  liqToCarrierMap_(pcm.liqToCarrierMap_),
118  liqToLiqMap_(pcm.liqToLiqMap_)
119 {}
120 
121 
122 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
123 
124 template<class CloudType>
126 {}
127 
128 
129 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
130 
131 template<class CloudType>
133 (
134  const scalar dt,
135  const label celli,
136  const scalar Re,
137  const scalar Pr,
138  const scalar d,
139  const scalar nu,
140  const scalar T,
141  const scalar Ts,
142  const scalar pc,
143  const scalar Tc,
144  const scalarField& X,
145  scalarField& dMassPC
146 ) const
147 {
148  // immediately evaporate mass that has reached critical condition
149  if ((liquids_.Tc(X) - T) < small)
150  {
151  if (debug)
152  {
154  << "Parcel reached critical conditions: "
155  << "evaporating all available mass" << endl;
156  }
157 
158  forAll(activeLiquids_, i)
159  {
160  const label lid = liqToLiqMap_[i];
161  dMassPC[lid] = great;
162  }
163 
164  return;
165  }
166 
167  // construct carrier phase species volume fractions for cell, celli
168  const scalarField Xc(calcXc(celli));
169 
170  // calculate mass transfer of each specie in liquid
171  forAll(activeLiquids_, i)
172  {
173  const label gid = liqToCarrierMap_[i];
174  const label lid = liqToLiqMap_[i];
175 
176  // vapour diffusivity [m2/s]
177  const scalar Dab = liquids_.properties()[lid].D(pc, Ts);
178 
179  // saturation pressure for species i [pa]
180  // - carrier phase pressure assumed equal to the liquid vapour pressure
181  // close to the surface
182  // NOTE: if pSat > pc then particle is superheated
183  // calculated evaporation rate will be greater than that of a particle
184  // at boiling point, but this is not a boiling model
185  const scalar pSat = liquids_.properties()[lid].pv(pc, T);
186 
187  // Schmidt number
188  const scalar Sc = nu/(Dab + rootVSmall);
189 
190  // Sherwood number
191  const scalar Sh = this->Sh(Re, Sc);
192 
193  // mass transfer coefficient [m/s]
194  const scalar kc = Sh*Dab/(d + rootVSmall);
195 
196  // vapour concentration at surface [kmol/m3] at film temperature
197  const scalar Cs = pSat/(RR*Ts);
198 
199  // vapour concentration in bulk gas [kmol/m3] at film temperature
200  const scalar Cinf = Xc[gid]*pc/(RR*Ts);
201 
202  // molar flux of vapour [kmol/m2/s]
203  const scalar Ni = max(kc*(Cs - Cinf), 0.0);
204 
205  // mass transfer [kg]
206  dMassPC[lid] += Ni*pi*sqr(d)*liquids_.properties()[lid].W()*dt;
207  }
208 }
209 
210 
211 template<class CloudType>
213 (
214  const label idc,
215  const label idl,
216  const scalar p,
217  const scalar T
218 ) const
219 {
220  scalar dh = 0;
221 
222  typedef PhaseChangeModel<CloudType> parent;
223  switch (parent::enthalpyTransfer_)
224  {
225  case (parent::etLatentHeat):
226  {
227  dh = liquids_.properties()[idl].hl(p, T);
228  break;
229  }
230  case (parent::etEnthalpyDifference):
231  {
232  scalar hc = this->owner().composition().carrier().Ha(idc, p, T);
233  scalar hp = liquids_.properties()[idl].h(p, T);
234 
235  dh = hc - hp;
236  break;
237  }
238  default:
239  {
241  << "Unknown enthalpyTransfer type" << abort(FatalError);
242  }
243  }
244 
245  return dh;
246 }
247 
248 
249 template<class CloudType>
251 (
252  const scalarField& X
253 ) const
254 {
255  return liquids_.Tpt(X);
256 }
257 
258 
259 template<class CloudType>
261 (
262  const scalar p,
263  const scalarField& X
264 ) const
265 {
266  return liquids_.pvInvert(p, X);
267 }
268 
269 
270 // ************************************************************************* //
const dictionary & properties() const
Return const access to the properties dictionary.
Definition: subModelBase.C:134
dictionary dict
Liquid evaporation model.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
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:59
Templated phase change model class.
Definition: ReactingCloud.H:55
List< word > activeLiquids_
List of active liquid names.
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.
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
virtual scalar TMax(const scalar p, const scalarField &X) const
Return maximum/limiting temperature.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
dimensionedSymmTensor sqr(const dimensionedVector &dv)
virtual scalar Tvap(const scalarField &X) const
Return vapourisation temperature.
dimensionedScalar sqrt(const dimensionedScalar &ds)
List< label > liqToCarrierMap_
Mapping between liquid and carrier species.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
List< label > liqToLiqMap_
Mapping between local and global liquid species.
rhoReactionThermo & thermo
Definition: createFields.H:28
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
const CloudType & owner() const
Return const access to the owner cloud.
mathematical constants.
dimensionedScalar cbrt(const dimensionedScalar &ds)
virtual scalar dh(const label idc, const label idl, const scalar p, const scalar T) const
Return the enthalpy per unit mass.
errorManip< error > abort(error &err)
Definition: errorManip.H:131
static const char nl
Definition: Ostream.H:265
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
virtual ~LiquidEvaporation()
Destructor.
LiquidEvaporation(const dictionary &dict, CloudType &cloud)
Construct from dictionary.
#define WarningInFunction
Report a warning using Foam::Warning.
PtrList< volScalarField > & Y
messageStream Info
const scalar RR
Universal gas constant (default in [J/(kmol K)])
A class for managing temporary objects.
Definition: PtrList.H:53
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:69
scalar Sh() const
Sherwood number.
tmp< scalarField > calcXc(const label celli) const
Calculate the carrier phase component volume fractions at celli.