Go to the documentation of this file.
34 template<
class CloudType>
45 this->owner().thermo().carrier().Y()[i][cellI]
46 /this->owner().thermo().carrier().W(i);
53 template<
class CloudType>
66 template<
class CloudType>
74 liquids_(owner.thermo().liquids()),
75 activeLiquids_(this->coeffDict().
lookup(
"activeLiquids")),
76 liqToCarrierMap_(activeLiquids_.size(), -1),
77 liqToLiqMap_(activeLiquids_.size(), -1)
79 if (activeLiquids_.size() == 0)
82 <<
"Evaporation model selected, but no active liquids defined"
87 Info<<
"Participating liquid species:" <<
endl;
92 Info<<
" " << activeLiquids_[i] <<
endl;
94 owner.composition().carrierId(activeLiquids_[i]);
98 const label idLiquid = owner.composition().idLiquid();
102 owner.composition().localId(idLiquid, activeLiquids_[i]);
108 template<
class CloudType>
115 liquids_(pcm.owner().thermo().liquids()),
124 template<
class CloudType>
131 template<
class CloudType>
149 if ((liquids_.Tc(X) -
T) < SMALL)
154 <<
"Parcel reached critical conditions: "
155 <<
"evaporating all avaliable mass" <<
endl;
160 const label lid = liqToLiqMap_[i];
161 dMassPC[lid] = GREAT;
168 scalar ps = liquids_.pv(pc, Ts, X);
171 scalar rhos = ps*liquids_.W(X)/(
RR*Ts);
183 scalar Yc = this->owner().thermo().carrier().Y()[i][cellI];
184 Hc += Yc*this->owner().thermo().carrier().Ha(i, pc, Tc);
185 Hsc += Yc*this->owner().thermo().carrier().Ha(i, ps, Ts);
186 Cpc += Yc*this->owner().thermo().carrier().Cp(i, ps, Ts);
187 kappac += Yc*this->owner().thermo().carrier().kappa(i, ps, Ts);
193 const label gid = liqToCarrierMap_[i];
194 const label lid = liqToLiqMap_[i];
197 const scalar TBoil = liquids_.properties()[lid].pvInvert(pc);
200 const scalar Td =
min(
T, 0.999*TBoil);
203 const scalar
pSat = liquids_.properties()[lid].pv(pc, Td);
206 const scalar Xc = XcMix[gid];
216 const scalar Dab = liquids_.properties()[lid].D(ps, Ts);
219 const scalar Sc =
nu/(Dab + ROOTVSMALL);
222 const scalar
Sh = this->
Sh(Re, Sc);
229 const scalar deltaT =
max(
T - TBoil, 0.5);
232 const scalar hv = liquids_.properties()[lid].hl(pc, Td);
238 alphaS = 760.0*
pow(deltaT, 0.26);
240 else if (deltaT < 25.0)
242 alphaS = 27.0*
pow(deltaT, 2.33);
246 alphaS = 13800.0*
pow(deltaT, 0.39);
250 const scalar Gf = alphaS*deltaT*
pi*
sqr(d)/hv;
254 const scalar
A = (Hc - Hsc)/hv;
255 const scalar B =
pi*kappac/Cpc*d*
Sh;
264 for (
label i=0; i<50; i++)
268 G = B/(1.0 + Gr)*
log(1.0 +
A*(1.0 + Gr));
271 if (
mag(Gr - GrDash)/GrDash < 1
e-3)
278 dMassPC[lid] += (
G + Gf)*dt;
285 const scalar Xs = X[lid]*
pSat/pc;
288 const scalar Xr = (Xs - Xc)/
max(SMALL, 1.0 - Xs);
293 dMassPC[lid] +=
pi*d*
Sh*Dab*rhos*
log(1.0 + Xr)*dt;
301 template<
class CloudType>
313 if (liquids_.properties()[idl].pv(
p,
T) >= 0.999*
p)
315 TDash = liquids_.properties()[idl].pvInvert(
p);
319 switch (parent::enthalpyTransfer_)
321 case (parent::etLatentHeat):
323 dh = liquids_.properties()[idl].hl(
p, TDash);
326 case (parent::etEnthalpyDifference):
328 scalar hc = this->owner().composition().carrier().Ha(idc,
p, TDash);
329 scalar hp = liquids_.properties()[idl].h(
p, TDash);
345 template<
class CloudType>
351 return liquids_.Tpt(X);
355 template<
class CloudType>
362 return liquids_.pvInvert(
p, X);
const scalar RR
Universal gas constant (default in [J/(kmol K)])
#define forAll(list, i)
Loop across all elements in list.
const dimensionedScalar G
Newtonian constant of gravitation.
A class for managing temporary objects.
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.
List< label > liqToLiqMap_
Mapping between local and global liquid species.
scalarField Re(const UList< complex > &cf)
Ostream & endl(Ostream &os)
Add newline and flush stream.
virtual scalar dh(const label idc, const label idl, const scalar p, const scalar T) const
Return the enthalpy per unit mass.
dimensioned< scalar > mag(const dimensioned< Type > &)
simpleMatrix< scalar > A(Nc)
virtual scalar Tvap(const scalarField &X) const
Return vapourisation temperature.
virtual ~LiquidEvaporationBoil()
Destructor.
virtual scalar TMax(const scalar p, const scalarField &X) const
Return maximum/limiting temperature.
List< label > liqToCarrierMap_
Mapping between liquid and carrier species.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
List< word > activeLiquids_
List of active liquid names.
Pre-declare SubField and related Field type.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
LiquidEvaporationBoil(const dictionary &dict, CloudType &cloud)
Construct from dictionary.
Liquid evaporation model.
dimensionedScalar Pr("Pr", dimless, laminarTransport)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Templated base class for dsmc cloud.
const dimensionedScalar & pSat
A list of keyword definitions, which are a keyword followed by any number of values (e....
dimensionedScalar log(const dimensionedScalar &ds)
errorManip< error > abort(error &err)
scalar Sh() const
Sherwood number.
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
tmp< scalarField > calcXc(const label cellI) const
Calculate the carrier phase component volume fractions at cellI.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
dimensionedScalar cbrt(const dimensionedScalar &ds)
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
#define WarningInFunction
Report a warning using Foam::Warning.
PtrList< volScalarField > & Y
stressControl lookup("compactNormalStress") >> compactNormalStress