ThermoParcel< ParcelType > Class Template Reference

Thermodynamic parcel class with one/two-way coupling with the continuous phase. Includes Kinematic parcel sub-models, plus: More...

Inheritance diagram for ThermoParcel< ParcelType >:
[legend]
Collaboration diagram for ThermoParcel< ParcelType >:
[legend]

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

 TypeName ("ThermoParcel")
 Runtime type information. More...
 
 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, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti, const label typeId, const scalar nParticle0, const scalar d0, const scalar dTarget0, const vector &U0, const vector &f0, const vector &angularMomentum0, const vector &torque0, const constantProperties &constProps)
 Construct from components. More...
 
 ThermoParcel (const polyMesh &mesh, Istream &is, bool readFields=true, bool newFormat=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< particleclone () const
 Construct and return a (basic particle) clone. More...
 
virtual autoPtr< particleclone (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 hs () const
 Return the parcel sensible enthalpy. 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...
 
void writeProperties (Ostream &os, const wordRes &filters, const word &delim, const bool namesOnly=false) const
 Write individual parcel properties to stream. 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 >
static void readObjects (CloudType &c, const objectRegistry &obr)
 Read particle fields as objects from the obr registry. More...
 
template<class CloudType >
static void writeObjects (const CloudType &c, objectRegistry &obr)
 Write particle fields as objects into the obr registry. More...
 

Static Public Attributes

static const std::size_t sizeofFields
 Size in bytes of the fields. 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

Ostreamoperator (Ostream &, const ThermoParcel< ParcelType > &)
 

Detailed Description

template<class ParcelType>
class Foam::ThermoParcel< ParcelType >

Thermodynamic parcel class with one/two-way coupling with the continuous phase. Includes Kinematic parcel sub-models, plus:

- heat transfer

Source files

Definition at line 71 of file ThermoParcel.H.

Constructor & Destructor Documentation

◆ ThermoParcel() [1/6]

ThermoParcel ( const polyMesh mesh,
const barycentric coordinates,
const label  celli,
const label  tetFacei,
const label  tetPti 
)
inline

Construct from mesh, coordinates and topology.

Other properties initialised as null

Definition at line 76 of file ThermoParcelI.H.

◆ ThermoParcel() [2/6]

ThermoParcel ( const polyMesh mesh,
const vector position,
const label  celli 
)
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 92 of file ThermoParcelI.H.

◆ ThermoParcel() [3/6]

ThermoParcel ( const polyMesh mesh,
const barycentric coordinates,
const label  celli,
const label  tetFacei,
const label  tetPti,
const label  typeId,
const scalar  nParticle0,
const scalar  d0,
const scalar  dTarget0,
const vector U0,
const vector f0,
const vector angularMomentum0,
const vector torque0,
const constantProperties constProps 
)
inline

Construct from components.

Definition at line 106 of file ThermoParcelI.H.

◆ ThermoParcel() [4/6]

ThermoParcel ( const polyMesh mesh,
Istream is,
bool  readFields = true,
bool  newFormat = true 
)

◆ ThermoParcel() [5/6]

ThermoParcel ( const ThermoParcel< ParcelType > &  p)

Construct as a copy.

◆ ThermoParcel() [6/6]

ThermoParcel ( const ThermoParcel< ParcelType > &  p,
const polyMesh mesh 
)

Construct as a copy.

Member Function Documentation

◆ calcHeatTransfer() [1/2]

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 
)
protected

Calculate new particle temperature.

◆ TypeName()

TypeName ( "ThermoParcel< ParcelType >"  )

Runtime type information.

◆ AddToPropertyList()

AddToPropertyList ( ParcelType  ,
" T"+" Cp  
)

String representation of properties.

◆ clone() [1/2]

virtual autoPtr< particle > clone ( ) const
inlinevirtual

Construct and return a (basic particle) clone.

Definition at line 351 of file ThermoParcel.H.

◆ clone() [2/2]

virtual autoPtr< particle > clone ( const polyMesh mesh) const
inlinevirtual

Construct and return a (basic particle) clone.

Definition at line 357 of file ThermoParcel.H.

References mesh.

◆ T() [1/2]

Foam::scalar T
inline

Return const access to temperature.

Definition at line 207 of file ThermoParcelI.H.

◆ Cp() [1/2]

Foam::scalar Cp
inline

Return const access to specific heat capacity.

Definition at line 214 of file ThermoParcelI.H.

◆ hs()

Foam::scalar hs
inline

Return the parcel sensible enthalpy.

Definition at line 221 of file ThermoParcelI.H.

◆ T() [2/2]

Foam::scalar & T
inline

Return access to temperature.

Definition at line 228 of file ThermoParcelI.H.

◆ Cp() [2/2]

Foam::scalar & Cp
inline

Return access to specific heat capacity.

Definition at line 235 of file ThermoParcelI.H.

◆ setCellValues()

void setCellValues ( TrackCloudType &  cloud,
trackingData td 
)

◆ cellValueSourceCorrection()

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 >::trackingData::CpInterp(), Foam::endl(), Foam::nl, interpolation< Type >::psi(), ThermoParcel< ParcelType >::trackingData::Tc(), and WarningInFunction.

Here is the call graph for this function:

◆ calcSurfaceValues()

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 101 of file ThermoParcel.C.

References coordinates(), ThermoParcel< ParcelType >::trackingData::Cpc(), Foam::endl(), interpolation< Type >::interpolate(), ThermoParcel< ParcelType >::trackingData::kappaInterp(), max(), Foam::nl, Pr(), T, ThermoParcel< ParcelType >::trackingData::Tc(), and WarningInFunction.

Here is the call graph for this function:

◆ calc()

void calc ( TrackCloudType &  cloud,
trackingData td,
const scalar  dt 
)

Update parcel properties over the time interval.

Definition at line 144 of file ThermoParcel.C.

References Foam::pow4(), Pr(), Foam::Re(), Su, and T0.

Here is the call graph for this function:

◆ readFields()

void readFields ( CloudType c)
static

Read.

Definition at line 91 of file ThermoParcelIO.C.

References Cp, IOobject::MUST_READ, p, and T.

Referenced by ThermoParcel< ParcelType >::ThermoParcel().

Here is the caller graph for this function:

◆ writeFields()

void writeFields ( const CloudType c)
static

Write.

Definition at line 117 of file ThermoParcelIO.C.

References Cp, IOobject::NO_READ, p, T, regIOobject::write(), and fieldAverage::writeFields().

Here is the call graph for this function:

◆ writeProperties()

void writeProperties ( Ostream os,
const wordRes filters,
const word delim,
const bool  namesOnly = false 
) const

Write individual parcel properties to stream.

Definition at line 142 of file ThermoParcelIO.C.

References os(), writeProp, and injectedParticle::writeProperties().

Here is the call graph for this function:

◆ readObjects()

void readObjects ( CloudType c,
const objectRegistry obr 
)
static

Read particle fields as objects from the obr registry.

Definition at line 165 of file ThermoParcelIO.C.

References Cp, p, and T.

◆ writeObjects()

void writeObjects ( const CloudType c,
objectRegistry obr 
)
static

Write particle fields as objects into the obr registry.

Definition at line 191 of file ThermoParcelIO.C.

References Cp, p, and T.

◆ calcHeatTransfer() [2/2]

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 253 of file ThermoParcel.C.

References coordinates(), epsilon, ThermoParcel< ParcelType >::trackingData::GInterp(), interpolation< Type >::interpolate(), max(), Foam::min(), Foam::pow4(), Pr(), Foam::Re(), rho, Foam::constant::physicoChemical::sigma, ThermoParcel< ParcelType >::trackingData::Tc(), and dimensioned< Type >::value().

Here is the call graph for this function:

Friends And Related Function Documentation

◆ operator

Ostream & operator ( Ostream ,
const ThermoParcel< ParcelType > &   
)
friend

Member Data Documentation

◆ sizeofFields

const std::size_t sizeofFields
static

Size in bytes of the fields.

Definition at line 78 of file ThermoParcel.H.

Referenced by ThermoParcel< ParcelType >::ThermoParcel().

◆ T_

scalar T_
protected

Temperature [K].

Definition at line 253 of file ThermoParcel.H.

Referenced by ThermoParcel< ParcelType >::ThermoParcel().

◆ Cp_

scalar Cp_
protected

Specific heat capacity [J/(kg.K)].

Definition at line 256 of file ThermoParcel.H.

Referenced by ThermoParcel< ParcelType >::trackingData::Cp(), and ThermoParcel< ParcelType >::ThermoParcel().


The documentation for this class was generated from the following files: