trackedParticle Class Reference

Particle class that marks cells it passes through. Used to mark cells visited by feature edges. More...

Inheritance diagram for trackedParticle:
[legend]
Collaboration diagram for trackedParticle:
[legend]

Classes

class  iNew
 Factory class to read-construct particles used for. More...
 
class  trackingData
 Class used to pass tracking data to the trackToFace function. More...
 

Public Member Functions

 trackedParticle (const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPtI, const point &end, const label level, const label i, const label j, const label k)
 Construct from components. More...
 
 trackedParticle (const polyMesh &mesh, const vector &position, const label celli, const point &end, const label level, const label i, const label j, const label k)
 Construct from a position and a cell, searching for the rest of the. More...
 
 trackedParticle (const polyMesh &mesh, Istream &is, bool readFields=true, bool newFormat=true)
 Construct from Istream. More...
 
autoPtr< particleclone () const
 Construct and return a clone. More...
 
pointstart ()
 Point to track from. More...
 
pointend ()
 Point to track to. More...
 
label i () const
 Transported label. More...
 
label & i ()
 Transported label. More...
 
label j () const
 Transported label. More...
 
label & j ()
 Transported label. More...
 
label k () const
 Transported label. More...
 
label & k ()
 Transported label. More...
 
bool move (Cloud< trackedParticle > &, trackingData &, const scalar)
 Track all particles to their end point. More...
 
bool hitPatch (Cloud< trackedParticle > &, trackingData &)
 Overridable function to handle the particle hitting a patch. More...
 
void hitWedgePatch (Cloud< trackedParticle > &, trackingData &)
 Overridable function to handle the particle hitting a wedge. More...
 
void hitSymmetryPlanePatch (Cloud< trackedParticle > &, trackingData &)
 Overridable function to handle the particle hitting a. More...
 
void hitSymmetryPatch (Cloud< trackedParticle > &, trackingData &)
 Overridable function to handle the particle hitting a. More...
 
void hitCyclicPatch (Cloud< trackedParticle > &, trackingData &)
 Overridable function to handle the particle hitting a cyclic. More...
 
void hitCyclicAMIPatch (Cloud< trackedParticle > &, trackingData &, const vector &)
 Overridable function to handle the particle hitting a cyclicAMI. More...
 
void hitCyclicACMIPatch (Cloud< trackedParticle > &, trackingData &, const vector &)
 Overridable function to handle the particle hitting a cyclicACMI. More...
 
void hitProcessorPatch (Cloud< trackedParticle > &, trackingData &)
 
void hitWallPatch (Cloud< trackedParticle > &, trackingData &)
 Overridable function to handle the particle hitting a wallPatch. More...
 
void correctAfterParallelTransfer (const label, trackingData &)
 Convert processor patch addressing to the global equivalents. More...
 
- Public Member Functions inherited from particle< Type >
 TypeName ("particle")
 Runtime type information. More...
 
 particle (const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti)
 Construct from components. More...
 
 particle (const polyMesh &mesh, const vector &position, const label celli=-1)
 Construct from a position and a cell. More...
 
 particle (const polyMesh &mesh, const vector &position, const label celli, const label tetFacei, const label tetPti, const bool doLocate=true)
 Construct from position components. More...
 
 particle (const polyMesh &mesh, Istream &, bool readFields=true, bool newFormat=true)
 Construct from Istream. More...
 
 particle (const particle &p)
 Construct as a copy. More...
 
 particle (const particle &p, const polyMesh &mesh)
 Construct as a copy with reference to a new mesh. More...
 
virtual ~particle ()=default
 Destructor. More...
 
label getNewParticleID () const
 Get unique particle creation id. More...
 
const polyMeshmesh () const
 Return the mesh database. More...
 
const barycentriccoordinates () const
 Return current particle coordinates. More...
 
label cell () const
 Return current cell particle is in. More...
 
label & cell ()
 Return current cell particle is in for manipulation. More...
 
label tetFace () const
 Return current tet face particle is in. More...
 
label & tetFace ()
 Return current tet face particle is in for manipulation. More...
 
label tetPt () const
 Return current tet face particle is in. More...
 
label & tetPt ()
 Return current tet face particle is in for manipulation. More...
 
label face () const
 Return current face particle is on otherwise -1. More...
 
label & face ()
 Return current face particle is on for manipulation. More...
 
scalar stepFraction () const
 Return the fraction of time-step completed. More...
 
