Peano
Loading...
Searching...
No Matches
toolbox::particles Namespace Reference

Namespaces

namespace  assignmentchecks
 Correctness checks of particle-mesh assignment.
 
namespace  internal
 
namespace  memorypool
 
namespace  potentials
 
namespace  tests
 

Data Structures

class  FileReader
 A file reader. More...
 
class  FileReaderHDF5
 An HDF5 file reader This class works if and only if you have compiled Peano using –with-hdf5. More...
 
struct  ParticleAccessorIdentity
 
class  ParticleSet
 Abstract base class for all flavours of particle sets. More...
 
class  SieveParticles
 Utility class for global sorting for all flavours of particle sets. More...
 
class  TrajectoryDatabase
 A simple particle database. More...
 

Concepts

concept  ParticleTimeSteppingAccessor
 

Typedefs

typedef int ParticleReassociationInstruction
 

Functions

template<typename Particle >
ParticleReassociationInstruction liftParticleAssociatedWithVertex (const Particle &p, const peano4::datamanagement::VertexMarker marker)
 Take particle from current vertex and lift it or keep it local.
 
ParticleReassociationInstruction liftParticleAssociatedWithVertex (bool isLocal, double searchRadius, const tarch::la::Vector< Dimensions, double > &x, const peano4::datamanagement::VertexMarker marker)
 Implementation of liftParticleAssociatedWithVertex() without the template.
 
template<typename Particle >
bool sieveParticle (const Particle &particle, const peano4::datamanagement::VertexMarker &marker)
 A particle is to be sieved into a vertex if.
 
bool sieveParticle (double searchRadius, const tarch::la::Vector< Dimensions, double > &x, const peano4::datamanagement::VertexMarker &marker)
 
template<typename Particle >
bool dropParticle (const Particle &particle, const peano4::datamanagement::VertexMarker &marker)
 
bool dropParticle (double searchRadius, const tarch::la::Vector< Dimensions, double > &x, const peano4::datamanagement::VertexMarker &marker)
 Alternative version of dropParticle() without templates.
 
template<typename Particle >
bool particleWillBeDroppedFurther (const Particle &particle, const peano4::datamanagement::CellMarker &marker)
 Will the particle be dropped further throughout the traversal.
 
bool particleWillBeDroppedFurther (double searchRadius, const peano4::datamanagement::CellMarker &marker)
 
template<typename Particle >
bool particleWillBeDroppedFurther (const Particle &particle, const peano4::datamanagement::VertexMarker &marker)
 Will a particle be dropped further.
 
bool particleWillBeDroppedFurther (double searchRadius, const peano4::datamanagement::VertexMarker &marker)
 
double relativeSpatialOwnershipTolerance (const ::peano4::datamanagement::CellMarker &marker)
 Ownership tolerance.
 
template<typename Particle >
ParticleReassociationInstruction getParticleReassociationInstructionWithinCellWithIntraCellReassignment (const Particle &p, const peano4::datamanagement::CellMarker &marker, int numberOfVertexWithinCell)
 If you use this operation, you can be sure that every particle is associated to the right vertex after we've finished the traversal, or it is dumped into the set of particles that have to be sieved.
 
ParticleReassociationInstruction getParticleReassociationInstructionWithinCellWithIntraCellReassignment (bool isLocal, double searchRadius, const tarch::la::Vector< Dimensions, double > &x, const peano4::datamanagement::CellMarker &marker, int numberOfVertexWithinCell)
 
std::bitset< TwoPowerDgetAdjacentCellsOwningParticle (const tarch::la::Vector< Dimensions, double > &x, const peano4::datamanagement::VertexMarker &marker)
 Find out which adjacent cell will hold a particle.
 
template<class T >
void init (T &particle, const tarch::la::Vector< Dimensions, double > &x, double searchRadius)
 Init particle.
 
template<class T >
std::vector< T * > createEquallySpacedParticles (double spacingH, const tarch::la::Vector< Dimensions, double > &voxelX, const tarch::la::Vector< Dimensions, double > &voxelH, bool addNoise)
 Create equally spaced particles for one cell.
 
template<class T >
std::vector< T * > createParticlesAlignedWithGlobalCartesianMesh (double spacingH, const tarch::la::Vector< Dimensions, double > &voxelX, const tarch::la::Vector< Dimensions, double > &voxelH, const tarch::la::Vector< Dimensions, double > &domainOffset, const tarch::la::Vector< Dimensions, double > &domainSize, bool addNoise)
 Insert particles that are aligned along a Cartesian grid globally.
 
template<typename T >
void ensureAllParticleListsAreGatheredOrEmpty (peano4::datamanagement::VertexEnumerator< T > vertices)
 Check if all particles are either gathered or host the empt set.
 
template<typename T >
void ensureThatParticleListsEncodeContinuousMemoryChunks (const std::list< T * > activeParticles, const std::vector< int > &numberOfParticlesPerVertexAdded)
 Validate encoding of particle lists.
 
template<typename Iterator >
Iterator particleIsDuplicate (const typename Iterator::value_type &particle, Iterator begin, const Iterator &end)
 Check if there's already a particle at p's position in the set.
 
template<typename Particle >
bool particleWillStayWithinComputationalDomain (const Particle &particle, const tarch::la::Vector< Dimensions, double > &domainOffset, const tarch::la::Vector< Dimensions, double > &domainWidth)
 Return true if and only if.
 
