Peano
Loading...
Searching...
No Matches
swift2::kernels Namespace Reference

Namespaces

namespace  coalesced
 
namespace  gravity
 
namespace  internal
 
namespace  legacy
 Legacy SPH implementation.
 
namespace  ompoffloading
 

Concepts

concept  ParticleUnaryOperatorOnVertex
 Definition of particle update (unary operation)
 
concept  ParticleUnaryOperatorOnCell
 
concept  ParticleBinaryOperator
 Definition of particle-particle interaction.
 

Typedefs

template<typename ParticleContainer >
using PCParticle = typename std::remove_pointer<typename ParticleContainer::value_type>::type
 
template<typename Particle >
using UpdateParticlePairWithinCellPredicate = std::function<bool(const peano4::datamanagement::CellMarker& marker, const Particle& localParticle, const Particle& activeParticle)>
 
template<typename Particle >
using UpdateParticleAssignedToCellPredicate = std::function<bool(const peano4::datamanagement::CellMarker& marker, const Particle& localParticle)>
 
template<typename Particle >
using UpdateParticleAssignedToVertexPredicate = std::function<bool(const peano4::datamanagement::VertexMarker& marker, const Particle& localParticle)>
 

Functions

template<typename Particle >
void adoptInteractionRadiusAndTriggerRerun (const std::list< Particle * > &localParticles, const std::list< Particle * > &activeParticles, int targetNumberOfNeighbourParticles, double maxGrowthPerSweep=2.0, double shrinkingFactor=0.8)
 This routine runs over all the local particles and tries to ensure that the number of particles equals roughly targetNumberOfNeighbourParticles.
 
template<typename ParticleContainer >
void flagBoundaryParticles (const ParticleContainer &localParticles, const double nparts)
 Flag boundary particles.
 
template<typename ParticleContainer >
void flagBoundaryParticles (const ParticleContainer &localParticles, const double nparts, const tarch::la::Vector< Dimensions, double > &domainSize, const tarch::la::Vector< Dimensions, double > &domainOffset)
 
template<typename ParticleContainer >
void forAllParticles (const peano4::datamanagement::VertexMarker &marker, ParticleContainer &assignedParticles, int numberOfCoalescedAssignedParticles, ParticleUnaryOperatorOnVertex< typename ParticleContainer::DoFType > auto f, UpdateParticleAssignedToVertexPredicate< typename ParticleContainer::DoFType > predicate)
 Run over all particles and update them independent of each other.
 
template<typename ParticleContainer >
void forAllParticles (const peano4::datamanagement::CellMarker &marker, ParticleContainer &localParticles, const std::vector< int > &numberOfCoalescedLocalParticlesPerVertex, ParticleUnaryOperatorOnCell< PCParticle< ParticleContainer > > auto f, UpdateParticleAssignedToCellPredicate< PCParticle< ParticleContainer > > predicate)
 Loop over all particles within a cell.
 
template<typename LocalParticleContainer , typename ActiveParticleContainer >
void forAllParticlePairs (const peano4::datamanagement::CellMarker &marker, LocalParticleContainer &localParticles, ActiveParticleContainer &activeParticles, const std::vector< int > &numberOfCoalescedLocalParticlesPerVertex, const std::vector< int > &numberOfCoalescedActiveParticlesPerVertex, ParticleBinaryOperator< PCParticle< LocalParticleContainer >, PCParticle< ActiveParticleContainer > > auto f, UpdateParticleAssignedToCellPredicate< PCParticle< LocalParticleContainer > > localParticlePredicate=::swift2::kernels::localParticleCanBeUpdatedInCellKernel< PCParticle< LocalParticleContainer > >, UpdateParticlePairWithinCellPredicate< PCParticle< LocalParticleContainer > > particlePairPredicate=::swift2::kernels::localParticleCanBeUpdatedInCellKernelFromAnyOtherParticle< PCParticle< LocalParticleContainer > >)
 Run over all local particle-active particle combinations.
 
template<typename Particle >
bool alwaysUpdateInVertexKernel (const peano4::datamanagement::VertexMarker &marker, const Particle &localParticle)
 Degenerated predicate which always allows for an update.
 
