33 template<
class ParcelType>
34 template<
class TrackData>
42 ParcelType::setCellValues(td, dt, cellI);
46 Cpc_ = td.CpInterp().interpolate(this->position(), tetIs);
48 Tc_ = td.TInterp().interpolate(this->position(), tetIs);
50 if (Tc_ < td.cloud().constProps().TMin())
56 "void Foam::ThermoParcel<ParcelType>::setCellValues" 62 ) <<
"Limiting observed temperature in cell " << cellI <<
" to " 63 << td.cloud().constProps().TMin() <<
nl <<
endl;
66 Tc_ = td.cloud().constProps().TMin();
71 template<
class ParcelType>
72 template<
class TrackData>
80 this->Uc_ += td.cloud().UTrans()[cellI]/this->massCell(cellI);
82 const scalar CpMean = td.CpInterp().psi()[cellI];
83 Tc_ += td.cloud().hsTrans()[cellI]/(CpMean*this->massCell(cellI));
85 if (Tc_ < td.cloud().constProps().TMin())
91 "void Foam::ThermoParcel<ParcelType>::cellValueSourceCorrection" 97 ) <<
"Limiting observed temperature in cell " << cellI <<
" to " 98 << td.cloud().constProps().TMin() <<
nl <<
endl;
101 Tc_ = td.cloud().constProps().TMin();
106 template<
class ParcelType>
107 template<
class TrackData>
121 Ts = (2.0*T + Tc_)/3.0;
123 if (Ts < td.cloud().constProps().TMin())
129 "void Foam::ThermoParcel<ParcelType>::calcSurfaceValues" 140 ) <<
"Limiting parcel surface temperature to " 141 << td.cloud().constProps().TMin() <<
nl <<
endl;
144 Ts = td.cloud().constProps().TMin();
148 const scalar TRatio = Tc_/Ts;
150 rhos = this->rhoc_*TRatio;
153 mus = td.muInterp().interpolate(this->position(), tetIs)/TRatio;
154 kappas = td.kappaInterp().interpolate(this->position(), tetIs)/TRatio;
156 Pr = Cpc_*mus/kappas;
157 Pr =
max(ROOTVSMALL, Pr);
161 template<
class ParcelType>
162 template<
class TrackData>
172 const scalar np0 = this->nParticle_;
173 const scalar mass0 = this->mass();
176 const scalar T0 = this->T_;
181 scalar Ts, rhos, mus,
Pr, kappas;
182 calcSurfaceValues(td, cellI, this->T_, Ts, rhos, mus, Pr, kappas);
185 scalar Re = this->
Re(this->U_, this->d_, rhos, mus);
207 scalar dhsTrans = 0.0;
218 this->calcHeatTransfer
238 this->calcVelocity(td, dt, cellI, Re, mus, mass0, Su, dUTrans, Spu);
243 if (td.cloud().solution().coupled())
246 td.cloud().UTrans()[cellI] += np0*dUTrans;
249 td.cloud().UCoeff()[cellI] += np0*Spu;
252 td.cloud().hsTrans()[cellI] += np0*dhsTrans;
255 td.cloud().hsCoeff()[cellI] += np0*Sph;
258 if (td.cloud().radiation())
260 const scalar ap = this->areaP();
261 const scalar T4 =
pow4(T0);
262 td.cloud().radAreaP()[cellI] += dt*np0*ap;
263 td.cloud().radT4()[cellI] += dt*np0*T4;
264 td.cloud().radAreaPT4()[cellI] += dt*np0*ap*T4;
270 template<
class ParcelType>
271 template<
class TrackData>
286 if (!td.cloud().heatTransfer().active())
291 const scalar d = this->d();
292 const scalar
rho = this->
rho();
295 scalar htc = td.cloud().heatTransfer().htc(d, Re, Pr, kappa, NCpW);
297 if (
mag(htc) < ROOTVSMALL && !td.cloud().radiation())
302 T_ + dt*Sh/(this->volume(d)*rho*Cp_),
303 td.cloud().constProps().TMin()
307 htc =
max(htc, ROOTVSMALL);
308 const scalar As = this->areaS(d);
310 scalar ap = Tc_ + Sh/(As*htc);
311 scalar bp = 6.0*(Sh/As + htc*(Tc_ - T_));
312 if (td.cloud().radiation())
315 const scalar Gc = td.GInterp().interpolate(this->position(), tetIs);
317 const scalar
epsilon = td.cloud().constProps().epsilon0();
320 scalar
s = epsilon*(Gc/4.0 - sigma*
pow4(T_));
325 bp /= rho*d*Cp_*(ap - T_) + ROOTVSMALL;
329 td.cloud().TIntegrator().integrate(T_, dt, ap*bp, bp);
337 td.cloud().constProps().TMin()
339 td.cloud().constProps().TMax()
344 dhsTrans += Sph*(Tres.average() - Tc_);
352 template<
class ParcelType>
366 template<
class ParcelType>
dimensionedScalar Pr("Pr", dimless, laminarTransport)
scalar Tc_
Temperature [K].
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){gmvFile<< iter().position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject( name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE ))
dimensioned< scalar > mag(const dimensioned< Type > &)
scalar Cp_
Specific heat capacity [J/(kg.K)].
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m2/K4].
Helper class to supply results of integration.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
void calc(TrackData &td, const scalar dt, const label cellI)
Update parcel properties over the time interval.
tmp< GeometricField< Type, fvPatchField, volMesh > > Su(const GeometricField< Type, fvPatchField, volMesh > &su, const GeometricField< Type, fvPatchField, volMesh > &vf)
Ostream & endl(Ostream &os)
Add newline and flush stream.
void cellValueSourceCorrection(TrackData &td, const scalar dt, const label cellI)
Correct cell values using latest transfer information.
#define WarningIn(functionName)
Report a warning using Foam::Warning.
Thermodynamic parcel class with one/two-way coupling with the continuous phase. Includes Kinematic pa...
scalar calcHeatTransfer(TrackData &td, const scalar dt, const label cellI, const scalar Re, const scalar Pr, const scalar kappa, const scalar NCpW, const scalar Sh, scalar &dhsTrans, scalar &Sph)
Calculate new particle temperature.
void setCellValues(TrackData &td, const scalar dt, const label cellI)
Set cell values.
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
scalar T_
Temperature [K].
Mesh consisting of general polyhedral cells.
dimensionedScalar pow4(const dimensionedScalar &ds)
void calcSurfaceValues(TrackData &td, const label cellI, const scalar T, scalar &Ts, scalar &rhos, scalar &mus, scalar &Pr, scalar &kappas) const
Calculate surface thermo properties.
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
Storage and named access for the indices of a tet which is part of the decomposition of a cell...
ThermoParcel(const polyMesh &mesh, const vector &position, const label cellI, const label tetFaceI, const label tetPtI)
Construct from owner, position, and cloud owner.
scalar Cpc_
Specific heat capacity [J/(kg.K)].
const Type & value() const
Return const reference to value.
scalarField Re(const UList< complex > &cf)