particle Class Reference

Base particle class. More...

Inheritance diagram for particle:
Collaboration diagram for particle:

Classes

class  trackingData
 

Public Member Functions

 TypeName ("particle")
 Runtime type information. More...
 
 particle (const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti, const label facei)
 Construct from components. More...
 
 particle (const polyMesh &mesh, const vector &position, const label celli)
 Construct from a position and a cell, searching for the rest of the. More...
 
 particle (Istream &, bool readFields=true)
 Construct from Istream. More...
 
 particle (const particle &p)
 Construct as a copy. More...
 
virtual autoPtr< particleclone () const
 Construct and return a clone. More...
 
virtual ~particle ()
 Destructor. More...
 
label getNewParticleID () const
 Get unique particle creation id. More...
 
const barycentriccoordinates () const
 Return current particle coordinates. More...
 
label cell () const
 Return current cell particle is in. More...
 
label tetFace () const
 Return current tet face particle is in. More...
 
label tetPt () const
 Return current tet face particle is in. More...
 
label face () const
 Return current face particle is on otherwise -1. More...
 
labelface ()
 Return current face particle is on otherwise -1. 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...
 
labelorigProc ()
 Return the originating processor ID. More...
 
label origId () const
 Return the particle ID on the originating processor. More...
 
labelorigId ()
 Return the particle ID on the originating processor. More...
 
tetIndices currentTetIndices (const polyMesh &mesh) const
 Return the indices of the current tet that the. More...
 
barycentricTensor currentTetTransform (const polyMesh &mesh) const
 Return the current tet transformation tensor. More...
 
vector normal (const polyMesh &mesh) const
 Return the normal of the tri on tetFacei_ for the. More...
 
bool onFace () const
 Is the particle on a face? More...
 
bool onInternalFace (const polyMesh &mesh) const
 Is the particle on an internal face? More...
 
bool onBoundaryFace (const polyMesh &mesh) const
 Is the particle on a boundary face? More...
 
label patch (const polyMesh &mesh) const
 Return the index of patch that the particle is on. More...
 
vector position (const polyMesh &mesh) const
 Return current particle position. More...
 
void reset (const scalar stepFraction)
 Set the step fraction and clear the behind data in preparation. More...
 
scalar track (const polyMesh &mesh, const vector &displacement, const scalar fraction)
 Track along the displacement for a given fraction of the overall. More...
 
scalar trackToCell (const polyMesh &mesh, const vector &displacement, const scalar fraction)
 As particle::track, but stops when a new cell is reached. More...
 
scalar trackToFace (const polyMesh &mesh, const vector &displacement, const scalar fraction)
 As particle::track, but stops when a face is hit. More...
 
scalar trackToTri (const polyMesh &mesh, const vector &displacement, const scalar fraction, label &tetTriI)
 As particle::trackToFace, but stops when a tet triangle is hit. More...
 
scalar trackToStationaryTri (const polyMesh &mesh, const vector &displacement, const scalar fraction, label &tetTriI)
 As particle::trackToTri, but for stationary meshes. More...
 
scalar trackToMovingTri (const polyMesh &mesh, const vector &displacement, const scalar fraction, label &tetTriI)
 As particle::trackToTri, but for moving meshes. More...
 
template<class TrackCloudType >
void hitFace (const vector &displacement, const scalar fraction, TrackCloudType &cloud, trackingData &td)
 Hit the current face. If the current face is internal than this. More...
 
template<class TrackCloudType >
scalar trackToAndHitFace (const vector &displacement, const scalar fraction, TrackCloudType &cloud, trackingData &td)
 Convenience function. Combines trackToFace and hitFace. More...
 
vector deviationFromMeshCentre (const polyMesh &mesh) const
 Get the displacement from the mesh centre. Used to correct the. More...
 
void patchData (const polyMesh &mesh, vector &normal, vector &displacement) const
 Get the normal and displacement of the current patch location. More...
 
virtual void transformProperties (const transformer &)
 Transform the physical properties of the particle. More...
 
template<class TrackCloudType >
void prepareForParallelTransfer (TrackCloudType &, trackingData &)
 Make changes prior to a parallel transfer. Runs either. More...
 
template<class TrackCloudType >
void correctAfterParallelTransfer (TrackCloudType &, trackingData &)
 Make changes following a parallel transfer. Runs either. More...
 
void prepareForProcessorTransfer (trackingData &td)
 Make changes prior to a transfer across a processor boundary. More...
 
void correctAfterProcessorTransfer (trackingData &td)
 Make changes following a transfer across a processor boundary. More...
 
void prepareForNonConformalCyclicTransfer (const polyMesh &mesh, const label sendToPatch, const label sendToPatchFace)
 Make changes prior to a transfer across a non conformal cyclic. More...
 
void correctAfterNonConformalCyclicTransfer (const polyMesh &mesh, const label sendToPatch)
 Make changes following a transfer across a non conformal cyclic. More...
 
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. More...
 
template<class TrackCloudType >
void hitCyclicPatch (TrackCloudType &, trackingData &)
 Overridable function to handle the particle hitting a. More...
 
template<class TrackCloudType >
bool hitNonConformalCyclicPatch (const vector &displacement, const scalar fraction, const label patchi, TrackCloudType &cloud, trackingData &td)
 Overridable function to handle the particle hitting an. More...
 
template<class TrackCloudType >
void hitProcessorPatch (TrackCloudType &, trackingData &)
 Overridable function to handle the particle hitting a. More...
 
template<class TrackCloudType >
void hitWallPatch (TrackCloudType &, trackingData &)
 Overridable function to handle the particle hitting a wallPatch. More...
 
template<class TrackCloudType >
void hitBasicPatch (TrackCloudType &, trackingData &)
 Overridable function to handle the particle hitting a basic. More...
 
void prepareForInteractionListReferral (const polyMesh &mesh, const transformer &transform)
 Break the topology and store the cartesian position so that the. More...
 
void correctAfterInteractionListReferral (const polyMesh &mesh, const label celli)
 Correct the topology after referral. Locates the particle. More...
 
label procTetPt (const polyMesh &mesh, const polyMesh &procMesh, const label procCell, const label procTetFace) const
 Return the tet point appropriate for decomposition or. More...
 
void map (const polyMesh &mesh, const point &position, const label celli)
 Map after a mesh change. More...
 
void writePosition (Ostream &) const
 Write the particle position and cell. More...
 
template<class TrackCloudType >
Foam::scalar trackToAndHitFace (const vector &displacement, const scalar fraction, TrackCloudType &cloud, trackingData &td)
 

Static Public Member Functions

static string propertyList ()
 
static autoPtr< particleNew (Istream &is)
 Construct from Istream and return. More...
 
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...
 

Static Public Attributes

static string propertyList_ = Foam::particle::propertyList()
 String representation of properties. More...
 
static label particleCount_ = 0
 Cumulative particle counter - used to provide unique ID. More...
 

Friends

Ostreamoperator<< (Ostream &, const particle &)
 
bool operator== (const particle &pA, const particle &pB)
 
bool operator!= (const particle &pA, const particle &pB)
 

Detailed Description

Base particle class.

Definition at line 80 of file particle.H.

Constructor & Destructor Documentation

◆ particle() [1/4]

particle ( const polyMesh mesh,
const barycentric coordinates,
const label  celli,
const label  tetFacei,
const label  tetPti,
const label  facei 
)

Construct from components.

Definition at line 477 of file particle.C.

Referenced by particle::clone(), and particle::New().

Here is the caller graph for this function:

◆ particle() [2/4]

particle ( const polyMesh mesh,
const vector position,
const label  celli 
)

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

required topology

Definition at line 500 of file particle.C.

References particle::position().

Here is the call graph for this function:

◆ particle() [3/4]

particle ( Istream is,
bool  readFields = true 
)

Construct from Istream.

Definition at line 46 of file particleIO.C.

References IOstream::ASCII, IOstream::check(), IOstream::format(), Istream::read(), and particle::readFields().

Here is the call graph for this function:

◆ particle() [4/4]

particle ( const particle p)

Construct as a copy.

Definition at line 529 of file particle.C.

◆ ~particle()

virtual ~particle ( )
inlinevirtual

Destructor.

Definition at line 374 of file particle.H.

Member Function Documentation

◆ TypeName()

TypeName ( "particle"  )

Runtime type information.

◆ propertyList()

static string propertyList ( )
inlinestatic

Definition at line 326 of file particle.H.

◆ clone()

virtual autoPtr<particle> clone ( ) const
inlinevirtual

Construct and return a clone.

Reimplemented in sampledSetParticle, trackedParticle, solidParticle, molecule, passiveParticle, indexedParticle, streamlinesParticle, and findCellParticle.

Definition at line 361 of file particle.H.

References particle::particle().

Here is the call graph for this function:

◆ New()

static autoPtr<particle> New ( Istream is)
inlinestatic

Construct from Istream and return.

Definition at line 367 of file particle.H.

References particle::particle().

Here is the call graph for this function:

◆ getNewParticleID()

Foam::label getNewParticleID ( ) const
inline

Get unique particle creation id.

Definition at line 114 of file particleI.H.

References Foam::endl(), Foam::labelMax, and WarningInFunction.

Here is the call graph for this function:

◆ coordinates()

const Foam::barycentric & coordinates ( ) const
inline

Return current particle coordinates.

Definition at line 128 of file particleI.H.

Referenced by solidParticle::move().

Here is the caller graph for this function:

◆ cell()

Foam::label cell ( ) const
inline

Return current cell particle is in.

Definition at line 134 of file particleI.H.

◆ tetFace()

Foam::label tetFace ( ) const
inline

Return current tet face particle is in.

Definition at line 140 of file particleI.H.

◆ tetPt()

Foam::label tetPt ( ) const
inline

Return current tet face particle is in.

Definition at line 146 of file particleI.H.

◆ face() [1/2]

Foam::label face ( ) const
inline

Return current face particle is on otherwise -1.

Definition at line 152 of file particleI.H.

◆ face() [2/2]

Foam::label & face ( )
inline

Return current face particle is on otherwise -1.

Definition at line 158 of file particleI.H.

◆ stepFraction() [1/2]

Foam::scalar stepFraction ( ) const
inline

Return the fraction of time-step completed.

Definition at line 164 of file particleI.H.

Referenced by solidParticle::move(), and streamlinesParticle::move().

Here is the caller graph for this function:

◆ stepFraction() [2/2]

Foam::scalar & stepFraction ( )
inline

Return the fraction of time-step completed.

Definition at line 170 of file particleI.H.

◆ origProc() [1/2]

Foam::label origProc ( ) const
inline

Return the originating processor ID.

Definition at line 176 of file particleI.H.

Referenced by Foam::operator==().

Here is the caller graph for this function:

◆ origProc() [2/2]

Foam::label & origProc ( )
inline

Return the originating processor ID.

Definition at line 182 of file particleI.H.

◆ origId() [1/2]

Foam::label origId ( ) const
inline

Return the particle ID on the originating processor.

Definition at line 188 of file particleI.H.

Referenced by Foam::operator==(), and particle::readFields().

Here is the caller graph for this function:

◆ origId() [2/2]

Foam::label & origId ( )
inline

Return the particle ID on the originating processor.

Definition at line 194 of file particleI.H.

◆ currentTetIndices()

Foam::tetIndices currentTetIndices ( const polyMesh mesh) const
inline

Return the indices of the current tet that the.

particle occupies.

Definition at line 200 of file particleI.H.

Referenced by solidParticle::move().

Here is the caller graph for this function:

◆ currentTetTransform()

Foam::barycentricTensor currentTetTransform ( const polyMesh mesh) const
inline

Return the current tet transformation tensor.

Definition at line 209 of file particleI.H.

References polyMesh::moving().

Here is the call graph for this function:

◆ normal()

Foam::vector normal ( const polyMesh mesh) const
inline

Return the normal of the tri on tetFacei_ for the.

current tet.

Definition at line 225 of file particleI.H.

◆ onFace()

bool onFace ( ) const
inline

Is the particle on a face?

Definition at line 231 of file particleI.H.

Referenced by lineFace::calcSamples().

Here is the caller graph for this function:

◆ onInternalFace()

bool onInternalFace ( const polyMesh mesh) const
inline

Is the particle on an internal face?

Definition at line 237 of file particleI.H.

References primitiveMesh::isInternalFace().

Here is the call graph for this function:

◆ onBoundaryFace()

bool onBoundaryFace ( const polyMesh mesh) const
inline

Is the particle on a boundary face?

Definition at line 243 of file particleI.H.

References primitiveMesh::isInternalFace().

Here is the call graph for this function:

◆ patch()

Foam::label patch ( const polyMesh mesh) const
inline

Return the index of patch that the particle is on.

Definition at line 249 of file particleI.H.

References polyMesh::boundaryMesh(), and polyBoundaryMesh::whichPatch().

Here is the call graph for this function:

◆ position()

Foam::vector position ( const polyMesh mesh) const
inline

Return current particle position.

Definition at line 255 of file particleI.H.

Referenced by nearWallFields::calcAddressing(), lineFace::calcSamples(), and particle::particle().

Here is the caller graph for this function:

◆ reset()

void reset ( const scalar  stepFraction)
inline

Set the step fraction and clear the behind data in preparation.

for a new track

Definition at line 261 of file particleI.H.

◆ track()

Foam::scalar track ( const polyMesh mesh,
const vector displacement,
const scalar  fraction 
)

Track along the displacement for a given fraction of the overall.

step. End when the track is complete, or when a boundary is hit. On exit, stepFraction_ will have been incremented to the current position, and facei_ will be set to the index of the boundary face that was hit, or -1 if the track completed within a cell. The proportion of the displacement still to be completed is returned.

Definition at line 546 of file particle.C.

References Foam::endl(), f(), FUNCTION_NAME, Foam::Info, and Foam::nl.

Here is the call graph for this function:

◆ trackToCell()

Foam::scalar trackToCell ( const polyMesh mesh,
const vector displacement,
const scalar  fraction 
)

As particle::track, but stops when a new cell is reached.

Definition at line 571 of file particle.C.

References Foam::endl(), f(), FUNCTION_NAME, Foam::Info, and Foam::nl.

Here is the call graph for this function:

◆ trackToFace()

Foam::scalar trackToFace ( const polyMesh mesh,
const vector displacement,
const scalar  fraction 
)

As particle::track, but stops when a face is hit.

Definition at line 594 of file particle.C.

References Foam::endl(), f(), FUNCTION_NAME, Foam::Info, Foam::nl, and WarningInFunction.

Referenced by lineFace::calcSamples(), and streamlinesParticle::move().

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

◆ trackToTri()

Foam::scalar trackToTri ( const polyMesh mesh,
const vector displacement,
const scalar  fraction,
label tetTriI 
)

As particle::trackToFace, but stops when a tet triangle is hit.

On exit, tetTriI is set to the index of the tet triangle that was hit, or -1 if the end position was reached.

Definition at line 983 of file particle.C.

References polyMesh::moving().

Here is the call graph for this function:

◆ trackToStationaryTri()

Foam::scalar trackToStationaryTri ( const polyMesh mesh,
const vector displacement,
const scalar  fraction,
label tetTriI 
)

As particle::trackToTri, but for stationary meshes.

Definition at line 649 of file particle.C.

References b, Foam::cmptSum(), Foam::endl(), Foam::Info, Foam::mag(), Foam::max(), Foam::constant::physicoChemical::mu, VectorSpace< Form, Cmpt, Ncmpts >::replace(), Foam::T(), and Foam::y0().

Here is the call graph for this function:

◆ trackToMovingTri()

Foam::scalar trackToMovingTri ( const polyMesh mesh,
const vector displacement,
const scalar  fraction,
label tetTriI 
)

◆ hitFace()

void hitFace ( const vector displacement,
const scalar  fraction,
TrackCloudType &  cloud,
trackingData td 
)

Hit the current face. If the current face is internal than this.

crosses into the next cell. If it is a boundary face then this will interact the particle with the relevant patch.

Definition at line 155 of file particleTemplates.C.

References polyMesh::boundaryMesh(), Foam::endl(), forAll, FUNCTION_NAME, Foam::Info, particle::trackingData::mesh, Foam::nl, and p.

Here is the call graph for this function:

◆ trackToAndHitFace() [1/2]

scalar trackToAndHitFace ( const vector displacement,
const scalar  fraction,
TrackCloudType &  cloud,
trackingData td 
)

Convenience function. Combines trackToFace and hitFace.

Referenced by solidParticle::move().

Here is the caller graph for this function:

◆ deviationFromMeshCentre()

Foam::vector deviationFromMeshCentre ( const polyMesh mesh) const

Get the displacement from the mesh centre. Used to correct the.

particle position in cases with reduced dimensionality. Returns a zero vector for three-dimensional cases.

Definition at line 1002 of file particle.C.

References Foam::cmptMin(), Foam::meshTools::constrainToMeshCentre(), polyMesh::geometricD(), Foam::pos(), and VectorSpace< Form, Cmpt, Ncmpts >::zero.

Here is the call graph for this function:

◆ patchData()

void patchData ( const polyMesh mesh,
vector normal,
vector displacement 
) const
inline

Get the normal and displacement of the current patch location.

Definition at line 269 of file particleI.H.

References Foam::exit(), Foam::FatalError, FatalErrorInFunction, polyMesh::moving(), and Foam::Zero.

Here is the call graph for this function:

◆ transformProperties()

void transformProperties ( const transformer )
virtual

Transform the physical properties of the particle.

according to the given transformation tensor

Reimplemented in solidParticle, molecule, and streamlinesParticle.

Definition at line 1020 of file particle.C.

Referenced by molecule::transformProperties(), and solidParticle::transformProperties().

Here is the caller graph for this function:

◆ prepareForParallelTransfer()

void prepareForParallelTransfer ( TrackCloudType &  cloud,
trackingData td 
)

Make changes prior to a parallel transfer. Runs either.

processor or nonConformalCyclic variant below.

Definition at line 106 of file particleTemplates.C.

References particle::trackingData::mesh, particle::trackingData::sendFromPatch, and particle::trackingData::sendToPatchFace.

◆ correctAfterParallelTransfer()

void correctAfterParallelTransfer ( TrackCloudType &  cloud,
trackingData td 
)

Make changes following a parallel transfer. Runs either.

processor or nonConformalCyclic variant below.

Definition at line 129 of file particleTemplates.C.

References polyMesh::boundaryMesh(), Foam::exit(), Foam::FatalError, FatalErrorInFunction, particle::trackingData::mesh, and particle::trackingData::sendToPatch.

Referenced by trackedParticle::correctAfterParallelTransfer(), and sampledSetParticle::correctAfterParallelTransfer().

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

◆ prepareForProcessorTransfer()

void prepareForProcessorTransfer ( trackingData td)

Make changes prior to a transfer across a processor boundary.

Stores the local patch face index (in facei_) so that the mesh face index can be determined on the other side.

Definition at line 1024 of file particle.C.

References particle::trackingData::sendToPatchFace.

◆ correctAfterProcessorTransfer()

void correctAfterProcessorTransfer ( trackingData td)

Make changes following a transfer across a processor boundary.

Converts the stored patch index to a mesh index. Accounts for the receiving face being reversed relative to the sending face.

Definition at line 1031 of file particle.C.

References polyMesh::boundaryMesh(), polyPatch::faceCells(), polyMesh::faces(), particle::trackingData::mesh, particle::trackingData::sendToPatch, List< T >::size(), polyPatch::start(), processorPolyPatch::transform(), and transformer::transformsPosition().

Here is the call graph for this function:

◆ prepareForNonConformalCyclicTransfer()

void prepareForNonConformalCyclicTransfer ( const polyMesh mesh,
const label  sendToPatch,
const label  sendToPatchFace 
)

Make changes prior to a transfer across a non conformal cyclic.

boundary. Stores the receiving patch face (in facei_). Breaks the topology and stores the cartesian position.

Definition at line 1065 of file particle.C.

References polyMesh::boundaryMesh(), Foam::cmptSum(), transformer::invTransformPosition(), nonConformalCyclicPolyPatch::nbrPatch(), Foam::pos(), nonConformalCyclicPolyPatch::transform(), and transformer::transformsPosition().

Here is the call graph for this function:

◆ correctAfterNonConformalCyclicTransfer()

void correctAfterNonConformalCyclicTransfer ( const polyMesh mesh,
const label  sendToPatch 
)

Make changes following a transfer across a non conformal cyclic.

boundary. Locates the particle using the stored face index and cartesian position.

Definition at line 1101 of file particle.C.

References polyMesh::boundaryMesh(), polyMesh::faceOwner(), patchIdentifier::name(), nonConformalCyclicPolyPatch::nbrPatch(), nonConformalPolyPatch::origPatch(), and polyPatch::start().

Here is the call graph for this function:

◆ hitPatch()

bool hitPatch ( TrackCloudType &  ,
trackingData  
)

Overridable function to handle the particle hitting a patch.

Executed before other patch-hitting functions.

Definition at line 261 of file particleTemplates.C.

◆ hitWedgePatch()

void hitWedgePatch ( TrackCloudType &  cloud,
trackingData td 
)

Overridable function to handle the particle hitting a wedgePatch.

Definition at line 268 of file particleTemplates.C.

References Foam::abort(), Foam::FatalError, and FatalErrorInFunction.

Here is the call graph for this function:

◆ hitSymmetryPlanePatch()

void hitSymmetryPlanePatch ( TrackCloudType &  cloud,
trackingData td 
)

Overridable function to handle the particle hitting a.

symmetryPlanePatch

Definition at line 279 of file particleTemplates.C.

◆ hitSymmetryPatch()

void hitSymmetryPatch ( TrackCloudType &  ,
trackingData td 
)

Overridable function to handle the particle hitting a.

symmetryPatch

Definition at line 290 of file particleTemplates.C.

References Foam::I, particle::trackingData::mesh, and transformer::rotation().

Here is the call graph for this function:

◆ hitCyclicPatch()

void hitCyclicPatch ( TrackCloudType &  ,
trackingData td 
)

Overridable function to handle the particle hitting a.

cyclicPatch

Definition at line 298 of file particleTemplates.C.

References polyMesh::boundaryMesh(), polyMesh::faceOwner(), polyMesh::faces(), particle::trackingData::mesh, cyclicPolyPatch::nbrPatch(), List< T >::size(), cyclicPolyPatch::transform(), cyclicPolyPatch::transformGlobalFace(), and transformer::transformsPosition().

Referenced by streamlinesParticle::hitCyclicPatch().

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

◆ hitNonConformalCyclicPatch()

bool hitNonConformalCyclicPatch ( const vector displacement,
const scalar  fraction,
const label  patchi,
TrackCloudType &  cloud,
trackingData td 
)

◆ hitProcessorPatch()

void hitProcessorPatch ( TrackCloudType &  cloud,
trackingData td 
)

Overridable function to handle the particle hitting a.

processorPatch

Definition at line 392 of file particleTemplates.C.

References polyMesh::boundaryMesh(), particle::trackingData::mesh, particle::trackingData::sendFromPatch, particle::trackingData::sendToPatch, particle::trackingData::sendToPatchFace, and particle::trackingData::sendToProc.

Referenced by findCellParticle::hitProcessorPatch(), sampledSetParticle::hitProcessorPatch(), and streamlinesParticle::hitProcessorPatch().

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

◆ hitWallPatch()

void hitWallPatch ( TrackCloudType &  ,
trackingData  
)

Overridable function to handle the particle hitting a wallPatch.

Definition at line 403 of file particleTemplates.C.

◆ hitBasicPatch()

void hitBasicPatch ( TrackCloudType &  ,
trackingData td 
)

Overridable function to handle the particle hitting a basic.

patch. Fall-through for the above.

Definition at line 408 of file particleTemplates.C.

References particle::trackingData::keepParticle.

◆ prepareForInteractionListReferral()

void prepareForInteractionListReferral ( const polyMesh mesh,
const transformer transform 
)

Break the topology and store the cartesian position so that the.

particle can be referred.

Definition at line 1139 of file particle.C.

References Foam::cmptSum(), Foam::inv(), Foam::pos(), and Foam::transform().

Here is the call graph for this function:

◆ correctAfterInteractionListReferral()

void correctAfterInteractionListReferral ( const polyMesh mesh,
const label  celli 
)

Correct the topology after referral. Locates the particle.

relative to a nearby cell/tet. The particle may end up outside this cell/tet and cannot therefore be tracked.

Definition at line 1164 of file particle.C.

References primitiveMesh::cells(), polyMesh::moving(), Foam::pos(), and Foam::T().

Here is the call graph for this function:

◆ procTetPt()

Foam::label procTetPt ( const polyMesh mesh,
const polyMesh procMesh,
const label  procCell,
const label  procTetFace 
) const

Return the tet point appropriate for decomposition or.

reconstruction to or from the given mesh.

Definition at line 1204 of file particle.C.

References polyMesh::faceOwner(), polyMesh::faces(), and List< T >::size().

Here is the call graph for this function:

◆ map()

void map ( const polyMesh mesh,
const point position,
const label  celli 
)

Map after a mesh change.

Definition at line 1231 of file particle.C.

◆ readFields()

◆ writeFields()

void writeFields ( const TrackCloudType &  c)
static

Write the fields associated with the owner cloud.

Definition at line 74 of file particleTemplates.C.

References Foam::constant::universal::c, forAllConstIter, IOobject::NO_READ, IOPosition< CloudType >::write(), and regIOobject::write().

Referenced by molecule::writeFields(), solidParticle::writeFields(), and streamlinesParticle::writeFields().

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

◆ writePosition()

void writePosition ( Ostream os) const

Write the particle position and cell.

Definition at line 86 of file particleIO.C.

References IOstream::ASCII, IOstream::check(), IOstream::format(), token::SPACE, and Ostream::write().

Here is the call graph for this function:

◆ trackToAndHitFace() [2/2]

Foam::scalar trackToAndHitFace ( const vector displacement,
const scalar  fraction,
TrackCloudType &  cloud,
trackingData td 
)

Definition at line 239 of file particleTemplates.C.

References Foam::endl(), f(), FUNCTION_NAME, Foam::Info, particle::trackingData::mesh, and Foam::nl.

Here is the call graph for this function:

Friends And Related Function Documentation

◆ operator<<

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

◆ operator==

bool operator== ( const particle pA,
const particle pB 
)
friend

◆ operator!=

bool operator!= ( const particle pA,
const particle pB 
)
friend

Member Data Documentation

◆ propertyList_

Foam::string propertyList_ = Foam::particle::propertyList()
static

String representation of properties.

Definition at line 326 of file particle.H.

◆ particleCount_

Foam::label particleCount_ = 0
static

Cumulative particle counter - used to provide unique ID.

Definition at line 329 of file particle.H.


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