particle< Type > Class Template Reference

Base particle class. More...

Inheritance diagram for particle< Type >:
[legend]
Collaboration diagram for particle< Type >:
[legend]

Classes

class  iNew
 Factory class to read-construct particles (for parallel transfer) More...
 
struct  positionsCompat1706
 Old particle positions content for OpenFOAM-1706 and earlier. More...
 
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)
 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 &, const bool readFields=true, const bool newFormat=true, const bool doLocate=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 autoPtr< particleclone () const
 Construct a clone. More...
 
virtual ~particle ()=default
 Destructor. More...
 
label getNewParticleID () const
 Get unique particle creation id. More...
 
const polyMeshmesh () const noexcept
 Return the mesh database. More...
 
const barycentriccoordinates () const noexcept
 Return current particle coordinates. More...
 
label cell () const noexcept
 Return current cell particle is in. More...
 
label & cell () noexcept
 Return current cell particle is in for manipulation. More...
 
label tetFace () const noexcept
 Return current tet face particle is in. More...
 
label & tetFace () noexcept
 Return current tet face particle is in for manipulation. More...
 
label tetPt () const noexcept
 Return current tet face particle is in. More...
 
label & tetPt () noexcept
 Return current tet face particle is in for manipulation. More...
 
label face () const noexcept
 Return current face particle is on otherwise -1. More...
 
label & face () noexcept
 Return current face particle is on for manipulation. More...
 
scalar stepFraction () const noexcept
 Return the fraction of time-step completed. More...
 
scalar & stepFraction () noexcept
 Return the fraction of time-step completed. More...
 
label origProc () const noexcept
 Return the originating processor ID. More...
 
label & origProc () noexcept
 Return the originating processor ID. More...
 
label origId () const noexcept
 Return the particle ID on the originating processor. More...
 
label & origId () noexcept
 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 noexcept
 Return indices of the current tet that the particle occupies. 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 noexcept
 Is the particle on a face? More...
 
bool onInternalFace () const noexcept
 Is the particle on an internal face? More...
 
bool onBoundaryFace () const noexcept
 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 Member Functions

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...
 

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...
 
static bool writeLagrangianCoordinates = true
 
static bool writeLagrangianPositions
 

Protected Member Functions

void readData (Istream &is, point &position, const bool readFields, const bool newFormat, const bool doLocate)
 Read particle from stream. Optionally (for old format) return. 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 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...
 
template<class TrackCloudType >
void hitBoundaryFace (const vector &direction, TrackCloudType &cloud, trackingData &td)
 Dispatch function for boundary face interaction. Calls one of. More...
 

Friends

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

Additional Inherited Members

Detailed Description

template<class Type>
class Foam::particle< Type >

Base particle class.

Definition at line 76 of file particle.H.

Constructor & Destructor Documentation

◆ particle() [1/6]

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

Construct from components.

Definition at line 515 of file particle.C.

◆ particle() [2/6]

particle ( const polyMesh mesh,
const vector position,
const label  celli = -1 
)

Construct from a position and a cell.

Searches for the rest of the required topology.

Definition at line 538 of file particle.C.

References particle< Type >::position().

Here is the call graph for this function:

◆ particle() [3/6]

particle ( const polyMesh mesh,
const vector position,
const label  celli,
const label  tetFacei,
const label  tetPti,
const bool  doLocate = true 
)

Construct from position components.

Definition at line 568 of file particle.C.

References particle< Type >::position().

Here is the call graph for this function:

◆ particle() [4/6]

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

Construct from Istream.

Definition at line 49 of file particleIO.C.

References particle< Type >::position(), particle< Type >::readData(), and particle< Type >::readFields().

Here is the call graph for this function:

◆ particle() [5/6]

particle ( const particle< Type > &  p)

Construct as a copy.

Definition at line 604 of file particle.C.

◆ particle() [6/6]

particle ( const particle< Type > &  p,
const polyMesh mesh 
)

Construct as a copy with reference to a new mesh.