template<typename Particle >
bool alwaysUpdateInCellKernel (const peano4::datamanagement::CellMarker &marker, const Particle &localParticle)
 
template<typename Particle >
bool alwaysUpdateParticlePairs (const peano4::datamanagement::CellMarker &marker, const Particle &localParticle, const Particle &activeParticle)
 
template<typename Particle >
bool particleIsLocal (const peano4::datamanagement::VertexMarker &marker, const Particle &localParticle)
 is this localParticle a local particle (in the ParallelState sense)?
 
template<typename Particle >
bool localParticleCanBeUpdatedAndMovedInVertexKernel (const peano4::datamanagement::VertexMarker &marker, const Particle &localParticle)
 Can we move (drift) this particle?
 
template<typename Particle >
bool localParticleCanBeUpdatedInCellKernelFromAnyOtherParticle (const peano4::datamanagement::CellMarker &marker, const Particle &localParticle, const Particle &activeParticle)
 Can we do work on this particle during a cell kernel sweep stage?
 
template<typename Particle >
bool localParticleCanBeUpdatedInCellKernelFromAnyOtherParticleWithinIterationRange (const peano4::datamanagement::CellMarker &marker, const Particle &localParticle, const Particle &activeParticle)
 
template<typename Particle >
bool localParticleCanBeUpdatedInCellKernel (const peano4::datamanagement::CellMarker &marker, const Particle &localParticle)
 
template<typename Particle >
bool localParticleCanBeUpdatedInVertexKernel (const peano4::datamanagement::VertexMarker &marker, const Particle &localParticle)
 Can we do work on this particle during a vertex kernel sweep stage?
 

Typedef Documentation

◆ PCParticle

template<typename ParticleContainer >
using swift2::kernels::PCParticle = typename std::remove_pointer<typename ParticleContainer::value_type>::type

Definition at line 18 of file ParticleSetIterators.h.

◆ UpdateParticleAssignedToCellPredicate

template<typename Particle >
using swift2::kernels::UpdateParticleAssignedToCellPredicate = std::function<bool(const peano4::datamanagement::CellMarker& marker, const Particle& localParticle)>

Definition at line 13 of file ParticleUpdatePredicates.h.

◆ UpdateParticleAssignedToVertexPredicate

template<typename Particle >
using swift2::kernels::UpdateParticleAssignedToVertexPredicate = std::function<bool(const peano4::datamanagement::VertexMarker& marker, const Particle& localParticle)>

Definition at line 16 of file ParticleUpdatePredicates.h.

◆ UpdateParticlePairWithinCellPredicate

template<typename Particle >
using swift2::kernels::UpdateParticlePairWithinCellPredicate = std::function<bool(const peano4::datamanagement::CellMarker& marker, const Particle& localParticle, const Particle& activeParticle)>

Definition at line 10 of file ParticleUpdatePredicates.h.

Function Documentation

◆ adoptInteractionRadiusAndTriggerRerun()

template<typename Particle >
void swift2::kernels::adoptInteractionRadiusAndTriggerRerun ( const std::list< Particle * > & localParticles,
const std::list< Particle * > & activeParticles,
int targetNumberOfNeighbourParticles,
double maxGrowthPerSweep = 2.0,
double shrinkingFactor = 0.8 )

This routine runs over all the local particles and tries to ensure that the number of particles equals roughly targetNumberOfNeighbourParticles.

The algorithm is fairly simple:

  • Find out how many other particles are currently within the search radius.
  • If this number is bigger than the target value, decrease the search radius slightly. Count again. If we are still above the target value, accept the slightly reduced search. Otherwise, undo the alteration.
  • If the number of interaction partners is too small, increase it and study the effect:
    • If we meet the target neighbour count now, we are happy.
    • If we have increased the number of neighbours but are not there yet, we keep this increase, but ask Peano/Swift to run this analysis again.
    • If the increase has not paid off, we keep the old search radius. It means that likely the interaction radii overall all are too small or the particle density is just not there.

