33 template<
class CloudType>
39 scalarField Xc(this->owner().composition().carrier().
Y().size());
44 this->owner().composition().carrier().Y()[i][celli]
45 /this->owner().composition().carrier().WiValue(i);
52 template<
class CloudType>
65 template<
class CloudType>
73 liquids_(owner.
thermo().liquids()),
76 this->coeffDict().template lookupOrDefault<
Switch>
82 activeLiquids_(this->coeffDict().lookup(
"activeLiquids")),
83 liqToCarrierMap_(activeLiquids_.size(), -1),
84 liqToLiqMap_(activeLiquids_.size(), -1)
89 <<
"Evaporation model selected, but no active liquids defined"
94 Info<<
"Participating liquid species:" <<
endl;
105 const label idLiquid =
owner.composition().idLiquid();
115 template<
class CloudType>
122 liquids_(pcm.owner().
thermo().liquids()),
123 activeLiquids_(pcm.activeLiquids_),
124 liqToCarrierMap_(pcm.liqToCarrierMap_),
125 liqToLiqMap_(pcm.liqToLiqMap_)
131 template<
class CloudType>
138 template<
class CloudType>
142 const typename CloudType::parcelType::trackingData& td,
157 if ((liquids_.Tc(X) -
T) < small)
162 <<
"Parcel reached critical conditions: "
163 <<
"evaporating all available mass" <<
endl;
168 const label lid = liqToLiqMap_[i];
169 dMassPC[lid] = great;
181 const label gid = liqToCarrierMap_[i];
182 const label lid = liqToLiqMap_[i];
185 const scalar Dab = liquids_.properties()[lid].D(pc, Ts);
193 const scalar pSat = liquids_.properties()[lid].pv(pc,
T);
196 const scalar Sc = nu/(Dab + rootVSmall);
199 const scalar Sh = this->Sh(
Re, Sc);
202 const scalar kc = Sh*Dab/(d + rootVSmall);
205 const scalar Cs = X[lid]*pSat/(
RR*Ts);
208 const scalar Cinf = Xc[gid]*pc/(
RR*Ts);
211 scalar Ni = kc*(Cs - Cinf);
220 dMassPC[lid] += Ni*
pi*
sqr(d)*liquids_.properties()[lid].W()*dt;
225 template<
class CloudType>
237 switch (parent::enthalpyTransfer_)
239 case (parent::etLatentHeat):
244 case (parent::etEnthalpyDifference):
246 scalar hc = this->owner().composition().carrier().hai(idc,
p,
T);
247 scalar hp = liquids_.properties()[idl].ha(
p,
T);
263 template<
class CloudType>
269 return liquids_.Tpt(X);
273 template<
class CloudType>
280 return liquids_.pvInvert(
p, X);
#define forAll(list, i)
Loop across all elements in list.
const CloudType & owner() const
Return const access to the owner cloud.
Templated base class for dsmc cloud.
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Liquid evaporation model.
List< word > activeLiquids_
List of active liquid names.
virtual scalar dh(const label idc, const label idl, const scalar p, const scalar T) const
Return the enthalpy per unit mass.
virtual ~LiquidEvaporation()
Destructor.
List< label > liqToLiqMap_
Mapping between local and global liquid species.
List< label > liqToCarrierMap_
Mapping between liquid and carrier species.
LiquidEvaporation(const dictionary &dict, CloudType &cloud)
Construct from dictionary.
virtual scalar Tvap(const scalarField &X) const
Return vapourisation temperature.
virtual void calculate(const typename CloudType::parcelType &p, const typename CloudType::parcelType::trackingData &td, const scalar dt, 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.
tmp< scalarField > calcXc(const label celli) const
Calculate the carrier phase component volume fractions at celli.
virtual scalar TMax(const scalar p, const scalarField &X) const
Return maximum/limiting temperature.
void size(const label)
Override size to be inconsistent with allocated storage.
Templated phase change model class.
scalar Sh() const
Sherwood number.
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
A list of keyword definitions, which are a keyword followed by any number of values (e....
const dictionary & properties() const
Return const access to the properties dictionary.
A class for managing temporary objects.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
#define WarningInFunction
Report a warning using Foam::Warning.
const dimensionedScalar RR
Universal gas constant: default SI units: [J/kmol/K].
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Ostream & endl(Ostream &os)
Add newline and flush stream.
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
errorManip< error > abort(error &err)
dimensionedScalar sqrt(const dimensionedScalar &ds)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
dimensionedScalar cbrt(const dimensionedScalar &ds)
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
scalarField Re(const UList< complex > &cf)
PtrList< volScalarField > & Y
fluidMulticomponentThermo & thermo