scalar & stepFraction ()
 Return the fraction of time-step completed. More...
 
label origProc () const
 Return the originating processor ID. More...
 
label & origProc ()
 Return the originating processor ID. More...
 
label origId () const
 Return the particle ID on the originating processor. More...
 
label & origId ()
 Return the particle ID on the originating processor. More...
 
Pair< scalar > stepFractionSpan () const
 Return the step fraction change within the overall time-step. More...
 
scalar currentTimeFraction () const
 Return the current fraction within the timestep. This differs. More...
 
tetIndices currentTetIndices () const
 Return the indices of the current tet that the. More...
 
barycentricTensor currentTetTransform () const
 Return the current tet transformation tensor. More...
 
vector normal () const
 The (unit) normal of the tri on tetFacei_ for the current tet. More...
 
bool onFace () const
 Is the particle on a face? More...
 
bool onInternalFace () const
 Is the particle on an internal face? More...
 
bool onBoundaryFace () const
 Is the particle on a boundary face? More...
 
label patch () const
 Return the index of patch that the particle is on. More...
 
vector position () const
 Return current particle position. More...
 
void reset ()
 Reset particle data. More...
 
scalar track (const vector &displacement, const scalar fraction)
 Track along the displacement for a given fraction of the overall. More...
 
scalar trackToFace (const vector &displacement, const scalar fraction)
 As particle::track, but also stops on internal faces. More...
 
scalar trackToTri (const vector &displacement, const scalar fraction, label &tetTriI)
 As particle::trackToFace, but also stops on tet triangles. On. More...
 
scalar trackToStationaryTri (const vector &displacement, const scalar fraction, label &tetTriI)
 As particle::trackToTri, but for stationary meshes. More...
 
scalar trackToMovingTri (const vector &displacement, const scalar fraction, label &tetTriI)
 As particle::trackToTri, but for moving meshes. More...
 
template<class TrackCloudType >
void hitFace (const vector &direction, TrackCloudType &cloud, trackingData &td)
 Hit the current face. If the current face is internal than this. More...
 
template<class TrackCloudType >
void trackToAndHitFace (const vector &direction, const scalar fraction, TrackCloudType &cloud, trackingData &td)
 Convenience function. Combines trackToFace and hitFace. More...
 
vector deviationFromMeshCentre () const
 Get the displacement from the mesh centre. Used to correct the. More...
 
void patchData (vector &n, vector &U) const
 Get the normal and velocity of the current patch location. 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...
 
void prepareForParallelTransfer ()
 Convert global addressing to the processor patch local equivalents. More...
 
void correctAfterParallelTransfer (const label patchi, trackingData &td)
 Convert processor patch addressing to the global equivalents. More...
 
void prepareForInteractionListReferral (const vectorTensorTransform &transform)
 Break the topology and store the particle position so that the. More...
 
void correctAfterInteractionListReferral (const label celli)
 Correct the topology after referral. The particle may still be. More...
 
label procTetPt (const polyMesh &procMesh, const label procCell, const label procTetFace) const
 Return the tet point appropriate for decomposition or reconstruction. More...
 
void autoMap (const vector &position, const mapPolyMesh &mapper)
 Map after a topology change. More...
 
void relocate (const point &position, const label celli=-1)
 Set the addressing based on the provided position. More...
 
void writeProperties (Ostream &os, const wordRes &filters, const word &delim, const bool namesOnly) const
 Write individual particle properties to stream. More...
 
void writeCoordinates (Ostream &os) const
 Write the particle barycentric coordinates and cell info. More...
 
virtual void writePosition (Ostream &os) const
 Write the particle position and cell id. More...
 

Static Public Attributes

static const std::size_t sizeofFields_
 Size in bytes of the fields. More...
 
- Static Public Attributes inherited from particle< Type >
static string propertyList_ = Foam::particle::propertyList()
 String representation of properties. More...
 
static label particleCount_ = 0
 Cumulative particle counter - used to provide unique ID. More...
 
static bool writeLagrangianCoordinates = true
 
static bool writeLagrangianPositions
 

Friends

class Cloud< trackedParticle >
 
Ostreamoperator<< (Ostream &, const trackedParticle &)
 

Additional Inherited Members

- Static Public Member Functions inherited from particle< Type >
static string propertyList ()
 
template<class Type >
static void writePropertyName (Ostream &os, const word &name, const word &delim)
 Write the name representation to stream. More...
 
template<class Type >
static void writeProperty (Ostream &os, const word &name, const Type &value, const bool nameOnly, const word &delim, const wordRes &filters=wordRes::null())
 
template<class Type >
static void writeProperty (Ostream &os, const word &name, const Field< Type > &values, const bool nameOnly, const word &delim, const wordRes &filters=wordRes::null())
 
template<class TrackCloudType >
static void readFields (TrackCloudType &c)
 Read the fields associated with the owner cloud. More...
 
template<class TrackCloudType >
static void writeFields (const TrackCloudType &c)
 Write the fields associated with the owner cloud. More...
 
template<class CloudType >
static void readObjects (CloudType &c, const objectRegistry &obr)
 Read particle fields as objects from the obr registry. More...
 
template<class CloudType >
static void writeObjects (const CloudType &c, objectRegistry &obr)
 Write particle fields as objects into the obr registry. More...
 
- Protected Member Functions inherited from particle< Type >
template<class TrackCloudType >
bool hitPatch (TrackCloudType &, trackingData &)
 Overridable function to handle the particle hitting a patch. More...
 
template<class TrackCloudType >
void hitWedgePatch (TrackCloudType &, trackingData &)
 Overridable function to handle the particle hitting a wedgePatch. More...
 
template<class TrackCloudType >
void hitSymmetryPlanePatch (TrackCloudType &, trackingData &)
 Overridable function to handle the particle hitting a. More...
 
template<class TrackCloudType >
void hitSymmetryPatch (TrackCloudType &, trackingData &)
 Overridable function to handle the particle hitting a symmetryPatch. More...
 
template<class TrackCloudType >
void hitCyclicPatch (TrackCloudType &, trackingData &)
 Overridable function to handle the particle hitting a cyclicPatch. More...
 
template<class TrackCloudType >
void hitCyclicAMIPatch (TrackCloudType &, trackingData &, const vector &)
 Overridable function to handle the particle hitting a cyclicAMIPatch. More...
 
template<class TrackCloudType >
void hitCyclicACMIPatch (TrackCloudType &, trackingData &, const vector &)
 Overridable function to handle the particle hitting a. More...
 
template<class TrackCloudType >
void hitProcessorPatch (TrackCloudType &, trackingData &)
 Overridable function to handle the particle hitting a processorPatch. More...
 
template<class TrackCloudType >
void hitWallPatch (TrackCloudType &, trackingData &)
 Overridable function to handle the particle hitting a wallPatch. More...
 

Detailed Description

Particle class that marks cells it passes through. Used to mark cells visited by feature edges.

Source files

Definition at line 62 of file trackedParticle.H.

Constructor & Destructor Documentation

◆ trackedParticle() [1/3]

trackedParticle ( const polyMesh mesh,
const barycentric coordinates,
const label  celli,
const label  tetFacei,
const label  tetPtI,
const point end,
const label  level,
const label  i,
const label  j,
const label  k 
)

Construct from components.

Definition at line 42 of file trackedParticle.C.

◆ trackedParticle() [2/3]

trackedParticle ( const polyMesh mesh,
const vector position,
const label  celli,
const point end,
const label  level,
const label  i,
const label  j,
const label  k 
)

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

required topology

Definition at line 66 of file trackedParticle.C.

◆ trackedParticle() [3/3]

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

Construct from Istream.

Definition at line 88 of file trackedParticle.C.

References Istream::beginRawRead(), IOstream::check(), IOstream::checkLabelSize(), IOstream::checkScalarSize(), Istream::endRawRead(), IOstreamOption::format(), FUNCTION_NAME, Istream::read(), Foam::readFields(), and Foam::readRawLabel().

Here is the call graph for this function:

Member Function Documentation

◆ clone()

autoPtr<particle> clone ( ) const
inlinevirtual

Construct and return a clone.

Reimplemented from particle< Type >.

Definition at line 167 of file trackedParticle.H.

◆ start()

point& start ( )
inline

Point to track from.

Definition at line 195 of file trackedParticle.H.

◆ end()

point& end ( )
inline

Point to track to.

Definition at line 201 of file trackedParticle.H.

◆ i() [1/2]

label i ( ) const
inline

Transported label.

Definition at line 207 of file trackedParticle.H.

◆ i() [2/2]

label& i ( )
inline

Transported label.

Definition at line 213 of file trackedParticle.H.

