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  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 "LiquidEvaporation.H"
30 #include "specie.H"
31 #include "mathematicalConstants.H"
32 
33 using namespace Foam::constant::mathematical;
34 
35 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
36 
37 template<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 
56 template<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 
69 template<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;
96  liqToCarrierMap_[i] =
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 
111 template<class CloudType>
113 (
114  const LiquidEvaporation<CloudType>& pcm
115 )
116 :
117  PhaseChangeModel<CloudType>(pcm),
118  liquids_(pcm.owner().thermo().liquids()),
119  activeLiquids_(pcm.activeLiquids_),
120  liqToCarrierMap_(pcm.liqToCarrierMap_),
121  liqToLiqMap_(pcm.liqToLiqMap_)
122 {}
123 
124 
125 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
126 
127 template<class CloudType>
129 {}
130 
131 
132 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
133 
134 template<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  // construct carrier phase species volume fractions for cell, celli
174  const scalarField Xc(calcXc(celli));
175 
176  // calculate mass transfer of each specie in liquid
177  forAll(activeLiquids_, i)
178  {
179  const label gid = liqToCarrierMap_[i];
180  const label lid = liqToLiqMap_[i];
181 
182  // vapour diffusivity [m2/s]
183  const scalar Dab = liquids_.properties()[lid].D(pc, Ts);
184 
185  // saturation pressure for species i [pa]
186  // - carrier phase pressure assumed equal to the liquid vapour pressure
187  // close to the surface
188  // NOTE: if pSat > pc then particle is superheated
189  // calculated evaporation rate will be greater than that of a particle
190  // at boiling point, but this is not a boiling model
191  const scalar pSat = liquids_.properties()[lid].pv(pc, T);
192 
193  // Schmidt number
194  const scalar Sc = nu/(Dab + ROOTVSMALL);
195 
196  // Sherwood number
197  const scalar Sh = this->Sh(Re, Sc);
198 
199  // mass transfer coefficient [m/s]
200  const scalar kc = Sh*Dab/(d + ROOTVSMALL);
201 
202  // vapour concentration at surface [kmol/m3] at film temperature
203  const scalar Cs = pSat/(RR*Ts);
204 
205  // vapour concentration in bulk gas [kmol/m3] at film temperature
206  const scalar Cinf = Xc[gid]*pc/(RR*Ts);
207 
208  // molar flux of vapour [kmol/m2/s]
209  const scalar Ni = max(kc*(Cs - Cinf), 0.0);
210 
211  // mass transfer [kg]
212  dMassPC[lid] += Ni*pi*sqr(d)*liquids_.properties()[lid].W()*dt;
213  }
214 }
215 
216 
217 template<class CloudType>
219 (
220  const label idc,
221  const label idl,
222  const scalar p,
223  const scalar T
224 ) const
225 {
226  scalar dh = 0;
227 
228  typedef PhaseChangeModel<CloudType> parent;
229  switch (parent::enthalpyTransfer_)
230  {
231  case (parent::etLatentHeat):
232  {
233  dh = liquids_.properties()[idl].hl(p, T);
234  break;
235  }
236  case (parent::etEnthalpyDifference):
237  {
238  scalar hc = this->owner().composition().carrier().Ha(idc, p, T);
239  scalar hp = liquids_.properties()[idl].h(p, T);
240 
241  dh = hc - hp;
242  break;
243  }
244  default:
245  {
247  << "Unknown enthalpyTransfer type" << abort(FatalError);
248  }
249  }
251  return dh;
252 }
253 
254 
255 template<class CloudType>
257 (
258  const scalarField& X
259 ) const
260 {
261  return liquids_.Tpt(X);
262 }
263 
264 
265 template<class CloudType>
267 (
268  const scalar p,
269  const scalarField& X
270 ) const
271 {
272  return liquids_.pvInvert(p, X);
273 }
274 
275 
276 // ************************************************************************* //
dimensionedScalar Pr("Pr", dimless, laminarTransport)
dictionary dict
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &f1)
DSMCCloud< dsmcParcel > CloudType
Templated phase change model class.
Definition: ReactingCloud.H:57
List< word > activeLiquids_
List of active liquid names.
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
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:598
virtual void calculate(const scalar dt, const label celli, const scalar Re, const scalar Pr, const scalar d, const scalar nu, const scalar rho, const scalar T, const scalar Ts, const scalar pc, const scalar Tc, const scalarField &X, const scalarField &solMass, const scalarField &liqMass, scalarField &dMassPC) const
Update model.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
dimensionedSymmTensor sqr(const dimensionedVector &dv)
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
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:531
List< label > liqToLiqMap_
Mapping between local and global liquid species.
Lookup type of boundary radiation properties.
Definition: lookup.H:57
const CloudType & owner() const
Return const access to the owner cloud.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
psiReactionThermo & thermo
Definition: createFields.H:28
Mathematical constants.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
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.
constexpr scalar pi(M_PI)
errorManip< error > abort(error &err)
Definition: errorManip.H:139
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
int debug
Static debugging option.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
virtual ~LiquidEvaporation()
Destructor.
scalarField Re(const UList< complex > &cmplx)
Extract real component.
Definition: complexField.C:207
LiquidEvaporation(const dictionary &dict, CloudType &cloud)
Construct from dictionary.
#define WarningInFunction
Report a warning using Foam::Warning.
PtrList< volScalarField > & Y
messageStream Info
Information stream (stdout output on master, null elsewhere)
const scalar RR
Universal gas constant: default in [J/(kmol K)].
volScalarField & p
A class for managing temporary objects.
Definition: HashPtrTable.H:50
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:67
volScalarField & nu
const scalarField & Cs
scalar Sh() const
Sherwood number.
tmp< scalarField > calcXc(const label celli) const
Calculate the carrier phase component volume fractions at celli.