template<typename Particle , typename ParticleSet >
void insertParticlesIntoCell (const peano4::datamanagement::CellMarker &marker, const std::vector< Particle * > &newParticles, peano4::datamanagement::VertexEnumerator< ParticleSet > &fineGridVertices, int spacetreeId)
 Insert particle into a cell.
 
template<typename Particle , typename ParticleSet >
void insertParticleIntoCell (const peano4::datamanagement::CellMarker &marker, Particle *newParticles, peano4::datamanagement::VertexEnumerator< ParticleSet > &fineGridVertices, int spacetreeId)
 Insert particle into cell.
 
void updateContainedInAdjacentCellProperty (const peano4::datamanagement::CellMarker &marker, const std::bitset< Dimensions > &assignedVertex, auto &particle)
 Update helper flag per particle for cell.
 
void updateContainedInAdjacentCellProperty (const peano4::datamanagement::VertexMarker &marker, auto &particle)
 We are.
 
void updateContainedInAdjacentCellProperty (const tarch::la::Vector< Dimensions, double > &x, const tarch::la::Vector< Dimensions, double > &h, auto &particle)
 Update internal marker of particle.
 
bool particleAssignedToVertexWillBeLocal (const tarch::la::Vector< Dimensions, double > &x, const peano4::datamanagement::VertexMarker &marker)
 
tarch::la::Vector< Dimensions, doublemirrorParticleAlongPeriodicDomains (const tarch::la::Vector< Dimensions, double > &x, const peano4::datamanagement::VertexMarker &marker, const tarch::la::Vector< Dimensions, double > domainOffset, const tarch::la::Vector< Dimensions, double > domainSize, const std::bitset< Dimensions > periodicBC)
 Mirror a particle along the periodic domain boundaries.
 
tarch::la::Vector< Dimensions, doubleapplyPeriodicBoundaryConditions (const tarch::la::Vector< Dimensions, double > &x, const tarch::la::Vector< Dimensions, double > domainOffset, const tarch::la::Vector< Dimensions, double > domainSize, const std::bitset< Dimensions > periodicBC)
 Applies the periodic boundary condition on particle coordinates x.
 
template<typename Particle >
void staticPosition (const Particle &particle, double timeStepSize)
 Leave the tracer particles where they are.
 
template<typename Particle >
void explicitEuler (Particle &particle, double timeStepSize)
 Simple explicit Euler.
 
tarch::tests::TestCasegetUnitTests ()
 Please destroy after usage.
 

Variables

constexpr double ReleaseOwnershipSpatialSortingTolerance = 1.01
 Sorting tolerance.
 
constexpr double GrabOwnershipSpatialSortingTolerance = 1.001
 
constexpr int ParticleReassociationInstruction_Keep = -1
 
constexpr int ParticleReassociationInstruction_SieveGlobally = -2
 
constexpr double SpatialDuplicateTolerance = 1e-6
 Default tolerance to decide when two particles are the same.
 

Typedef Documentation

◆ ParticleReassociationInstruction

Function Documentation

◆ applyPeriodicBoundaryConditions()

tarch::la::Vector< Dimensions, double > toolbox::particles::applyPeriodicBoundaryConditions ( const tarch::la::Vector< Dimensions, double > & x,
const tarch::la::Vector< Dimensions, double > domainOffset,
const tarch::la::Vector< Dimensions, double > domainSize,
const std::bitset< Dimensions > periodicBC )

Applies the periodic boundary condition on particle coordinates x.

In mirrorParticleAlongPeriodicDomains(), we check whether a vertex is located on the boundary, and if it is, apply the boundary conditions such that the vertex would "see" the particles. In contrast, here we blindly apply the periodic boundary conditions to a particle. This function is intended to be run on global particle lists, such as is used e.g. in exchangeSieveListsGlobally().

Definition at line 44 of file particles.cpp.

References tarch::la::greater(), periodicBC, and tarch::la::smaller().

Referenced by toolbox::particles::SieveParticles< T >::exchangeSieveListsGlobally().

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

◆ createEquallySpacedParticles()

template<class T >
std::vector< T * > toolbox::particles::createEquallySpacedParticles ( double spacingH,
const tarch::la::Vector< Dimensions, double > & voxelX,
const tarch::la::Vector< Dimensions, double > & voxelH,
bool addNoise )

Create equally spaced particles for one cell.

This routine creates equally spaced particles within one cell. It is a solely cell-local factory method, i.e. it does not take any adjacent cell into account. As a consequence, you will get a locally equidistant layout of particles, but you will not get a globally equidistant layout: The particles will not be spaced the same way across a cell face as they are spaced out inside a cell. If you need a globally equidistant layout, you will have to use createParticlesAlignedWithGlobalCartesianMesh().

The routine creates all the resulting particles on the heap via a plain new(). That is, the particles will be scattered in memory.

Parameters
spacingHSpacing of particles. Pick it smaller or equal to 0 and we'll always create one particle in the center.
voxelXCentre of voxel in space
voxelHSize of voxel, i.e. subcell, into which the particles shall be embedded

Definition at line 28 of file ParticleFactory.cpph.

References assertion, assertion4, dfor, tarch::la::min(), and tarch::la::smallerEquals().

Here is the call graph for this function:

◆ createParticlesAlignedWithGlobalCartesianMesh()

template<class T >
std::vector< T * > toolbox::particles::createParticlesAlignedWithGlobalCartesianMesh ( double spacingH,
const tarch::la::Vector< Dimensions, double > & voxelX,
const tarch::la::Vector< Dimensions, double > & voxelH,
const tarch::la::Vector< Dimensions, double > & domainOffset,
const tarch::la::Vector< Dimensions, double > & domainSize,
bool addNoise )

Insert particles that are aligned along a Cartesian grid globally.

Biased insertion

We compute the left number of the particle and the right one. If the left particle ends up exactly on the face, then we add it. If the right one ends up exactly on the face, we skip it. The rationale is that the right neighbour will insert that guy later on, and we don't want to insert particles twice.

Definition at line 71 of file ParticleFactory.cpph.

References _log, dfor, tarch::la::greaterEquals(), logTraceInWith4Arguments, logTraceOutWith1Argument, tarch::la::min(), tarch::la::remainder(), and tarch::la::smaller().

Here is the call graph for this function:

◆ dropParticle() [1/2]

◆ dropParticle() [2/2]

bool toolbox::particles::dropParticle ( double searchRadius,
const tarch::la::Vector< Dimensions, double > & x,
const peano4::datamanagement::VertexMarker & marker )

Alternative version of dropParticle() without templates.

This is the actual implementation, i.e. the template delegates to this one. Through the method, we can unittest the drop logic.

Definition at line 92 of file MultiscaleTransitions.cpp.

References assertion3, liftParticleAssociatedWithVertex(), tarch::la::NUMERICAL_ZERO_DIFFERENCE, and ParticleReassociationInstruction_SieveGlobally.

Here is the call graph for this function:

◆ ensureAllParticleListsAreGatheredOrEmpty()

template<typename T >
void toolbox::particles::ensureAllParticleListsAreGatheredOrEmpty ( peano4::datamanagement::VertexEnumerator< T > vertices)

Check if all particles are either gathered or host the empt set.

Definition at line 2 of file ParticleListValidation.cpph.

References assertion2, toString(), and TwoPowerD.

Here is the call graph for this function:

◆ ensureThatParticleListsEncodeContinuousMemoryChunks()

template<typename T >
void toolbox::particles::ensureThatParticleListsEncodeContinuousMemoryChunks ( const std::list< T * > activeParticles,
const std::vector< int > & numberOfParticlesPerVertexAdded )

Validate encoding of particle lists.

The flag numberOfParticlesPerVertexAdded encodes how many particles each vertex contributes to a list. If these particles are continuous in memory (per vertex), then we can iterate over activeParticles and validate that the resulting particles match the particles "reconstructed" via pointer arithmetic. If this is not the case, we dump an assertion.

Definition at line 16 of file ParticleListValidation.cpph.

References assertion, and assertionEquals2.

◆ explicitEuler()

template<typename Particle >
void toolbox::particles::explicitEuler ( Particle & particle,
double timeStepSize )

Simple explicit Euler.

We realise the explicit Euler through a policy template technique: The routine is generic as it accepts a template type as particle. However, it does not expect the template argument to have a particular velocity attribute. Instead, it takes a second template argument which takes the particle and provides setters and getters.

While the construct seems to be complicated, it allows us to write the Euler once, but to use it for various particles. Some of them might have a field called velocity, others might not have such a field and get their data from somewhere else.

The term policy is not uniquely determined in C++ and thus this name might not be optimal.

Template argument

ParticleAccesor has to be a class (or struct) which defines a typedef called Particle.

It has to to provide the routines

  • v()
  • v() const
  • x()
  • x() const

which return a copy or provide access to the velocity or position, respectively.

Definition at line 2 of file Tracer.cpph.

◆ getAdjacentCellsOwningParticle()

std::bitset< TwoPowerD > toolbox::particles::getAdjacentCellsOwningParticle ( const tarch::la::Vector< Dimensions, double > & x,
const peano4::datamanagement::VertexMarker & marker )

Find out which adjacent cell will hold a particle.

Let the routine find out which of the neighbouring cells will try to claim ownership of a particle. In exact arithmetics, this routine would be the same as internal::getParticleAssociationWithinCell(). However, we work with floating point precision and cells with a slight overlap. Therefore, multiple adjacent cells might try to claim ownership of the particle potentially, and this is important to find out a priori whether a particle will be labelled as inside or not.

For the floating point reason, i.e. the fact that the particle ownership is not unique, VertexMarker::getRelativeAdjacentCell() is a more strict version of this routine and cannot be used in some cases.

Parameters
particleThe particle of interest to be studied.
xVertex centre.

Definition at line 129 of file MultiscaleTransitions.cpp.

References dfor2, enddforx, peano4::datamanagement::CellMarker::isContained(), tarch::la::multiplyComponents(), and tarch::la::NUMERICAL_ZERO_DIFFERENCE.

Referenced by particleAssignedToVertexWillBeLocal().

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

◆ getParticleReassociationInstructionWithinCellWithIntraCellReassignment() [1/2]

toolbox::particles::ParticleReassociationInstruction toolbox::particles::getParticleReassociationInstructionWithinCellWithIntraCellReassignment ( bool isLocal,
double searchRadius,
const tarch::la::Vector< Dimensions, double > & x,
const peano4::datamanagement::CellMarker & marker,
int numberOfVertexWithinCell )