The reduction starts from the assumption that we should be careful with reducing the search radii. So we only slightly decrease the search radius. If we could have reduced it more aggressively, we accept that and hope that subsequent time steps or sweeps will eventually bring the search radius down. But we do not enforce it here, which might mean that we work with too big interaction sets.

The increase is different: If we

If an interaction radius has to be increased, the code sets the flag setRerunPreviousGridSweep() on the underlying particle species.

In prin

Todo
Mladen This docu is messed up, and we have to take into account here that we have to accommodate multiscale effects.

Definition at line 5 of file ParticleSearchRadiusCalculation.cpph.

References _log, assertion, logDebug, and tarch::la::norm2Squared().

Here is the call graph for this function:

◆ alwaysUpdateInCellKernel()

template<typename Particle >
bool swift2::kernels::alwaysUpdateInCellKernel ( const peano4::datamanagement::CellMarker & marker,
const Particle & localParticle )

Definition at line 29 of file ParticleUpdatePredicates.h.

◆ alwaysUpdateInVertexKernel()

template<typename Particle >
bool swift2::kernels::alwaysUpdateInVertexKernel ( const peano4::datamanagement::VertexMarker & marker,
const Particle & localParticle )

Degenerated predicate which always allows for an update.

Definition at line 23 of file ParticleUpdatePredicates.h.

◆ alwaysUpdateParticlePairs()

template<typename Particle >
bool swift2::kernels::alwaysUpdateParticlePairs ( const peano4::datamanagement::CellMarker & marker,
const Particle & localParticle,
const Particle & activeParticle )

Definition at line 35 of file ParticleUpdatePredicates.h.

◆ flagBoundaryParticles() [1/2]

template<typename ParticleContainer >
void swift2::kernels::flagBoundaryParticles ( const ParticleContainer & localParticles,
const double nparts )

Flag boundary particles.

These particles we will not be updated by any algorithmic step.

Parameters
ParticleContainerTypicaly either std::list<Particle *> or std::unordered_set<Particle *>

◆ flagBoundaryParticles() [2/2]

template<typename ParticleContainer >
void swift2::kernels::flagBoundaryParticles ( const ParticleContainer & localParticles,
const double nparts,
const tarch::la::Vector< Dimensions, double > & domainSize,
const tarch::la::Vector< Dimensions, double > & domainOffset )
Todo
Das ist falsch

Definition at line 10 of file ParticleSelfInteraction.cpph.

References tarch::la::oneGreater(), and tarch::la::oneSmaller().

Here is the call graph for this function:

◆ forAllParticlePairs()

template<typename LocalParticleContainer , typename ActiveParticleContainer >
void swift2::kernels::forAllParticlePairs ( const peano4::datamanagement::CellMarker & marker,
LocalParticleContainer & localParticles,
ActiveParticleContainer & activeParticles,
const std::vector< int > & numberOfCoalescedLocalParticlesPerVertex,
const std::vector< int > & numberOfCoalescedActiveParticlesPerVertex,
ParticleBinaryOperator< PCParticle< LocalParticleContainer >, PCParticle< ActiveParticleContainer > > auto f,
UpdateParticleAssignedToCellPredicate< PCParticle< LocalParticleContainer > > localParticlePredicate = ::swift2::kernels::localParticleCanBeUpdatedInCellKernel<PCParticle<LocalParticleContainer>>,
UpdateParticlePairWithinCellPredicate< PCParticle< LocalParticleContainer > > particlePairPredicate = ::swift2::kernels::localParticleCanBeUpdatedInCellKernelFromAnyOtherParticle<PCParticle<LocalParticleContainer>> )

Run over all local particle-active particle combinations.

The generic implementation realises a nested loops: We call interaction per pair of local and active particles. The interaction functor may modify the local particle. It is the responsibility of the update predicate to ensure that no particle is updated twice, even if it is sitting right at the face in-between two cell. The predicate also should check if there's a self-interaction and mask it out if appropriate.

The routine accepts a container over particle pointers, but the functor f actually accepts references. It is the responsibility of this routine to map pointers onto references. We use

for (auto& localParticle : localParticles) {

here, but it is important to note that this reference still becomes a reference to a pointer as localParticle is a container over pointers.

In Swift's terminology, this is a non-symmetric XXX calculation where XXX is the compute kernel (such as density): Our implementations loops over all local particles (outer loop) and, per local particle then over all active particles (inner loop). Inside this second loop, we compute a force acting from the active particle onto the local particle.

To avoid self-interaction, users have to ensure that \( f(p,p) :=0 \): If a particle acts on itself, the arising force equals zero. This should always be the case, as the predicate passed in here is an optimisation, i.e. users should not rely on it to mask out self-interaction a priori.

Equally important is a second case distinction: If a local and an active particle reside on the same level, we know that the symmetric force component will be computed later on. Let L be the local one and A the active one, we have

\( f(L,A) = -f(A,L) \)

but this antisymmetry is not exploited in this kernel. We know that there will be another loop iteration or cell interaction kernel call which takes care of the negative complementary force.

However, if A resides on a coarser level than L, this assumption is wrong. We need to add the force explicitly as discussed in Mesh traversal.

Multiscale interactions

This routine follows Peano's generic multiscale discussion, where we point out that we have correctly take particles into account which reside on coarser levels. This is necessary if some particles have larger cut-off radii than the others or if we work with adaptive meshes.

Predicates

This would be a generic predicate as most kernels require it:

localParticlePredicate = ::swift2::kernels::alwaysUpdateInCellKernel<PCParticle<LocalParticleContainer>>,
particlePairPredicate = ::swift2::kernels::localParticleCanBeUpdatedInCellKernelFromAnyOtherParticle<PCParticle<LocalParticleContainer>>
std::function< bool(const peano4::datamanagement::CellMarker &marker, const Particle &localParticle)> UpdateParticleAssignedToCellPredicate
std::function< bool(const peano4::datamanagement::CellMarker &marker, const Particle &localParticle, const Particle &activeParticle)> UpdateParticlePairWithinCellPredicate

However, you can use a more bespoke version as follows:

localParticlePredicate = ::swift2::kernels::localParticleCanBeUpdatedInCellKernel<PCParticle<LocalParticleContainer>>,
particlePairPredicate = ::swift2::kernels::localParticleCanBeUpdatedInCellKernelFromAnyOtherParticle<PCParticle<LocalParticleContainer>>

Hydro Force calculation particle-particle kernel

Todo
This documentation is likely outdated

Non-symmetric force interaction between two particles (as per Swift terminology). In this, only the Local particle is updated (in the symmetric case, both are updated at the same time). Time derivative of the internal energy, du/dt, is also updated here. Subscripts a and b represent indices of Local and Active particles, respectively:

h_a: Smoothing length of Local particle h_b: Smoothign length of Active particle dx : Vector separating the two particles (x_loc - x_act) r : Distance between the two particles (|| dx ||)

To inject the hydro force, we typically use the code snippet similar to

cell_kernel_for_force_calculation = "::swift2::forAllParticlePairs(marker, localParticles, activeParticles,::swift2::kernels::forceKernel<globaldata::HydroPart>);"

As we have to know the particle type for the injected functor, the SPH particle for example replaces this similar to

cell_kernel_for_force_calculation = "::swift2::forAllParticlePairs(marker, localParticles, activeParticles,::swift2::kernels::forceKernel<globaldata::{}>);".format(self.name)
Parameters
LocalParticleContainerA subtype of std::list<Particle *>
ActiveParticleContainerA subtype of std::list<Particle *>

Definition at line 43 of file ParticleSetIterators.cpph.

References f(), and tarch::la::oneGreater().

Referenced by swift2::kernels::ompoffloading::forAllParticlePairs(), and runBenchmark().

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

◆ forAllParticles() [1/2]

template<typename ParticleContainer >
void swift2::kernels::forAllParticles ( const peano4::datamanagement::CellMarker & marker,
ParticleContainer & localParticles,
const std::vector< int > & numberOfCoalescedLocalParticlesPerVertex,
ParticleUnaryOperatorOnCell< PCParticle< ParticleContainer > > auto f,
UpdateParticleAssignedToCellPredicate< PCParticle< ParticleContainer > > predicate )

Loop over all particles within a cell.

Predicates

The most popular predicates in use are

::swift2::kernels::alwaysUpdateInCellKernel<typename ParticleContainer::DoFType>

and

::swift2::kernels::localParticleCanBeUpdatedInCellKernel<PCParticle<ParticleContainer>>

Definition at line 26 of file ParticleSetIterators.cpph.

References f().

Here is the call graph for this function:

◆ forAllParticles() [2/2]

template<typename ParticleContainer >
void swift2::kernels::forAllParticles ( const peano4::datamanagement::VertexMarker & marker,
ParticleContainer & assignedParticles,
int numberOfCoalescedAssignedParticles,
ParticleUnaryOperatorOnVertex< typename ParticleContainer::DoFType > auto f,
UpdateParticleAssignedToVertexPredicate< typename ParticleContainer::DoFType > predicate )

Run over all particles and update them independent of each other.

The routine accepts a container over particle pointers, but the functor f actually accepts references. It is the responsibility of this routine to map pointers onto references.

Predicate

The predicate can be used to mask out certain updates.

We distinguish two use cases for the particle self-interactions:

  1. We want to execute a stand-alone computation over the particle, e.g. the kick step, which is does not require any particle-particle contribution.
  2. We want to execute computations after and before the particle-particle interactions kernels, e.g. to 'prepare' or 'end' a given calculation in SPH loops (e.g. density).

The most popular predicates therefore are:

::swift2::kernels::localParticleCanBeUpdatedInVertexKernel<typename ParticleContainer::DoFType>

and

::swift2::kernels::alwaysUpdateInVertexKernel<typename ParticleContainer::DoFType>

Consistency remarks

In order to ensure the self-interaction kernels execute consistently during the mesh traversals for these type of operations, the user should bear in mind the difference between these two cases:

  • The 'prepare' case should be mapped into touchVertexFirstTime(). Peano internally ensures that CellHasUpdatedParticle is false in this situation.
  • The 'end' case should be mapped into touchVertexLastTime. Peano internally ensures that CellHasUpdatedParticle is true in this situation.

This routine does not vectorise over the particles. If any vectorisation is used, you'll see the vector instructions arise from the actual compute kernel.

Parameters
ParticleContainerA subtype of std::list<Particle *>

Definition at line 10 of file ParticleSetIterators.cpph.

References f().

Referenced by swift2::kernels::ompoffloading::forAllParticles(), swift2::kernels::ompoffloading::forAllParticles(), prepareForceBenchmark(), and runBenchmark().

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

◆ localParticleCanBeUpdatedAndMovedInVertexKernel()

template<typename Particle >
bool swift2::kernels::localParticleCanBeUpdatedAndMovedInVertexKernel ( const peano4::datamanagement::VertexMarker & marker,
const Particle & localParticle )

Can we move (drift) this particle?

This is a more restrictive version compared to localParticleCanBeUpdatedInVertexKernel(), as it allows the underlying kernel to move a particle, too. Hence, the particle has to be local, and we have to check if it has not been moved yet. It is important that we distinguish this more restrictive version from its counterpart, as not each and every mesh traversal might reset the moved marker.

Parameters
localParticleParticle to check for
markerIdentifier for this vertex

Definition at line 63 of file ParticleUpdatePredicates.h.

References particleIsLocal().

Referenced by swift2::timestepping::computeExplicitEulerCromerWithGlobalTimeStepSize(), swift2::timestepping::computeExplicitEulerWithGlobalTimeStepSize(), swift2::kernels::legacy::leapfrogDriftWithGlobalTimeStepSize(), and swift2::timestepping::leapfrogDriftWithGlobalTimeStepSize().

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

◆ localParticleCanBeUpdatedInCellKernel()

template<typename Particle >
bool swift2::kernels::localParticleCanBeUpdatedInCellKernel ( const peano4::datamanagement::CellMarker & marker,
const Particle & localParticle )

◆ localParticleCanBeUpdatedInCellKernelFromAnyOtherParticle()

template<typename Particle >
bool swift2::kernels::localParticleCanBeUpdatedInCellKernelFromAnyOtherParticle ( const peano4::datamanagement::CellMarker & marker,
const Particle & localParticle,
const Particle & activeParticle )

Can we do work on this particle during a cell kernel sweep stage?

A particle is to be updated if and only if

  • it is located within the cell of interest;
  • it has not yet been updated by a cell;
  • the particle interaction triggering the update is not a self-interaction.

The second point is important. A particle might be located right at the face in-between two cells. In this case, it is not clear to which cell is actually belong to. So we are fine if either cell updates it, but it should be only one cel at a time.

Further remarks

I originally thought that this predicate should resemble

return (
not localParticle->getCellHasUpdatedParticle()
and
marker.isContained(localParticle->getX())
and
);
bool particleWillBeDroppedFurther(const Particle &particle, const peano4::datamanagement::CellMarker &marker)
Will the particle be dropped further throughout the traversal.

However, that's a poor idea, as it does not work along AMR boundaries for particles which reside in a refined cell yet would be dropped into a hanging vertex (which we don't). The file peano4.toolbox.particles.api.AbstractParticleGridAssociation provides some examples on this.

Parameters
localParticleparticle to check for
markerthe cell's CellMarker

Definition at line 109 of file ParticleUpdatePredicates.h.

References tarch::la::NUMERICAL_ZERO_DIFFERENCE, and toolbox::particles::internal::relativeGrabOwnershipSpatialSortingTolerance().

Referenced by swift2::kernels::legacy::forceKernel().

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

◆ localParticleCanBeUpdatedInCellKernelFromAnyOtherParticleWithinIterationRange()

template<typename Particle >
bool swift2::kernels::localParticleCanBeUpdatedInCellKernelFromAnyOtherParticleWithinIterationRange ( const peano4::datamanagement::CellMarker & marker,
const Particle & localParticle,
const Particle & activeParticle )

◆ localParticleCanBeUpdatedInVertexKernel()

template<typename Particle >
bool swift2::kernels::localParticleCanBeUpdatedInVertexKernel ( const peano4::datamanagement::VertexMarker & marker,
const Particle & localParticle )

Can we do work on this particle during a vertex kernel sweep stage?

This predicate filters out all halo (virtual) particle. It implicitly assumes that the particle-vertex association is correct. Therefore, we really only have to mask out virtual particles. The predicate breaks down if the association is not correct, which means it does not work if particles move.

Parameters
localParticleParticle to check for
markerIdentifier for this vertex

Definition at line 162 of file ParticleUpdatePredicates.h.

References particleIsLocal().

Referenced by swift2::kernels::legacy::computeNextTimestep(), swift2::kernels::legacy::endDensityCalculation(), swift2::kernels::legacy::endHydroForceCalculation(), swift2::kernels::legacy::firstInitParticle(), swift2::kernels::legacy::hydroPredictExtra(), swift2::kernels::legacy::leapfrogKickExtraWithGlobalTimeStepSize(), swift2::kernels::legacy::leapfrogKickWithGlobalTimeStepSize(), swift2::timestepping::leapfrogKickWithGlobalTimeStepSize(), swift2::kernels::legacy::predictHydro(), swift2::kernels::legacy::prepareDensity(), swift2::kernels::legacy::prepareHydroForce(), swift2::statistics::reduceVelocityAndSearchRadius(), swift2::kernels::legacy::resetAcceleration(), swift2::timestepping::resetMovedParticleMarker(), swift2::kernels::legacy::resetPredictedValues(), swift2::kernels::legacy::resetSmoothingLengthIterationCounter(), and swift2::kernels::legacy::updateSmoothingLengthAndRerunIfRequired().

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

◆ particleIsLocal()

template<typename Particle >
bool swift2::kernels::particleIsLocal ( const peano4::datamanagement::VertexMarker & marker,
const Particle & localParticle )

is this localParticle a local particle (in the ParallelState sense)?

Definition at line 44 of file ParticleUpdatePredicates.h.

Referenced by localParticleCanBeUpdatedAndMovedInVertexKernel(), and localParticleCanBeUpdatedInVertexKernel().

Here is the caller graph for this function: