Classes | Public Member Functions | Static Public Member Functions | Protected Attributes | Friends | List of all members
DSMCParcel< ParcelType > Class Template Reference

DSMC parcel class. More...

Inheritance diagram for DSMCParcel< ParcelType >:
Inheritance graph
[legend]
Collaboration diagram for DSMCParcel< ParcelType >:
Collaboration graph
[legend]

Classes

class  constantProperties
 Class to hold DSMC particle constant properties. More...
 
class  iNew
 Factory class to read-construct particles used for. More...
 
class  trackingData
 Class used to pass kinematic tracking data to the trackToFace function. More...
 

Public Member Functions

 TypeName ("DSMCParcel")
 Runtime type information. More...
 
 DSMCParcel (const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti, const vector &U, const scalar Ei, const label typeId)
 Construct from components. 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 (const polyMesh &mesh, 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 TrackData >
bool move (TrackData &td, const scalar trackTime)
 Move the parcel. More...
 
template<class TrackData >
bool hitPatch (const polyPatch &, TrackData &td, const label patchi, const scalar trackFraction, const tetIndices &tetIs)
 Overridable function to handle the particle hitting a patch. More...
 
template<class TrackData >
void hitProcessorPatch (const processorPolyPatch &, TrackData &td)
 Overridable function to handle the particle hitting a. More...
 
template<class TrackData >
void hitWallPatch (const wallPolyPatch &, TrackData &td, const tetIndices &)
 Overridable function to handle the particle hitting a wallPatch. More...
 
template<class TrackData >
void hitPatch (const polyPatch &, TrackData &td)
 Overridable function to handle the particle hitting a polyPatch. More...
 
virtual void transformProperties (const tensor &T)
 Transform the physical properties of the particle. More...
 
virtual void transformProperties (const vector &separation)
 Transform the physical properties of the particle. More...
 

Static Public Member Functions

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 51 of file DSMCParcel.H.

Constructor & Destructor Documentation

◆ DSMCParcel() [1/3]

DSMCParcel ( const polyMesh mesh,
const barycentric coordinates,
const label  celli,
const label  tetFacei,
const label  tetPti,
const vector U,
const scalar  Ei,
const label  typeId 
)
inline

Construct from components.

Definition at line 56 of file DSMCParcelI.H.

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

Here is the caller graph for this function:

◆ DSMCParcel() [2/3]

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 76 of file DSMCParcelI.H.

◆ DSMCParcel() [3/3]

DSMCParcel ( const polyMesh mesh,
Istream is,
bool  readFields = true 
)

Construct from Istream.

Definition at line 44 of file DSMCParcelIO.C.

References IOstream::check(), IOstream::format(), Istream::read(), Foam::readLabel(), and readScalar.

Here is the call graph for this function:

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 210 of file DSMCParcel.H.

◆ typeId()

Foam::label typeId ( ) const
inline

Return type id.

Definition at line 137 of file DSMCParcelI.H.

References DSMCParcel< ParcelType >::typeId_.

Referenced by DSMCParcel< ParcelType >::iNew::operator()(), and DSMCParcel< ParcelType >::writeFields().

Here is the caller graph for this function:

◆ U() [1/2]

const Foam::vector & U ( ) const
inline

Return const access to velocity.

Definition at line 144 of file DSMCParcelI.H.

References DSMCParcel< ParcelType >::U_.

Referenced by DSMCParcel< ParcelType >::iNew::operator()(), and DSMCParcel< ParcelType >::writeFields().

Here is the caller graph for this function:

◆ Ei() [1/2]

Foam::scalar Ei ( ) const
inline

Return const access to internal energy.

Definition at line 151 of file DSMCParcelI.H.

References DSMCParcel< ParcelType >::Ei_.

Referenced by DSMCParcel< ParcelType >::iNew::operator()(), and DSMCParcel< ParcelType >::writeFields().

Here is the caller graph for this function:

◆ U() [2/2]

Foam::vector & U ( )
inline

Return access to velocity.

Definition at line 158 of file DSMCParcelI.H.

References DSMCParcel< ParcelType >::U_.

◆ Ei() [2/2]

Foam::scalar & Ei ( )
inline

Return access to internal energy.

Definition at line 165 of file DSMCParcelI.H.

References DSMCParcel< ParcelType >::Ei_.

◆ move()

bool move ( TrackData &  td,
const scalar  trackTime 
)

Move the parcel.

Definition at line 33 of file DSMCParcel.C.

References polyMesh::boundaryMesh(), Foam::meshTools::constrainDirection(), f(), DSMCParcel< ParcelType >::hitPatch(), mesh, p, and polyMesh::solutionD().

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

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

◆ hitPatch() [1/2]

bool hitPatch ( const polyPatch ,
TrackData &  td,
const label  patchi,
const scalar  trackFraction,
const tetIndices tetIs 
)

Overridable function to handle the particle hitting a patch.

Executed before other patch-hitting functions

Definition at line 80 of file DSMCParcel.C.

References DSMCParcel< ParcelType >::hitProcessorPatch().

Referenced by DSMCParcel< ParcelType >::move(), and DSMCParcel< ParcelType >::iNew::operator()().

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

◆ hitProcessorPatch()

void hitProcessorPatch ( const processorPolyPatch ,
TrackData &  td 
)

Overridable function to handle the particle hitting a.

processorPatch

Definition at line 95 of file DSMCParcel.C.

References DSMCParcel< ParcelType >::hitWallPatch().

Referenced by DSMCParcel< ParcelType >::hitPatch(), and DSMCParcel< ParcelType >::iNew::operator()().

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

◆ hitWallPatch()

void hitWallPatch ( const wallPolyPatch wpp,
TrackData &  td,
const tetIndices tetIs 
)

Overridable function to handle the particle hitting a wallPatch.

Definition at line 107 of file DSMCParcel.C.

References polyPatch::faceAreas(), patchIdentifier::index(), Foam::mag(), DSMCParcel< ParcelType >::constantProperties::mass(), Foam::max(), nw, and polyPatch::whichFace().

Referenced by DSMCParcel< ParcelType >::hitProcessorPatch(), and DSMCParcel< ParcelType >::iNew::operator()().

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

◆ hitPatch() [2/2]

void hitPatch ( const polyPatch ,
TrackData &  td 
)

Overridable function to handle the particle hitting a polyPatch.

Definition at line 199 of file DSMCParcel.C.

◆ transformProperties() [1/2]

void transformProperties ( const tensor T)
virtual

Transform the physical properties of the particle.

according to the given transformation tensor

Definition at line 206 of file DSMCParcel.C.

References Foam::transform().

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

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

◆ transformProperties() [2/2]

void transformProperties ( const vector separation)
virtual

Transform the physical properties of the particle.

according to the given separation vector

Definition at line 215 of file DSMCParcel.C.

◆ readFields()

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

Definition at line 79 of file DSMCParcelIO.C.

References Foam::constant::universal::c, DSMCParcel< ParcelType >::Ei_, forAllIter, p, Foam::readFields(), DSMCParcel< ParcelType >::typeId_, U, DSMCParcel< ParcelType >::U_, and DSMCParcel< ParcelType >::writeFields().

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

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 113 of file DSMCParcelIO.C.

References Foam::constant::universal::c, DSMCParcel< ParcelType >::Ei(), forAllConstIter(), p, DSMCParcel< ParcelType >::typeId(), U, and DSMCParcel< ParcelType >::U().

Referenced by DSMCParcel< ParcelType >::iNew::operator()(), and DSMCParcel< ParcelType >::readFields().

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

Friends And Related Function Documentation

◆ Cloud< ParcelType >

friend class Cloud< ParcelType >
friend

Definition at line 171 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 156 of file DSMCParcel.H.

Referenced by DSMCParcel< ParcelType >::readFields(), and DSMCParcel< ParcelType >::U().

◆ Ei_

scalar Ei_
protected

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

degrees of freedom [J]

Definition at line 160 of file DSMCParcel.H.

Referenced by DSMCParcel< ParcelType >::Ei(), and DSMCParcel< ParcelType >::readFields().

◆ typeId_

label typeId_
protected

Parcel type id.

Definition at line 163 of file DSMCParcel.H.

Referenced by DSMCParcel< ParcelType >::readFields(), and DSMCParcel< ParcelType >::typeId().


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