◆ getParticleReassociationInstructionWithinCellWithIntraCellReassignment() [2/2]

template<typename Particle >
ParticleReassociationInstruction toolbox::particles::getParticleReassociationInstructionWithinCellWithIntraCellReassignment ( const Particle & p,
const peano4::datamanagement::CellMarker & marker,
int numberOfVertexWithinCell )

If you use this operation, you can be sure that every particle is associated to the right vertex after we've finished the traversal, or it is dumped into the set of particles that have to be sieved.

Therefore, you don't have to check any particle anymore in touchVertexLastTime(). All should be in place.

We return

  • ParticleReassociationInstruction_Keep if the particle should remain where it is.
  • A value between 0 and 2^d-1 if we should assign it to another vertex of this cell (on the same level).

Algorithm

  • If the particle is not local or actually not covered by the current cell, we return keep, as we don't dare to make a decision.
  • When we identify the destination cell, we violate our halo discussion above, i.e. we move the association within a cell given harsh boundaries not subject to MultilevelSortingTolerance.

◆ getUnitTests()

tarch::tests::TestCase * toolbox::particles::getUnitTests ( )

Please destroy after usage.

Definition at line 12 of file UnitTests.cpp.

References tarch::tests::TreeTestCaseCollection::addTestCase().

Referenced by runTests().

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

◆ init()

template<class T >
void toolbox::particles::init ( T & particle,
const tarch::la::Vector< Dimensions, double > & x,
double searchRadius )

Init particle.

This is the most generic, minimalist initialisation of a particle that one can think of. The minimal data we need per particle is

  • a position x which defines where the particle is located,
  • a search radius which defines on which resolution level to hold the particle.

The further attribute that is set is the parallel state. At this point, we set it to Virtual. This highlights that fact that you should not add particles directly to the mesh, but instead use insertParticleIntoCell() which then in turn will toggle to Local.

Definition at line 14 of file ParticleFactory.cpph.

◆ insertParticleIntoCell()

template<typename Particle , typename ParticleSet >
void toolbox::particles::insertParticleIntoCell ( const peano4::datamanagement::CellMarker & marker,
Particle * newParticles,
peano4::datamanagement::VertexEnumerator< ParticleSet > & fineGridVertices,
int spacetreeId )

Insert particle into cell.

This routine assumes that the particle indeed fits into a cell. Otherwise, I cannot assign in properly. If it falls into a cell, then I can in return rightfully assume that the particle is local. Setting the flag to local is important, as the locality analysis of Peano usually sets particles to local in touchCellFirstTime(). You don't know if the particle insertion is called before or after this analysis step, so we might just have missed the local insertion and therefore miss out on particles.

I recommend never to insert a particle directly, but always to run through this routine and to combine it with a particle sorting step. On the one hand, thisensures that particles are inserted in a consistent way, i.e. stored next to their closest neighbour. On the other hand, this ensures that a particle is never assigned to a hanging vertex and then lost. If you assign to a hanging vertex in a step where no grid re-sorting is active, you will effectively loose the particle, as it will disappear with the elimination of the hanging vertex.

Realisation

  • Find which vertex in the cell is the closest vertex.
  • Set the particle state to local. As we know that we insert only into local cells, we can do this without any fear.

Definition at line 126 of file particles.cpph.

References _log, assertion3, assertion4, toolbox::particles::assignmentchecks::assignParticleToVertex(), tarch::la::greaterEquals(), logDebug, tarch::la::min(), tarch::la::multiplyComponents(), tarch::la::smallerEquals(), and updateContainedInAdjacentCellProperty().

Referenced by insertParticlesIntoCell().

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

◆ insertParticlesIntoCell()

template<typename Particle , typename ParticleSet >
void toolbox::particles::insertParticlesIntoCell ( const peano4::datamanagement::CellMarker & marker,
const std::vector< Particle * > & newParticles,
peano4::datamanagement::VertexEnumerator< ParticleSet > & fineGridVertices,
int spacetreeId )

Insert particle into a cell.

Shallow copy, i.e. you can remove the input container newParticles, but you may not delete the particles it is pointing to.

This routine basically runs through the vector passed and calls the other insertParticleIntoCell() operation per particle.

Definition at line 53 of file particles.cpph.

References insertParticleIntoCell().

Here is the call graph for this function:

◆ liftParticleAssociatedWithVertex() [1/2]

toolbox::particles::ParticleReassociationInstruction toolbox::particles::liftParticleAssociatedWithVertex ( bool isLocal,
double searchRadius,
const tarch::la::Vector< Dimensions, double > & x,
const peano4::datamanagement::VertexMarker marker )

◆ liftParticleAssociatedWithVertex() [2/2]

template<typename Particle >
ParticleReassociationInstruction toolbox::particles::liftParticleAssociatedWithVertex ( const Particle & p,
const peano4::datamanagement::VertexMarker marker )

Take particle from current vertex and lift it or keep it local.

This routine should be called by touchVertexLastTime(). Different to its cousin getParticleReassociationInstructionWithinCellWithIntraCellReassignment(), it should not be called within the cell. Usually, users do not invoke the routine directly. Instead the action set peano4.toolbox.particles.api.UpdateParticleGridAssociation_ListDrop injects calls to this routine into the generated source code.

