LiquidEvaporation.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 "LiquidEvaporation.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  // construct carrier phase species volume fractions for cell, celli
170  const scalarField Xc(calcXc(celli));
171 
172  // calculate mass transfer of each specie in liquid
173  forAll(activeLiquids_, i)
174  {
175  const label gid = liqToCarrierMap_[i];
176  const label lid = liqToLiqMap_[i];
177 
178  // vapour diffusivity [m2/s]
179  const scalar Dab = liquids_.properties()[lid].D(pc, Ts);
180 
181  // saturation pressure for species i [pa]
182  // - carrier phase pressure assumed equal to the liquid vapour pressure
183  // close to the surface
184  // NOTE: if pSat > pc then particle is superheated
185  // calculated evaporation rate will be greater than that of a particle
186  // at boiling point, but this is not a boiling model
187  const scalar pSat = liquids_.properties()[lid].pv(pc, T);
188 
189  // Schmidt number
190  const scalar Sc = nu/(Dab + ROOTVSMALL);
191 
192  // Sherwood number
193  const scalar Sh = this->Sh(Re, Sc);
194 
195  // mass transfer coefficient [m/s]
196  const scalar kc = Sh*Dab/(d + ROOTVSMALL);
197 
198  // vapour concentration at surface [kmol/m3] at film temperature
199  const scalar Cs = pSat/(RR*Ts);
200 
201  // vapour concentration in bulk gas [kmol/m3] at film temperature
202  const scalar Cinf = Xc[gid]*pc/(RR*Ts);
203 
204  // molar flux of vapour [kmol/m2/s]
205  const scalar Ni = max(kc*(Cs - Cinf), 0.0);
206 
207  // mass transfer [kg]
208  dMassPC[lid] += Ni*pi*sqr(d)*liquids_.properties()[lid].W()*dt;
209  }
210 }
211 
212 
213 template<class CloudType>
215 (
216  const label idc,
217  const label idl,
218  const scalar p,
219  const scalar T
220 ) const
221 {
222  scalar dh = 0;
223 
224  typedef PhaseChangeModel<CloudType> parent;
225  switch (parent::enthalpyTransfer_)
226  {
227  case (parent::etLatentHeat):
228  {
229  dh = liquids_.properties()[idl].hl(p, T);
230  break;
231  }
232  case (parent::etEnthalpyDifference):
233  {
234  scalar hc = this->owner().composition().carrier().Ha(idc, p, T);
235  scalar hp = liquids_.properties()[idl].h(p, T);
236 
237  dh = hc - hp;
238  break;
239  }
240  default:
241  {
243  << "Unknown enthalpyTransfer type" << abort(FatalError);
244  }
245  }
246 
247  return dh;
248 }
249 
250 
251 template<class CloudType>
253 (
254  const scalarField& X
255 ) const
256 {
257  return liquids_.Tpt(X);
258 }
259 
260 
261 template<class CloudType>
263 (
264  const scalar p,
265  const scalarField& X
266 ) const
267 {
268  return liquids_.pvInvert(p, X);
269 }
270 
271 
272 // ************************************************************************* //
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::tmp
A class for managing temporary objects.
Definition: PtrList.H:59
Foam::LiquidEvaporation::TMax
virtual scalar TMax(const scalar p, const scalarField &X) const
Return maximum/limiting temperature.
Definition: LiquidEvaporation.C:263
thermo
psiReactionThermo & thermo
Definition: createFields.H:28
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:350
Foam::LiquidEvaporation::liqToLiqMap_
List< label > liqToLiqMap_
Mapping between local and global liquid species.
Definition: LiquidEvaporation.H:70
specie.H
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
nu
volScalarField & nu
Definition: readMechanicalProperties.H:176
Foam::Field< scalar >
Foam::LiquidEvaporation::liqToCarrierMap_
List< label > liqToCarrierMap_
Mapping between liquid and carrier species.
Definition: LiquidEvaporation.H:67
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::LiquidEvaporation::activeLiquids_
List< word > activeLiquids_
List of active liquid names.
Definition: LiquidEvaporation.H:64
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::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:137
Y
PtrList< volScalarField > & Y
Definition: createFieldRefs.H:7
Foam::PhaseChangeModel::Sh
scalar Sh() const
Sherwood number.
Foam::LiquidEvaporation::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: LiquidEvaporation.C:135
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::LiquidEvaporation::dh
virtual scalar dh(const label idc, const label idl, const scalar p, const scalar T) const
Return the enthalpy per unit mass.
Definition: LiquidEvaporation.C:215
Foam::LiquidEvaporation
Liquid evaporation model.
Definition: LiquidEvaporation.H:52
Foam::LiquidEvaporation::LiquidEvaporation
LiquidEvaporation(const dictionary &dict, CloudType &cloud)
Construct from dictionary.
Definition: LiquidEvaporation.C:70
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
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
Cs
const scalarField & Cs
Definition: solveBulkSurfactant.H:10
Foam::LiquidEvaporation::Tvap
virtual scalar Tvap(const scalarField &X) const
Return vapourisation temperature.
Definition: LiquidEvaporation.C:253
LiquidEvaporation.H
Foam::cbrt
dimensionedScalar cbrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:155
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:298
Foam::LiquidEvaporation::~LiquidEvaporation
virtual ~LiquidEvaporation()
Destructor.
Definition: LiquidEvaporation.C:127
Foam::LiquidEvaporation::calcXc
tmp< scalarField > calcXc(const label celli) const
Calculate the carrier phase component volume fractions at celli.
Definition: LiquidEvaporation.C:38