Go to the documentation of this file.
37 template<
class CloudType>
48 this->owner().thermo().carrier().Y()[i][celli]
49 /this->owner().thermo().carrier().W(i);
56 template<
class CloudType>
69 template<
class CloudType>
77 liquids_(owner.
thermo().liquids()),
78 activeLiquids_(this->coeffDict().
lookup(
"activeLiquids")),
79 liqToCarrierMap_(activeLiquids_.size(), -1),
80 liqToLiqMap_(activeLiquids_.size(), -1)
85 <<
"Evaporation model selected, but no active liquids defined"
90 Info<<
"Participating liquid species:" <<
endl;
101 const label idLiquid = owner.composition().idLiquid();
111 template<
class CloudType>
114 const LiquidEvaporationBoil<CloudType>& pcm
118 liquids_(pcm.owner().
thermo().liquids()),
119 activeLiquids_(pcm.activeLiquids_),
120 liqToCarrierMap_(pcm.liqToCarrierMap_),
121 liqToLiqMap_(pcm.liqToLiqMap_)
127 template<
class CloudType>
134 template<
class CloudType>
155 if ((liquids_.Tc(X) -
T) < SMALL)
160 <<
"Parcel reached critical conditions: "
161 <<
"evaporating all available mass" <<
endl;
166 const label lid = liqToLiqMap_[i];
167 dMassPC[lid] = GREAT;
174 scalar ps = liquids_.pv(pc, Ts, X);
177 scalar rhos = ps*liquids_.W(X)/(
RR*Ts);
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);
199 const label gid = liqToCarrierMap_[i];
200 const label lid = liqToLiqMap_[i];
203 const scalar TBoil = liquids_.properties()[lid].pvInvert(pc);
206 const scalar Td =
min(
T, 0.999*TBoil);
209 const scalar pSat = liquids_.properties()[lid].pv(pc, Td);
212 const scalar Xc = XcMix[gid];
222 const scalar Dab = liquids_.properties()[lid].D(ps, Ts);
225 const scalar Sc =
nu/(Dab + ROOTVSMALL);
228 const scalar Sh = this->Sh(
Re, Sc);
235 const scalar deltaT =
max(
T - TBoil, 0.5);
238 const scalar hv = liquids_.properties()[lid].hl(pc, Td);
244 alphaS = 760.0*
pow(deltaT, 0.26);
246 else if (deltaT < 25.0)
248 alphaS = 27.0*
pow(deltaT, 2.33);
252 alphaS = 13800.0*
pow(deltaT, 0.39);
256 const scalar Gf = alphaS*deltaT*
pi*
sqr(d)/hv;
260 const scalar
A = (Hc - Hsc)/hv;
261 const scalar
B =
pi*kappac/Cpc*d*Sh;
270 for (label i=0; i<50; i++)
274 G =
B/(1.0 + Gr)*
log(1.0 +
A*(1.0 + Gr));
277 if (
mag(Gr - GrDash)/GrDash < 1
e-3)
284 dMassPC[lid] += (
G + Gf)*dt;
291 const scalar Xs = X[lid]*pSat/pc;
294 const scalar
Xr = (Xs - Xc)/
max(SMALL, 1.0 - Xs);
299 dMassPC[lid] +=
pi*d*Sh*Dab*rhos*
log(1.0 +
Xr)*dt;
307 template<
class CloudType>
319 if (liquids_.properties()[idl].pv(
p,
T) >= 0.999*
p)
321 TDash = liquids_.properties()[idl].pvInvert(
p);
324 typedef PhaseChangeModel<CloudType> parent;
325 switch (parent::enthalpyTransfer_)
327 case (parent::etLatentHeat):
329 dh = liquids_.properties()[idl].hl(
p, TDash);
332 case (parent::etEnthalpyDifference):
334 scalar hc = this->owner().composition().carrier().Ha(idc,
p, TDash);
335 scalar hp = liquids_.properties()[idl].h(
p, TDash);
351 template<
class CloudType>
357 return liquids_.Tpt(X);
361 template<
class CloudType>
368 return liquids_.pvInvert(
p, X);
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
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const dimensionedScalar G
A class for managing temporary objects.
List< label > liqToLiqMap_
psiReactionThermo & thermo
static const Foam::dimensionedScalar B("", Foam::dimless, 18.678)
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
static const Foam::dimensionedScalar A("", Foam::dimPressure, 611.21)
Ostream & endl(Ostream &os)
virtual scalar dh(const label idc, const label idl, const scalar p, const scalar T) const
label min(const labelHashSet &set, label minValue=labelMax)
virtual scalar Tvap(const scalarField &X) const
virtual ~LiquidEvaporationBoil()
virtual scalar TMax(const scalar p, const scalarField &X) const
List< label > liqToCarrierMap_
List< word > activeLiquids_
Generic templated field type.
dimensionedScalar Pr("Pr", dimless, laminarTransport)
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
LiquidEvaporationBoil(const dictionary &dict, CloudType &cloud)
Liquid evaporation model.
DSMCCloud< dsmcParcel > CloudType
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
label max(const labelHashSet &set, label maxValue=labelMin)
Lookup type of boundary radiation properties.
Templated base class for dsmc cloud.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
tmp< scalarField > calcXc(const label celli) const
dimensionedScalar log(const dimensionedScalar &ds)
errorManip< error > abort(error &err)
PtrList< volScalarField > & Y
spatialTransform Xr(const vector &a, const scalar omega)
#define FatalErrorInFunction
dimensionedSymmTensor sqr(const dimensionedVector &dv)
constexpr scalar pi(M_PI)
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
const dimensionedScalar e
scalarField Re(const UList< complex > &cf)
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
dimensionedScalar cbrt(const dimensionedScalar &ds)
#define WarningInFunction