Cabana 0.8.0-dev
 
Loading...
Searching...
No Matches
Cabana::Slice< DataType, MemorySpace, MemoryAccessType, VectorLength, Stride > Class Template Reference

A slice of an array-of-structs-of-arrays with data access to a single multidimensional member. More...

#include <Cabana_Slice.hpp>

Public Types

enum  { rank = std::rank<DataType>::value + 1 }
 
using slice_type
 Slice type.
 
using memory_space = typename MemorySpace::memory_space
 Memory space.
 
using execution_space = typename memory_space::execution_space
 Default execution space.
 
using memory_access_type = MemoryAccessType
 Memory access type.
 
using size_type = typename memory_space::size_type
 Memory space size type.
 
using index_type = Impl::Index<vector_length>
 Index type.
 
using view_wrapper
 Kokkos view wrapper.
 
using kokkos_view
 Kokkos view type.
 
using reference_type = typename kokkos_view::reference_type
 View reference type alias.
 
using value_type = typename kokkos_view::value_type
 View value type alias.
 
using pointer_type = typename kokkos_view::pointer_type
 View pointer type alias.
 
using view_layout = typename kokkos_view::array_layout
 View array layout type alias.
 
using default_access_slice
 Default memory access slice type.
 
using atomic_access_slice
 Atomic memory access slice type.
 
using random_access_slice
 Random memory access slice type.
 

Public Member Functions

 Slice ()
 Default constructor.
 
 Slice (const pointer_type data, const size_type size, const size_type num_soa, const std::string &label="")
 Constructor.
 
template<class MAT>
 Slice (const Slice< DataType, MemorySpace, MAT, VectorLength, Stride > &rhs)
 Shallow copy constructor for different memory spaces for assigning new memory access traits to the view.
 
template<class MAT>
Sliceoperator= (const Slice< DataType, MemorySpace, MAT, VectorLength, Stride > &rhs)
 Assignment operator for different memory spaces for assigning new memory access traits to the view.
 
std::string label () const
 Returns the data structure label.
 
KOKKOS_INLINE_FUNCTION size_type size () const
 Returns the total number tuples in the slice.
 
KOKKOS_INLINE_FUNCTION size_type numSoA () const
 Get the number of structs-of-arrays in the container.
 
KOKKOS_INLINE_FUNCTION size_type arraySize (const size_type s) const
 Get the size of the data array at a given struct member index.
 
template<typename U = DataType>
KOKKOS_FORCEINLINE_FUNCTION std::enable_if<(0==std::rank< U >::value &&std::is_same< U, DataType >::value), reference_type >::type access (const size_type s, const size_type a) const
 2d access for Rank 0
 
template<typename U = DataType>
KOKKOS_FORCEINLINE_FUNCTION std::enable_if<(1==std::rank< U >::value &&std::is_same< U, DataType >::value), reference_type >::type access (const size_type s, const size_type a, const size_type d0) const
 2d access for Rank 1
 
template<typename U = DataType>
KOKKOS_FORCEINLINE_FUNCTION std::enable_if<(2==std::rank< U >::value &&std::is_same< U, DataType >::value), reference_type >::type access (const size_type s, const size_type a, const size_type d0, const size_type d1) const
 2d access for Rank 2
 
template<typename U = DataType>
KOKKOS_FORCEINLINE_FUNCTION std::enable_if<(3==std::rank< U >::value &&std::is_same< U, DataType >::value), reference_type >::type access (const size_type s, const size_type a, const size_type d0, const size_type d1, const size_type d2) const
 2d access for Rank 3
 
template<typename U = DataType>
KOKKOS_FORCEINLINE_FUNCTION std::enable_if<(0==std::rank< U >::value &&std::is_same< U, DataType >::value), reference_type >::type operator() (const size_type i) const
 1d access for Rank 0
 
template<typename U = DataType>
KOKKOS_FORCEINLINE_FUNCTION std::enable_if<(1==std::rank< U >::value &&std::is_same< U, DataType >::value), reference_type >::type operator() (const size_type i, const size_type d0) const
 1d access for Rank 1
 
template<typename U = DataType>
KOKKOS_FORCEINLINE_FUNCTION std::enable_if<(2==std::rank< U >::value &&std::is_same< U, DataType >::value), reference_type >::type operator() (const size_type i, const size_type d0, const size_type d1) const
 1d access for Rank 2
 
template<typename U = DataType>
KOKKOS_FORCEINLINE_FUNCTION std::enable_if<(3==std::rank< U >::value &&std::is_same< U, DataType >::value), reference_type >::type operator() (const size_type i, const size_type d0, const size_type d1, const size_type d2) const
 1d access for Rank 3
 
KOKKOS_INLINE_FUNCTION pointer_type data () const
 Get a raw pointer to the data for this member.
 
KOKKOS_INLINE_FUNCTION constexpr size_type viewRank () const
 Get the rank of the raw data for this slice. This includes the struct dimension, array dimension, and all tuple slice dimensions.
 
KOKKOS_INLINE_FUNCTION size_type extent (const size_type d) const
 Get the extent of a given raw slice data dimension. This includes the struct dimension, array dimension, and all tuple slice dimensions.
 
KOKKOS_INLINE_FUNCTION size_type stride (const size_type d) const
 Get the stride of a given raw slice dimension. This includes the struct dimension, array dimension, and all tuple slice dimensions.
 
KOKKOS_INLINE_FUNCTION kokkos_view view () const
 Get the underlying Kokkos View managing the slice data.
 

Static Public Attributes

static constexpr int vector_length = VectorLength
 Vector length.
 
static constexpr int soa_stride = Stride
 SoA stride.
 
static constexpr std::size_t max_supported_rank = 3
 Maximum supported rank.
 
static constexpr std::size_t max_label_length = 128
 Maximum label length.
 

Friends

class Slice< DataType, MemorySpace, DefaultAccessMemory, VectorLength, Stride >
 
class Slice< DataType, MemorySpace, AtomicAccessMemory, VectorLength, Stride >
 
class Slice< DataType, MemorySpace, RandomAccessMemory, VectorLength, Stride >
 

Detailed Description

template<typename DataType, typename MemorySpace, typename MemoryAccessType, int VectorLength, int Stride>
class Cabana::Slice< DataType, MemorySpace, MemoryAccessType, VectorLength, Stride >

A slice of an array-of-structs-of-arrays with data access to a single multidimensional member.

Member Typedef Documentation

◆ atomic_access_slice

template<typename DataType, typename MemorySpace, typename MemoryAccessType, int VectorLength, int Stride>
using Cabana::Slice< DataType, MemorySpace, MemoryAccessType, VectorLength, Stride >::atomic_access_slice
Initial value:
Slice()
Default constructor.
Definition Cabana_Slice.hpp:598

Atomic memory access slice type.

◆ default_access_slice

template<typename DataType, typename MemorySpace, typename MemoryAccessType, int VectorLength, int Stride>
using Cabana::Slice< DataType, MemorySpace, MemoryAccessType, VectorLength, Stride >::default_access_slice
Initial value:

Default memory access slice type.

◆ kokkos_view

template<typename DataType, typename MemorySpace, typename MemoryAccessType, int VectorLength, int Stride>
using Cabana::Slice< DataType, MemorySpace, MemoryAccessType, VectorLength, Stride >::kokkos_view
Initial value:
Kokkos::View<typename view_wrapper::data_type,
typename view_wrapper::cabana_layout, MemorySpace,
typename MemoryAccessType::kokkos_memory_traits>

Kokkos view type.

◆ random_access_slice

template<typename DataType, typename MemorySpace, typename MemoryAccessType, int VectorLength, int Stride>
using Cabana::Slice< DataType, MemorySpace, MemoryAccessType, VectorLength, Stride >::random_access_slice
Initial value:

Random memory access slice type.

◆ slice_type

template<typename DataType, typename MemorySpace, typename MemoryAccessType, int VectorLength, int Stride>
using Cabana::Slice< DataType, MemorySpace, MemoryAccessType, VectorLength, Stride >::slice_type

◆ view_wrapper

template<typename DataType, typename MemorySpace, typename MemoryAccessType, int VectorLength, int Stride>
using Cabana::Slice< DataType, MemorySpace, MemoryAccessType, VectorLength, Stride >::view_wrapper
Initial value:
Impl::KokkosViewWrapper<DataType, vector_length, soa_stride>

Kokkos view wrapper.

Constructor & Destructor Documentation

◆ Slice() [1/2]

template<typename DataType, typename MemorySpace, typename MemoryAccessType, int VectorLength, int Stride>
Cabana::Slice< DataType, MemorySpace, MemoryAccessType, VectorLength, Stride >::Slice ( const pointer_type data,
const size_type size,
const size_type num_soa,
const std::string & label = "" )
inline

Constructor.

Parameters
dataPointer to the first member of the slice.
sizeThe number of tuples in the slice.
num_soaThe number of structs in the slice.
labelAn optional label for the slice.

◆ Slice() [2/2]

template<typename DataType, typename MemorySpace, typename MemoryAccessType, int VectorLength, int Stride>
template<class MAT>
Cabana::Slice< DataType, MemorySpace, MemoryAccessType, VectorLength, Stride >::Slice ( const Slice< DataType, MemorySpace, MAT, VectorLength, Stride > & rhs)
inline

Shallow copy constructor for different memory spaces for assigning new memory access traits to the view.

Template Parameters
MATMemory access type.
Parameters
rhsThe slice to shallow copy with a potentially different memory space.

Member Function Documentation

◆ arraySize()

template<typename DataType, typename MemorySpace, typename MemoryAccessType, int VectorLength, int Stride>
KOKKOS_INLINE_FUNCTION size_type Cabana::Slice< DataType, MemorySpace, MemoryAccessType, VectorLength, Stride >::arraySize ( const size_type s) const
inline

Get the size of the data array at a given struct member index.

Parameters
sThe struct index to get the array size for.
Returns
The size of the array at the given struct index.

◆ data()

template<typename DataType, typename MemorySpace, typename MemoryAccessType, int VectorLength, int Stride>
KOKKOS_INLINE_FUNCTION pointer_type Cabana::Slice< DataType, MemorySpace, MemoryAccessType, VectorLength, Stride >::data ( ) const
inline

Get a raw pointer to the data for this member.

Returns
A raw pointer to the data for this slice.

◆ extent()

template<typename DataType, typename MemorySpace, typename MemoryAccessType, int VectorLength, int Stride>
KOKKOS_INLINE_FUNCTION size_type Cabana::Slice< DataType, MemorySpace, MemoryAccessType, VectorLength, Stride >::extent ( const size_type d) const
inline

Get the extent of a given raw slice data dimension. This includes the struct dimension, array dimension, and all tuple slice dimensions.

Parameters
dThe member data dimension to get the extent for.
Returns
The extent of the given member data dimension.

◆ label()

template<typename DataType, typename MemorySpace, typename MemoryAccessType, int VectorLength, int Stride>
std::string Cabana::Slice< DataType, MemorySpace, MemoryAccessType, VectorLength, Stride >::label ( ) const
inline

Returns the data structure label.

Returns
A string identifying the data structure.

◆ numSoA()

template<typename DataType, typename MemorySpace, typename MemoryAccessType, int VectorLength, int Stride>
KOKKOS_INLINE_FUNCTION size_type Cabana::Slice< DataType, MemorySpace, MemoryAccessType, VectorLength, Stride >::numSoA ( ) const
inline

Get the number of structs-of-arrays in the container.

Returns
The number of structs-of-arrays in the container.

◆ operator=()

template<typename DataType, typename MemorySpace, typename MemoryAccessType, int VectorLength, int Stride>
template<class MAT>
Slice & Cabana::Slice< DataType, MemorySpace, MemoryAccessType, VectorLength, Stride >::operator= ( const Slice< DataType, MemorySpace, MAT, VectorLength, Stride > & rhs)
inline

Assignment operator for different memory spaces for assigning new memory access traits to the view.

Template Parameters
MATMemory access type
Parameters
rhsThe slice to shallow copy with a potentially different memory space.
Returns
A reference to a new slice with a potentially different memory space.

◆ size()

template<typename DataType, typename MemorySpace, typename MemoryAccessType, int VectorLength, int Stride>
KOKKOS_INLINE_FUNCTION size_type Cabana::Slice< DataType, MemorySpace, MemoryAccessType, VectorLength, Stride >::size ( ) const
inline

Returns the total number tuples in the slice.

Returns
The number of tuples in the slice.

◆ stride()

template<typename DataType, typename MemorySpace, typename MemoryAccessType, int VectorLength, int Stride>
KOKKOS_INLINE_FUNCTION size_type Cabana::Slice< DataType, MemorySpace, MemoryAccessType, VectorLength, Stride >::stride ( const size_type d) const
inline

Get the stride of a given raw slice dimension. This includes the struct dimension, array dimension, and all tuple slice dimensions.

Parameters
dThe member data dimension to get the stride for.
Returns
The stride of the given member data dimension.

◆ viewRank()

template<typename DataType, typename MemorySpace, typename MemoryAccessType, int VectorLength, int Stride>
KOKKOS_INLINE_FUNCTION constexpr size_type Cabana::Slice< DataType, MemorySpace, MemoryAccessType, VectorLength, Stride >::viewRank ( ) const
inlineconstexpr

Get the rank of the raw data for this slice. This includes the struct dimension, array dimension, and all tuple slice dimensions.

Returns
The rank of the data for this slice.

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