Definition at line 620 of file particle.C.

◆ ~particle()

virtual ~particle ( )
virtualdefault

Destructor.

Member Function Documentation

◆ readData()

void readData ( Istream is,
point position,
const bool  readFields,
const bool  newFormat,
const bool  doLocate 
)
protected

Read particle from stream. Optionally (for old format) return.

read position. Used by construct-from-Istream

Definition at line 75 of file particleIO.C.

References STLCore::ASCII, Istream::beginRawRead(), IOstream::check(), IOstream::checkLabelSize(), IOstream::checkScalarSize(), UList< T >::data(), Istream::endRawRead(), IOstreamOption::format(), FUNCTION_NAME, pTraits< bool >::nComponents, p, Istream::read(), Foam::readFields(), Foam::readRawLabel(), and s().

Referenced by particle< Type >::particle(), and passivePositionParticle::passivePositionParticle().

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

◆ hitPatch()

bool hitPatch ( TrackCloudType &  ,
trackingData  
)
protected

Overridable function to handle the particle hitting a patch.

Executed before other patch-hitting functions.

Definition at line 393 of file particleTemplates.C.

◆ hitWedgePatch()

void hitWedgePatch ( TrackCloudType &  cloud,
trackingData td 
)
protected

Overridable function to handle the particle hitting a wedgePatch.

Definition at line 400 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 
)
protected

Overridable function to handle the particle hitting a.

symmetryPlanePatch

Definition at line 411 of file particleTemplates.C.

◆ hitSymmetryPatch()

void hitSymmetryPatch ( TrackCloudType &  ,
trackingData  
)
protected

Overridable function to handle the particle hitting a symmetryPatch.

Definition at line 422 of file particleTemplates.C.

References Foam::I.

◆ hitCyclicPatch()

void hitCyclicPatch ( TrackCloudType &  ,
trackingData  
)
protected

Overridable function to handle the particle hitting a cyclicPatch.

Definition at line 431 of file particleTemplates.C.

References coupledPolyPatch::forwardT(), cyclicPolyPatch::neighbPatch(), coupledPolyPatch::parallel(), s(), coupledPolyPatch::separated(), coupledPolyPatch::separation(), UList< T >::size(), T, cyclicPolyPatch::transformGlobalFace(), and polyPatch::whichFace().

Here is the call graph for this function:

◆ hitCyclicAMIPatch()

void hitCyclicAMIPatch ( TrackCloudType &  ,
trackingData td,
const vector displacement 
)
protected

◆ hitCyclicACMIPatch()

void hitCyclicACMIPatch ( TrackCloudType &  cloud,
trackingData td,
const vector displacement 
)
protected

Overridable function to handle the particle hitting a.

cyclicACMIPatch

Definition at line 571 of file particleTemplates.C.

References cyclicACMIPolyPatch::mask(), cyclicACMIPolyPatch::nonOverlapPatch(), cyclicAMIPolyPatch::pointFace(), Foam::pos(), polyPatch::start(), cyclicACMIPolyPatch::tolerance(), and polyPatch::whichFace().

Here is the call graph for this function:

◆ hitProcessorPatch()

void hitProcessorPatch ( TrackCloudType &  ,
trackingData  
)
protected

Overridable function to handle the particle hitting a processorPatch.

Definition at line 616 of file particleTemplates.C.

◆ hitWallPatch()

void hitWallPatch ( TrackCloudType &  ,
trackingData  
)
protected

Overridable function to handle the particle hitting a wallPatch.

Definition at line 621 of file particleTemplates.C.

◆ hitBoundaryFace()

void hitBoundaryFace ( const vector direction,
TrackCloudType &  cloud,
trackingData td 
)
protected

Dispatch function for boundary face interaction. Calls one of.

the above (hitWedgePatch, hitCyclicPatch etc) depending on the patch type

Definition at line 294 of file particleTemplates.C.

References particle< Type >::trackingData::keepParticle, and p.

◆ TypeName()

