ReactingParcel< ParcelType > Class Template Reference

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

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

Classes

class  constantProperties
 Class to hold reacting parcel constant properties. More...
 

Public Types

typedef ParcelType::trackingData trackingData
 Use base tracking data. More...
 

Public Member Functions

 AddToPropertyList (ParcelType, " nPhases(Y1..YN)")
 String representation of properties. More...
 
 ReactingParcel (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...
 
 ReactingParcel (const polyMesh &mesh, const vector &position, const label celli)
 Construct from a position and a cell, searching for the rest of the. More...
 
 ReactingParcel (Istream &is, bool readFields=true)
 Construct from Istream. More...
 
 ReactingParcel (const ReactingParcel &p)
 Construct as a copy. More...
 
virtual autoPtr< particleclone () const
 Construct and return a clone. More...
 
const scalarFieldY () const
 Return const access to mass fractions of mixture []. More...
 
scalarFieldY ()
 Return access to mass fractions of mixture []. 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 correctSurfaceValues (TrackCloudType &cloud, trackingData &td, const scalar T, const scalarField &Cs, scalar &rhos, scalar &mus, scalar &Prs, scalar &kappas)
 Correct surface values due to emitted species. More...
 
template<class TrackCloudType >
void calc (TrackCloudType &cloud, trackingData &td, const scalar dt)
 Update parcel properties over the time interval. More...
 

Static Public Member Functions

static autoPtr< ReactingParcelNew (Istream &is)
 Construct from Istream and return. More...
 
template<class CloudType , class CompositionType >
static void readFields (CloudType &c, const CompositionType &compModel)
 Read. More...
 
template<class CloudType >
static void readFields (CloudType &c)
 Read - no composition. More...
 
template<class CloudType , class CompositionType >
static void writeFields (const CloudType &c, const CompositionType &compModel)
 Write. More...
 
template<class CloudType >
static void writeFields (const CloudType &c)
 Write - composition supplied. More...
 

Protected Member Functions

template<class TrackCloudType >
void calcPhaseChange (TrackCloudType &cloud, trackingData &td, const scalar dt, const scalar Re, const scalar Pr, const scalar Ts, const scalar nus, const scalar d, const scalar T, const scalar mass, const label idPhase, const scalar YPhase, const scalarField &YComponents, scalarField &dMassPC, scalar &Sh, scalar &N, scalar &NCpW, scalarField &Cs)
 Calculate Phase change. More...
 
scalar updateMassFraction (const scalar mass0, const scalarField &dMass, scalarField &Y) const
 Update mass fraction. More...
 

Protected Attributes

scalarField Y_
 Mass fractions of mixture []. More...
 

Friends

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

Detailed Description

template<class ParcelType>
class Foam::ReactingParcel< ParcelType >

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

Source files

Definition at line 73 of file ReactingParcel.H.

Member Typedef Documentation

◆ trackingData

typedef ParcelType::trackingData trackingData

Use base tracking data.

Definition at line 125 of file ReactingParcel.H.

Constructor & Destructor Documentation

◆ ReactingParcel() [1/4]

ReactingParcel ( 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 65 of file ReactingParcelI.H.

Referenced by ReactingParcel< ParcelType >::New().

Here is the caller graph for this function:

◆ ReactingParcel() [2/4]

ReactingParcel ( 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 81 of file ReactingParcelI.H.

◆ ReactingParcel() [3/4]

ReactingParcel ( Istream is,
bool  readFields = true 
)

Construct from Istream.

Definition at line 45 of file ReactingParcelIO.C.

References IOstream::check(), ReactingParcel< ParcelType >::readFields(), List< T >::transfer(), and ReactingParcel< ParcelType >::Y_.

Here is the call graph for this function:

◆ ReactingParcel() [4/4]

ReactingParcel ( const ReactingParcel< ParcelType > &  p)

Construct as a copy.

Member Function Documentation

◆ calcPhaseChange()

void calcPhaseChange ( TrackCloudType &  cloud,
trackingData td,
const scalar  dt,
const scalar  Re,
const scalar  Pr,
const scalar  Ts,
const scalar  nus,
const scalar  d,
const scalar  T,
const scalar  mass,
const label  idPhase,
const scalar  YPhase,
const scalarField YComponents,
scalarField dMassPC,
scalar &  Sh,
scalar &  N,
scalar &  NCpW,
scalarField Cs 
)
protected

◆ updateMassFraction()

Foam::scalar updateMassFraction ( const scalar  mass0,
const scalarField dMass,
scalarField Y 
) const
protected

Update mass fraction.

Definition at line 156 of file ReactingParcel.C.

References forAll, Foam::sum(), and Y.

Here is the call graph for this function:

◆ AddToPropertyList()

AddToPropertyList ( ParcelType  ,
" nPhases(Y1..YN)"   
)

String representation of properties.

◆ clone()

virtual autoPtr<particle> clone ( ) const
inlinevirtual

Construct and return a clone.

Definition at line 215 of file ReactingParcel.H.

◆ New()

static autoPtr<ReactingParcel> New ( Istream is)
inlinestatic

Construct from Istream and return.

Definition at line 221 of file ReactingParcel.H.

References ReactingParcel< ParcelType >::ReactingParcel().

Here is the call graph for this function:

◆ Y() [1/2]

const Foam::scalarField & Y
inline

Return const access to mass fractions of mixture [].

Definition at line 114 of file ReactingParcelI.H.

◆ Y() [2/2]

Foam::scalarField & Y
inline

Return access to mass fractions of mixture [].

Definition at line 121 of file ReactingParcelI.H.

◆ setCellValues()

void setCellValues ( TrackCloudType &  cloud,
trackingData td 
)

Set cell values.

Definition at line 195 of file ReactingParcel.C.

References Foam::endl(), Foam::nl, and WarningInFunction.

Here is the call graph for this function:

◆ cellValueSourceCorrection()

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

Correct cell values using latest transfer information.

Definition at line 225 of file ReactingParcel.C.

References Foam::endl(), forAll, Foam::mag(), Foam::max(), Foam::nl, WarningInFunction, and Y.

Here is the call graph for this function:

◆ correctSurfaceValues()

void correctSurfaceValues ( TrackCloudType &  cloud,
trackingData td,
const scalar  T,
const scalarField Cs,
scalar &  rhos,
scalar &  mus,
scalar &  Prs,
scalar &  kappas 
)

◆ calc()

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

Update parcel properties over the time interval.

Definition at line 373 of file ReactingParcel.C.

References Foam::cbrt(), composition, basicSpecieMixture::Cp(), forAll, basicSpecieMixture::Hs(), Foam::constant::mathematical::pi(), Foam::pow4(), Foam::Re(), List< T >::size(), basicSpecieMixture::species(), Foam::fvc::Su(), T0, and Foam::Zero.

Here is the call graph for this function:

◆ readFields() [1/2]

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

Read.

Definition at line 82 of file ReactingParcelIO.C.

References Foam::constant::universal::c, forAll, forAllIter, IOobject::MUST_READ, p, Foam::readFields(), List< T >::setSize(), List< T >::size(), Foam::blendedInterfacialModel::valid(), and Y.

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

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

◆ readFields() [2/2]

void readFields ( CloudType c)
static

Read - no composition.

Definition at line 74 of file ReactingParcelIO.C.

References Foam::constant::universal::c, and Foam::readFields().

Here is the call graph for this function:

◆ writeFields() [1/2]

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

Write.

Definition at line 141 of file ReactingParcelIO.C.

References Foam::constant::universal::c, forAll, forAllConstIter, IOobject::NO_READ, p, List< T >::size(), and Y.

Here is the call graph for this function:

◆ writeFields() [2/2]

void writeFields ( const CloudType c)
static

Write - composition supplied.

Definition at line 133 of file ReactingParcelIO.C.

References Foam::constant::universal::c.

Friends And Related Function Documentation

◆ operator

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

Member Data Documentation

◆ Y_

scalarField Y_
protected

Mass fractions of mixture [].

Definition at line 135 of file ReactingParcel.H.

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


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