The function returns one of the following values:

  • particles::internal::ParticleReassociationInstruction_Keep. All is fine, just analyse the next particles.
  • particles::internal::ParticleReassociationInstruction_SieveGlobally. Particle should not be associated with this vertex, but you can't just move it one level up. Remove it completely from mesh and push it into a global list of particles which shall be sieved in the next step.
  • A value between 0 and \( 2^d-1 \). In this case, the particle should not be associated to this vertex. Remove the association and associate to the next coarsest level. The vertex within this coarse cell is identified by the result number.

Algorithm

The algorithm starts from a simple observation:

  • If a particle is not local, we keep it, as we don't dare to make a comment to which mesh entity is should be moved to. Actually, non-local particles should never be moved. Hence, we could run the association analysis and then check if a non-local particle is kept where it should be, but we do it the other way round: We don't continue to check for virtual particles. We use the local flag as an optimisation opportunity.
  • If a particle is local, it might have moved and it might have to be assigned to another vertex.

Particles are always associated to their closest vertex. If the blue particle moves, it remains within that area around the vertex where it belongs to, and we can return a particles::internal::ParticleReassociationInstruction_Keep command. peano4::datatraversal::VertexMarker::isContainedInAdjacentCells() allows us to make this decision quickly: We check if a particle is contained within the h/2 environment around a vertex. If a particle has left this (greenish) area around a vertex, we have to associate it to another vertex. We have to lift. Lifts however can only be realised within the local spacetree, where we have local and thread-safe access to the coarse level data structures of the tree.

liftParticleAssociatedWithVertex() is called by touchVertexLastTime() which in turn is called from within a cell: We are about to leave a cell, find out that a vertex won't be needed anymore, and consequently trigger touchVertexLastTime(). Each cell has a unique parent within the spacetree.

If the parent cell is not local, we have to assume that we just did hit a vertical tree cut. In this case, we are on the safe side if we just add the particles to the global sieve particle list. We return ParticleReassociationInstruction_SieveGlobally.

If the parent cell is local, we can determine which of the \( 2^d \) coarse grid vertices would hold the particle. If this coarse grid is not adjacent to the local domain, we have to sieve. In the illustration below, we can see this behaviour on the right-hand side:

  • The particle (red solid dot) moves from left to right.
  • Due to the move, the particle's current vertex association (blurry blue) becomes invalid.
  • We have would like to lift it to the red blurry vertex.
  • The parent cell belongs to another rank (yellow instaed of blue).
  • Therefore, we have to add it to the global sieve list.

There is a second, more delicate case which is illustrated on the left-hand side of the sketch:

  • Let the part move from left to right on the fine grid owned by the blue rank.
  • Therefore, the particle changes its vertex association from blurry blue to blurry green.
  • We cannot realise this reassociation directly within touchVertexLastTime(). Therefore, we lift the particle to the red coarse grid vertex.
  • This vertex also is adjacent to the blue rank. Nevertheless, we have to assign is to the global sieve list, as we move it through the horizontal tree cut.

If we reassigned it locally, it would be labelled as remote on the coarser level in the next iteration. The yellow rank in return would receive a copy from the blue-yellow coarse grid boundary and then label it as local there. While this is correct, we now have to sieve it through a horizontal tree cut again. All not nice. We better get it straight and assign it to the global sieve list right from the start.

The need for special care becomes apparent once we take into account which cell triggers the touchVertexLastTime(): Should it be a cell above the two fine grid vertices of interest, we would always associate the particle to the global sieve set, as we spot the multilevel transition. If it is called from one of the two adjacent cells below, the parent cell is local and we cannot spot the tree cut.

Originally, I thought we could get away with an analysis of the coarse vertex and use particleAssignedToVertexWillBeLocal() to figure out whether we lift thorough a tree cut or not. However, we do not have the marker of the coarse grid vertices. So we cannot access the information we'd need. So the only thing that we can do is to lift if and only if we lift within the coarse cell and therefore can be sure that our analysis of the parent cell's owner provides sufficient information.

See also
toolbox::particles::internal::relativeReleaseOwnershipSpatialSortingTolerance()

Referenced by dropParticle(), toolbox::particles::tests::MultiscaleTransitionsTest::testLiftDropOfParticleAssociatedWithVertex01(), toolbox::particles::tests::MultiscaleTransitionsTest::testLiftDropOfParticleAssociatedWithVertex02(), toolbox::particles::tests::MultiscaleTransitionsTest::testLiftDropOfParticleAssociatedWithVertex03(), and toolbox::particles::tests::MultiscaleTransitionsTest::testLiftDropOfParticleAssociatedWithVertex04().

Here is the caller graph for this function:

◆ mirrorParticleAlongPeriodicDomains()

tarch::la::Vector< Dimensions, double > toolbox::particles::mirrorParticleAlongPeriodicDomains ( const tarch::la::Vector< Dimensions, double > & x,
const peano4::datamanagement::VertexMarker & marker,
const tarch::la::Vector< Dimensions, double > domainOffset,
const tarch::la::Vector< Dimensions, double > domainSize,
const std::bitset< Dimensions > periodicBC )

Mirror a particle along the periodic domain boundaries.

We check for every boundary whether it is periodic. Furthermore, we check whether the particle would be already contained within a vertex along each dimension. This is necessary for vertices in corners of the domain: In these cases, we need to mirror the particle only along the boundary that isn't contained in the vertex. Otherwise, the particle is displaced too far away. E.g. For a vertex at (0, 0) a particle at (0.95, 0.) would be displaced to (-0.05, -1.) instead of (-0.5, 0.).

Parameters
xPosition of the particle as received along a periodic boundary exchange.
markerMarker of vertex which triggers the receive, i.e. the one to which the incoming particle at x will be attached to.
Positionof particle as "seen" from receiving side

Definition at line 19 of file particles.cpp.

References tarch::la::equals(), and periodicBC.

Here is the call graph for this function:

◆ particleAssignedToVertexWillBeLocal()

bool toolbox::particles::particleAssignedToVertexWillBeLocal ( const tarch::la::Vector< Dimensions, double > & x,
const peano4::datamanagement::VertexMarker & marker )

◆ particleIsDuplicate()

template<typename Iterator >
Iterator toolbox::particles::particleIsDuplicate ( const typename Iterator::value_type & particle,
Iterator begin,
const Iterator & end )

Check if there's already a particle at p's position in the set.

As we return a non-const iterator, the operation may not be labelled as const. It is essential

Parameters
particleYou want to search the container identified by begin and end for this particle.
Returns
container's end() if there is no duplicate. Otherwise, return iterator pointing to existing particle.

Definition at line 13 of file particles.cpph.

References _log, tarch::la::equals(), logDebug, and SpatialDuplicateTolerance.

Here is the call graph for this function:

◆ particleWillBeDroppedFurther() [1/4]

template<typename Particle >
bool toolbox::particles::particleWillBeDroppedFurther ( const Particle & particle,
const peano4::datamanagement::CellMarker & marker )

Will the particle be dropped further throughout the traversal.

The sorting in Peano usually realises pull semantics. That is, the fine grid vertex pulls particles down from coarser levels into its own domain. Particles recursively drop down the resolution cascade. Throughout the top-down tree traversal, it is often important to know if a particle will be pulled down further later on throughout the mesh traversal. This routine tells the user if this will happen. It is somehow the top-down counterpart of dropParticle() and sieveParticle().

Consult the generic mesh traversal discussion for some further "use cases" where this operation becomes absolutely mandatory.

Realisation

There are two major checks which determine if a particle will be dropped further or not: The mesh size vs the cut-off radius of a particle, and the mesh cell state. If a mesh cell is not refined and will not be refined throughout this very traversal, then the particle will also not be dropped. It will however be dropped if it resides within a cell that's already refined or will be refined in this mesh sweep, and the cut-off radius will also fit into the next finer level.

A cell is refined or will be refined

We note that hasBeenRefined() can be wrong and willBeRefined() might hold. In this case, the new children in the tree are added in this very mesh traversal. It can also happen that hasBeenRefined() holds and willBeRefined() is wrong. In this case, this very mesh traversal will eliminate the children, i.e. coarsen. Nevertheless, the particles still are to be dropped into their respective mesh level - even though we might then lift them straightaway again as we coarsen.

Usage in combination with particle-particle interactions

Deciding whether to update a particle within a mesh sweep that also introduces drops is delicate: Quickly, we run into situations where particles are updated multiple times. The documentation of peano4.toolbox.particles.api.AbstractUpdateParticleGridAssociation provides examples for this.

Definition at line 30 of file MultiscaleTransitions.cpph.

References particleWillBeDroppedFurther().

Referenced by particleWillBeDroppedFurther(), and particleWillBeDroppedFurther().

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

◆ particleWillBeDroppedFurther() [2/4]

template<typename Particle >
bool toolbox::particles::particleWillBeDroppedFurther ( const Particle & particle,
const peano4::datamanagement::VertexMarker & marker )

Will a particle be dropped further.

A a first glance, this routine seems to be a replica of dropParticle(). However, this is not the case. dropParticle() is asked from a child for a parent vertex if it shall drop a particle from the parent into this particular child vertex. This routine is called for a plain vertex without any multiscale relation and tells the callee if a vertex holding this particular particle will later on be stipped off the particle.

Definition at line 36 of file MultiscaleTransitions.cpph.

References particleWillBeDroppedFurther().

Here is the call graph for this function:

◆ particleWillBeDroppedFurther() [3/4]

bool toolbox::particles::particleWillBeDroppedFurther ( double searchRadius,
const peano4::datamanagement::CellMarker & marker )

Definition at line 119 of file MultiscaleTransitions.cpp.

◆ particleWillBeDroppedFurther() [4/4]

bool toolbox::particles::particleWillBeDroppedFurther ( double searchRadius,
const peano4::datamanagement::VertexMarker & marker )

Definition at line 124 of file MultiscaleTransitions.cpp.

◆ particleWillStayWithinComputationalDomain()

template<typename Particle >
bool toolbox::particles::particleWillStayWithinComputationalDomain ( const Particle & particle,
const tarch::la::Vector< Dimensions, double > & domainOffset,
const tarch::la::Vector< Dimensions, double > & domainWidth )

Return true if and only if.

  • we deal with a local particle (virtual ones are copies, to we are not really entitled to say something about them), and
  • the particle resides within the computational domain.

Definition at line 40 of file particles.cpph.

References tarch::la::allGreaterEquals(), and tarch::la::allSmallerEquals().

Here is the call graph for this function:

◆ relativeSpatialOwnershipTolerance()

double toolbox::particles::relativeSpatialOwnershipTolerance ( const ::peano4::datamanagement::CellMarker & marker)

Ownership tolerance.

A particle is "owned" by a cell (or vertex, respectively), if it resides within its area. However, we slightly augment this area, and this tolerance reflects the augmentation. It basically invokes internal::relativeReleaseOwnershipSpatialSortingTolerance() and multiply it with the numerical zero tarch::la::NUMERICAL_ZERO_DIFFERENCE.

Definition at line 13 of file MultiscaleTransitions.cpp.

References tarch::la::NUMERICAL_ZERO_DIFFERENCE.

◆ sieveParticle() [1/2]

template<typename Particle >
bool toolbox::particles::sieveParticle ( const Particle & particle,
const peano4::datamanagement::VertexMarker & marker )

A particle is to be sieved into a vertex if.

  • this vertex is closer to the particle than half the mesh size of an adjacent mesh element, as we always store particles next to their closest vertex; AND
  • this vertex's search radius would not fit into the next finer mesh level OR there is no next finer mesh level; AND
  • the particle's search radius is not too big for the current mesh level.

This routine does not check in any way if the particle will be local.

Rationale and usage

I assume that the tree sieving particles removes a particle from the sieve set as soon as it is inserted into the domain. The association of particles to vertices is not unique. In line with this, the sieving is not unique either. So it is important that you remove the particle from the sieve set once it is inserted.

In theory, we could insert particles on the top level of a segment and then rely on the drops to move them down. This all works fine in a serial code. Once you parallelise, we however encounter an issue when particles drop through vertical tree cuts, i.e. if the fine grid child of a cell is remote. In this case, we might have inserted a particle into a rather coarse resolution level and rely on the drop. The remote child will insert the particle as well. The father will not be able to drop the particle all the way down into a remote tree.

When we sieve a particle into the next level, we have to find out to which particle we do so. To find this out, we run over the vertices as we encounter touchVertexFirstTime(). As particles are not uniquely assigned to cells or vertices, we have some flexibility here. There are two things to weight up:

  • We have to ensure that the particles actually go down;
  • We have to ensure that we don't immediately lift the particle again.

The latter property is essential: Usually, drops, sieves and lifts are done within one action set. However, we always need a second sweep after each lift. If the user runs two sweeps with the action set, the users relies upon the fact that all is sorted completely. If you lift again immediately due to an invalid drop, we mess up this logic. Therefore, it is important that this sorting tolerance is always smaller than the lift tolerance.

We always try to sieve a particle directly into the right level. Once more, it is absolutely essential that we don't sieve into a cell or vertex where we then have to lift stuff again.

Lifting will always be subject to a comparison similar to

  marker.h() < 2.0 * particle.getSearchRadius()

This threshold is hard. However, we can sieve into a coarser level obviously to ensure that we don't violate it.

Optimisation

I started off with this nice simple code

const bool fitsIntoThisGridResolution = tarch::la::min( marker.h() ) * MultilevelSortingTolerance *
AsymmetricSortingDamping >= 2.0 * particle.getSearchRadius(); const bool fitsIntoNextFinerGridResolution =
tarch::la::min( marker.h() )/3.0 > 2.0 * particle.getSearchRadius();
return marker.isContainedInAdjacentCells( particle.getX(), 0.5, relativeReleaseOwnershipSpatialSortingTolerance(marker) *
AsymmetricSortingDamping * tarch::la::NUMERICAL_ZERO_DIFFERENCE) and fitsIntoThisGridResolution and ( not
fitsIntoNextFinerGridResolution or not marker.hasBeenRefined()
);
constexpr double NUMERICAL_ZERO_DIFFERENCE
Definition Scalar.h:17
Scalar min(const Vector< Size, Scalar > &vector)
Returns the element with minimal value (NOT absolute value).

However, it turns out a lot of runtime is spent in getSearchRadius(). It is not clear if the compiler manages to exploit the non-strict evaluation if we first define the two const bools. So I moved them now into the boolean expression, i.e. work without the explicit temporary variables.

See also
sieveParticle()

Referenced by toolbox::particles::SieveParticles< T >::getParticlesToBeSievedIntoVertex(), toolbox::particles::tests::MultiscaleTransitionsTest::testLiftDropOfParticleAssociatedWithVertex03(), toolbox::particles::tests::MultiscaleTransitionsTest::testLiftDropOfParticleAssociatedWithVertex04(), and toolbox::particles::tests::MultiscaleTransitionsTest::testSievePredicate().

Here is the caller graph for this function:

◆ sieveParticle() [2/2]

bool toolbox::particles::sieveParticle ( double searchRadius,
const tarch::la::Vector< Dimensions, double > & x,
const peano4::datamanagement::VertexMarker & marker )

Definition at line 107 of file MultiscaleTransitions.cpp.

References tarch::la::NUMERICAL_ZERO_DIFFERENCE.

◆ staticPosition()

template<typename Particle >
void toolbox::particles::staticPosition ( const Particle & particle,
double timeStepSize )

Leave the tracer particles where they are.

