Loading...
Searching...
No Matches
DSMCParcel< ParcelType > Class Template Reference

DSMC parcel class. More...

#include <DSMCParcel.H>

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

Classes

class  constantProperties
 Class to hold DSMC particle constant properties. More...
class  iNew
 Factory class to read-construct particles (for parallel transfer). More...

Public Types

typedef ParcelType::trackingData trackingData
 Use base tracking data.

Public Member Functions

 TypeName ("DSMCParcel")
 Runtime type information.
 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.
 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.
 DSMCParcel (const polyMesh &mesh, Istream &is, bool readFields=true, bool newFormat=true)
 Construct from Istream.
virtual autoPtr< particleclone () const
 Return a clone.
label typeId () const
 Return type id.
const vectorU () const
 Return const access to velocity.
scalar Ei () const
 Return const access to internal energy.
vectorU ()
 Return access to velocity.
scalar & Ei ()
 Return access to internal energy.
template<class TrackCloudType>
bool move (TrackCloudType &cloud, trackingData &td, const scalar trackTime)
 Move the parcel.
template<class TrackCloudType>
bool hitPatch (TrackCloudType &, trackingData &)
 Overridable function to handle the particle hitting a patch.
template<class TrackCloudType>
void hitProcessorPatch (TrackCloudType &, trackingData &)
 Overridable function to handle the particle hitting a.
template<class TrackCloudType>
void hitWallPatch (TrackCloudType &, trackingData &)
 Overridable function to handle the particle hitting a wallPatch.
virtual void transformProperties (const tensor &T)
 Transform the physical properties of the particle.
virtual void transformProperties (const vector &separation)
 Transform the physical properties of the particle.

Static Public Member Functions

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

Static Public Attributes

static const std::size_t sizeofFields
 Size in bytes of the fields.

Protected Attributes

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

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

Member Typedef Documentation

◆ trackingData

template<class ParcelType>
typedef ParcelType::trackingData trackingData

Use base tracking data.

Definition at line 155 of file DSMCParcel.H.

Constructor & Destructor Documentation

◆ DSMCParcel() [1/3]

template<class ParcelType>
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 47 of file DSMCParcelI.H.

References coordinates(), Ei(), Ei_, mesh, typeId(), typeId_, U(), and U_.

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

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

◆ DSMCParcel() [2/3]

template<class ParcelType>
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 67 of file DSMCParcelI.H.

References Ei(), Ei_, mesh, typeId(), typeId_, U(), and U_.

Here is the call graph for this function:

◆ DSMCParcel() [3/3]

template<class ParcelType>
DSMCParcel ( const polyMesh & mesh,
Istream & is,
bool readFields = true,
bool newFormat = true )

Construct from Istream.

Definition at line 39 of file DSMCParcelIO.C.

References IOstreamOption::ASCII, IOstream::check(), Ei_, IOstream::fatalCheckNativeSizes(), IOstreamOption::format(), FUNCTION_NAME, mesh, Istream::read(), readFields(), sizeofFields, typeId_, U_, and Foam::Zero.

Here is the call graph for this function:

Member Function Documentation

◆ TypeName()

template<class ParcelType>
TypeName ( "DSMCParcel< ParcelType >" )

Runtime type information.

◆ clone()

template<class ParcelType>
virtual autoPtr< particle > clone ( ) const
inlinevirtual

Return a clone.

Definition at line 238 of file DSMCParcel.H.

◆ typeId()

◆ U() [1/2]

template<class ParcelType>
const Foam::vector & U ( ) const
inline

Return const access to velocity.

Definition at line 136 of file DSMCParcelI.H.

Referenced by LarsenBorgnakkeVariableHardSphere< CloudType >::collide(), VariableHardSphere< CloudType >::collide(), ORourkeCollision< CloudType >::collideParcels(), TrajectoryCollision< CloudType >::collideParcels(), ORourkeCollision< CloudType >::collideSorted(), DSMCParcel(), DSMCParcel(), PairSpringSliderDashpot< CloudType >::evaluatePair(), InjectionModel< CloudType >::injectSteadyState(), InflationInjection< CloudType >::parcelsToInject(), RecycleInteraction< CloudType >::postEvolve(), readFields(), CellZoneInjection< CloudType >::setProperties(), ConeInjection< CloudType >::setProperties(), ConeNozzleInjection< CloudType >::setProperties(), FieldActivatedInjection< CloudType >::setProperties(), InflationInjection< CloudType >::setProperties(), InjectedParticleDistributionInjection< CloudType >::setProperties(), InjectedParticleInjection< CloudType >::setProperties(), KinematicLookupTableInjection< CloudType >::setProperties(), ManualInjection< CloudType >::setProperties(), NoInjection< CloudType >::setProperties(), PatchFlowRateInjection< CloudType >::setProperties(), PatchInjection< CloudType >::setProperties(), ReactingLookupTableInjection< CloudType >::setProperties(), ReactingMultiphaseLookupTableInjection< CloudType >::setProperties(), LarsenBorgnakkeVariableHardSphere< CloudType >::sigmaTcR(), VariableHardSphere< CloudType >::sigmaTcR(), KinematicSurfaceFilm< CloudType >::splashInteraction(), and writeFields().

