Wrapper around kinematic parcel types to add MPPIC modelling. More...
Classes | |
class | iNew |
Factory class to read-construct particles used for. More... | |
class | trackingData |
Public Member Functions | |
TypeName ("MPPICParcel") | |
Runtime type information. More... | |
AddToPropertyList (ParcelType, "(UCorrectx UCorrecty UCorrectz)") | |
String representation of properties. More... | |
MPPICParcel (const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti) | |
Construct from mesh, coordinates and topology. More... | |
MPPICParcel (const polyMesh &mesh, const vector &position, const label celli) | |
Construct from a position and a cell, searching for the rest of the. More... | |
MPPICParcel (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 &UCorrect0, const typename ParcelType::constantProperties &constProps) | |
Construct from components. More... | |
MPPICParcel (const polyMesh &mesh, Istream &is, bool readFields=true) | |
Construct from Istream. More... | |
MPPICParcel (const MPPICParcel &p) | |
Construct as a copy. More... | |
MPPICParcel (const MPPICParcel &p, const polyMesh &mesh) | |
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... | |
const vector & | UCorrect () const |
Return const access to correction velocity. More... | |
vector & | UCorrect () |
Return access to correction velocity. More... | |
template<class TrackCloudType > | |
bool | move (TrackCloudType &cloud, trackingData &td, const scalar trackTime) |
Move the parcel. More... | |
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... | |
Protected Attributes | |
vector | UCorrect_ |
Velocity correction due to collisions [m/s]. More... | |
Friends | |
Ostream & | operator (Ostream &, const MPPICParcel< ParcelType > &) |
Wrapper around kinematic parcel types to add MPPIC modelling.
Definition at line 52 of file MPPICParcel.H.
|
inline |
Construct from mesh, coordinates and topology.
Other properties initialised as null
Definition at line 30 of file MPPICParcelI.H.
Referenced by MPPICParcel< ParcelType >::clone(), and MPPICParcel< ParcelType >::MPPICParcel().
|
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 45 of file MPPICParcelI.H.
References MPPICParcel< ParcelType >::MPPICParcel().
|
inline |
Construct from components.
Definition at line 58 of file MPPICParcelI.H.
MPPICParcel | ( | const polyMesh & | mesh, |
Istream & | is, | ||
bool | readFields = true |
||
) |
Construct from Istream.
Definition at line 47 of file MPPICParcelIO.C.
References IOstream::check(), IOstream::format(), and Istream::read().
MPPICParcel | ( | const MPPICParcel< ParcelType > & | p | ) |
Construct as a copy.
MPPICParcel | ( | const MPPICParcel< ParcelType > & | p, |
const polyMesh & | mesh | ||
) |
Construct as a copy.
TypeName | ( | "MPPICParcel< ParcelType >" | ) |
Runtime type information.
AddToPropertyList | ( | ParcelType | , |
"(UCorrectx UCorrecty UCorrectz)" | |||
) |
String representation of properties.
Construct and return a (basic particle) clone.
Definition at line 232 of file MPPICParcel.H.
References MPPICParcel< ParcelType >::MPPICParcel().
Construct and return a (basic particle) clone.
Definition at line 238 of file MPPICParcel.H.
References MPPICParcel< ParcelType >::MPPICParcel().
|
inline |
Return const access to correction velocity.
Definition at line 94 of file MPPICParcelI.H.
Referenced by MPPICParcel< ParcelType >::iNew::operator()(), and MPPICParcel< ParcelType >::writeFields().
|
inline |
Return access to correction velocity.
Definition at line 101 of file MPPICParcelI.H.
bool move | ( | TrackCloudType & | cloud, |
trackingData & | td, | ||
const scalar | trackTime | ||
) |
Move the parcel.
Definition at line 58 of file MPPICParcel.C.
References f(), p, MPPICParcel< ParcelType >::trackingData::part(), and Foam::Swap().
Referenced by MPPICParcel< ParcelType >::iNew::operator()().
|
static |
Read.
Definition at line 78 of file MPPICParcelIO.C.
References Cloud< ParticleType >::checkFieldIOobject(), Cloud< ParticleType >::fieldIOobject(), forAllIter, p, Foam::readFields(), Cloud< ParticleType >::size(), and MPPICParcel< ParcelType >::UCorrect_.
Referenced by MPPICParcel< ParcelType >::iNew::operator()().
|
static |
Write.
Definition at line 106 of file MPPICParcelIO.C.
References Cloud< ParticleType >::fieldIOobject(), forAllConstIter(), p, Cloud< ParticleType >::size(), MPPICParcel< ParcelType >::UCorrect(), and regIOobject::write().
Referenced by MPPICParcel< ParcelType >::iNew::operator()().
|
friend |
|
protected |
Velocity correction due to collisions [m/s].
Definition at line 160 of file MPPICParcel.H.
Referenced by MPPICParcel< ParcelType >::readFields().