ThermoParcel< ParcelType > Class Template Reference

Thermodynamic parcel class with one/two-way coupling with the continuous phase. More...

Inheritance diagram for ThermoParcel< ParcelType >:
Collaboration diagram for ThermoParcel< ParcelType >:

Classes

class  constantProperties
 Class to hold thermo particle constant properties. 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, const label facei)
 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 (Istream &is, bool readFields=true)
 Construct from Istream. More...
 
 ThermoParcel (const ThermoParcel &p)
 Construct as a copy. More...
 
virtual autoPtr< particleclone () const
 Construct and return a 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 (const TrackCloudType &cloud, const 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

static autoPtr< ThermoParcelNew (Istream &is)
 Construct from Istream and return. More...
 
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

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.

Source files

Definition at line 73 of file ThermoParcel.H.

Constructor & Destructor Documentation

◆ ThermoParcel() [1/4]

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

Construct from mesh, coordinates and topology.

Other properties initialised as null

Definition at line 76 of file ThermoParcelI.H.

Referenced by ThermoParcel< ParcelType >::clone(), and ThermoParcel< ParcelType >::New().

Here is the caller graph for this function:

◆ ThermoParcel() [2/4]

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 93 of file ThermoParcelI.H.

◆ ThermoParcel() [3/4]

ThermoParcel ( Istream is,
bool  readFields = true 
)

◆ ThermoParcel() [4/4]

ThermoParcel ( const ThermoParcel< ParcelType > &  p)

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.

◆ AddToPropertyList()

AddToPropertyList ( ParcelType  ,
" T"+" Cp  
)

String representation of properties.

◆ clone()

virtual autoPtr<particle> clone ( ) const
inlinevirtual

Construct and return a clone.

Definition at line 336 of file ThermoParcel.H.

References ThermoParcel< ParcelType >::ThermoParcel().

Here is the call graph for this function:

◆ New()

static autoPtr<ThermoParcel> New ( Istream is)
inlinestatic

Construct from Istream and return.

Definition at line 342 of file ThermoParcel.H.

References ThermoParcel< ParcelType >::ThermoParcel().

Here is the call graph for this function:

◆ T() [1/2]

Foam::scalar T
inline

Return const access to temperature.

Definition at line 167 of file ThermoParcelI.H.

◆ Cp() [1/2]

Foam::scalar Cp
inline

Return const access to specific heat capacity.

Definition at line 174 of file ThermoParcelI.H.

◆ T() [2/2]

Foam::scalar & T
inline

Return access to temperature.

Definition at line 181 of file ThermoParcelI.H.

◆ Cp() [2/2]

Foam::scalar & Cp
inline

Return access to specific heat capacity.

Definition at line 188 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 66 of file ThermoParcel.C.

References ThermoParcel< ParcelType >::trackingData::CpInterp(), Foam::endl(), Foam::nl, ThermoParcel< ParcelType >::trackingData::Tc(), and WarningInFunction.

Here is the call graph for this function:

◆ calcSurfaceValues()

void calcSurfaceValues ( const TrackCloudType &  cloud,
const trackingData td,
const scalar  T,
scalar &  Ts,
scalar &  rhos,
scalar &  mus,
scalar &  Pr,
scalar &  kappas 
) const

◆ calc()

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

Update parcel properties over the time interval.

Definition at line 137 of file ThermoParcel.C.

References ThermoParcel< ParcelType >::trackingData::pc(), Foam::pow4(), Foam::Re(), Foam::fvc::Su(), T0, Y, and Foam::Zero.

Here is the call graph for this function:

◆ readFields()

void readFields ( CloudType c)
static

Read.

Definition at line 75 of file ThermoParcelIO.C.

References Foam::constant::universal::c, Cp(), forAllIter, IOobject::MUST_READ, p, Foam::readFields(), Foam::T(), and Foam::blendedInterfacialModel::valid().

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

Here is the call graph for this function:
Here is the caller graph for this function:

◆ writeFields() [1/2]

void writeFields ( const CloudType c)
static

Write.

Definition at line 102 of file ThermoParcelIO.C.

References Foam::constant::universal::c, Cp(), forAllConstIter, IOobject::NO_READ, p, and Foam::T().

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

Here is the call graph for this function:
Here is the caller graph for this function:

◆ writeFields() [2/2]

void writeFields ( const CloudType c,
const CompositionType &  compModel 
)
static

Write.

Definition at line 128 of file ThermoParcelIO.C.

References Foam::constant::universal::c, and ThermoParcel< ParcelType >::writeFields().

Here is the call graph for this function:

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

Friends And Related Function Documentation

◆ operator

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

Member Data Documentation

◆ T_

scalar T_
protected

Temperature [K].

Definition at line 268 of file ThermoParcel.H.

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

◆ Cp_

scalar Cp_
protected

Specific heat capacity [J/kg/K].

Definition at line 271 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: