KinematicParcel< ParcelType > Class Template Reference

Kinematic parcel class with rotational motion (as spherical particles only) and one/two-way coupling with the continuous phase. More...

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

Classes

class  constantProperties
 Class to hold kinematic particle constant properties. More...
 
class  iNew
 Factory class to read-construct particles used for. More...
 
class  trackingData
 

Public Member Functions

 TypeName ("KinematicParcel")
 Runtime type information. More...
 
 AddToPropertyList (ParcelType, " active"+" typeId"+" nParticle"+" d"+" dTarget "+" (Ux Uy Uz)"+" rho"+" age"+" tTurb"+" (UTurbx UTurby UTurbz)")
 String representation of properties. More...
 
 KinematicParcel (const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti)
 Construct from mesh, coordinates and topology. More...
 
 KinematicParcel (const polyMesh &mesh, const vector &position, const label celli)
 Construct from a position and a cell, searching for the rest of the. More...
 
 KinematicParcel (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 constantProperties &constProps)
 Construct from components. More...
 
 KinematicParcel (const polyMesh &mesh, Istream &is, bool readFields=true)
 Construct from Istream. More...
 
 KinematicParcel (const KinematicParcel &p)
 Construct as a copy. More...
 
 KinematicParcel (const KinematicParcel &p, const polyMesh &mesh)
 Construct as a copy. More...
 
virtual autoPtr< particleclone () const
 Construct and return a (basic particle) clone. More...
 
virtual autoPtr< particleclone (const polyMesh &mesh) const
 Construct and return a (basic particle) clone. More...
 
bool active () const
 Return const access to active flag. More...
 
label typeId () const
 Return const access to type id. More...
 
scalar nParticle () const
 Return const access to number of particles. More...
 
scalar d () const
 Return const access to diameter. More...
 
scalar dTarget () const
 Return const access to target diameter. More...
 
const vectorU () const
 Return const access to velocity. More...
 
scalar rho () const
 Return const access to density. More...
 
scalar age () const
 Return const access to the age. More...
 
scalar tTurb () const
 Return const access to time spent in turbulent eddy. More...
 
const vectorUTurb () const
 Return const access to turbulent velocity fluctuation. More...
 
bool & active ()
 Return const access to active flag. More...
 
labeltypeId ()
 Return access to type id. More...
 
scalar & nParticle ()
 Return access to number of particles. More...
 
scalar & d ()
 Return access to diameter. More...
 
scalar & dTarget ()
 Return access to target diameter. More...
 
vectorU ()
 Return access to velocity. More...
 
scalar & rho ()
 Return access to density. More...
 
scalar & age ()
 Return access to the age. More...
 
scalar & tTurb ()
 Return access to time spent in turbulent eddy. More...
 
vectorUTurb ()
 Return access to turbulent velocity fluctuation. More...
 
scalar massCell (const trackingData &td) const
 Cell owner mass. More...
 
scalar mass () const
 Particle mass. More...
 
scalar momentOfInertia () const
 Particle moment of inertia around diameter axis. More...
 
scalar volume () const
 Particle volume. More...
 
scalar areaP () const
 Particle projected area. More...
 
scalar areaS () const
 Particle surface area. More...
 
scalar Re (const trackingData &td) const
 Reynolds number. More...
 
scalar We (const trackingData &td, const scalar sigma) const
 Weber number. More...
 
scalar Eo (const trackingData &td, const scalar sigma) const
 Eotvos number. More...
 
template<class TrackCloudType >
void setCellValues (TrackCloudType &cloud, trackingData &td)
 Set cell values. More...
 
template<class TrackCloudType >
void calcDispersion (TrackCloudType &cloud, trackingData &td, const scalar dt)
 Apply dispersion to the carrier phase velocity and update. More...
 
template<class TrackCloudType >
void cellValueSourceCorrection (TrackCloudType &cloud, trackingData &td, const scalar dt)
 Correct cell values using latest transfer information. More...
 
template<class TrackCloudType >
void calc (TrackCloudType &cloud, trackingData &td, const scalar dt)
 Update parcel properties over the time interval. More...
 
template<class TrackCloudType >
bool move (TrackCloudType &cloud, trackingData &td, const scalar trackTime)
 Move the parcel. More...
 
template<class TrackCloudType >
bool hitPatch (TrackCloudType &cloud, trackingData &td)
 Overridable function to handle the particle hitting a patch. More...
 
template<class TrackCloudType >
void hitProcessorPatch (TrackCloudType &cloud, trackingData &td)
 Overridable function to handle the particle hitting a. More...
 
template<class TrackCloudType >
void hitWallPatch (TrackCloudType &cloud, trackingData &td)
 Overridable function to handle the particle hitting a wallPatch. 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...
 
template<class TrackCloudType >
const Foam::vector calcVelocity (TrackCloudType &cloud, trackingData &td, const scalar dt, const scalar Re, const scalar mu, const scalar mass, const vector &Su, vector &dUTrans, scalar &Spu) const
 
template<class CloudType >
void readFields (CloudType &c)
 
template<class CloudType >
void writeFields (const CloudType &c)
 

Static Public Member Functions

static scalar volume (const scalar d)
 Particle volume for a given diameter. More...
 
static scalar areaP (const scalar d)
 Projected area for given diameter. More...
 
static scalar areaS (const scalar d)
 Surface area for given diameter. More...
 
static scalar Re (const scalar rhoc, const vector &U, const vector &Uc, const scalar d, const scalar muc)
 Reynolds number for given conditions. More...
 
static scalar We (const scalar rhoc, const vector &U, const vector &Uc, const scalar d, const scalar sigma)
 Weber number for given conditions. More...
 
static scalar Eo (const vector &g, const scalar rho, const scalar rhoc, const vector &U, const scalar d, const scalar sigma)
 Eotvos number for given conditions. More...
 
template<class TrackCloudType >
static void readFields (TrackCloudType &c)
 Read. More...
 
template<class TrackCloudType >
static void writeFields (const TrackCloudType &c)
 Write. More...
 

Protected Member Functions

template<class TrackCloudType >
const vector calcVelocity (TrackCloudType &cloud, trackingData &td, const scalar dt, const scalar Re, const scalar mu, const scalar mass, const vector &Su, vector &dUTrans, scalar &Spu) const
 Calculate new particle velocity. More...
 

Protected Attributes

bool active_
 Active flag - tracking inactive when active = false. More...
 
label typeId_
 Parcel type id. More...
 
scalar nParticle_
 Number of particles in Parcel. More...
 
scalar d_
 Diameter [m]. More...
 
scalar dTarget_
 Target diameter [m]. More...
 
vector U_
 Velocity of Parcel [m/s]. More...
 
scalar rho_
 Density [kg/m^3]. More...
 
scalar age_
 Age [s]. More...
 
scalar tTurb_
 Time spent in turbulent eddy [s]. More...
 
vector UTurb_
 Turbulent velocity fluctuation [m/s]. More...
 

Friends

Ostreamoperator (Ostream &, const KinematicParcel< ParcelType > &)
 

Detailed Description

template<class ParcelType>
class Foam::KinematicParcel< ParcelType >

Kinematic parcel class with rotational motion (as spherical particles only) and one/two-way coupling with the continuous phase.

Sub-models include:

  • drag
  • turbulent dispersion
  • wall interactions
Source files

Definition at line 60 of file KinematicParcel.H.

Constructor & Destructor Documentation

◆ KinematicParcel() [1/6]

KinematicParcel ( const polyMesh mesh,
const barycentric coordinates,
const label  celli,
const label  tetFacei,
const label  tetPti 
)
inline

Construct from mesh, coordinates and topology.

Other properties initialised as null

Definition at line 74 of file KinematicParcelI.H.

