IPPL (Independent Parallel Particle Layer)
IPPL
Loading...
Searching...
No Matches
ippl::ParticleAttrib< T, Properties > Class Template Reference

#include <ParticleAttrib.h>

Inheritance diagram for ippl::ParticleAttrib< T, Properties >:
Collaboration diagram for ippl::ParticleAttrib< T, Properties >:

Public Types

typedef T value_type
using Base = typename detail::ParticleAttribBase<>::with_properties<Properties...>
using hash_type = typename Base::hash_type
using view_type = typename detail::ViewType<T, 1, Properties...>::view_type
using HostMirror = typename view_type::host_mirror_type
using memory_space = typename view_type::memory_space
using execution_space = typename view_type::execution_space
using size_type = detail::size_type

Public Member Functions

void create (size_type) override
void destroy (const hash_type &deleteIndex, const hash_type &keepIndex, size_type invalidCount) override
void pack (const hash_type &) override
void unpack (size_type) override
void serialize (detail::Archive< memory_space > &ar, size_type nsends) override
void deserialize (detail::Archive< memory_space > &ar, size_type nrecvs) override
virtual ~ParticleAttrib ()=default
size_type size () const override
size_type packedSize (const size_type count) const override
void resize (size_type n)
void realloc (size_type n)
void print ()
KOKKOS_INLINE_FUNCTION Toperator() (const size_t i) const
view_typegetView ()
const view_typegetView () const
HostMirror getHostMirror () const
void set_name (const std::string &name_) override
std::string get_name () const override
ParticleAttrib< T, Properties... > & operator= (T x)
template<typename E, size_t N>
ParticleAttrib< T, Properties... > & operator= (detail::Expression< E, N > const &expr)
template<typename Field, typename P2, typename policy_type>
void scatter (Field &f, const ParticleAttrib< Vector< P2, Field::dim >, Properties... > &pp, policy_type iteration_policy, hash_type hash_array={}) const
 Scatter particle attribute data onto a field.
template<typename Field, typename P2>
void gather (Field &f, const ParticleAttrib< Vector< P2, Field::dim >, Properties... > &pp, const bool addToAttribute=false)
 Gather field data into the particle attribute.
T sum ()
T max ()
T min ()
T prod ()
void applyPermutation (const hash_type &permutation) override
 Sort the attribute according to a permutation.
void internalCopy (const hash_type &indices) override
 Copy and create values of given indices.
template<typename Field, class PT, typename policy_type>
void scatter (Field &f, const ParticleAttrib< Vector< PT, Field::dim >, Properties... > &pp, policy_type iteration_policy, hash_type hash_array) const
KOKKOS_INLINE_FUNCTION auto operator[] (size_t i) const

Static Public Attributes

static constexpr unsigned dim = 1

Private Attributes

view_type dview_m
view_type buf_m

Detailed Description

template<typename T, class... Properties>
class ippl::ParticleAttrib< T, Properties >

Definition at line 28 of file ParticleAttrib.h.

Member Typedef Documentation

◆ Base

template<typename T, class... Properties>
using ippl::ParticleAttrib< T, Properties >::Base = typename detail::ParticleAttribBase<>::with_properties<Properties...>

Definition at line 36 of file ParticleAttrib.h.

◆ execution_space

template<typename T, class... Properties>
using ippl::ParticleAttrib< T, Properties >::execution_space = typename view_type::execution_space

Definition at line 45 of file ParticleAttrib.h.

◆ hash_type

template<typename T, class... Properties>
using ippl::ParticleAttrib< T, Properties >::hash_type = typename Base::hash_type

Definition at line 38 of file ParticleAttrib.h.

◆ HostMirror

template<typename T, class... Properties>
using ippl::ParticleAttrib< T, Properties >::HostMirror = typename view_type::host_mirror_type

Definition at line 42 of file ParticleAttrib.h.

◆ memory_space

template<typename T, class... Properties>
using ippl::ParticleAttrib< T, Properties >::memory_space = typename view_type::memory_space

Definition at line 44 of file ParticleAttrib.h.

◆ size_type

template<typename T, class... Properties>
using ippl::ParticleAttrib< T, Properties >::size_type = detail::size_type

Definition at line 47 of file ParticleAttrib.h.

◆ value_type

template<typename T, class... Properties>
typedef T ippl::ParticleAttrib< T, Properties >::value_type

Definition at line 33 of file ParticleAttrib.h.

◆ view_type

template<typename T, class... Properties>
using ippl::ParticleAttrib< T, Properties >::view_type = typename detail::ViewType<T, 1, Properties...>::view_type

Definition at line 40 of file ParticleAttrib.h.

Constructor & Destructor Documentation

◆ ~ParticleAttrib()

template<typename T, class... Properties>
virtual ippl::ParticleAttrib< T, Properties >::~ParticleAttrib ( )
virtualdefault

Member Function Documentation

◆ applyPermutation()

template<typename T, class... Properties>
void ippl::ParticleAttrib< T, Properties >::applyPermutation ( const hash_type & permutation)
override

Sort the attribute according to a permutation.

This function sorts the attribute according to a given permutation such that afterward attr(permutation(i)) = attr(i).

Note
this cannot be done inplace and is not very cache efficient
Parameters
permutationThe permutation to apply

Definition at line 220 of file ParticleAttrib.hpp.

References getView(), and size().

Here is the call graph for this function:

◆ create()

template<typename T, class... Properties>
void ippl::ParticleAttrib< T, Properties >::create ( size_type n)
override

Definition at line 25 of file ParticleAttrib.hpp.

References ippl::Comm, realloc(), and size().

Referenced by internalCopy().

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

◆ deserialize()

template<typename T, class... Properties>
void ippl::ParticleAttrib< T, Properties >::deserialize ( detail::Archive< memory_space > & ar,
size_type nrecvs )
inlineoverride

Definition at line 71 of file ParticleAttrib.h.

◆ destroy()

template<typename T, class... Properties>
void ippl::ParticleAttrib< T, Properties >::destroy ( const hash_type & deleteIndex,
const hash_type & keepIndex,
size_type invalidCount )
override

Particle deletion function. Partition the particles into a valid region and an invalid region.

Parameters
deleteIndexList of indices of invalid particles in the valid region
keepIndexList of indices of valid particles in the invalid region
invalidCountNumber of invalid particles in the valid region

Definition at line 34 of file ParticleAttrib.hpp.

References dview_m.

◆ gather()

template<typename T, class... Properties>
template<typename Field, typename P2>
void ippl::ParticleAttrib< T, Properties >::gather ( Field & f,
const ParticleAttrib< Vector< P2, Field::dim >, Properties... > & pp,
const bool addToAttribute = false )

Gather field data into the particle attribute.

This function gathers data from the given field into the particle attribute by iterating over all particles. Depending on the parameter addToAttribute, the gathered field value is either added to the existing attribute value (using "+=") or used to overwrite the attribute value.

Note
This behavior exists to give the OPAL-X field solver the ablity to gather field data per "energy bin".
Template Parameters
FieldThe type of the field.
P2The particle type for the position attribute.
Parameters
fThe field from which data is gathered.
ppThe ParticleAttrib representing particle positions.
addToAttributeIf true, the gathered value is added to the current attribute value; otherwise, the attribute value is overwritten.

Definition at line 167 of file ParticleAttrib.hpp.

References Dim, ippl::Field< T, Dim, Mesh_t< Dim >, Centering_t< Dim >, ViewArgs... >::dim, and IpplTimings::getTimer().

Here is the call graph for this function:

◆ get_name()

template<typename T, class... Properties>
std::string ippl::ParticleAttrib< T, Properties >::get_name ( ) const
inlineoverride

Definition at line 105 of file ParticleAttrib.h.

◆ getHostMirror()

template<typename T, class... Properties>
HostMirror ippl::ParticleAttrib< T, Properties >::getHostMirror ( ) const
inline

Definition at line 101 of file ParticleAttrib.h.

◆ getView() [1/2]

template<typename T, class... Properties>
view_type & ippl::ParticleAttrib< T, Properties >::getView ( )
inline

Definition at line 97 of file ParticleAttrib.h.

Referenced by ippl::detail::ParticleLayout< T, Dim, PositionProperties >::applyBC(), applyPermutation(), and internalCopy().

Here is the caller graph for this function:

◆ getView() [2/2]

template<typename T, class... Properties>
const view_type & ippl::ParticleAttrib< T, Properties >::getView ( ) const
inline

Definition at line 99 of file ParticleAttrib.h.

◆ internalCopy()

template<typename T, class... Properties>
void ippl::ParticleAttrib< T, Properties >::internalCopy ( const hash_type & indices)
override

Copy and create values of given indices.

This functions creates new particles with the values copied from the given indices. The values will be at the highest new indices.

Parameters
indicesThe indices to copy.

Definition at line 239 of file ParticleAttrib.hpp.

References create(), getView(), and size().

Here is the call graph for this function:

◆ max()

template<typename T, class... Properties>
T ippl::ParticleAttrib< T, Properties >::max ( )

◆ min()

template<typename T, class... Properties>
T ippl::ParticleAttrib< T, Properties >::min ( )

◆ operator()()

template<typename T, class... Properties>
KOKKOS_INLINE_FUNCTION T & ippl::ParticleAttrib< T, Properties >::operator() ( const size_t i) const
inline

Definition at line 95 of file ParticleAttrib.h.

◆ operator=() [1/2]

template<typename T, class... Properties>
template<typename E, size_t N>
ParticleAttrib< T, Properties... > & ippl::ParticleAttrib< T, Properties >::operator= ( detail::Expression< E, N > const & expr)

Assign an arbitrary particle attribute expression

Template Parameters
Eexpression type
Nsize of the expression, this is necessary for running on the device since otherwise it does not allocate enough memory
Parameters
expris the expression

Definition at line 92 of file ParticleAttrib.hpp.

References dview_m.

◆ operator=() [2/2]

template<typename T, class... Properties>
ParticleAttrib< T, Properties... > & ippl::ParticleAttrib< T, Properties >::operator= ( T x)

Assign the same value to the whole attribute.

Definition at line 81 of file ParticleAttrib.hpp.

References dview_m.

◆ operator[]()

KOKKOS_INLINE_FUNCTION auto ippl::detail::Expression< ParticleAttrib< T, Properties... >, N >::operator[] ( size_t i) const
inlineinherited

Access single element of the expression

Definition at line 32 of file IpplExpressions.h.

◆ pack()

template<typename T, class... Properties>
void ippl::ParticleAttrib< T, Properties >::pack ( const hash_type & hash)
override

Definition at line 48 of file ParticleAttrib.hpp.

References buf_m, and size().

Here is the call graph for this function:

◆ packedSize()

template<typename T, class... Properties>
size_type ippl::ParticleAttrib< T, Properties >::packedSize ( const size_type count) const
inlineoverride

Definition at line 79 of file ParticleAttrib.h.

◆ print()

template<typename T, class... Properties>
void ippl::ParticleAttrib< T, Properties >::print ( )
inline

Definition at line 87 of file ParticleAttrib.h.

◆ prod()

template<typename T, class... Properties>
T ippl::ParticleAttrib< T, Properties >::prod ( )

◆ realloc()

template<typename T, class... Properties>
void ippl::ParticleAttrib< T, Properties >::realloc ( size_type n)
inline

Definition at line 85 of file ParticleAttrib.h.

Referenced by create().

Here is the caller graph for this function:

◆ resize()

template<typename T, class... Properties>
void ippl::ParticleAttrib< T, Properties >::resize ( size_type n)
inline

Definition at line 83 of file ParticleAttrib.h.

◆ scatter() [1/2]

template<typename T, class... Properties>
template<typename Field, typename P2, typename policy_type>
void ippl::ParticleAttrib< T, Properties >::scatter ( Field & f,
const ParticleAttrib< Vector< P2, Field::dim >, Properties... > & pp,
policy_type iteration_policy,
hash_type hash_array = {} ) const

Scatter particle attribute data onto a field.

This function scatters data from this attribute onto the given field, using the given position attribute. The function can be used together with a custom iteration policy to iterate over a specified range and, optionally, an ippl::hash_type array to remap iteration indices.

When a non-empty hash_array is provided, the function:

  • Checks that the iteration policy's range does not exceed the size of hash_array.

  • Maps the current index to the appropriate index using the hash_array.

  • Careful: access pattern optimization might be lost when using hash_array.

Note
This custom iteration functionality is needed to support energy binning in the field solver of OPAL-X, allowing only particles within a specific bin to be scattered.
Template Parameters
FieldThe type of the field.
P2The type for the position attribute.
policy_typeThe type of the Kokkos iteration policy when using hash_array.
Parameters
fThe field onto which the particle data is scattered.
ppThe ParticleAttrib representing particle positions.
iteration_policyA custom Kokkos::range_policy defining the iteration range.
hash_arrayAn optional ippl::hash_type array for index mapping. If empty, no map is used.

◆ scatter() [2/2]

template<typename T, class... Properties>
template<typename Field, class PT, typename policy_type>
void ippl::ParticleAttrib< T, Properties >::scatter ( Field & f,
const ParticleAttrib< Vector< PT, Field::dim >, Properties... > & pp,
policy_type iteration_policy,
hash_type hash_array ) const

Definition at line 106 of file ParticleAttrib.hpp.

◆ serialize()

template<typename T, class... Properties>
void ippl::ParticleAttrib< T, Properties >::serialize ( detail::Archive< memory_space > & ar,
size_type nsends )
inlineoverride

Definition at line 67 of file ParticleAttrib.h.

◆ set_name()

template<typename T, class... Properties>
void ippl::ParticleAttrib< T, Properties >::set_name ( const std::string & name_)
inlineoverride

Definition at line 103 of file ParticleAttrib.h.

◆ size()

template<typename T, class... Properties>
size_type ippl::ParticleAttrib< T, Properties >::size ( ) const
inlineoverride

Definition at line 77 of file ParticleAttrib.h.

Referenced by applyPermutation(), create(), internalCopy(), and pack().

Here is the caller graph for this function:

◆ sum()

template<typename T, class... Properties>
T ippl::ParticleAttrib< T, Properties >::sum ( )

◆ unpack()

template<typename T, class... Properties>
void ippl::ParticleAttrib< T, Properties >::unpack ( size_type nrecvs)
override

Definition at line 63 of file ParticleAttrib.hpp.

Member Data Documentation

◆ buf_m

template<typename T, class... Properties>
view_type ippl::ParticleAttrib< T, Properties >::buf_m
private

Definition at line 206 of file ParticleAttrib.h.

Referenced by pack().

◆ dim

template<typename T, class... Properties>
unsigned ippl::ParticleAttrib< T, Properties >::dim = 1
staticconstexpr

Definition at line 34 of file ParticleAttrib.h.

◆ dview_m

template<typename T, class... Properties>
view_type ippl::ParticleAttrib< T, Properties >::dview_m
private

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