37template<
class CloudType>
48 this->
owner().thermo().carrier().Y()[i][celli]
49 /this->
owner().thermo().carrier().W(i);
56template<
class CloudType>
69template<
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;
93 forAll(activeLiquids_, i)
95 Info<<
" " << activeLiquids_[i] << endl;
97 owner.composition().carrierId(activeLiquids_[i]);
101 const label idLiquid =
owner.composition().idLiquid();
111template<
class CloudType>
118 liquids_(pcm.owner().
thermo().liquids()),
119 activeLiquids_(pcm.activeLiquids_),
120 liqToCarrierMap_(pcm.liqToCarrierMap_),
127template<
class CloudType>
134template<
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;
179 const label gid = liqToCarrierMap_[i];
180 const label lid = liqToLiqMap_[i];
183 const scalar Dab = liquids_.properties()[lid].D(pc, Ts);
191 const scalar pSat = liquids_.properties()[lid].pv(pc,
T);
194 const scalar Sc =
nu/(Dab + ROOTVSMALL);
197 const scalar Sh = this->Sh(
Re, Sc);
200 const scalar kc = Sh*Dab/(d + ROOTVSMALL);
203 const scalar
Cs = pSat/(
RR*Ts);
206 const scalar Cinf = Xc[gid]*pc/(
RR*Ts);
209 const scalar Ni =
max(kc*(
Cs - Cinf), 0.0);
212 dMassPC[lid] += Ni*
pi*
sqr(d)*
liquids_.properties()[lid].W()*dt;
217template<
class CloudType>
229 switch (parent::enthalpyTransfer_)
231 case (parent::etLatentHeat):
233 dh = liquids_.properties()[idl].hl(
p,
T);
236 case (parent::etEnthalpyDifference):
238 scalar hc = this->owner().composition().carrier().Ha(idc,
p,
T);
239 scalar hp = liquids_.properties()[idl].h(
p,
T);
255template<
class CloudType>
265template<
class CloudType>
272 return liquids_.pvInvert(
p, X);
const CloudType & owner() const
Return const access to the owner cloud.
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.
scalar Sh(const scalar Re, const scalar Sc) const
Sherwood number as a function of Reynolds and Schmidt numbers.
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.
const liquidMixtureProperties & liquids_
Global liquid properties data.
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.
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.
Templated phase change model class.
scalar Sh() const
Sherwood number.
PhaseChangeModel(CloudType &owner)
Construct null from owner.
void size(const label n)
Older name for setAddressableSize.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Lookup type of boundary radiation properties.
const dictionary & coeffDict() const
Return const access to the coefficients dictionary.
const dictionary & dict() const
Return const access to the cloud dictionary.
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
A class for managing temporary objects.
PtrList< volScalarField > & Y
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
#define WarningInFunction
Report a warning using Foam::Warning.
constexpr scalar pi(M_PI)
const scalar RR
Universal gas constant: default in [J/(kmol K)].
Namespace for handling debugging switches.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
DSMCCloud< dsmcParcel > CloudType
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
Ostream & endl(Ostream &os)
Add newline and flush stream.
dimensionedScalar sqrt(const dimensionedScalar &ds)
errorManip< error > abort(error &err)
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &f1, const label comm)
scalarField Re(const UList< complex > &cmplx)
Extract real component.
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
dimensionedScalar cbrt(const dimensionedScalar &ds)
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
psiReactionThermo & thermo
dimensionedScalar Pr("Pr", dimless, laminarTransport)
#define forAll(list, i)
Loop across all elements in list.