Referenced by KinematicParcel< ParcelType >::calcVelocity(), KinematicParcel< ParcelType >::clone(), KinematicParcel< ParcelType >::constantProperties::constantProperties(), and KinematicParcel< ParcelType >::KinematicParcel().

Here is the caller graph for this function:

◆ KinematicParcel() [2/6]

KinematicParcel ( const polyMesh mesh,
const vector position,
const label  celli 
)
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 98 of file KinematicParcelI.H.

References KinematicParcel< ParcelType >::KinematicParcel().

Here is the call graph for this function:

◆ KinematicParcel() [3/6]

KinematicParcel ( 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 constantProperties constProps 
)
inline

Construct from components.

Definition at line 120 of file KinematicParcelI.H.

◆ KinematicParcel() [4/6]

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

Construct from Istream.

Definition at line 49 of file KinematicParcelIO.C.

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

Here is the call graph for this function:

◆ KinematicParcel() [5/6]

KinematicParcel ( const KinematicParcel< ParcelType > &  p)

Construct as a copy.

◆ KinematicParcel() [6/6]

KinematicParcel ( const KinematicParcel< ParcelType > &  p,
const polyMesh mesh 
)

Construct as a copy.

Member Function Documentation

◆ calcVelocity() [1/2]

const vector calcVelocity ( TrackCloudType &  cloud,
trackingData td,
const scalar  dt,
const scalar  Re,
const scalar  mu,
const scalar  mass,
const vector Su,
vector dUTrans,
scalar &  Spu 
) const
protected

Calculate new particle velocity.

Referenced by KinematicParcel< ParcelType >::calc().

Here is the caller graph for this function:

◆ TypeName()

TypeName ( "KinematicParcel< ParcelType >"  )

Runtime type information.

◆ AddToPropertyList()

AddToPropertyList ( ParcelType  ,
" active"+" typeId"+" nParticle"+" d"+" dTarget "+" (Ux Uy Uz)"+" rho"+" age"+" tTurb"+" (UTurbx UTurby UTurbz)"   
)

String representation of properties.

◆ clone() [1/2]

virtual autoPtr<particle> clone ( ) const
inlinevirtual

Construct and return a (basic particle) clone.

Definition at line 391 of file KinematicParcel.H.

References KinematicParcel< ParcelType >::KinematicParcel().

Here is the call graph for this function:

◆ clone() [2/2]

virtual autoPtr<particle> clone ( const polyMesh mesh) const
inlinevirtual

Construct and return a (basic particle) clone.

Definition at line 397 of file KinematicParcel.H.

References KinematicParcel< ParcelType >::KinematicParcel().

Here is the call graph for this function:

◆ active() [1/2]

bool active ( ) const
inline

Return const access to active flag.

Definition at line 193 of file KinematicParcelI.H.

References KinematicParcel< ParcelType >::active_.

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

Here is the caller graph for this function:

◆ typeId() [1/2]

Foam::label typeId ( ) const
inline

Return const access to type id.

Definition at line 200 of file KinematicParcelI.H.

References KinematicParcel< ParcelType >::typeId_.

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

Here is the caller graph for this function:

◆ nParticle() [1/2]

Foam::scalar nParticle ( ) const
inline

Return const access to number of particles.

Definition at line 207 of file KinematicParcelI.H.

References KinematicParcel< ParcelType >::nParticle_.

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

Here is the caller graph for this function:

◆ d() [1/2]

Foam::scalar d ( ) const
inline

Return const access to diameter.

Definition at line 214 of file KinematicParcelI.H.

References KinematicParcel< ParcelType >::d_.

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

Here is the caller graph for this function:

◆ dTarget() [1/2]

Foam::scalar dTarget ( ) const
inline

Return const access to target diameter.

Definition at line 221 of file KinematicParcelI.H.

References KinematicParcel< ParcelType >::dTarget_.

Referenced by KinematicParcel< ParcelType >::iNew::operator()(), and KinematicParcel< 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 228 of file KinematicParcelI.H.

References KinematicParcel< ParcelType >::U_.

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

Here is the caller graph for this function:

◆ rho() [1/2]

Foam::scalar rho ( ) const
inline

Return const access to density.

Definition at line 235 of file KinematicParcelI.H.

References KinematicParcel< ParcelType >::rho_.

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

Here is the caller graph for this function:

◆ age() [1/2]

Foam::scalar age ( ) const
inline

Return const access to the age.

Definition at line 242 of file KinematicParcelI.H.

References KinematicParcel< ParcelType >::age_.

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

Here is the caller graph for this function:

◆ tTurb() [1/2]

Foam::scalar tTurb ( ) const
inline

Return const access to time spent in turbulent eddy.

Definition at line 249 of file KinematicParcelI.H.

References KinematicParcel< ParcelType >::tTurb_.

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

Here is the caller graph for this function:

◆ UTurb() [1/2]

const Foam::vector & UTurb ( ) const
inline

Return const access to turbulent velocity fluctuation.

Definition at line 256 of file KinematicParcelI.H.

References KinematicParcel< ParcelType >::UTurb_.

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

Here is the caller graph for this function:

◆ active() [2/2]

bool & active ( )
inline

Return const access to active flag.

Definition at line 263 of file KinematicParcelI.H.

References KinematicParcel< ParcelType >::active_.

◆ typeId() [2/2]

Foam::label & typeId ( )
inline

Return access to type id.

Definition at line 270 of file KinematicParcelI.H.

References KinematicParcel< ParcelType >::typeId_.

◆ nParticle() [2/2]

Foam::scalar & nParticle ( )
inline

Return access to number of particles.

Definition at line 277 of file KinematicParcelI.H.

References KinematicParcel< ParcelType >::nParticle_.

◆ d() [2/2]

Foam::scalar & d ( )
inline

Return access to diameter.

Definition at line 284 of file KinematicParcelI.H.

References KinematicParcel< ParcelType >::d_.

◆ dTarget() [2/2]

Foam::scalar & dTarget ( )
inline

Return access to target diameter.

Definition at line 291 of file KinematicParcelI.H.

References KinematicParcel< ParcelType >::dTarget_.

◆ U() [2/2]

Foam::vector & U ( )
inline

Return access to velocity.

Definition at line 298 of file KinematicParcelI.H.

References KinematicParcel< ParcelType >::U_.

◆ rho() [2/2]

Foam::scalar & rho ( )
inline

Return access to density.

Definition at line 305 of file KinematicParcelI.H.

References KinematicParcel< ParcelType >::rho_.

◆ age() [2/2]

Foam::scalar & age ( )
inline

Return access to the age.

Definition at line 312 of file KinematicParcelI.H.

References KinematicParcel< ParcelType >::age_.

◆ tTurb() [2/2]

Foam::scalar & tTurb ( )
inline

Return access to time spent in turbulent eddy.

Definition at line 319 of file KinematicParcelI.H.

References KinematicParcel< ParcelType >::tTurb_.

◆ UTurb() [2/2]

Foam::vector & UTurb ( )
inline

Return access to turbulent velocity fluctuation.

Definition at line 326 of file KinematicParcelI.H.

References KinematicParcel< ParcelType >::massCell(), and KinematicParcel< ParcelType >::UTurb_.

Here is the call graph for this function:

◆ massCell()

Foam::scalar massCell ( const trackingData td) const
inline

Cell owner mass.

Definition at line 334 of file KinematicParcelI.H.

References mesh, and KinematicParcel< ParcelType >::trackingData::rhoc().

Referenced by KinematicParcel< ParcelType >::iNew::operator()(), and KinematicParcel< ParcelType >::UTurb().

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

◆ mass()

Foam::scalar mass ( ) const
inline

Particle mass.

Definition at line 343 of file KinematicParcelI.H.

References KinematicParcel< ParcelType >::rho_, and KinematicParcel< ParcelType >::volume().

Referenced by KinematicParcel< ParcelType >::momentOfInertia(), and KinematicParcel< ParcelType >::iNew::operator()().

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