Essentially, this is a complicated way to write down a `‘not an operation’' (nope). The routine simply does nothing with the particle passed.

Definition at line 20 of file Tracer.h.

◆ updateContainedInAdjacentCellProperty() [1/3]

void toolbox::particles::updateContainedInAdjacentCellProperty ( const peano4::datamanagement::CellMarker & marker,
const std::bitset< Dimensions > & assignedVertex,
auto & particle )

Update helper flag per particle for cell.

Each particle is associated with a vertex of a cell. This vertex is identified by assignedVertex. We now compute a bitset that tells us, from the vertex' point of view, which of the adjacent cells will host this particle. This might not be unique, as we work with a fuzzy geometric ownership.

See also
the explanation in Particle.py

Definition at line 65 of file particles.cpph.

References tarch::la::multiplyComponents(), and updateContainedInAdjacentCellProperty().

Referenced by insertParticleIntoCell(), updateContainedInAdjacentCellProperty(), and updateContainedInAdjacentCellProperty().

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

◆ updateContainedInAdjacentCellProperty() [2/3]

void toolbox::particles::updateContainedInAdjacentCellProperty ( const peano4::datamanagement::VertexMarker & marker,
auto & particle )

We are.

Definition at line 85 of file particles.cpph.

References updateContainedInAdjacentCellProperty().

Here is the call graph for this function:

◆ updateContainedInAdjacentCellProperty() [3/3]

void toolbox::particles::updateContainedInAdjacentCellProperty ( const tarch::la::Vector< Dimensions, double > & x,
const tarch::la::Vector< Dimensions, double > & h,
auto & particle )

Update internal marker of particle.

Each particle is associated to a vertex. The vertex is identified by its position x and the surrounding cells with size h. Each vertex has \( 2^d \) adjacent cells. The particle's ContainedInAdjacentCellProperty bitset flags which of these adjacent cells consider the particle to fall in their area.

As we work with a weak ownership, there can be more than two bits (if the particle sits exactly on the face between two adjacent cells).

Particle movement

We have to call updateContainedInAdjacentCellProperty() after each particle movement. Even if we don't change the vertex association as a particle drifts, it might stay with the same vertex and still move from one cell into another one and hence require an update of the flags.

See also
the explanation in toolbox.particles.Particle

Definition at line 97 of file particles.cpph.

References dfor2, enddforx, peano4::datamanagement::CellMarker::isContained(), tarch::la::multiplyComponents(), tarch::la::NUMERICAL_ZERO_DIFFERENCE, and toolbox::particles::internal::relativeGrabOwnershipSpatialSortingTolerance().

Here is the call graph for this function:

Variable Documentation

◆ GrabOwnershipSpatialSortingTolerance

constexpr double toolbox::particles::GrabOwnershipSpatialSortingTolerance = 1.001
constexpr
See also
SpatialReassignmentTolerance

Definition at line 63 of file internal.h.

Referenced by toolbox::particles::internal::relativeGrabOwnershipSpatialSortingTolerance().

◆ ParticleReassociationInstruction_Keep

◆ ParticleReassociationInstruction_SieveGlobally

◆ ReleaseOwnershipSpatialSortingTolerance

constexpr double toolbox::particles::ReleaseOwnershipSpatialSortingTolerance = 1.01
constexpr

Sorting tolerance.

As the particle association is not unique and not sharp, as we work with floating point numbers. The correlation

\( centre_{left} + h/2 < centre_{right} - h/2 \)

might not hold for two neighbouring voxels, i.e. there might be a floating-point gap in-between two voxels. If particles fall into this gap, they would not belong to any voxel. We therefore give each voxel a small overlap. The relative size of this overlap is quantified by ReleaseOwnershipSpatialSortingTolerance. This is out first principle: particles are always assigned to a cell/volume relative to a relative threshold. We issue a resort if and only if the particle is outside of its respective control volume subject to the additional tolerance. If a particle leaves a cell or-vertex environment subject to the ReleaseOwnershipSpatialSortingTolerance, we either assign it to a neighbouring voxel or lift it into a coarser level or assign it to a global sorting list.

There is a second threshold: Our sorting should be idempotent. That is, once we assign a particle to a vertex or cell, we should not immediately require yet another resorting step. Two rules formalise this:

  1. If a vertex releases a stationary particle, it will never grab it back.
  2. If a vertex grabs a stationary particle, it will never release it.

The two principles imply that we work with two thresholds: We trigger a resort if a particle leaves its control volume. The region in which we grab a particle however is slightly smaller. We work with two thresholds, which all span areas which are slighter than the control volumes. However, releases are way more conservative/careful than grabs or drops.

Related action sets

The action set peano4.toolbox.particles.api.UpdateParallelState employs the threshold to decide if a particle is local or not. It has to use the most generic environment, i.e. the release radius.

See also
AsymmetricSortingDamping

Definition at line 58 of file internal.h.

Referenced by getParticleReassociationInstructionWithinCellWithIntraCellReassignment(), and toolbox::particles::internal::relativeReleaseOwnershipSpatialSortingTolerance().

◆ SpatialDuplicateTolerance

constexpr double toolbox::particles::SpatialDuplicateTolerance = 1e-6
constexpr

Default tolerance to decide when two particles are the same.

See also
toolbox::particles::particleIsDuplicate()

Definition at line 18 of file particles.h.

Referenced by particleIsDuplicate().