Reacting parcel class with one/two-way coupling with the continuous phase. More...
Classes | |
class | constantProperties |
Class to hold reacting parcel constant properties. More... | |
class | iNew |
Factory class to read-construct particles used for. More... | |
class | TrackingData |
Public Member Functions | |
TypeName ("ReactingParcel") | |
Runtime type information. More... | |
AddToPropertyList (ParcelType, " mass0"+" nPhases(Y1..YN)") | |
String representation of properties. More... | |
ReactingParcel (const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti) | |
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 (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 scalarField &Y0, const constantProperties &constProps) | |
Construct from components. More... | |
ReactingParcel (const polyMesh &mesh, Istream &is, bool readFields=true) | |
Construct from Istream. More... | |
ReactingParcel (const ReactingParcel &p, const polyMesh &mesh) | |
Construct as a copy. More... | |
ReactingParcel (const ReactingParcel &p) | |
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 | mass0 () const |
Return const access to initial mass [kg]. More... | |
const scalarField & | Y () const |
Return const access to mass fractions of mixture []. More... | |
scalar | pc () const |
Return the owner cell pressure [Pa]. More... | |
scalar & | pc () |
Return reference to the owner cell pressure [Pa]. More... | |
scalar & | mass0 () |
Return access to initial mass [kg]. More... | |
scalarField & | Y () |
Return access to mass fractions of mixture []. More... | |
template<class TrackData > | |
void | setCellValues (TrackData &td, const scalar dt, const label celli) |
Set cell values. More... | |
template<class TrackData > | |
void | cellValueSourceCorrection (TrackData &td, const scalar dt, const label celli) |
Correct cell values using latest transfer information. More... | |
template<class TrackData > | |
void | correctSurfaceValues (TrackData &td, const label celli, const scalar T, const scalarField &Cs, scalar &rhos, scalar &mus, scalar &Prs, scalar &kappas) |
Correct surface values due to emitted species. More... | |
template<class TrackData > | |
void | calc (TrackData &td, const scalar dt, const label celli) |
Update parcel properties over the time interval. More... | |
Static Public Member Functions | |
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 TrackData > | |
void | calcPhaseChange (TrackData &td, const scalar dt, const label celli, 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 | |
scalar | mass0_ |
Initial mass [kg]. More... | |
scalarField | Y_ |
Mass fractions of mixture []. More... | |
scalar | pc_ |
Pressure [Pa]. More... | |
Friends | |
Ostream & | operator (Ostream &, const ReactingParcel< ParcelType > &) |
Reacting parcel class with one/two-way coupling with the continuous phase.
Definition at line 50 of file ReactingParcel.H.
|
inline |
Construct from mesh, coordinates and topology.
Other properties initialised as null
Definition at line 64 of file ReactingParcelI.H.
Referenced by ReactingParcel< ParcelType >::constantProperties::constantProperties(), ReactingParcel< ParcelType >::ReactingParcel(), and ReactingParcel< ParcelType >::updateMassFraction().
|
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.
References ReactingParcel< ParcelType >::ReactingParcel().
|
inline |
Construct from components.
Definition at line 96 of file ReactingParcelI.H.
References ReactingParcel< ParcelType >::mass0_.
ReactingParcel | ( | const polyMesh & | mesh, |
Istream & | is, | ||
bool | readFields = true |
||
) |
Construct from Istream.
Definition at line 46 of file ReactingParcelIO.C.
References IOstream::check(), IOstream::format(), Istream::read(), and DynamicList< T, SizeInc, SizeMult, SizeDiv >::transfer().
ReactingParcel | ( | const ReactingParcel< ParcelType > & | p, |
const polyMesh & | mesh | ||
) |
Construct as a copy.
ReactingParcel | ( | const ReactingParcel< ParcelType > & | p | ) |
Construct as a copy.
|
protected |
Calculate Phase change.
Definition at line 39 of file ReactingParcel.C.
References subModelBase::active(), PhaseChangeModel< CloudType >::addToPhaseChangeMass(), PhaseChangeModel< CloudType >::calculate(), CompositionModel< CloudType >::carrier(), composition, PhaseChangeModel< CloudType >::dh(), forAll, CompositionModel< CloudType >::liquids(), CompositionModel< CloudType >::localToCarrierId(), Foam::min(), Foam::constant::thermodynamic::RR, Foam::sum(), PhaseChangeModel< CloudType >::TMax(), PhaseChangeModel< CloudType >::Tvap(), and ReactingParcel< ParcelType >::updateMassFraction().
|
protected |
Update mass fraction.
Definition at line 151 of file ReactingParcel.C.
References forAll, ReactingParcel< ParcelType >::mass0_, mesh, p, ReactingParcel< ParcelType >::pc_, ReactingParcel< ParcelType >::ReactingParcel(), ReactingParcel< ParcelType >::setCellValues(), Foam::sum(), and ReactingParcel< ParcelType >::Y_.
Referenced by ReactingParcel< ParcelType >::calcPhaseChange().
TypeName | ( | "ReactingParcel< ParcelType >" | ) |
Runtime type information.
AddToPropertyList | ( | ParcelType | , |
" mass0"+" nPhases(Y1..YN)" | |||
) |
String representation of properties.
Construct and return a (basic particle) clone.
Definition at line 285 of file ReactingParcel.H.
Construct and return a (basic particle) clone.
Definition at line 291 of file ReactingParcel.H.
References mesh.
|
inline |
Return const access to initial mass [kg].
Definition at line 161 of file ReactingParcelI.H.
References ReactingParcel< ParcelType >::mass0_.
Referenced by ReactingParcel< ParcelType >::iNew::operator()().
|
inline |
Return const access to mass fractions of mixture [].
Definition at line 168 of file ReactingParcelI.H.
References ReactingParcel< ParcelType >::Y_.
Referenced by ReactingParcel< ParcelType >::iNew::operator()(), and ReactingParcel< ParcelType >::writeFields().
|
inline |
Return the owner cell pressure [Pa].
Definition at line 175 of file ReactingParcelI.H.
References ReactingParcel< ParcelType >::pc_.
Referenced by ReactingParcel< ParcelType >::iNew::operator()().
|
inline |
Return reference to the owner cell pressure [Pa].
Definition at line 182 of file ReactingParcelI.H.
References ReactingParcel< ParcelType >::pc_.
|
inline |
Return access to initial mass [kg].
Definition at line 189 of file ReactingParcelI.H.
References ReactingParcel< ParcelType >::mass0_.
|
inline |
Return access to mass fractions of mixture [].
Definition at line 196 of file ReactingParcelI.H.
References ReactingParcel< ParcelType >::Y_.
void setCellValues | ( | TrackData & | td, |
const scalar | dt, | ||
const label | celli | ||
) |
Set cell values.
Definition at line 206 of file ReactingParcel.C.
References ReactingParcel< ParcelType >::cellValueSourceCorrection(), coordinates(), Foam::endl(), Foam::nl, and WarningInFunction.
Referenced by ReactingParcel< ParcelType >::iNew::operator()(), and ReactingParcel< ParcelType >::updateMassFraction().
void cellValueSourceCorrection | ( | TrackData & | td, |
const scalar | dt, | ||
const label | celli | ||
) |
Correct cell values using latest transfer information.
Definition at line 237 of file ReactingParcel.C.
References ReactingParcel< ParcelType >::correctSurfaceValues(), Foam::endl(), forAll, Foam::mag(), Foam::max(), Foam::nl, and WarningInFunction.
Referenced by ReactingParcel< ParcelType >::iNew::operator()(), and ReactingParcel< ParcelType >::setCellValues().
void correctSurfaceValues | ( | TrackData & | td, |
const label | celli, | ||
const scalar | T, | ||
const scalarField & | Cs, | ||
scalar & | rhos, | ||
scalar & | mus, | ||
scalar & | Prs, | ||
scalar & | kappas | ||
) |
Correct surface values due to emitted species.
Definition at line 298 of file ReactingParcel.C.
References ReactingParcel< ParcelType >::calc(), SLGThermo::carrier(), Foam::cbrt(), basicSpecieMixture::Cp(), forAll, basicSpecieMixture::kappa(), Foam::max(), Foam::min(), basicSpecieMixture::mu(), Foam::constant::thermodynamic::RR, List< T >::size(), basicMultiComponentMixture::species(), Foam::sqrt(), Foam::sum(), Foam::T(), basicSpecieMixture::W(), and basicMultiComponentMixture::Y().
Referenced by ReactingParcel< ParcelType >::cellValueSourceCorrection(), and ReactingParcel< ParcelType >::iNew::operator()().
void calc | ( | TrackData & | td, |
const scalar | dt, | ||
const label | celli | ||
) |
Update parcel properties over the time interval.
Definition at line 390 of file ReactingParcel.C.
References CompositionModel< CloudType >::carrier(), Foam::cbrt(), composition, CompositionModel< CloudType >::Cp(), forAll, CompositionModel< CloudType >::localToCarrierId(), Foam::constant::mathematical::pi(), Foam::pow4(), Foam::Re(), Su, and Foam::Zero.
Referenced by ReactingParcel< ParcelType >::correctSurfaceValues(), and ReactingParcel< ParcelType >::iNew::operator()().
|
static |
Read.
Definition at line 98 of file ReactingParcelIO.C.
References Cloud< ParticleType >::checkFieldIOobject(), Cloud< ParticleType >::fieldIOobject(), forAll, forAllIter, ReactingParcel< ParcelType >::mass0_, p, Foam::readFields(), List< T >::size(), Cloud< ParticleType >::size(), Y, and ReactingParcel< ParcelType >::Y_.
Referenced by ReactingParcel< ParcelType >::iNew::operator()(), and ReactingParcel< ParcelType >::readFields().
|
static |
Read - no composition.
Definition at line 89 of file ReactingParcelIO.C.
References Foam::readFields(), and ReactingParcel< ParcelType >::readFields().
|
static |
Write.
Definition at line 172 of file ReactingParcelIO.C.
References Cloud< ParticleType >::fieldIOobject(), forAll, forAllConstIter(), ReactingParcel< ParcelType >::mass0_, p, List< T >::size(), Cloud< ParticleType >::size(), Y, and ReactingParcel< ParcelType >::Y().
Referenced by ReactingParcel< ParcelType >::iNew::operator()(), and ReactingParcel< ParcelType >::writeFields().
|
static |
Write - composition supplied.
Definition at line 163 of file ReactingParcelIO.C.
References ReactingParcel< ParcelType >::writeFields().
|
friend |
|
protected |
Initial mass [kg].
Definition at line 161 of file ReactingParcel.H.
Referenced by ReactingParcel< ParcelType >::mass0(), ReactingParcel< ParcelType >::ReactingParcel(), ReactingParcel< ParcelType >::readFields(), ReactingParcel< ParcelType >::updateMassFraction(), and ReactingParcel< ParcelType >::writeFields().
|
protected |
Mass fractions of mixture [].
Definition at line 164 of file ReactingParcel.H.
Referenced by ReactingParcel< ParcelType >::readFields(), ReactingParcel< ParcelType >::updateMassFraction(), and ReactingParcel< ParcelType >::Y().
|
protected |
Pressure [Pa].
Definition at line 170 of file ReactingParcel.H.
Referenced by ReactingParcel< ParcelType >::pc(), and ReactingParcel< ParcelType >::updateMassFraction().