◆ momentOfInertia()

Foam::scalar momentOfInertia ( ) const
inline

Particle moment of inertia around diameter axis.

Definition at line 350 of file KinematicParcelI.H.

References KinematicParcel< ParcelType >::d_, KinematicParcel< ParcelType >::mass(), and Foam::sqr().

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

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

◆ volume() [1/2]

Foam::scalar volume ( ) const
inline

Particle volume.

Definition at line 357 of file KinematicParcelI.H.

References KinematicParcel< ParcelType >::d_.

Referenced by KinematicParcel< ParcelType >::mass(), and KinematicParcel< ParcelType >::iNew::operator()().

Here is the caller graph for this function:

◆ volume() [2/2]

Foam::scalar volume ( const scalar  d)
inlinestatic

Particle volume for a given diameter.

Definition at line 364 of file KinematicParcelI.H.

References Foam::constant::mathematical::pi(), and Foam::pow3().

Here is the call graph for this function:

◆ areaP() [1/2]

Foam::scalar areaP ( ) const
inline

Particle projected area.

Definition at line 371 of file KinematicParcelI.H.

References KinematicParcel< ParcelType >::d_.

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

Here is the caller graph for this function:

◆ areaP() [2/2]

Foam::scalar areaP ( const scalar  d)
inlinestatic

Projected area for given diameter.

Definition at line 378 of file KinematicParcelI.H.

References KinematicParcel< ParcelType >::areaS().

Here is the call graph for this function:

◆ areaS() [1/2]

Foam::scalar areaS ( ) const
inline

Particle surface area.

Definition at line 385 of file KinematicParcelI.H.

References KinematicParcel< ParcelType >::d_.

Referenced by KinematicParcel< ParcelType >::areaP(), and KinematicParcel< ParcelType >::iNew::operator()().

Here is the caller graph for this function:

◆ areaS() [2/2]

Foam::scalar areaS ( const scalar  d)
inlinestatic

Surface area for given diameter.

Definition at line 392 of file KinematicParcelI.H.

References KinematicParcel< ParcelType >::d(), Foam::constant::mathematical::pi(), and KinematicParcel< ParcelType >::Re().

Here is the call graph for this function:

◆ Re() [1/2]

Foam::scalar Re ( const trackingData td) const
inline

◆ Re() [2/2]

Foam::scalar Re ( const scalar  rhoc,
const vector U,
const vector Uc,
const scalar  d,
const scalar  muc 
)
inlinestatic

Reynolds number for given conditions.

Definition at line 410 of file KinematicParcelI.H.

References Foam::mag(), Foam::max(), and KinematicParcel< ParcelType >::We().

Here is the call graph for this function:

◆ We() [1/2]

Foam::scalar We ( const trackingData td,
const scalar  sigma 
) const
inline

◆ We() [2/2]

Foam::scalar We ( const scalar  rhoc,
const vector U,
const vector Uc,
const scalar  d,
const scalar  sigma 
)
inlinestatic

Weber number for given conditions.

Definition at line 435 of file KinematicParcelI.H.

References KinematicParcel< ParcelType >::Eo(), Foam::magSqr(), and Foam::max().

Here is the call graph for this function:

◆ Eo() [1/2]

Foam::scalar Eo ( const trackingData td,
const scalar  sigma 
) const
inline

◆ Eo() [2/2]

Foam::scalar Eo ( const vector g,
const scalar  rho,
const scalar  rhoc,
const vector U,
const scalar  d,
const scalar  sigma 
)
inlinestatic

Eotvos number for given conditions.

Definition at line 460 of file KinematicParcelI.H.

References Foam::mag(), Foam::max(), and Foam::sqr().

Here is the call graph for this function:

◆ setCellValues()

◆ calcDispersion()

void calcDispersion ( TrackCloudType &  cloud,
trackingData td,
const scalar  dt 
)

Apply dispersion to the carrier phase velocity and update.

parcel turbulence parameters