◆ j() [1/2]

label j ( ) const
inline

Transported label.

Definition at line 219 of file trackedParticle.H.

◆ j() [2/2]

label& j ( )
inline

Transported label.

Definition at line 225 of file trackedParticle.H.

◆ k() [1/2]

label k ( ) const
inline

Transported label.

Definition at line 231 of file trackedParticle.H.

◆ k() [2/2]

label& k ( )
inline

Transported label.

Definition at line 237 of file trackedParticle.H.

◆ move()

bool move ( Cloud< trackedParticle > &  cloud,
trackingData td,
const scalar  trackTime 
)

Track all particles to their end point.

Definition at line 134 of file trackedParticle.C.

References f(), particle< Type >::trackingData::keepParticle, Foam::max(), trackedParticle::trackingData::maxLevel_, s, and particle< Type >::trackingData::switchProcessor.

Here is the call graph for this function:

◆ hitPatch()

bool hitPatch ( Cloud< trackedParticle > &  ,
trackingData  
)

Overridable function to handle the particle hitting a patch.

Executed before other patch-hitting functions

Definition at line 172 of file trackedParticle.C.

◆ hitWedgePatch()

void hitWedgePatch ( Cloud< trackedParticle > &  ,
trackingData td 
)

Overridable function to handle the particle hitting a wedge.

Definition at line 179 of file trackedParticle.C.

References particle< Type >::trackingData::keepParticle.

◆ hitSymmetryPlanePatch()

void hitSymmetryPlanePatch ( Cloud< trackedParticle > &  ,
trackingData td 
)

Overridable function to handle the particle hitting a.

symmetry plane

Definition at line 190 of file trackedParticle.C.

References particle< Type >::trackingData::keepParticle.

◆ hitSymmetryPatch()

void hitSymmetryPatch ( Cloud< trackedParticle > &  ,
trackingData td 
)

Overridable function to handle the particle hitting a.

symmetry patch

Definition at line 201 of file trackedParticle.C.

References particle< Type >::trackingData::keepParticle.

◆ hitCyclicPatch()

void hitCyclicPatch ( Cloud< trackedParticle > &  ,
trackingData td 
)

Overridable function to handle the particle hitting a cyclic.

Definition at line 212 of file trackedParticle.C.

References particle< Type >::trackingData::keepParticle.

◆ hitCyclicAMIPatch()

void hitCyclicAMIPatch ( Cloud< trackedParticle > &  ,
trackingData td,
const vector direction 
)

Overridable function to handle the particle hitting a cyclicAMI.

Definition at line 223 of file trackedParticle.C.

References particle< Type >::trackingData::keepParticle.

◆ hitCyclicACMIPatch()

void hitCyclicACMIPatch ( Cloud< trackedParticle > &  ,
trackingData td,
const vector  
)

Overridable function to handle the particle hitting a cyclicACMI.

Definition at line 235 of file trackedParticle.C.

References particle< Type >::trackingData::keepParticle.

◆ hitProcessorPatch()

void hitProcessorPatch ( Cloud< trackedParticle > &  ,
trackingData td 
)

Overridable function to handle the particle hitting a processorPatch

Definition at line 247 of file trackedParticle.C.

References particle< Type >::trackingData::switchProcessor.

◆ hitWallPatch()

void hitWallPatch ( Cloud< trackedParticle > &  ,
trackingData td 
)

Overridable function to handle the particle hitting a wallPatch.

Definition at line 258 of file trackedParticle.C.

References particle< Type >::trackingData::keepParticle.

◆ correctAfterParallelTransfer()

void correctAfterParallelTransfer ( const label  patchi,
trackingData td 
)

Convert processor patch addressing to the global equivalents.

and set the celli to the face-neighbour

Definition at line 269 of file trackedParticle.C.

References particle< Type >::correctAfterParallelTransfer(), trackedParticle::trackingData::featureEdgeVisited_, and k.

Here is the call graph for this function:

Friends And Related Function Documentation

◆ Cloud< trackedParticle >

friend class Cloud< trackedParticle >
friend

Definition at line 89 of file trackedParticle.H.

◆ operator<<

Ostream& operator<< ( Ostream ,
const trackedParticle  
)
friend

Member Data Documentation

◆ sizeofFields_

const std::size_t sizeofFields_
static

Size in bytes of the fields.

Definition at line 123 of file trackedParticle.H.

Referenced by Foam::operator<<().


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