Cabana 0.8.0-dev
 
Loading...
Searching...
No Matches
Cabana::Grid::Experimental::SparseHalo< MemorySpace, DataTypes, EntityType, NumSpaceDim, cellBitsPerTileDim, Value, Key > Class Template Reference

#include <Cabana_Grid_SparseHalo.hpp>

Public Types

enum  KeyValue { invalid_key = ~static_cast<key_type>( 0 ) }
 invalid key in sparse map
 
enum  Index { own = 0 , ghost = 1 , total = 2 }
 index (own or ghost)
 
using memory_space = MemorySpace
 memory space to store the sparse grid
 
using entity_type = EntityType
 entity type on sparse grid
 
using halo_pattern_type = NodeHaloPattern<NumSpaceDim>
 sparse grid halo pattern (TODO currently reusing Halo's Node pattern)
 
using value_type = Value
 value type of entities on sparse grid
 
using key_type = Key
 key type in sparse map
 
using aosoa_member_types = DataTypes
 data members in AoSoA structure
 
using tuple_type = Cabana::Tuple<aosoa_member_types>
 AoSoA tuple type.
 
template<std::size_t M>
using member_data_type
 AoSoA member data type.
 
using buffer_view = Kokkos::View<tuple_type*, memory_space>
 communication data buffer view type
 
using steering_view = Kokkos::View<key_type*, memory_space>
 
using tile_index_space
 tile index space type TODO
 
using counting_view
 

Public Member Functions

template<class SparseArrayType>
 SparseHalo (const halo_pattern_type pattern, const std::shared_ptr< SparseArrayType > &sparse_array)
 constructor
 
template<class SparseArrayType>
MPI_Comm getComm (const SparseArrayType sparse_array) const
 Get the communicator.
 
template<class DecompositionTag, class LocalGridType>
void buildCommData (DecompositionTag decomposition_tag, const std::shared_ptr< LocalGridType > &local_grid, const std::array< int, num_space_dim > &nid, std::vector< buffer_view > &buffers, std::vector< steering_view > &steering, std::vector< tile_index_space > &spaces)
 Build communication data.
 
template<class LocalGridType>
void updateTileSpace (const std::shared_ptr< LocalGridType > &local_grid)
 update tile index space according to current partition
 
template<class ExecSpace, class SparseMapType>
void register_halo (SparseMapType &map)
 register valid halos (according to grid activation status in sparse map) in the steerings
 
void clear (MPI_Comm comm)
 clear guiding information in sparse halo,
 
void collectNeighborCounting (MPI_Comm comm, const bool is_neighbor_counting_collected=false) const
 neighbor tile counting, communication needed only if the counting is non-zero
 
void scatterValidSendAndRecvRanks (MPI_Comm comm, std::vector< int > &valid_sends, std::vector< int > &valid_recvs, const bool is_neighbor_counting_collected=false) const
 collect all valid ranks for sparse grid scatter operations
 
void gatherValidSendAndRecvRanks (MPI_Comm comm, std::vector< int > &valid_sends, std::vector< int > &valid_recvs, const bool is_neighbor_counting_collected=false) const
 collect all valid ranks for sparse grid gather operations
 
template<class ExecSpace, class SparseArrayType>
void gather (const ExecSpace &exec_space, SparseArrayType &sparse_array, const bool is_neighbor_counting_collected=false) const
 Gather data into our ghosted share space from their owners.
 
template<class ExecSpace, class ReduceOp, class SparseArrayType>
void scatter (const ExecSpace &exec_space, const ReduceOp &reduce_op, SparseArrayType &sparse_array, const bool is_neighbor_counting_collected=false) const
 Scatter data from our ghosts to their owners using the given type of reduce operation.
 
template<class ExecSpace, class SparseArrayType>
void packBuffer (const ExecSpace &exec_space, const buffer_view &buffer, const steering_view &tile_steering, SparseArrayType &sparse_array, const int count) const
 Pack sparse arrays at halo regions into a buffer.
 
template<class ReduceOp, class ExecSpace, class SparseArrayType, class SparseMapType>
void unpackBuffer (const ReduceOp &reduce_op, const ExecSpace &exec_space, const buffer_view &buffer, const steering_view &tile_steering, const SparseArrayType &sparse_array, SparseMapType &map, const int count) const
 Unpack a sparse array communication buffer.
 

Static Public Member Functions

template<class T>
static KOKKOS_INLINE_FUNCTION void unpackOp (ScatterReduce::Sum, const T &buffer_val, T &array_val)
 Reduce an element into the buffer. Sum reduction.
 
template<class T>
static KOKKOS_INLINE_FUNCTION void unpackOp (ScatterReduce::Min, const T &buffer_val, T &array_val)
 Reduce an element into the buffer. Min reduction.
 
template<class T>
static KOKKOS_INLINE_FUNCTION void unpackOp (ScatterReduce::Max, const T &buffer_val, T &array_val)
 Reduce an element into the buffer. Max reduction.
 
template<class T>
static KOKKOS_INLINE_FUNCTION void unpackOp (ScatterReduce::Replace, const T &buffer_val, T &array_val)
 Reduce an element into the buffer. Replace reduction.
 
template<class ReduceOp, std::size_t N, std::size_t M, class SoAType>
static KOKKOS_FORCEINLINE_FUNCTION std::enable_if_t< 3==M, void > unpackTupleMember (const ReduceOp &reduce_op, const tuple_type &src_tuple, SoAType &dst_soa, const int soa_idx, const Kokkos::Array< std::size_t, M > &extents, const std::integral_constant< std::size_t, N >, const std::integral_constant< std::size_t, M >)
 Unpack a sparse arrays element (a tuple) in a buffer (for case tuple members with rank == 3)
 
template<class ReduceOp, std::size_t N, std::size_t M, class SoAType>
static KOKKOS_FORCEINLINE_FUNCTION std::enable_if_t< 2==M, void > unpackTupleMember (const ReduceOp &reduce_op, const tuple_type &src_tuple, SoAType &dst_soa, const int soa_idx, const Kokkos::Array< std::size_t, M > &extents, const std::integral_constant< std::size_t, N >, const std::integral_constant< std::size_t, M >)
 Unpack a sparse arrays element (a tuple) in a buffer (for case tuple members with rank == 2)
 
template<class ReduceOp, std::size_t N, std::size_t M, class SoAType>
static KOKKOS_FORCEINLINE_FUNCTION std::enable_if_t< 1==M, void > unpackTupleMember (const ReduceOp &reduce_op, const tuple_type &src_tuple, SoAType &dst_soa, const int soa_idx, const Kokkos::Array< std::size_t, M > &extents, const std::integral_constant< std::size_t, N >, const std::integral_constant< std::size_t, M >)
 Unpack a sparse arrays element (a tuple) in a buffer (for case tuple members with rank == 1)
 
template<class ReduceOp, std::size_t N, std::size_t M, class SoAType>
static KOKKOS_FORCEINLINE_FUNCTION std::enable_if_t< 0==M, void > unpackTupleMember (const ReduceOp &reduce_op, const tuple_type &src_tuple, SoAType &dst_soa, const int soa_idx, const Kokkos::Array< std::size_t, M > &, const std::integral_constant< std::size_t, N >, const std::integral_constant< std::size_t, M >)
 Unpack a sparse arrays element (a tuple) in a buffer (for case tuple members with rank == 0)
 
template<class ReduceOp, class SoAType>
static KOKKOS_FORCEINLINE_FUNCTION void unpackTuple (const ReduceOp &reduce_op, const tuple_type &src_tuple, SoAType &dst_soa, const int soa_idx, const std::integral_constant< std::size_t, 0 >)
 Unpack a sparse arrays tuple for it's member with index 0.
 
template<class ReduceOp, std::size_t N, class SoAType>
static KOKKOS_FORCEINLINE_FUNCTION void unpackTuple (const ReduceOp &reduce_op, const tuple_type &src_tuple, SoAType &dst_soa, const int soa_idx, const std::integral_constant< std::size_t, N >)
 Unpack a sparse arrays tuple for all members when element ID!=0.
 

Static Public Attributes

static constexpr std::size_t num_space_dim = NumSpaceDim
 sparse array dimension number
 
static constexpr std::size_t member_num = aosoa_member_types::size
 AoSoA member #.
 
static constexpr unsigned long long cell_bits_per_tile_dim
 sparse grid hierarchy: cell id bit# per dimension
 
static constexpr unsigned long long cell_num_per_tile
 sparse grid hierarchy: cell # per dimension
 

Detailed Description

template<class MemorySpace, class DataTypes, class EntityType, std::size_t NumSpaceDim, unsigned long long cellBitsPerTileDim, typename Value = int, typename Key = uint64_t>
class Cabana::Grid::Experimental::SparseHalo< MemorySpace, DataTypes, EntityType, NumSpaceDim, cellBitsPerTileDim, Value, Key >

General multiple array halo communication plan for migrating shared data between sparse grid blocks. Arrays may be defined on different entity types and have different data types.

The halo operates on an arbitrary set of arrays. Each of these arrays must be defined on the same local grid meaning they that share the same communicator and halo size. The arrays must also reside in the same memory space. These requirements are checked at construction.

Member Typedef Documentation

◆ counting_view

template<class MemorySpace, class DataTypes, class EntityType, std::size_t NumSpaceDim, unsigned long long cellBitsPerTileDim, typename Value = int, typename Key = uint64_t>
using Cabana::Grid::Experimental::SparseHalo< MemorySpace, DataTypes, EntityType, NumSpaceDim, cellBitsPerTileDim, Value, Key >::counting_view
Initial value:
Kokkos::View<int[2], memory_space,
Kokkos::MemoryTraits<Kokkos::Atomic>>
MemorySpace memory_space
memory space to store the sparse grid
Definition Cabana_Grid_SparseHalo.hpp:58

view type used to count common sparse grid # with neighbors [0] shared_owned_num [1] shared_ghost_num

◆ member_data_type

template<class MemorySpace, class DataTypes, class EntityType, std::size_t NumSpaceDim, unsigned long long cellBitsPerTileDim, typename Value = int, typename Key = uint64_t>
template<std::size_t M>
using Cabana::Grid::Experimental::SparseHalo< MemorySpace, DataTypes, EntityType, NumSpaceDim, cellBitsPerTileDim, Value, Key >::member_data_type
Initial value:
Get the type of the member at a given index.
Definition Cabana_MemberTypes.hpp:75

AoSoA member data type.

◆ steering_view

template<class MemorySpace, class DataTypes, class EntityType, std::size_t NumSpaceDim, unsigned long long cellBitsPerTileDim, typename Value = int, typename Key = uint64_t>
using Cabana::Grid::Experimental::SparseHalo< MemorySpace, DataTypes, EntityType, NumSpaceDim, cellBitsPerTileDim, Value, Key >::steering_view = Kokkos::View<key_type*, memory_space>

communication steering view type, used to check whether there are common sparse grids in the halo region

◆ tile_index_space

template<class MemorySpace, class DataTypes, class EntityType, std::size_t NumSpaceDim, unsigned long long cellBitsPerTileDim, typename Value = int, typename Key = uint64_t>
using Cabana::Grid::Experimental::SparseHalo< MemorySpace, DataTypes, EntityType, NumSpaceDim, cellBitsPerTileDim, Value, Key >::tile_index_space
Initial value:
Index space with tile as unit; _min and _max forms the tile range. Note this is for sparse grid only,...
Definition Cabana_Grid_SparseIndexSpace.hpp:1137

tile index space type TODO

Constructor & Destructor Documentation

◆ SparseHalo()

template<class MemorySpace, class DataTypes, class EntityType, std::size_t NumSpaceDim, unsigned long long cellBitsPerTileDim, typename Value = int, typename Key = uint64_t>
template<class SparseArrayType>
Cabana::Grid::Experimental::SparseHalo< MemorySpace, DataTypes, EntityType, NumSpaceDim, cellBitsPerTileDim, Value, Key >::SparseHalo ( const halo_pattern_type pattern,
const std::shared_ptr< SparseArrayType > & sparse_array )
inline

constructor

Template Parameters
LocalGridTypelocal grid type
Parameters
patternThe halo pattern to use for halo communication
sparse_arraySparse array to communicate

Member Function Documentation

◆ buildCommData()

template<class MemorySpace, class DataTypes, class EntityType, std::size_t NumSpaceDim, unsigned long long cellBitsPerTileDim, typename Value = int, typename Key = uint64_t>
template<class DecompositionTag, class LocalGridType>
void Cabana::Grid::Experimental::SparseHalo< MemorySpace, DataTypes, EntityType, NumSpaceDim, cellBitsPerTileDim, Value, Key >::buildCommData ( DecompositionTag decomposition_tag,
const std::shared_ptr< LocalGridType > & local_grid,
const std::array< int, num_space_dim > & nid,
std::vector< buffer_view > & buffers,
std::vector< steering_view > & steering,
std::vector< tile_index_space > & spaces )
inline

Build communication data.

Template Parameters
DecompositionTagdecomposition tag type
LocalGridTypesparse local grid type
Parameters
decomposition_tagtag to indicate if it's owned or ghosted halo
local_gridsparse local grid shared pointer
nidneighbor local id (ijk in pattern)
buffersbuffer to be used to store communicated data
steeringsteering to be used to guide communications
spacessparse tile index spaces

◆ clear()

template<class MemorySpace, class DataTypes, class EntityType, std::size_t NumSpaceDim, unsigned long long cellBitsPerTileDim, typename Value = int, typename Key = uint64_t>
void Cabana::Grid::Experimental::SparseHalo< MemorySpace, DataTypes, EntityType, NumSpaceDim, cellBitsPerTileDim, Value, Key >::clear ( MPI_Comm comm)
inline

clear guiding information in sparse halo,

Parameters
commMPI communicator

This information needs to be cleared before steering and counting, then recollected before halo communication in each step.

◆ collectNeighborCounting()

template<class MemorySpace, class DataTypes, class EntityType, std::size_t NumSpaceDim, unsigned long long cellBitsPerTileDim, typename Value = int, typename Key = uint64_t>
void Cabana::Grid::Experimental::SparseHalo< MemorySpace, DataTypes, EntityType, NumSpaceDim, cellBitsPerTileDim, Value, Key >::collectNeighborCounting ( MPI_Comm comm,
const bool is_neighbor_counting_collected = false ) const
inline

neighbor tile counting, communication needed only if the counting is non-zero

Parameters
commMPI communicator
is_neighbor_counting_collectedlabel if the neighbor has already been collected; if true, it means all neighbor counting information is up-to-date and there's no need for recollection

◆ gather()

template<class MemorySpace, class DataTypes, class EntityType, std::size_t NumSpaceDim, unsigned long long cellBitsPerTileDim, typename Value = int, typename Key = uint64_t>
template<class ExecSpace, class SparseArrayType>
void Cabana::Grid::Experimental::SparseHalo< MemorySpace, DataTypes, EntityType, NumSpaceDim, cellBitsPerTileDim, Value, Key >::gather ( const ExecSpace & exec_space,
SparseArrayType & sparse_array,
const bool is_neighbor_counting_collected = false ) const
inline

Gather data into our ghosted share space from their owners.

Template Parameters
ExecSpaceexecution space
SparseArrayTypesparse array type
SparseMapTypesparse map type
Parameters
exec_spaceexecution space
sparse_arraysparse AoSoA array used to store grid data
is_neighbor_counting_collectedlabel if the neighbor has already been collected; if true, it means all neighbor counting information is up-to-date and there's no need for recollection

◆ gatherValidSendAndRecvRanks()

template<class MemorySpace, class DataTypes, class EntityType, std::size_t NumSpaceDim, unsigned long long cellBitsPerTileDim, typename Value = int, typename Key = uint64_t>
void Cabana::Grid::Experimental::SparseHalo< MemorySpace, DataTypes, EntityType, NumSpaceDim, cellBitsPerTileDim, Value, Key >::gatherValidSendAndRecvRanks ( MPI_Comm comm,
std::vector< int > & valid_sends,
std::vector< int > & valid_recvs,
const bool is_neighbor_counting_collected = false ) const
inline

collect all valid ranks for sparse grid gather operations

Parameters
commMPI communicator
valid_sendsneighbor array id that requires data from current rank
valid_recvsneighbor array id the current ranks requires data from
is_neighbor_counting_collectedlabel if the neighbor has already been collected; if true, it means all neighbor counting information is up-to-date and there's no need for recollection

◆ packBuffer()

template<class MemorySpace, class DataTypes, class EntityType, std::size_t NumSpaceDim, unsigned long long cellBitsPerTileDim, typename Value = int, typename Key = uint64_t>
template<class ExecSpace, class SparseArrayType>
void Cabana::Grid::Experimental::SparseHalo< MemorySpace, DataTypes, EntityType, NumSpaceDim, cellBitsPerTileDim, Value, Key >::packBuffer ( const ExecSpace & exec_space,
const buffer_view & buffer,
const steering_view & tile_steering,
SparseArrayType & sparse_array,
const int count ) const
inline

Pack sparse arrays at halo regions into a buffer.

Template Parameters
ExecSpaceexecution space type
SparseArrayTypesparse array type
Parameters
exec_spaceexecution space
bufferbuffer to store sparse array data and to communicate
tile_steeringKokkos view to store halo tile keys
sparse_arraysparse array (all sparse grids on current rank)
countnumber of halo grids to pack

◆ register_halo()

template<class MemorySpace, class DataTypes, class EntityType, std::size_t NumSpaceDim, unsigned long long cellBitsPerTileDim, typename Value = int, typename Key = uint64_t>
template<class ExecSpace, class SparseMapType>
void Cabana::Grid::Experimental::SparseHalo< MemorySpace, DataTypes, EntityType, NumSpaceDim, cellBitsPerTileDim, Value, Key >::register_halo ( SparseMapType & map)
inline

register valid halos (according to grid activation status in sparse map) in the steerings

Template Parameters
ExecSpaceexecution space
SparseMapTypesparse map type
scalar_typescalar type (type of dx)
Parameters
mapsparse map

◆ scatter()

template<class MemorySpace, class DataTypes, class EntityType, std::size_t NumSpaceDim, unsigned long long cellBitsPerTileDim, typename Value = int, typename Key = uint64_t>
template<class ExecSpace, class ReduceOp, class SparseArrayType>
void Cabana::Grid::Experimental::SparseHalo< MemorySpace, DataTypes, EntityType, NumSpaceDim, cellBitsPerTileDim, Value, Key >::scatter ( const ExecSpace & exec_space,
const ReduceOp & reduce_op,
SparseArrayType & sparse_array,
const bool is_neighbor_counting_collected = false ) const
inline

Scatter data from our ghosts to their owners using the given type of reduce operation.

Template Parameters
ExecSpaceexecution space
ReduceOpThe type of reduction functor
SparseArrayTypesparse array type
SparseMapTypesparse map type
Parameters
exec_spaceexecution space
reduce_opThe functor used to reduce the results
sparse_arraysparse AoSoA array used to store grid data
is_neighbor_counting_collectedlabel if the neighbor has already been collected; if true, it means all neighbor counting information is up-to-date and there's no need for recollection

◆ scatterValidSendAndRecvRanks()

template<class MemorySpace, class DataTypes, class EntityType, std::size_t NumSpaceDim, unsigned long long cellBitsPerTileDim, typename Value = int, typename Key = uint64_t>
void Cabana::Grid::Experimental::SparseHalo< MemorySpace, DataTypes, EntityType, NumSpaceDim, cellBitsPerTileDim, Value, Key >::scatterValidSendAndRecvRanks ( MPI_Comm comm,
std::vector< int > & valid_sends,
std::vector< int > & valid_recvs,
const bool is_neighbor_counting_collected = false ) const
inline

collect all valid ranks for sparse grid scatter operations

Parameters
commMPI communicator
valid_sendsneighbor array id that requires data from current rank
valid_recvsneighbor array id the current ranks requires data from
is_neighbor_counting_collectedlabel if the neighbor has already been collected; if true, it means all neighbor counting information is up-to-date and there's no need for recollection

◆ unpackBuffer()

template<class MemorySpace, class DataTypes, class EntityType, std::size_t NumSpaceDim, unsigned long long cellBitsPerTileDim, typename Value = int, typename Key = uint64_t>
template<class ReduceOp, class ExecSpace, class SparseArrayType, class SparseMapType>
void Cabana::Grid::Experimental::SparseHalo< MemorySpace, DataTypes, EntityType, NumSpaceDim, cellBitsPerTileDim, Value, Key >::unpackBuffer ( const ReduceOp & reduce_op,
const ExecSpace & exec_space,
const buffer_view & buffer,
const steering_view & tile_steering,
const SparseArrayType & sparse_array,
SparseMapType & map,
const int count ) const
inline

Unpack a sparse array communication buffer.

Template Parameters
ReduceOpreduce functor type
ExecSpaceexecution space type
SparseArrayTypesparse array type
SparseMapTypesparse map type
Parameters
reduce_opreduce operation
exec_spaceexecution space
bufferbuffer to store sparse array data and to communicate
tile_steeringKokkos view to store halo tile keys
sparse_arraysparse array (all sparse grids on current rank)
mapsparse map that has valid grids registered
countnumber of halo grids to unpack

◆ unpackTuple() [1/2]

template<class MemorySpace, class DataTypes, class EntityType, std::size_t NumSpaceDim, unsigned long long cellBitsPerTileDim, typename Value = int, typename Key = uint64_t>
template<class ReduceOp, class SoAType>
static KOKKOS_FORCEINLINE_FUNCTION void Cabana::Grid::Experimental::SparseHalo< MemorySpace, DataTypes, EntityType, NumSpaceDim, cellBitsPerTileDim, Value, Key >::unpackTuple ( const ReduceOp & reduce_op,
const tuple_type & src_tuple,
SoAType & dst_soa,
const int soa_idx,
const std::integral_constant< std::size_t, 0 >  )
inlinestatic

Unpack a sparse arrays tuple for it's member with index 0.

Template Parameters
ReduceOpreduce functor type
SoATypeSoA type in sparse array (which is an AoSoA)
Parameters
reduce_opreduce operation
src_tuplesource tuple
dst_soadestination SoA to store copied data
soa_idxtuple index inside the destination SoA

◆ unpackTuple() [2/2]

template<class MemorySpace, class DataTypes, class EntityType, std::size_t NumSpaceDim, unsigned long long cellBitsPerTileDim, typename Value = int, typename Key = uint64_t>
template<class ReduceOp, std::size_t N, class SoAType>
static KOKKOS_FORCEINLINE_FUNCTION void Cabana::Grid::Experimental::SparseHalo< MemorySpace, DataTypes, EntityType, NumSpaceDim, cellBitsPerTileDim, Value, Key >::unpackTuple ( const ReduceOp & reduce_op,
const tuple_type & src_tuple,
SoAType & dst_soa,
const int soa_idx,
const std::integral_constant< std::size_t, N >  )
inlinestatic

Unpack a sparse arrays tuple for all members when element ID!=0.

Template Parameters
ReduceOpreduce functor type
SoATypeSoA type in sparse array (which is an AoSoA)
NUnpack N-th data member in this call
Parameters
reduce_opreduce operation
src_tuplesource tuple
dst_soadestination SoA to store copied data
soa_idxtuple index inside the destination SoA

◆ unpackTupleMember() [1/4]

template<class MemorySpace, class DataTypes, class EntityType, std::size_t NumSpaceDim, unsigned long long cellBitsPerTileDim, typename Value = int, typename Key = uint64_t>
template<class ReduceOp, std::size_t N, std::size_t M, class SoAType>
static KOKKOS_FORCEINLINE_FUNCTION std::enable_if_t< 0==M, void > Cabana::Grid::Experimental::SparseHalo< MemorySpace, DataTypes, EntityType, NumSpaceDim, cellBitsPerTileDim, Value, Key >::unpackTupleMember ( const ReduceOp & reduce_op,
const tuple_type & src_tuple,
SoAType & dst_soa,
const int soa_idx,
const Kokkos::Array< std::size_t, M > & ,
const std::integral_constant< std::size_t, N > ,
const std::integral_constant< std::size_t, M >  )
inlinestatic

Unpack a sparse arrays element (a tuple) in a buffer (for case tuple members with rank == 0)

Template Parameters
ReduceOpreduce functor type
Nelement ID inside a SoA tuple (N-th data member)
Mrank number of the current element (N-th data member)
SoATypeSoA type in sparse array (which is an AoSoA)
Parameters
reduce_opreduce operation
src_tuplesource tuple
dst_soadestination SoA to store copied data
soa_idxtuple index inside the destination SoA

◆ unpackTupleMember() [2/4]

template<class MemorySpace, class DataTypes, class EntityType, std::size_t NumSpaceDim, unsigned long long cellBitsPerTileDim, typename Value = int, typename Key = uint64_t>
template<class ReduceOp, std::size_t N, std::size_t M, class SoAType>
static KOKKOS_FORCEINLINE_FUNCTION std::enable_if_t< 3==M, void > Cabana::Grid::Experimental::SparseHalo< MemorySpace, DataTypes, EntityType, NumSpaceDim, cellBitsPerTileDim, Value, Key >::unpackTupleMember ( const ReduceOp & reduce_op,
const tuple_type & src_tuple,
SoAType & dst_soa,
const int soa_idx,
const Kokkos::Array< std::size_t, M > & extents,
const std::integral_constant< std::size_t, N > ,
const std::integral_constant< std::size_t, M >  )
inlinestatic

Unpack a sparse arrays element (a tuple) in a buffer (for case tuple members with rank == 3)

Template Parameters
ReduceOpreduce functor type
Nelement ID inside a SoA tuple (N-th data member)
Mrank number of the current element (N-th data member)
SoATypeSoA type in sparse array (which is an AoSoA)
Parameters
reduce_opreduce operation
src_tuplesource tuple
dst_soadestination SoA to store copied data
soa_idxtuple index inside the destination SoA
extentselement member extents in all ranks

◆ unpackTupleMember() [3/4]

template<class MemorySpace, class DataTypes, class EntityType, std::size_t NumSpaceDim, unsigned long long cellBitsPerTileDim, typename Value = int, typename Key = uint64_t>
template<class ReduceOp, std::size_t N, std::size_t M, class SoAType>
static KOKKOS_FORCEINLINE_FUNCTION std::enable_if_t< 2==M, void > Cabana::Grid::Experimental::SparseHalo< MemorySpace, DataTypes, EntityType, NumSpaceDim, cellBitsPerTileDim, Value, Key >::unpackTupleMember ( const ReduceOp & reduce_op,
const tuple_type & src_tuple,
SoAType & dst_soa,
const int soa_idx,
const Kokkos::Array< std::size_t, M > & extents,
const std::integral_constant< std::size_t, N > ,
const std::integral_constant< std::size_t, M >  )
inlinestatic

Unpack a sparse arrays element (a tuple) in a buffer (for case tuple members with rank == 2)

Template Parameters
ReduceOpreduce functor type
Nelement ID inside a SoA tuple (N-th data member)
Mrank number of the current element (N-th data member)
SoATypeSoA type in sparse array (which is an AoSoA)
Parameters
reduce_opreduce operation
src_tuplesource tuple
dst_soadestination SoA to store copied data
soa_idxtuple index inside the destination SoA
extentselement member extents in all ranks

◆ unpackTupleMember() [4/4]

template<class MemorySpace, class DataTypes, class EntityType, std::size_t NumSpaceDim, unsigned long long cellBitsPerTileDim, typename Value = int, typename Key = uint64_t>
template<class ReduceOp, std::size_t N, std::size_t M, class SoAType>
static KOKKOS_FORCEINLINE_FUNCTION std::enable_if_t< 1==M, void > Cabana::Grid::Experimental::SparseHalo< MemorySpace, DataTypes, EntityType, NumSpaceDim, cellBitsPerTileDim, Value, Key >::unpackTupleMember ( const ReduceOp & reduce_op,
const tuple_type & src_tuple,
SoAType & dst_soa,
const int soa_idx,
const Kokkos::Array< std::size_t, M > & extents,
const std::integral_constant< std::size_t, N > ,
const std::integral_constant< std::size_t, M >  )
inlinestatic

Unpack a sparse arrays element (a tuple) in a buffer (for case tuple members with rank == 1)

Template Parameters
ReduceOpreduce functor type
Nelement ID inside a SoA tuple (N-th data member)
Mrank number of the current element (N-th data member)
SoATypeSoA type in sparse array (which is an AoSoA)
Parameters
reduce_opreduce operation
src_tuplesource tuple
dst_soadestination SoA to store copied data
soa_idxtuple index inside the destination SoA
extentselement member extents in all ranks

◆ updateTileSpace()

template<class MemorySpace, class DataTypes, class EntityType, std::size_t NumSpaceDim, unsigned long long cellBitsPerTileDim, typename Value = int, typename Key = uint64_t>
template<class LocalGridType>
void Cabana::Grid::Experimental::SparseHalo< MemorySpace, DataTypes, EntityType, NumSpaceDim, cellBitsPerTileDim, Value, Key >::updateTileSpace ( const std::shared_ptr< LocalGridType > & local_grid)
inline

update tile index space according to current partition

Template Parameters
LocalGridTypesparse local grid type
Parameters
local_gridsparse local grid pointer

Member Data Documentation

◆ cell_bits_per_tile_dim

template<class MemorySpace, class DataTypes, class EntityType, std::size_t NumSpaceDim, unsigned long long cellBitsPerTileDim, typename Value = int, typename Key = uint64_t>
unsigned long long Cabana::Grid::Experimental::SparseHalo< MemorySpace, DataTypes, EntityType, NumSpaceDim, cellBitsPerTileDim, Value, Key >::cell_bits_per_tile_dim
staticconstexpr
Initial value:
=
cellBitsPerTileDim

sparse grid hierarchy: cell id bit# per dimension

◆ cell_num_per_tile

template<class MemorySpace, class DataTypes, class EntityType, std::size_t NumSpaceDim, unsigned long long cellBitsPerTileDim, typename Value = int, typename Key = uint64_t>
unsigned long long Cabana::Grid::Experimental::SparseHalo< MemorySpace, DataTypes, EntityType, NumSpaceDim, cellBitsPerTileDim, Value, Key >::cell_num_per_tile
staticconstexpr
Initial value:
=
static constexpr unsigned long long cell_bits_per_tile_dim
Definition Cabana_Grid_SparseIndexSpace.hpp:1143

sparse grid hierarchy: cell # per dimension


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