36 template<
class ParcelType>
37 template<
class TrackCloudType>
40 TrackCloudType& cloud,
60 typedef typename TrackCloudType::thermoCloudType thermoCloudType;
64 typedef typename TrackCloudType::reactingCloudType reactingCloudType;
74 scalar Tvap = phaseChange.
Tvap(X);
81 const scalar TMax = phaseChange.
TMax(td.pc(), X);
82 const scalar Tdash =
min(T, TMax);
83 const scalar Tsdash =
min(Ts, TMax);
103 dMassPC =
min(mass*YPhase*Y, dMassPC);
105 const scalar dMassTot =
sum(dMassPC);
114 const scalar dh = phaseChange.
dh(cid, i, td.pc(), Tdash);
115 Sh -= dMassPC[i]*dh/dt;
120 if (cloud.heatTransfer().BirdCorrection())
123 const scalar Wc = td.rhoc()*
RR*td.Tc()/td.pc();
129 const scalar
Cp = composition.
carrier().Cp(cid, td.pc(), Tsdash);
130 const scalar
W = composition.
carrier().Wi(cid);
131 const scalar Ni = dMassPC[i]/(this->areaS(d)*dt*
W);
134 composition.
liquids().properties()[i].D(td.pc(), Tsdash, Wc);
143 Cs[cid] += Ni*d/(2.0*Dab);
149 template<
class ParcelType>
157 scalar mass1 = mass0 -
sum(dMass);
160 if (mass1 > rootVSmall)
164 Y[i] = (Y[i]*mass0 - dMass[i])/mass1;
174 template<
class ParcelType>
186 template<
class ParcelType>
201 template<
class ParcelType>
202 template<
class TrackCloudType>
205 TrackCloudType& cloud,
209 ParcelType::setCellValues(cloud, td);
211 td.pc() = td.pInterp().interpolate
214 this->currentTetIndices()
217 if (td.pc() < cloud.constProps().pMin())
222 <<
"Limiting observed pressure in cell " << this->
cell()
223 <<
" to " << cloud.constProps().pMin() <<
nl <<
endl;
226 td.pc() = cloud.constProps().pMin();
231 template<
class ParcelType>
232 template<
class TrackCloudType>
235 TrackCloudType& cloud,
240 scalar addedMass = 0.0;
241 scalar maxMassI = 0.0;
242 forAll(cloud.rhoTrans(), i)
244 scalar dm = cloud.rhoTrans(i)[this->
cell()];
245 maxMassI =
max(maxMassI,
mag(dm));
249 if (maxMassI < rootVSmall)
254 const scalar massCell = this->massCell(td);
256 td.rhoc() += addedMass/cloud.pMesh().cellVolumes()[this->
cell()];
258 const scalar massCellNew = massCell + addedMass;
259 td.Uc() = (td.Uc()*massCell + cloud.UTransRef()[this->
cell()])/massCellNew;
262 forAll(cloud.rhoTrans(), i)
264 scalar Y = cloud.rhoTrans(i)[this->
cell()]/addedMass;
265 CpEff += Y*cloud.composition().carrier().Cp(i, td.pc(), td.Tc());
268 const scalar Cpc = td.CpInterp().psi()[this->
cell()];
269 td.Cpc() = (massCell*Cpc + addedMass*CpEff)/massCellNew;
271 td.Tc() += cloud.hsTransRef()[this->
cell()]/(td.Cpc()*massCellNew);
273 if (td.Tc() < cloud.constProps().TMin())
278 <<
"Limiting observed temperature in cell " << this->
cell()
279 <<
" to " << cloud.constProps().TMin() <<
nl <<
endl;
282 td.Tc() = cloud.constProps().TMin();
287 template<
class ParcelType>
288 template<
class TrackCloudType>
291 TrackCloudType& cloud,
302 if (!cloud.heatTransfer().BirdCorrection() || (
sum(Cs) < small))
314 Xinf[i] = carrier.
Y(i)[this->
cell()]/carrier.
Wi(i);
319 const scalar Xsff = 1.0 -
min(
sum(Cs)*
RR*this->T_/td.pc(), 1.0);
322 const scalar CsTot = td.pc()/(
RR*this->T_);
333 const scalar Csi = Cs[i] + Xsff*Xinf[i]*CsTot;
335 Xs[i] = (2.0*Csi + Xinf[i]*CsTot)/3.0;
336 Ys[i] = Xs[i]*carrier.
Wi(i);
346 scalar sumYiSqrtW = 0;
347 scalar sumYiCbrtW = 0;
351 const scalar
W = carrier.
Wi(i);
352 const scalar sqrtW =
sqrt(W);
353 const scalar cbrtW =
cbrt(W);
356 mus += Ys[i]*sqrtW*carrier.
mu(i, td.pc(),
T);
357 kappas += Ys[i]*cbrtW*carrier.
kappa(i, td.pc(),
T);
358 Cps += Xs[i]*carrier.
Cp(i, td.pc(),
T);
360 sumYiSqrtW += Ys[i]*sqrtW;
361 sumYiCbrtW += Ys[i]*cbrtW;
364 Cps =
max(Cps, rootVSmall);
366 rhos *= td.pc()/(
RR*
T);
367 rhos =
max(rhos, rootVSmall);
370 mus =
max(mus, rootVSmall);
372 kappas /= sumYiCbrtW;
373 kappas =
max(kappas, rootVSmall);
375 Prs = Cps*mus/kappas;
379 template<
class ParcelType>
380 template<
class TrackCloudType>
383 TrackCloudType& cloud,
388 typedef typename TrackCloudType::thermoCloudType thermoCloudType;
396 const scalar np0 = this->nParticle_;
397 const scalar d0 = this->d_;
398 const vector& U0 = this->U_;
399 const scalar
T0 = this->T_;
400 const scalar mass0 = this->mass();
404 scalar Ts, rhos, mus, Prs, kappas;
405 this->calcSurfaceValues(cloud, td, T0, Ts, rhos, mus, Prs, kappas);
406 scalar Res = this->
Re(rhos, U0, td.Uc(), d0, mus);
428 scalar dhsTrans = 0.0;
477 scalar mass1 = updateMassFraction(mass0, dMass, Y_);
479 this->Cp_ = composition.
Cp(0, Y_, td.pc(),
T0);
482 if (cloud.constProps().constantVolume())
484 this->rho_ = mass1/this->volume();
488 this->d_ =
cbrt(mass1/this->rho_*6.0/
pi);
492 if (np0*mass1 < cloud.constProps().minParcelMass())
494 td.keepParticle =
false;
496 if (cloud.solution().coupled())
498 scalar dm = np0*mass0;
503 scalar dmi = dm*Y_[i];
505 scalar hs = composition.
carrier().Hs(gid, td.pc(),
T0);
507 cloud.rhoTrans(gid)[this->
cell()] += dmi;
508 cloud.hsTransRef()[this->
cell()] += dmi*hs;
510 cloud.UTransRef()[this->
cell()] += dm*U0;
512 cloud.phaseChange().addToPhaseChangeMass(np0*mass1);
519 correctSurfaceValues(cloud, td, Ts, Cs, rhos, mus, Prs, kappas);
520 Res = this->
Re(rhos, U0, td.Uc(), this->d(), mus);
531 this->calcHeatTransfer
545 this->Cp_ = composition.
Cp(0, Y_, td.pc(),
T0);
553 this->calcVelocity(cloud, td, dt, Res, mus, mass1, Su, dUTrans, Spu);
559 if (cloud.solution().coupled())
564 scalar dm = np0*dMass[i];
566 scalar hs = composition.
carrier().Hs(gid, td.pc(),
T0);
568 cloud.rhoTrans(gid)[this->
cell()] += dm;
569 cloud.UTransRef()[this->
cell()] += dm*U0;
570 cloud.hsTransRef()[this->
cell()] += dm*hs;
574 cloud.UTransRef()[this->
cell()] += np0*dUTrans;
575 cloud.UCoeffRef()[this->
cell()] += np0*Spu;
578 cloud.hsTransRef()[this->
cell()] += np0*dhsTrans;
579 cloud.hsCoeffRef()[this->
cell()] += np0*Sph;
582 if (cloud.radiation())
584 const scalar ap = this->areaP();
585 const scalar T4 =
pow4(T0);
586 cloud.radAreaP()[this->
cell()] += dt*np0*ap;
587 cloud.radT4()[this->
cell()] += dt*np0*T4;
588 cloud.radAreaPT4()[this->
cell()] += dt*np0*ap*T4;
virtual scalar TMax(const scalar p, const scalarField &X) const
Return maximum/limiting temperature.
#define forAll(list, i)
Loop across all elements in list.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Templated phase change model class.
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
virtual scalar kappa(const label speciei, const scalar p, const scalar T) const =0
Thermal conductivity [W/m/K].
void calcPhaseChange(TrackCloudType &cloud, trackingData &td, const scalar dt, const scalar Re, const scalar Pr, const scalar Ts, const scalar nus, const scalar d, const scalar T, const scalar mass, const label idPhase, const scalar YPhase, const scalarField &YComponents, scalarField &dMassPC, scalar &Sh, scalar &N, scalar &NCpW, scalarField &Cs)
Calculate Phase change.
basicSpecieMixture & composition
scalar updateMassFraction(const scalar mass0, const scalarField &dMass, scalarField &Y) const
Update mass fraction.
void size(const label)
Override size to be inconsistent with allocated storage.
scalarField Y_
Mass fractions of mixture [].
dimensionedScalar sqrt(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
ReactingParcel(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti)
Construct from mesh, coordinates and topology.
twoPhaseChangeModel & phaseChange
virtual scalar Wi(const label speciei) const =0
Molecular weight of the given specie [kg/kmol].
label localToCarrierId(const label phaseI, const label id, const bool allowNotFound=false) const
Return carrier id of component given local id.
void addToPhaseChangeMass(const scalar dMass)
Add to phase change mass.
virtual scalar Cp(const label phaseI, const scalarField &Y, const scalar p, const scalar T) const
Return specific heat capacity for the phase phaseI.
Specialisation of basicMixture for a mixture consisting of a number for molecular species...
void cellValueSourceCorrection(TrackCloudType &cloud, trackingData &td, const scalar dt)
Correct cell values using latest transfer information.
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
ParcelType::trackingData trackingData
Use base tracking data.
const dimensionedScalar RR
Universal gas constant: default SI units: [J/kmol/K].
virtual scalar Cp(const label speciei, const scalar p, const scalar T) const =0
Heat capacity at constant pressure [J/kg/K].
virtual scalar mu(const label speciei, const scalar p, const scalar T) const =0
Dynamic viscosity [kg/m/s].
dimensionedScalar cbrt(const dimensionedScalar &ds)
void correctSurfaceValues(TrackCloudType &cloud, trackingData &td, const scalar T, const scalarField &Cs, scalar &rhos, scalar &mus, scalar &Prs, scalar &kappas)
Correct surface values due to emitted species.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
virtual scalar dh(const label idc, const label idl, const scalar p, const scalar T) const
Return the enthalpy per unit mass.
#define WarningInFunction
Report a warning using Foam::Warning.
A cell is defined as a list of faces with extra functionality.
const scalarList W(::W(thermo))
dimensionedScalar pow4(const dimensionedScalar &ds)
scalar Cp(const scalar p, const scalar T) const
void setCellValues(TrackCloudType &cloud, trackingData &td)
Set cell values.
dimensioned< scalar > mag(const dimensioned< Type > &)
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 =0
Update model.
Mesh consisting of general polyhedral cells.
const tmp< volScalarField::Internal > & Su
Reacting parcel class with one/two-way coupling with the continuous phase.
scalar mass0_
Initial mass [kg].
PtrList< volScalarField > & Y()
Return the mass-fraction fields.
Templated reacting parcel composition model class Consists of carrier species (via thermo package)...
virtual scalar Tvap(const scalarField &X) const
Return vapourisation temperature.
scalarField Re(const UList< complex > &cf)
const basicSpecieMixture & carrier() const
Return the carrier components (wrapper function)
void calc(TrackCloudType &cloud, trackingData &td, const scalar dt)
Update parcel properties over the time interval.
const speciesTable & species() const
Return the table of species.
const liquidMixtureProperties & liquids() const
Return the global (additional) liquids.