Thermodynamic parcel class with one/two-way coupling with the continuous phase. More...
Classes | |
class | constantProperties |
Class to hold thermo particle constant properties. More... | |
class | iNew |
Factory class to read-construct particles used for. More... | |
class | trackingData |
Public Member Functions | |
AddToPropertyList (ParcelType, " T"+" Cp") | |
String representation of properties. More... | |
ThermoParcel (const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti) | |
Construct from mesh, coordinates and topology. More... | |
ThermoParcel (const polyMesh &mesh, const vector &position, const label celli) | |
Construct from a position and a cell, searching for the rest of the. More... | |
ThermoParcel (const polyMesh &mesh, Istream &is, bool readFields=true) | |
Construct from Istream. More... | |
ThermoParcel (const ThermoParcel &p) | |
Construct as a copy. More... | |
ThermoParcel (const ThermoParcel &p, const polyMesh &mesh) | |
Construct as a copy. More... | |
virtual autoPtr< particle > | clone () const |
Construct and return a (basic particle) clone. More... | |
virtual autoPtr< particle > | clone (const polyMesh &mesh) const |
Construct and return a (basic particle) clone. More... | |
scalar | T () const |
Return const access to temperature. More... | |
scalar | Cp () const |
Return const access to specific heat capacity. More... | |
scalar & | T () |
Return access to temperature. More... | |
scalar & | Cp () |
Return access to specific heat capacity. More... | |
template<class TrackCloudType > | |
void | setCellValues (TrackCloudType &cloud, trackingData &td) |
Set cell values. More... | |
template<class TrackCloudType > | |
void | cellValueSourceCorrection (TrackCloudType &cloud, trackingData &td, const scalar dt) |
Correct cell values using latest transfer information. More... | |
template<class TrackCloudType > | |
void | calcSurfaceValues (TrackCloudType &cloud, trackingData &td, const scalar T, scalar &Ts, scalar &rhos, scalar &mus, scalar &Pr, scalar &kappas) const |
Calculate surface thermo properties. More... | |
template<class TrackCloudType > | |
void | calc (TrackCloudType &cloud, trackingData &td, const scalar dt) |
Update parcel properties over the time interval. More... | |
template<class TrackCloudType > | |
Foam::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) |
Static Public Member Functions | |
template<class CloudType > | |
static void | readFields (CloudType &c) |
Read. More... | |
template<class CloudType > | |
static void | writeFields (const CloudType &c) |
Write. More... | |
template<class CloudType , class CompositionType > | |
static void | writeFields (const CloudType &c, const CompositionType &compModel) |
Write. More... | |
Protected Member Functions | |
template<class TrackCloudType > | |
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. More... | |
Protected Attributes | |
scalar | T_ |
Temperature [K]. More... | |
scalar | Cp_ |
Specific heat capacity [J/kg/K]. More... | |
Friends | |
Ostream & | operator (Ostream &, const ThermoParcel< ParcelType > &) |
Thermodynamic parcel class with one/two-way coupling with the continuous phase.
Definition at line 51 of file ThermoParcel.H.
|
inline |
Construct from mesh, coordinates and topology.
Other properties initialised as null
Definition at line 75 of file ThermoParcelI.H.
Referenced by ThermoParcel< ParcelType >::calcHeatTransfer(), ThermoParcel< ParcelType >::clone(), and ThermoParcel< ParcelType >::constantProperties::constantProperties().
|
inline |
Construct from a position and a cell, searching for the rest of the.
required topology. Other properties are initialised as null.
Definition at line 91 of file ThermoParcelI.H.
ThermoParcel | ( | const polyMesh & | mesh, |
Istream & | is, | ||
bool | readFields = true |
||
) |
Construct from Istream.
Definition at line 47 of file ThermoParcelIO.C.
References IOstream::check(), IOstream::format(), Istream::read(), and readScalar.
ThermoParcel | ( | const ThermoParcel< ParcelType > & | p | ) |
Construct as a copy.
ThermoParcel | ( | const ThermoParcel< ParcelType > & | p, |
const polyMesh & | mesh | ||
) |
Construct as a copy.
|
protected |
Calculate new particle temperature.
Referenced by ThermoParcel< ParcelType >::calc().
AddToPropertyList | ( | ParcelType | , |
" T"+" Cp" | |||
) |
String representation of properties.
Construct and return a (basic particle) clone.
Definition at line 343 of file ThermoParcel.H.
References ThermoParcel< ParcelType >::ThermoParcel().
Construct and return a (basic particle) clone.
Definition at line 349 of file ThermoParcel.H.
References ThermoParcel< ParcelType >::ThermoParcel().
|
inline |
Return const access to temperature.
Definition at line 164 of file ThermoParcelI.H.
References ThermoParcel< ParcelType >::T_.
Referenced by ThermoParcel< ParcelType >::iNew::operator()().
|
inline |
Return const access to specific heat capacity.
Definition at line 171 of file ThermoParcelI.H.
References ThermoParcel< ParcelType >::Cp_.
Referenced by ThermoParcel< ParcelType >::iNew::operator()().
|
inline |
Return access to temperature.
Definition at line 178 of file ThermoParcelI.H.
References ThermoParcel< ParcelType >::T_.
|
inline |
Return access to specific heat capacity.
Definition at line 185 of file ThermoParcelI.H.
References ThermoParcel< ParcelType >::Cp_.
void setCellValues | ( | TrackCloudType & | cloud, |
trackingData & | td | ||
) |
Set cell values.
Definition at line 37 of file ThermoParcel.C.
References ThermoParcel< ParcelType >::cellValueSourceCorrection(), ThermoParcel< ParcelType >::trackingData::Cpc(), ThermoParcel< ParcelType >::trackingData::CpInterp(), Foam::endl(), interpolation< Type >::interpolate(), Foam::nl, ThermoParcel< ParcelType >::trackingData::Tc(), ThermoParcel< ParcelType >::trackingData::TInterp(), and WarningInFunction.
Referenced by ThermoParcel< ParcelType >::iNew::operator()().
void cellValueSourceCorrection | ( | TrackCloudType & | cloud, |
trackingData & | td, | ||
const scalar | dt | ||
) |
Correct cell values using latest transfer information.
Definition at line 67 of file ThermoParcel.C.
References ThermoParcel< ParcelType >::calcSurfaceValues(), ThermoParcel< ParcelType >::trackingData::CpInterp(), Foam::endl(), Foam::nl, interpolation< Type >::psi(), ThermoParcel< ParcelType >::trackingData::Tc(), and WarningInFunction.
Referenced by ThermoParcel< ParcelType >::iNew::operator()(), and ThermoParcel< ParcelType >::setCellValues().
void calcSurfaceValues | ( | TrackCloudType & | cloud, |
trackingData & | td, | ||
const scalar | T, | ||
scalar & | Ts, | ||
scalar & | rhos, | ||
scalar & | mus, | ||
scalar & | Pr, | ||
scalar & | kappas | ||
) | const |
Calculate surface thermo properties.
Definition at line 95 of file ThermoParcel.C.
References ThermoParcel< ParcelType >::calc(), ThermoParcel< ParcelType >::trackingData::Cpc(), Foam::endl(), interpolation< Type >::interpolate(), ThermoParcel< ParcelType >::trackingData::kappaInterp(), Foam::max(), Foam::nl, ThermoParcel< ParcelType >::trackingData::Tc(), and WarningInFunction.
Referenced by ThermoParcel< ParcelType >::cellValueSourceCorrection(), and ThermoParcel< ParcelType >::iNew::operator()().
void calc | ( | TrackCloudType & | cloud, |
trackingData & | td, | ||
const scalar | dt | ||
) |
Update parcel properties over the time interval.
Definition at line 138 of file ThermoParcel.C.
References ThermoParcel< ParcelType >::calcHeatTransfer(), ThermoParcel< ParcelType >::trackingData::pc(), Foam::pow4(), Foam::Re(), Su, T0, Y, and Foam::Zero.
Referenced by ThermoParcel< ParcelType >::calcSurfaceValues(), and ThermoParcel< ParcelType >::iNew::operator()().
|
static |
Read.
Definition at line 80 of file ThermoParcelIO.C.
References Cloud< ParticleType >::checkFieldIOobject(), Cp(), ThermoParcel< ParcelType >::Cp_, Cloud< ParticleType >::fieldIOobject(), forAllIter, p, Foam::readFields(), Cloud< ParticleType >::size(), T, and ThermoParcel< ParcelType >::T_.
Referenced by ThermoParcel< ParcelType >::iNew::operator()().
|
static |
Write.
Definition at line 107 of file ThermoParcelIO.C.
References Cp(), ThermoParcel< ParcelType >::Cp_, Cloud< ParticleType >::fieldIOobject(), forAllConstIter(), p, Cloud< ParticleType >::size(), T, and ThermoParcel< ParcelType >::T_.
Referenced by ThermoParcel< ParcelType >::iNew::operator()().
|
static |
Foam::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 | ||
) |
Definition at line 251 of file ThermoParcel.C.
References ThermoParcel< ParcelType >::Cp_, epsilon, ThermoParcel< ParcelType >::trackingData::GInterp(), interpolation< Type >::interpolate(), Foam::isType(), Foam::max(), mesh, Foam::min(), p, Foam::pow4(), rho, Foam::constant::physicoChemical::sigma, ThermoParcel< ParcelType >::T_, ThermoParcel< ParcelType >::trackingData::Tc(), ThermoParcel< ParcelType >::ThermoParcel(), and dimensioned< Type >::value().
|
friend |
|
protected |
Temperature [K].
Definition at line 268 of file ThermoParcel.H.
Referenced by ThermoParcel< ParcelType >::calcHeatTransfer(), ThermoParcel< ParcelType >::readFields(), ThermoParcel< ParcelType >::T(), and ThermoParcel< ParcelType >::writeFields().
|
protected |
Specific heat capacity [J/kg/K].
Definition at line 271 of file ThermoParcel.H.
Referenced by ThermoParcel< ParcelType >::calcHeatTransfer(), ThermoParcel< ParcelType >::Cp(), ThermoParcel< ParcelType >::readFields(), and ThermoParcel< ParcelType >::writeFields().