Definition at line 72 of file KinematicParcel.C.

References KinematicParcel< ParcelType >::cellValueSourceCorrection(), and KinematicParcel< ParcelType >::trackingData::Uc().

Referenced by KinematicParcel< ParcelType >::iNew::operator()(), and KinematicParcel< ParcelType >::setCellValues().

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

◆ cellValueSourceCorrection()

void cellValueSourceCorrection ( TrackCloudType &  cloud,
trackingData td,
const scalar  dt 
)

Correct cell values using latest transfer information.

Definition at line 93 of file KinematicParcel.C.

References KinematicParcel< ParcelType >::calc(), and KinematicParcel< ParcelType >::trackingData::Uc().

Referenced by KinematicParcel< ParcelType >::calcDispersion(), and KinematicParcel< ParcelType >::iNew::operator()().

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

◆ calc()

void calc ( TrackCloudType &  cloud,
trackingData td,
const scalar  dt 
)

Update parcel properties over the time interval.

Definition at line 106 of file KinematicParcel.C.

References KinematicParcel< ParcelType >::calcVelocity(), KinematicParcel< ParcelType >::trackingData::muc(), Foam::Re(), Su, and Foam::Zero.

Referenced by KinematicParcel< ParcelType >::cellValueSourceCorrection(), and KinematicParcel< ParcelType >::iNew::operator()().

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

◆ move()

bool move ( TrackCloudType &  cloud,
trackingData td,
const scalar  trackTime 
)

Move the parcel.

Definition at line 279 of file KinematicParcel.C.

References f(), KinematicParcel< ParcelType >::hitPatch(), Foam::mag(), Foam::max(), maxCo, Foam::min(), and s().

Referenced by KinematicParcel< ParcelType >::calcVelocity(), and KinematicParcel< ParcelType >::iNew::operator()().

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

◆ hitPatch()

bool hitPatch ( TrackCloudType &  cloud,
trackingData td 
)

Overridable function to handle the particle hitting a patch.

Executed before other patch-hitting functions

Definition at line 372 of file KinematicParcel.C.

References polyPatch::boundaryMesh(), polyPatch::coupled(), and KinematicParcel< ParcelType >::hitProcessorPatch().

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

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

◆ hitProcessorPatch()

void hitProcessorPatch ( TrackCloudType &  cloud,
trackingData td 
)

Overridable function to handle the particle hitting a.

processorPatch

Definition at line 407 of file KinematicParcel.C.

References KinematicParcel< ParcelType >::hitWallPatch().

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

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

◆ hitWallPatch()

void hitWallPatch ( TrackCloudType &  cloud,
trackingData td 
)

Overridable function to handle the particle hitting a wallPatch.

Definition at line 419 of file KinematicParcel.C.

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

Here is the caller graph for this function:

◆ 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 429 of file KinematicParcel.C.

References Foam::transform().

Referenced by KinematicParcel< 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 439 of file KinematicParcel.C.

◆ readFields() [1/2]

static void readFields ( TrackCloudType &  c)
static

Read.

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

Here is the caller graph for this function:

◆ writeFields() [1/2]

static void writeFields ( const TrackCloudType &  c)
static

Write.

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

Here is the caller graph for this function:

◆ calcVelocity() [2/2]

◆ readFields() [2/2]

◆ writeFields() [2/2]

Friends And Related Function Documentation

◆ operator

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

Member Data Documentation

◆ active_

bool active_
protected

Active flag - tracking inactive when active = false.

Definition at line 266 of file KinematicParcel.H.

Referenced by KinematicParcel< ParcelType >::active(), KinematicParcel< ParcelType >::calcVelocity(), and KinematicParcel< ParcelType >::readFields().

◆ typeId_

◆ nParticle_

scalar nParticle_
protected

◆ d_

◆ dTarget_

scalar dTarget_
protected

◆ U_

◆ rho_

◆ age_

◆ tTurb_

scalar tTurb_
protected

◆ UTurb_

vector UTurb_
protected

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