Go to the documentation of this file.
39 #ifndef ThermoParcel_H
40 #define ThermoParcel_H
51 template<
class ParcelType>
54 template<
class ParcelType>
65 template<
class ParcelType>
81 public ParcelType::constantProperties
124 inline scalar
T0()
const;
127 inline scalar
TMin()
const;
130 inline scalar
TMax()
const;
137 inline scalar
Cp0()
const;
145 inline scalar
f0()
const;
149 template<
class CloudType>
249 template<
class TrackData>
290 const label tetFaceI,
300 const label tetFaceI,
303 const scalar nParticle0,
305 const scalar dTarget0,
308 const vector& angularMomentum0,
367 inline scalar
T()
const;
370 inline scalar
Cp()
const;
373 inline scalar
hs()
const;
376 inline scalar
Tc()
const;
379 inline scalar
Cpc()
const;
394 template<
class TrackData>
403 template<
class TrackData>
412 template<
class TrackData>
426 template<
class TrackData>
438 template<
class CloudType>
442 template<
class CloudType>
448 friend Ostream& operator<< <ParcelType>
const volScalarField kappa_
Local copy of carrier thermal conductivity field.
Class to hold thermo particle constant properties.
scalar epsilon0() const
Return const access to the particle emissivity [].
static void writeFields(const CloudType &c)
Write.
AddToPropertyList(ParcelType, " T"+" Cp")
String representation of properties.
scalarField Re(const UList< complex > &cf)
Factory class to read-construct particles used for.
void calc(TrackData &td, const scalar dt, const label cellI)
Update parcel properties over the time interval.
autoPtr< interpolation< scalar > > TInterp_
Temperature field interpolator.
TypeName("ThermoParcel")
Runtime type information.
bool cp(const fileName &src, const fileName &dst)
Copy, recursively if necessary, the source to the destination.
scalar Tc_
Temperature [K].
This function object reads fields from the time directories and adds them to the mesh database for fu...
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
Mesh consisting of general polyhedral cells.
void setCellValues(TrackData &td, const scalar dt, const label cellI)
Set cell values.
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.
autoPtr< ThermoParcel< ParcelType > > operator()(Istream &is) const
scalar f0() const
Return const access to the particle scattering factor [].
autoPtr< interpolation< scalar > > CpInterp_
Specific heat capacity field interpolator.
scalar Tc() const
Return const access to carrier temperature.
const volScalarField & Cp() const
Return access to the locally stored carrier Cp field.
void setTMax(const scalar TMax)
Set the maximum temperature [K].
scalar T0() const
Return const access to the particle initial temperature [K].
scalar Cpc() const
Return const access to carrier specific heat capacity.
const interpolation< scalar > & CpInterp() const
Return const access to the interpolator for continuous.
const volScalarField & kappa() const
Return access to the locally stored carrier kappa field.
iNew(const polyMesh &mesh)
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
virtual autoPtr< particle > clone() const
Construct and return a (basic particle) clone.
scalar T_
Temperature [K].
scalar Cp_
Specific heat capacity [J/(kg.K)].
scalar T() const
Return const access to temperature.
scalar hs() const
Return the parcel sensible enthalpy.
const interpolation< scalar > & kappaInterp() const
Return const access to the interpolator for continuous.
scalar Cp0() const
Return const access to the particle specific heat capacity.
dimensionedScalar Pr("Pr", dimless, laminarTransport)
Templated base class for dsmc cloud.
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.
demandDrivenEntry< scalar > f0_
Particle scattering factor [] (radiation)
A list of keyword definitions, which are a keyword followed by any number of values (e....
const volScalarField Cp_
Local copy of carrier specific heat field.
scalar TMax() const
Return const access to maximum temperature [K].
autoPtr< interpolation< scalar > > kappaInterp_
Thermal conductivity field interpolator.
const interpolation< scalar > & GInterp() const
Return const access to the interpolator for continuous.
demandDrivenEntry< scalar > T0_
Particle initial temperature [K].
demandDrivenEntry< scalar > TMax_
Maximum temperature [K].
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
A cloud is a collection of lagrangian particles.
constantProperties()
Null constructor.
demandDrivenEntry< scalar > TMin_
Minimum temperature [K].
ThermoParcel(const polyMesh &mesh, const vector &position, const label cellI, const label tetFaceI, const label tetPtI)
Construct from owner, position, and cloud owner.
Thermodynamic parcel class with one/two-way coupling with the continuous phase. Includes Kinematic pa...
autoPtr< interpolation< scalar > > GInterp_
Radiation field interpolator.
void cellValueSourceCorrection(TrackData &td, const scalar dt, const label cellI)
Correct cell values using latest transfer information.
scalar TMin() const
Return const access to minimum temperature [K].
const dimensionedScalar c
Speed of light in a vacuum.
demandDrivenEntry< scalar > Cp0_
Particle specific heat capacity [J/(kg.K)].
virtual autoPtr< particle > clone(const polyMesh &mesh) const
Construct and return a (basic particle) clone.
const interpolation< scalar > & TInterp() const
Return const access to the interpolator for continuous.
static const std::size_t sizeofFields_
Size in bytes of the fields.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Generic GeometricField class.
demandDrivenEntry< scalar > epsilon0_
Particle emissivity [] (radiation)
ParcelType::template TrackingData< CloudType >::trackPart trackPart
static void readFields(CloudType &c)
Read.
scalar Cp() const
Return const access to specific heat capacity.
scalar Cpc_
Specific heat capacity [J/(kg.K)].
TrackingData(CloudType &cloud, trackPart part=ParcelType::template TrackingData< CloudType >::tpLinearTrack)
Construct from components.