TypeName ( "particle< Type >"  )

Runtime type information.

◆ propertyList()

static string propertyList ( )
inlinestatic

Definition at line 372 of file particle.H.

◆ clone()

virtual autoPtr< particle > clone ( ) const
inlinevirtual

Construct a clone.

Reimplemented in passivePositionParticle, findCellParticle, streamLineParticle, wallBoundedParticle, wallBoundedStreamLineParticle, indexedParticle, injectedParticle, passiveParticle, molecule, solidParticle, trackedParticle, and passivePositionParticle.

Definition at line 436 of file particle.H.

References Time::New().

Here is the call graph for this function:

◆ getNewParticleID()

Foam::label getNewParticleID ( ) const
inline

Get unique particle creation id.

Definition at line 123 of file particleI.H.

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

Here is the call graph for this function:

◆ mesh()

const Foam::polyMesh & mesh ( ) const
inlinenoexcept

Return the mesh database.

Definition at line 137 of file particleI.H.

Referenced by injectedParticle::clone(), wallBoundedParticle::currentEdge(), solidParticle::move(), wallBoundedStreamLineParticle::move(), and wallBoundedParticle::patchInteraction().

Here is the caller graph for this function:

◆ coordinates()

const Foam::barycentric & coordinates ( ) const
inlinenoexcept

Return current particle coordinates.

Definition at line 143 of file particleI.H.

Referenced by lagrangianFieldDecomposer::lagrangianFieldDecomposer(), solidParticle::move(), and lagrangianReconstructor::reconstructPositions().

Here is the caller graph for this function:

◆ cell() [1/2]

Foam::label cell ( ) const
inlinenoexcept

Return current cell particle is in.

Definition at line 149 of file particleI.H.

Referenced by lagrangianReconstructor::reconstructPositions(), and passivePositionParticle::writePosition().

Here is the caller graph for this function:

◆ cell() [2/2]

Foam::label & cell ( )
inlinenoexcept

Return current cell particle is in for manipulation.

Definition at line 155 of file particleI.H.

◆ tetFace() [1/2]

Foam::label tetFace ( ) const
inlinenoexcept

Return current tet face particle is in.

Definition at line 161 of file particleI.H.

Referenced by wallBoundedParticle::currentEdge(), lagrangianFieldDecomposer::lagrangianFieldDecomposer(), and lagrangianReconstructor::reconstructPositions().

Here is the caller graph for this function:

◆ tetFace() [2/2]

Foam::label & tetFace ( )
inlinenoexcept

Return current tet face particle is in for manipulation.

Definition at line 167 of file particleI.H.

◆ tetPt() [1/2]

Foam::label tetPt ( ) const
inlinenoexcept

Return current tet face particle is in.

Definition at line 173 of file particleI.H.

◆ tetPt() [2/2]

Foam::label & tetPt ( )
inlinenoexcept

Return current tet face particle is in for manipulation.

Definition at line 179 of file particleI.H.

◆ face() [1/2]

Foam::label face ( ) const
inlinenoexcept

Return current face particle is on otherwise -1.

Definition at line 185 of file particleI.H.

Referenced by wallBoundedParticle::patchInteraction().

Here is the caller graph for this function:

◆ face() [2/2]

Foam::label & face ( )
inlinenoexcept

Return current face particle is on for manipulation.

Definition at line 191 of file particleI.H.

◆ stepFraction() [1/2]

Foam::scalar stepFraction ( ) const
inlinenoexcept

Return the fraction of time-step completed.

Definition at line 197 of file particleI.H.

Referenced by solidParticle::move(), streamLineParticle::move(), and wallBoundedStreamLineParticle::move().

Here is the caller graph for this function:

◆ stepFraction() [2/2]

Foam::scalar & stepFraction ( )
inlinenoexcept

Return the fraction of time-step completed.

Definition at line 203 of file particleI.H.

◆ origProc() [1/2]

Foam::label origProc ( ) const
inlinenoexcept

Return the originating processor ID.

Definition at line 209 of file particleI.H.

Referenced by Foam::operator==().

Here is the caller graph for this function:

◆ origProc() [2/2]

Foam::label & origProc ( )
inlinenoexcept

Return the originating processor ID.

Definition at line 215 of file particleI.H.

◆ origId() [1/2]

Foam::label origId ( ) const
inlinenoexcept

Return the particle ID on the originating processor.

Definition at line 221 of file particleI.H.

Referenced by Foam::operator==().

Here is the caller graph for this function:

◆ origId() [2/2]

Foam::label & origId ( )
inlinenoexcept

Return the particle ID on the originating processor.

Definition at line 227 of file particleI.H.

◆ stepFractionSpan()

Foam::Pair< Foam::scalar > stepFractionSpan ( ) const
inline

Return the step fraction change within the overall time-step.

Returns the start value and the change as a scalar pair. Always return Pair<scalar>(0, 1), unless sub-cycling is in effect, in which case the values will reflect the span of the sub-cycle within the time-step.

Definition at line 233 of file particleI.H.

References TimeState::deltaTValue(), and dimensioned< Type >::value().

Here is the call graph for this function:

◆ currentTimeFraction()

Foam::scalar currentTimeFraction ( ) const
inline

Return the current fraction within the timestep. This differs.

from the stored step fraction due to sub-cycling.

Definition at line 257 of file particleI.H.

References s().

Here is the call graph for this function:

◆ currentTetIndices()

Foam::tetIndices currentTetIndices ( ) const
inlinenoexcept

Return indices of the current tet that the particle occupies.

Definition at line 265 of file particleI.H.

Referenced by solidParticle::move().

Here is the caller graph for this function:

◆ currentTetTransform()

Foam::barycentricTensor currentTetTransform ( ) const
inline

Return the current tet transformation tensor.

Definition at line 271 of file particleI.H.

◆ normal()

Foam::vector normal ( ) const
inline

The (unit) normal of the tri on tetFacei_ for the current tet.

Definition at line 284 of file particleI.H.

◆ onFace()

bool onFace ( ) const
inlinenoexcept

Is the particle on a face?

Definition at line 290 of file particleI.H.

◆ onInternalFace()

bool onInternalFace ( ) const
inlinenoexcept

Is the particle on an internal face?

Definition at line 296 of file particleI.H.

◆ onBoundaryFace()

bool onBoundaryFace ( ) const
inlinenoexcept

Is the particle on a boundary face?

Definition at line 302 of file particleI.H.

◆ patch()

Foam::label patch ( ) const
inline

Return the index of patch that the particle is on.

Definition at line 308 of file particleI.H.

Referenced by wallBoundedParticle::patchInteraction().

Here is the caller graph for this function:

◆ position()

Foam::vector position ( ) const
inline

Return current particle position.

Definition at line 314 of file particleI.H.

Referenced by injectedParticle::injectedParticle(), wallBoundedStreamLineParticle::interpolateFields(), particle< Type >::particle(), and passivePositionParticle::writePosition().

Here is the caller graph for this function:

◆ reset()

void reset ( )
inline

Reset particle data.

Definition at line 320 of file particleI.H.

◆ track()