Here is the caller graph for this function:

◆ Ei() [1/2]

template<class ParcelType>
Foam::scalar Ei ( ) const
inline

Return const access to internal energy.

Definition at line 143 of file DSMCParcelI.H.

References Ei_.

Referenced by LarsenBorgnakkeVariableHardSphere< CloudType >::collide(), DSMCParcel(), DSMCParcel(), readFields(), and writeFields().

Here is the caller graph for this function:

◆ U() [2/2]

template<class ParcelType>
Foam::vector & U ( )
inline

Return access to velocity.

Definition at line 150 of file DSMCParcelI.H.

References U_.

◆ Ei() [2/2]

template<class ParcelType>
Foam::scalar & Ei ( )
inline

Return access to internal energy.

Definition at line 157 of file DSMCParcelI.H.

References Ei_.

◆ move()

template<class ParcelType>
template<class TrackCloudType>
bool move ( TrackCloudType & cloud,
trackingData & td,
const scalar trackTime )

Move the parcel.

Definition at line 28 of file DSMCParcel.C.

References Foam::meshTools::constrainDirection(), f(), mesh, p, td(), and U_.

Here is the call graph for this function:

◆ hitPatch()

template<class ParcelType>
template<class TrackCloudType>
bool hitPatch ( TrackCloudType & ,
trackingData &  )

Overridable function to handle the particle hitting a patch.

Executed before other patch-hitting functions

Definition at line 70 of file DSMCParcel.C.

◆ hitProcessorPatch()

template<class ParcelType>
template<class TrackCloudType>
void hitProcessorPatch ( TrackCloudType & ,
trackingData & td )

Overridable function to handle the particle hitting a.

processorPatch

Definition at line 78 of file DSMCParcel.C.

References td().

Here is the call graph for this function:

◆ hitWallPatch()

template<class ParcelType>
template<class TrackCloudType>
void hitWallPatch ( TrackCloudType & cloud,
trackingData &  )

◆ transformProperties() [1/2]

template<class ParcelType>
void transformProperties ( const tensor & T)
virtual

Transform the physical properties of the particle.

according to the given transformation tensor

Definition at line 181 of file DSMCParcel.C.

References Foam::T(), Foam::transform(), and U_.

Here is the call graph for this function:

◆ transformProperties() [2/2]

template<class ParcelType>
void transformProperties ( const vector & separation)
virtual

Transform the physical properties of the particle.

according to the given separation vector

Definition at line 189 of file DSMCParcel.C.

◆ readFields()

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

Definition at line 72 of file DSMCParcelIO.C.

References DSMCParcel(), Ei(), IOobjectOption::MUST_READ, p, typeId(), U(), and U.

Referenced by DSMCParcel().

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

◆ writeFields()

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

Definition at line 103 of file DSMCParcelIO.C.

References DSMCParcel(), Ei(), IOobjectOption::NO_READ, p, typeId(), U(), and U.

Here is the call graph for this function:

◆ Cloud< ParcelType >

template<class ParcelType>
friend class Cloud< ParcelType >
friend

Definition at line 187 of file DSMCParcel.H.

◆ operator<<

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

Member Data Documentation

◆ sizeofFields

template<class ParcelType>
const std::size_t sizeofFields
static

Size in bytes of the fields.

Definition at line 74 of file DSMCParcel.H.

Referenced by DSMCParcel().

◆ U_

template<class ParcelType>
vector U_
protected

Velocity of Parcel [m/s].

Definition at line 167 of file DSMCParcel.H.

Referenced by DSMCParcel(), DSMCParcel(), DSMCParcel(), hitWallPatch(), move(), transformProperties(), and U().

◆ Ei_

template<class ParcelType>
scalar Ei_
protected

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

degrees of freedom [J]

Definition at line 174 of file DSMCParcel.H.

Referenced by DSMCParcel(), DSMCParcel(), DSMCParcel(), Ei(), Ei(), and hitWallPatch().

◆ typeId_

template<class ParcelType>
label typeId_
protected

Parcel type id.

Definition at line 179 of file DSMCParcel.H.

Referenced by DSMCParcel(), DSMCParcel(), DSMCParcel(), hitWallPatch(), and typeId().


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