34 template<
class ParcelType>
35 template<
class TrackCloudType>
38 TrackCloudType& cloud,
42 ParcelType::setCellValues(cloud, td);
50 if (td.
Tc() < cloud.constProps().TMin())
55 <<
"Limiting observed temperature in cell " << this->
cell()
56 <<
" to " << cloud.constProps().TMin() <<
nl <<
endl;
59 td.
Tc() = cloud.constProps().TMin();
64 template<
class ParcelType>
65 template<
class TrackCloudType>
68 TrackCloudType& cloud,
73 td.Uc() += cloud.UTransRef()[this->
cell()]/this->massCell(td);
76 td.
Tc() += cloud.hsTransRef()[this->
cell()]/(CpMean*this->massCell(td));
78 if (td.
Tc() < cloud.constProps().TMin())
83 <<
"Limiting observed temperature in cell " << this->
cell()
84 <<
" to " << cloud.constProps().TMin() <<
nl <<
endl;
87 td.
Tc() = cloud.constProps().TMin();
92 template<
class ParcelType>
93 template<
class TrackCloudType>
96 const TrackCloudType& cloud,
107 Ts = (2.0*T + td.
Tc())/3.0;
109 if (Ts < cloud.constProps().TMin())
114 <<
"Limiting parcel surface temperature to " 115 << cloud.constProps().TMin() <<
nl <<
endl;
118 Ts = cloud.constProps().TMin();
122 const scalar TRatio = td.
Tc()/Ts;
124 rhos = td.rhoc()*TRatio;
127 mus = td.muInterp().interpolate(this->coordinates(), tetIs)/TRatio;
130 Pr = td.
Cpc()*mus/kappas;
131 Pr =
max(rootVSmall, Pr);
135 template<
class ParcelType>
136 template<
class TrackCloudType>
139 TrackCloudType& cloud,
146 const scalar np0 = this->nParticle_;
147 const scalar mass0 = this->mass();
150 const scalar
T0 = this->T_;
155 scalar Ts, rhos, mus, Pr, kappas;
156 calcSurfaceValues(cloud, td, this->T_, Ts, rhos, mus, Pr, kappas);
159 scalar Re = this->
Re(rhos, this->U_, td.Uc(), this->d_, mus);
181 scalar dhsTrans = 0.0;
192 this->calcHeatTransfer
208 Cp_ = cloud.composition().Cp(0, Y, td.
pc(),
T0);
216 this->calcVelocity(cloud, td, dt, Re, mus, mass0, Su, dUTrans, Spu);
221 if (cloud.solution().coupled())
224 cloud.UTransRef()[this->
cell()] += np0*dUTrans;
227 cloud.UCoeffRef()[this->
cell()] += np0*Spu;
230 cloud.hsTransRef()[this->
cell()] += np0*dhsTrans;
233 cloud.hsCoeffRef()[this->
cell()] += np0*Sph;
236 if (cloud.radiation())
238 const scalar ap = this->areaP();
239 const scalar T4 =
pow4(T0);
240 cloud.radAreaP()[this->
cell()] += dt*np0*ap;
241 cloud.radT4()[this->
cell()] += dt*np0*T4;
242 cloud.radAreaPT4()[this->
cell()] += dt*np0*ap*T4;
248 template<
class ParcelType>
249 template<
class TrackCloudType>
252 TrackCloudType& cloud,
275 const scalar d = this->d();
276 const scalar
rho = this->
rho();
277 const scalar As = this->areaS(d);
278 const scalar V = this->volume(d);
279 const scalar m = rho*V;
282 scalar htc = cloud.heatTransfer().htc(d, Re, Pr, kappa, NCpW);
285 const scalar bcp = htc*As/(m*Cp_);
286 const scalar acp = bcp*td.
Tc();
288 if (cloud.radiation())
290 const tetIndices tetIs = this->currentTetIndices();
293 const scalar
epsilon = cloud.constProps().epsilon0();
295 ancp += As*epsilon*(Gc/4.0 - sigma*
pow4(T_));
300 const scalar deltaT = cloud.TIntegrator().delta(T_, dt, acp + ancp, bcp);
301 const scalar deltaTncp = ancp*dt;
302 const scalar deltaTcp = deltaT - deltaTncp;
305 scalar Tnew = T_ + deltaT;
306 Tnew =
min(
max(Tnew, cloud.constProps().TMin()), cloud.constProps().TMax());
308 dhsTrans -= m*Cp_*deltaTcp;
318 template<
class ParcelType>
330 template<
class ParcelType>
const interpolation< scalar > & kappaInterp() const
Return const access to the interpolator for continuous.
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
scalar Cp_
Specific heat capacity [J/kg/K].
void calc(TrackCloudType &cloud, trackingData &td, const scalar dt)
Update parcel properties over the time interval.
Ostream & endl(Ostream &os)
Add newline and flush stream.
ThermoParcel(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti)
Construct from mesh, coordinates and topology.
const interpolation< scalar > & GInterp() const
Return const access to the interpolator for continuous.
const interpolation< scalar > & TInterp() const
Return const access to the interpolator for continuous.
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m^2/K^4].
virtual Type interpolate(const vector &position, const label celli, const label facei=-1) const =0
Interpolate field to the given point in the given cell.
const GeometricField< Type, fvPatchField, volMesh > & psi() const
Return the field to be interpolated.
scalar T_
Temperature [K].
bool isType(const Type &t)
Check the typeid.
void cellValueSourceCorrection(TrackCloudType &cloud, trackingData &td, const scalar dt)
Correct cell values using latest transfer information.
const Type & value() const
Return const reference to value.
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
Storage and named access for the indices of a tet which is part of the decomposition of a cell...
Dummy heat transfer model for 'none'.
scalar Tc() const
Return the continuous phase temperature.
void calcSurfaceValues(const TrackCloudType &cloud, const trackingData &td, const scalar T, scalar &Ts, scalar &rhos, scalar &mus, scalar &Pr, scalar &kappas) const
Calculate surface thermo properties.
const interpolation< scalar > & CpInterp() const
Return const access to the interpolator for continuous.
scalar calcHeatTransfer(TrackCloudType &cloud, trackingData &td, const scalar dt, const scalar Re, const scalar Pr, const scalar kappa, const scalar NCpW, const scalar Sh, scalar &dhsTrans, scalar &Sph)
Calculate new particle temperature.
#define WarningInFunction
Report a warning using Foam::Warning.
scalar Cpc() const
Return the continuous phase specific heat capacity.
scalar pc() const
Return the continuous phase pressure.
A cell is defined as a list of faces with extra functionality.
PtrList< volScalarField > & Y
dimensionedScalar pow4(const dimensionedScalar &ds)
void setCellValues(TrackCloudType &cloud, trackingData &td)
Set cell values.
Mesh consisting of general polyhedral cells.
const tmp< volScalarField::Internal > & Su
scalarField Re(const UList< complex > &cf)
Thermodynamic parcel class with one/two-way coupling with the continuous phase.