Foam::scalar track ( 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 638 of file particle.C.

References f().

Here is the call graph for this function:

◆ trackToFace()

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

As particle::track, but also stops on internal faces.

Definition at line 657 of file particle.C.

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

Referenced by streamLineParticle::move().

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

◆ trackToTri()

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

As particle::trackToFace, but also stops on tet triangles. 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 1039 of file particle.C.

◆ trackToStationaryTri()

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

As particle::trackToTri, but for stationary meshes.

Definition at line 713 of file particle.C.

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

Here is the call graph for this function:

◆ trackToMovingTri()

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

As particle::trackToTri, but for moving meshes.

Definition at line 845 of file particle.C.

References b, Foam::cmptSum(), Foam::endl(), Foam::exit(), Foam::FatalError, FatalErrorInFunction, forAll, Foam::mag(), Foam::max(), mu, Foam::name(), Foam::nl, Foam::Pout, Foam::roots::real, VectorSpace< Form, Cmpt, Ncmpts >::replace(), Foam::sqr(), T, cubicEqn::value(), and Foam::y0().

Here is the call graph for this function:

◆ hitFace()

void hitFace ( const vector direction,
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 351 of file particleTemplates.C.

◆ trackToAndHitFace()

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

Convenience function. Combines trackToFace and hitFace.

Definition at line 378 of file particleTemplates.C.

Referenced by solidParticle::move().

Here is the caller graph for this function:

◆ deviationFromMeshCentre()

Foam::vector deviationFromMeshCentre ( ) 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 1057 of file particle.C.

References Foam::cmptMin(), Foam::meshTools::constrainToMeshCentre(), and Foam::pos().

Here is the call graph for this function:

◆ patchData()

void patchData ( vector n,
vector U 
) const
inline

Get the normal and velocity of the current patch location.

Definition at line 328 of file particleI.H.

References Foam::exit(), Foam::FatalError, FatalErrorInFunction, n, U, and Foam::Zero.

Here is the call 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

Reimplemented in molecule, and solidParticle.

Definition at line 1072 of file particle.C.

Referenced by DSMCParcel< ParcelType >::transformProperties(), CollidingParcel< ParcelType >::transformProperties(), KinematicParcel< ParcelType >::transformProperties(), molecule::transformProperties(), and solidParticle::transformProperties().

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

Reimplemented in molecule, and solidParticle.

Definition at line 1076 of file particle.C.

◆ prepareForParallelTransfer()

void prepareForParallelTransfer ( )

Convert global addressing to the processor patch local equivalents.

Definition at line 1080 of file particle.C.

◆ 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 1087 of file particle.C.

References polyPatch::faceCells(), coupledPolyPatch::forwardT(), coupledPolyPatch::parallel(), s(), coupledPolyPatch::separated(), coupledPolyPatch::separation(), UList< T >::size(), polyPatch::start(), and T.

Referenced by trackedParticle::correctAfterParallelTransfer().

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

◆ prepareForInteractionListReferral()

void prepareForInteractionListReferral ( const vectorTensorTransform transform)

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

particle can be referred.

Definition at line 1136 of file particle.C.

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

Here is the call graph for this function:

◆ correctAfterInteractionListReferral()

void correctAfterInteractionListReferral ( const label  celli)

Correct the topology after referral. The particle may still be.

outside the stored tet and therefore not track-able.

Definition at line 1162 of file particle.C.

References Foam::pos(), and T.

Here is the call graph for this function:

◆ procTetPt()

Foam::label procTetPt ( 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 1199 of file particle.C.

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

Referenced by lagrangianFieldDecomposer::lagrangianFieldDecomposer(), and lagrangianReconstructor::reconstructPositions().

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

◆ autoMap()

void autoMap ( const vector position,
const mapPolyMesh mapper 
)

Map after a topology change.

Definition at line 1225 of file particle.C.

References mapPolyMesh::reverseCellMap().

Here is the call graph for this function:

◆ relocate()

void relocate ( const point position,
const label  celli = -1 
)

Set the addressing based on the provided position.

Definition at line 1242 of file particle.C.

◆ writePropertyName()

void writePropertyName ( Ostream os,
const word name,
const word delim 
)
static

Write the name representation to stream.

Definition at line 44 of file particleTemplates.C.

References Foam::name(), and os().

Here is the call graph for this function:

◆ writeProperty() [1/2]

void writeProperty ( Ostream os,
const word name,
const Type &  value,
const bool  nameOnly,
const word delim,
const wordRes filters = wordRes::null() 
)
static

Write a named particle property to stream, optionally filtered based on its name

Definition at line 70 of file particleTemplates.C.

References UList< T >::empty(), wordRes::match(), Foam::name(), and os().

Here is the call graph for this function:

◆ writeProperty() [2/2]

void writeProperty ( Ostream os,
const word name,
const Field< Type > &  values,
const bool  nameOnly,
const word delim,
const wordRes filters = wordRes::null() 
)
static

Write a named particle property list to stream, optionally filtered based on its name

Definition at line 98 of file particleTemplates.C.

References UList< T >::empty(), forAll, wordRes::match(), Foam::name(), and os().

Here is the call graph for this function:

◆ readFields()

void readFields ( TrackCloudType &  c)
static

Read the fields associated with the owner cloud.

Definition at line 140 of file particleTemplates.C.

References IOobject::MUST_READ, p, and IOobject::typeHeaderOk().

Referenced by findCellParticle::findCellParticle(), particle< Type >::particle(), passivePositionParticle::passivePositionParticle(), and trackedParticle::trackedParticle().

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

◆ writeFields()

void writeFields ( const TrackCloudType &  c)
static

Write the fields associated with the owner cloud.

Definition at line 170 of file particleTemplates.C.

References Foam::exit(), Foam::FatalError, FatalErrorInFunction, Foam::nl, IOobject::NO_READ, p, cloud::POSITIONS, IOPosition< CloudType >::write(), and regIOobject::write().

Here is the call graph for this function:

◆ writeProperties()

void writeProperties ( Ostream os,
const wordRes filters,
const word delim,
const bool  namesOnly 
) const

Write individual particle properties to stream.

Definition at line 219 of file particleIO.C.

References writeProp.

Referenced by injectedParticle::writeProperties().

Here is the caller graph for this function:

◆ readObjects()

void readObjects ( CloudType c,
const objectRegistry obr 
)
static

Read particle fields as objects from the obr registry.

Definition at line 224 of file particleTemplates.C.

References cloud::findIOPosition(), and p.

Referenced by injectedParticle::readObjects().

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

◆ writeObjects()

void writeObjects ( const CloudType c,
objectRegistry obr 
)
static

Write particle fields as objects into the obr registry.

Always writes "position", not "coordinate"

Definition at line 273 of file particleTemplates.C.

References p.

◆ writeCoordinates()

void writeCoordinates ( Ostream os) const

Write the particle barycentric coordinates and cell info.

Definition at line 245 of file particleIO.C.

References STLCore::ASCII, IOstream::check(), IOstreamOption::format(), FUNCTION_NAME, os(), token::SPACE, and OBJstream::write().

Here is the call graph for this function:

◆ writePosition()

void writePosition ( Ostream os) const
virtual

Write the particle position and cell id.

Reimplemented in injectedParticle, and passivePositionParticle.

Definition at line 264 of file particleIO.C.

References STLCore::ASCII, IOstream::check(), IOstreamOption::format(), FUNCTION_NAME, os(), p, s(), token::SPACE, and OBJstream::write().

Here is the call graph for this function:

Friends And Related Function Documentation

◆ operator<<

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

◆ operator==

bool operator== ( const particle< Type > &  pA,
const particle< Type > &  pB 
)
friend

◆ operator!=

bool operator!= ( const particle< Type > &  pA,
const particle< Type > &  pB 
)
friend

Member Data Documentation

◆ propertyList_

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

String representation of properties.

Definition at line 372 of file particle.H.

◆ particleCount_

Foam::label particleCount_ = 0
static

Cumulative particle counter - used to provide unique ID.

Definition at line 375 of file particle.H.

◆ writeLagrangianCoordinates

bool writeLagrangianCoordinates = true
static

Write particle coordinates file (v1712 and later) Default is true

Definition at line 379 of file particle.H.

Referenced by injectedParticle::writeFields().

◆ writeLagrangianPositions

bool writeLagrangianPositions
static

Write particle positions file (v1706 format and earlier) Default is true (disable in etc/controlDict)

Definition at line 383 of file particle.H.

Referenced by lagrangianReconstructor::reconstructPositions(), and injectedParticle::writeFields().


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