DSMCParcel< ParcelType > Class Template Reference

DSMC parcel class. More...

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

Classes

class  constantProperties
 Class to hold DSMC particle constant properties. More...
 

Public Types

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

Public Member Functions

 TypeName ("DSMCParcel")
 Runtime type information. More...
 
 DSMCParcel (const polyMesh &mesh, const vector &position, const label celli, const vector &U, const scalar Ei, const label typeId)
 Construct from a position and a cell, searching for the rest of the. More...
 
 DSMCParcel (Istream &is, bool readFields=true)
 Construct from Istream. More...
 
virtual autoPtr< particleclone () const
 Construct and return a clone. More...
 
label typeId () const
 Return type id. More...
 
const vectorU () const
 Return const access to velocity. More...
 
scalar Ei () const
 Return const access to internal energy. More...
 
vectorU ()
 Return access to velocity. More...
 
scalar & Ei ()
 Return access to internal energy. More...
 
template<class TrackCloudType >
bool move (TrackCloudType &cloud, trackingData &td)
 Move the parcel. More...
 
template<class TrackCloudType >
void hitWallPatch (TrackCloudType &, trackingData &)
 Overridable function to handle the particle hitting a wallPatch. More...
 
virtual void transformProperties (const transformer &)
 Transform the physical properties of the particle. More...
 

Static Public Member Functions

static autoPtr< DSMCParcelNew (Istream &is)
 Construct from Istream and return. More...
 
static void readFields (Cloud< DSMCParcel< ParcelType >> &c)
 
static void writeFields (const Cloud< DSMCParcel< ParcelType >> &c)
 

Protected Attributes

vector U_
 Velocity of Parcel [m/s]. More...
 
scalar Ei_
 Internal energy of the Parcel, covering all non-translational. More...
 
label typeId_
 Parcel type id. More...
 

Friends

class Cloud< ParcelType >
 
Ostreamoperator (Ostream &, const DSMCParcel< ParcelType > &)
 

Detailed Description

template<class ParcelType>
class Foam::DSMCParcel< ParcelType >

DSMC parcel class.

Source files

Definition at line 67 of file DSMCParcel.H.

Member Typedef Documentation

◆ trackingData

typedef ParcelType::trackingData trackingData

Use base tracking data.

Definition at line 130 of file DSMCParcel.H.

Constructor & Destructor Documentation

◆ DSMCParcel() [1/2]

DSMCParcel ( const polyMesh mesh,
const vector position,
const label  celli,
const vector U,
const scalar  Ei,
const label  typeId 
)
inline

Construct from a position and a cell, searching for the rest of the.

required topology

Definition at line 56 of file DSMCParcelI.H.

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

Here is the caller graph for this function:

◆ DSMCParcel() [2/2]

Member Function Documentation

◆ TypeName()

TypeName ( "DSMCParcel< ParcelType >"  )

Runtime type information.

◆ clone()

virtual autoPtr<particle> clone ( ) const
inlinevirtual

Construct and return a clone.

Definition at line 176 of file DSMCParcel.H.

◆ New()

static autoPtr<DSMCParcel> New ( Istream is)
inlinestatic

Construct from Istream and return.

Definition at line 182 of file DSMCParcel.H.

References DSMCParcel< ParcelType >::DSMCParcel().

Here is the call graph for this function:

◆ typeId()

Foam::label typeId
inline

Return type id.

Definition at line 118 of file DSMCParcelI.H.

◆ U() [1/2]

const Foam::vector & U
inline

Return const access to velocity.

Definition at line 125 of file DSMCParcelI.H.

◆ Ei() [1/2]

Foam::scalar Ei
inline

Return const access to internal energy.

Definition at line 132 of file DSMCParcelI.H.

◆ U() [2/2]

Foam::vector & U
inline

Return access to velocity.

Definition at line 139 of file DSMCParcelI.H.

◆ Ei() [2/2]

Foam::scalar & Ei
inline

Return access to internal energy.

Definition at line 146 of file DSMCParcelI.H.

◆ move()

bool move ( TrackCloudType &  cloud,
trackingData td 
)

Move the parcel.

Definition at line 33 of file DSMCParcel.C.

References Foam::meshTools::constrainDirection(), TimeState::deltaTValue(), f(), p, polyMesh::solutionD(), and objectRegistry::time().

Here is the call graph for this function:

◆ hitWallPatch()

void hitWallPatch ( TrackCloudType &  cloud,
trackingData  
)

◆ transformProperties()

void transformProperties ( const transformer transform)
virtual

Transform the physical properties of the particle.

according to the given transformation tensor

Definition at line 168 of file DSMCParcel.C.

References Foam::transform(), and dimensionSet::transform.

Here is the call graph for this function:

◆ readFields()

void readFields ( Cloud< DSMCParcel< ParcelType >> &  c)
static

Definition at line 74 of file DSMCParcelIO.C.

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

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

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

◆ writeFields()

void writeFields ( const Cloud< DSMCParcel< ParcelType >> &  c)
static

Definition at line 107 of file DSMCParcelIO.C.

References Foam::constant::universal::c, forAllConstIter, IOobject::NO_READ, p, U, and regIOobject::write().

Here is the call graph for this function:

Friends And Related Function Documentation

◆ Cloud< ParcelType >

friend class Cloud< ParcelType >
friend

Definition at line 153 of file DSMCParcel.H.

◆ operator

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

Member Data Documentation

◆ U_

vector U_
protected

Velocity of Parcel [m/s].

Definition at line 140 of file DSMCParcel.H.

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

◆ Ei_

scalar Ei_
protected

Internal energy of the Parcel, covering all non-translational.

degrees of freedom [J]

Definition at line 144 of file DSMCParcel.H.

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

◆ typeId_

label typeId_
protected

Parcel type id.

Definition at line 147 of file DSMCParcel.H.

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


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