Stokhos Package Browser (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Classes | Public Types | Public Member Functions | Static Public Member Functions | Static Public Attributes | Private Types | Private Attributes | Friends | List of all members
Stokhos::CrsProductTensor< ValueType, ExecutionSpace, Memory > Class Template Reference

Sparse product tensor with replicated entries to provide subsets with a given coordinate. More...

#include <Stokhos_CrsProductTensor.hpp>

Inheritance diagram for Stokhos::CrsProductTensor< ValueType, ExecutionSpace, Memory >:
Inheritance graph
[legend]

Classes

struct  CijkRowCount
 
struct  CompareCijkRowCount
 

Public Types

typedef ExecutionSpace execution_space
 
typedef int size_type
 
typedef ValueType value_type
 
typedef Memory memory_type
 
typedef Kokkos::ViewTraits
< size_type *, execution_space,
void, void >
::host_mirror_space 
host_mirror_space
 
typedef CrsProductTensor
< value_type,
host_mirror_space
HostMirror
 

Public Member Functions

KOKKOS_INLINE_FUNCTION ~CrsProductTensor ()
 
KOKKOS_INLINE_FUNCTION CrsProductTensor ()
 
template<class M >
KOKKOS_INLINE_FUNCTION CrsProductTensor (const CrsProductTensor< value_type, execution_space, M > &rhs)
 
template<class M >
KOKKOS_INLINE_FUNCTION
CrsProductTensor
operator= (const CrsProductTensor< value_type, execution_space, M > &rhs)
 
KOKKOS_INLINE_FUNCTION size_type dimension () const
 Dimension of the tensor. More...
 
KOKKOS_INLINE_FUNCTION bool is_empty () const
 Is the tensor empty. More...
 
KOKKOS_INLINE_FUNCTION size_type entry_count () const
 Number of sparse entries. More...
 
KOKKOS_INLINE_FUNCTION size_type entry_maximum () const
 Maximum sparse entries for any coordinate. More...
 
KOKKOS_INLINE_FUNCTION size_type entry_begin (size_type i) const
 Begin entries with a coordinate 'i'. More...
 
KOKKOS_INLINE_FUNCTION size_type entry_end (size_type i) const
 End entries with a coordinate 'i'. More...
 
KOKKOS_INLINE_FUNCTION size_type num_entry (size_type i) const
 Number of entries with a coordinate 'i'. More...
 
KOKKOS_INLINE_FUNCTION const
size_type
coord (const size_type entry, const size_type c) const
 Coordinates of an entry. More...
 
KOKKOS_INLINE_FUNCTION const
size_type
coord (const size_type entry) const
 Coordinates of an entry. More...
 
KOKKOS_INLINE_FUNCTION const
value_type
value (const size_type entry) const
 Value of an entry. More...
 
KOKKOS_INLINE_FUNCTION size_type num_non_zeros () const
 Number of non-zero's. More...
 
KOKKOS_INLINE_FUNCTION size_type num_flops () const
 Number flop's per multiply-add. More...
 
KOKKOS_INLINE_FUNCTION size_type avg_entries_per_row () const
 Number average number of entries per row. More...
 

Static Public Member Functions

template<typename OrdinalType >
static CrsProductTensor create (const Stokhos::ProductBasis< OrdinalType, ValueType > &basis, const Stokhos::Sparse3Tensor< OrdinalType, ValueType > &Cijk, const Teuchos::ParameterList &params=Teuchos::ParameterList())
 
static CrsProductTensor createMeanBased ()
 
static HostMirror create_mirror_view (const CrsProductTensor &tensor)
 
template<class DstDevice , class DstMemory >
static void deep_copy (const CrsProductTensor< ValueType, DstDevice, DstMemory > &dst, const CrsProductTensor &src)
 

Static Public Attributes

static const size_type host_vectorsize = 2
 
static const bool use_intrinsics = false
 
static const size_type num_entry_align = 1
 
static const size_type cuda_vectorsize = 32
 
static const bool is_cuda
 
static const size_type vectorsize = is_cuda ? cuda_vectorsize : host_vectorsize
 
static const size_type tensor_align = vectorsize
 

Private Types

typedef Kokkos::View
< value_type
*, Kokkos::LayoutLeft,
execution_space, memory_type
vec_type
 
typedef Kokkos::View
< size_type
*, Kokkos::LayoutLeft,
execution_space, memory_type
coord_array_type
 
typedef Kokkos::View
< size_type
*[2], Kokkos::LayoutLeft,
execution_space, memory_type
coord2_array_type
 
typedef Kokkos::View
< value_type
*, Kokkos::LayoutLeft,
execution_space, memory_type
value_array_type
 
typedef Kokkos::View
< size_type
*, Kokkos::LayoutLeft,
execution_space, memory_type
entry_array_type
 
typedef Kokkos::View
< size_type
*, Kokkos::LayoutLeft,
execution_space, memory_type
row_map_array_type
 

Private Attributes

coord_array_type m_coord
 
coord2_array_type m_coord2
 
value_array_type m_value
 
entry_array_type m_num_entry
 
row_map_array_type m_row_map
 
size_type m_dim
 
size_type m_entry_max
 
size_type m_nnz
 
size_type m_flops
 
size_type m_avg_entries_per_row
 

Friends

template<class , class , class >
class CrsProductTensor
 

Detailed Description

template<typename ValueType, class ExecutionSpace, class Memory = void>
class Stokhos::CrsProductTensor< ValueType, ExecutionSpace, Memory >

Sparse product tensor with replicated entries to provide subsets with a given coordinate.

This allows product tensor multiplication to be partitioned on a given coordinate values.

for ( size_type i = 0 ; i < p.dimension() ; ++i ) { y[i] = 0 ; for ( size_type e = p.entry_begin(i) ; e < p.entry_end(i) ; ++e ) { const size_type j = p.coord(e,0); const size_type k = p.coord(e,1); Scalar tmp = a[j] * x[k] ; if ( j != k ) tmp += a[k] * x[j] ; y[i] += p.value(e) * tmp ; } }

Definition at line 46 of file Stokhos_CrsProductTensor.hpp.

Member Typedef Documentation

template<typename ValueType, class ExecutionSpace, class Memory = void>
typedef ExecutionSpace Stokhos::CrsProductTensor< ValueType, ExecutionSpace, Memory >::execution_space

Definition at line 49 of file Stokhos_CrsProductTensor.hpp.

template<typename ValueType, class ExecutionSpace, class Memory = void>
typedef int Stokhos::CrsProductTensor< ValueType, ExecutionSpace, Memory >::size_type

Definition at line 50 of file Stokhos_CrsProductTensor.hpp.

template<typename ValueType, class ExecutionSpace, class Memory = void>
typedef ValueType Stokhos::CrsProductTensor< ValueType, ExecutionSpace, Memory >::value_type

Definition at line 51 of file Stokhos_CrsProductTensor.hpp.

template<typename ValueType, class ExecutionSpace, class Memory = void>
typedef Memory Stokhos::CrsProductTensor< ValueType, ExecutionSpace, Memory >::memory_type

Definition at line 52 of file Stokhos_CrsProductTensor.hpp.

template<typename ValueType, class ExecutionSpace, class Memory = void>
typedef Kokkos::ViewTraits< size_type*, execution_space,void,void >::host_mirror_space Stokhos::CrsProductTensor< ValueType, ExecutionSpace, Memory >::host_mirror_space

Definition at line 54 of file Stokhos_CrsProductTensor.hpp.

template<typename ValueType, class ExecutionSpace, class Memory = void>
typedef CrsProductTensor<value_type, host_mirror_space> Stokhos::CrsProductTensor< ValueType, ExecutionSpace, Memory >::HostMirror

Definition at line 55 of file Stokhos_CrsProductTensor.hpp.

template<typename ValueType, class ExecutionSpace, class Memory = void>
typedef Kokkos::View< value_type*, Kokkos::LayoutLeft, execution_space, memory_type > Stokhos::CrsProductTensor< ValueType, ExecutionSpace, Memory >::vec_type
private

Definition at line 87 of file Stokhos_CrsProductTensor.hpp.

template<typename ValueType, class ExecutionSpace, class Memory = void>
typedef Kokkos::View< size_type*, Kokkos::LayoutLeft, execution_space, memory_type > Stokhos::CrsProductTensor< ValueType, ExecutionSpace, Memory >::coord_array_type
private

Definition at line 88 of file Stokhos_CrsProductTensor.hpp.

template<typename ValueType, class ExecutionSpace, class Memory = void>
typedef Kokkos::View< size_type*[2], Kokkos::LayoutLeft, execution_space, memory_type > Stokhos::CrsProductTensor< ValueType, ExecutionSpace, Memory >::coord2_array_type
private

Definition at line 89 of file Stokhos_CrsProductTensor.hpp.

template<typename ValueType, class ExecutionSpace, class Memory = void>
typedef Kokkos::View< value_type*, Kokkos::LayoutLeft, execution_space, memory_type > Stokhos::CrsProductTensor< ValueType, ExecutionSpace, Memory >::value_array_type
private

Definition at line 90 of file Stokhos_CrsProductTensor.hpp.

template<typename ValueType, class ExecutionSpace, class Memory = void>
typedef Kokkos::View< size_type*, Kokkos::LayoutLeft, execution_space, memory_type > Stokhos::CrsProductTensor< ValueType, ExecutionSpace, Memory >::entry_array_type
private

Definition at line 91 of file Stokhos_CrsProductTensor.hpp.

template<typename ValueType, class ExecutionSpace, class Memory = void>
typedef Kokkos::View< size_type*, Kokkos::LayoutLeft, execution_space, memory_type > Stokhos::CrsProductTensor< ValueType, ExecutionSpace, Memory >::row_map_array_type
private

Definition at line 92 of file Stokhos_CrsProductTensor.hpp.

Constructor & Destructor Documentation

template<typename ValueType, class ExecutionSpace, class Memory = void>
KOKKOS_INLINE_FUNCTION Stokhos::CrsProductTensor< ValueType, ExecutionSpace, Memory >::~CrsProductTensor ( )
inline

Definition at line 124 of file Stokhos_CrsProductTensor.hpp.

template<typename ValueType, class ExecutionSpace, class Memory = void>
KOKKOS_INLINE_FUNCTION Stokhos::CrsProductTensor< ValueType, ExecutionSpace, Memory >::CrsProductTensor ( )
inline

Definition at line 127 of file Stokhos_CrsProductTensor.hpp.

template<typename ValueType, class ExecutionSpace, class Memory = void>
template<class M >
KOKKOS_INLINE_FUNCTION Stokhos::CrsProductTensor< ValueType, ExecutionSpace, Memory >::CrsProductTensor ( const CrsProductTensor< value_type, execution_space, M > &  rhs)
inline

Definition at line 141 of file Stokhos_CrsProductTensor.hpp.

Member Function Documentation

template<typename ValueType, class ExecutionSpace, class Memory = void>
template<class M >
KOKKOS_INLINE_FUNCTION CrsProductTensor& Stokhos::CrsProductTensor< ValueType, ExecutionSpace, Memory >::operator= ( const CrsProductTensor< value_type, execution_space, M > &  rhs)
inline

Definition at line 156 of file Stokhos_CrsProductTensor.hpp.

template<typename ValueType, class ExecutionSpace, class Memory = void>
KOKKOS_INLINE_FUNCTION size_type Stokhos::CrsProductTensor< ValueType, ExecutionSpace, Memory >::dimension ( ) const
inline

Dimension of the tensor.

Definition at line 173 of file Stokhos_CrsProductTensor.hpp.

template<typename ValueType, class ExecutionSpace, class Memory = void>
KOKKOS_INLINE_FUNCTION bool Stokhos::CrsProductTensor< ValueType, ExecutionSpace, Memory >::is_empty ( ) const
inline

Is the tensor empty.

Definition at line 177 of file Stokhos_CrsProductTensor.hpp.

template<typename ValueType, class ExecutionSpace, class Memory = void>
KOKKOS_INLINE_FUNCTION size_type Stokhos::CrsProductTensor< ValueType, ExecutionSpace, Memory >::entry_count ( ) const
inline

Number of sparse entries.

Definition at line 181 of file Stokhos_CrsProductTensor.hpp.

template<typename ValueType, class ExecutionSpace, class Memory = void>
KOKKOS_INLINE_FUNCTION size_type Stokhos::CrsProductTensor< ValueType, ExecutionSpace, Memory >::entry_maximum ( ) const
inline

Maximum sparse entries for any coordinate.

Definition at line 186 of file Stokhos_CrsProductTensor.hpp.

template<typename ValueType, class ExecutionSpace, class Memory = void>
KOKKOS_INLINE_FUNCTION size_type Stokhos::CrsProductTensor< ValueType, ExecutionSpace, Memory >::entry_begin ( size_type  i) const
inline

Begin entries with a coordinate 'i'.

Definition at line 191 of file Stokhos_CrsProductTensor.hpp.

template<typename ValueType, class ExecutionSpace, class Memory = void>
KOKKOS_INLINE_FUNCTION size_type Stokhos::CrsProductTensor< ValueType, ExecutionSpace, Memory >::entry_end ( size_type  i) const
inline

End entries with a coordinate 'i'.

Definition at line 196 of file Stokhos_CrsProductTensor.hpp.

template<typename ValueType, class ExecutionSpace, class Memory = void>
KOKKOS_INLINE_FUNCTION size_type Stokhos::CrsProductTensor< ValueType, ExecutionSpace, Memory >::num_entry ( size_type  i) const
inline

Number of entries with a coordinate 'i'.

Definition at line 201 of file Stokhos_CrsProductTensor.hpp.

template<typename ValueType, class ExecutionSpace, class Memory = void>
KOKKOS_INLINE_FUNCTION const size_type& Stokhos::CrsProductTensor< ValueType, ExecutionSpace, Memory >::coord ( const size_type  entry,
const size_type  c 
) const
inline

Coordinates of an entry.

Definition at line 206 of file Stokhos_CrsProductTensor.hpp.

template<typename ValueType, class ExecutionSpace, class Memory = void>
KOKKOS_INLINE_FUNCTION const size_type& Stokhos::CrsProductTensor< ValueType, ExecutionSpace, Memory >::coord ( const size_type  entry) const
inline

Coordinates of an entry.

Definition at line 211 of file Stokhos_CrsProductTensor.hpp.

template<typename ValueType, class ExecutionSpace, class Memory = void>
KOKKOS_INLINE_FUNCTION const value_type& Stokhos::CrsProductTensor< ValueType, ExecutionSpace, Memory >::value ( const size_type  entry) const
inline

Value of an entry.

Definition at line 216 of file Stokhos_CrsProductTensor.hpp.

template<typename ValueType, class ExecutionSpace, class Memory = void>
KOKKOS_INLINE_FUNCTION size_type Stokhos::CrsProductTensor< ValueType, ExecutionSpace, Memory >::num_non_zeros ( ) const
inline

Number of non-zero's.

Definition at line 221 of file Stokhos_CrsProductTensor.hpp.

template<typename ValueType, class ExecutionSpace, class Memory = void>
KOKKOS_INLINE_FUNCTION size_type Stokhos::CrsProductTensor< ValueType, ExecutionSpace, Memory >::num_flops ( ) const
inline

Number flop's per multiply-add.

Definition at line 226 of file Stokhos_CrsProductTensor.hpp.

template<typename ValueType, class ExecutionSpace, class Memory = void>
KOKKOS_INLINE_FUNCTION size_type Stokhos::CrsProductTensor< ValueType, ExecutionSpace, Memory >::avg_entries_per_row ( ) const
inline

Number average number of entries per row.

Definition at line 231 of file Stokhos_CrsProductTensor.hpp.

template<typename ValueType, class ExecutionSpace, class Memory = void>
template<typename OrdinalType >
static CrsProductTensor Stokhos::CrsProductTensor< ValueType, ExecutionSpace, Memory >::create ( const Stokhos::ProductBasis< OrdinalType, ValueType > &  basis,
const Stokhos::Sparse3Tensor< OrdinalType, ValueType > &  Cijk,
const Teuchos::ParameterList params = Teuchos::ParameterList() 
)
inlinestatic

Definition at line 236 of file Stokhos_CrsProductTensor.hpp.

template<typename ValueType, class ExecutionSpace, class Memory = void>
static CrsProductTensor Stokhos::CrsProductTensor< ValueType, ExecutionSpace, Memory >::createMeanBased ( )
inlinestatic

Definition at line 392 of file Stokhos_CrsProductTensor.hpp.

template<typename ValueType, class ExecutionSpace, class Memory = void>
static HostMirror Stokhos::CrsProductTensor< ValueType, ExecutionSpace, Memory >::create_mirror_view ( const CrsProductTensor< ValueType, ExecutionSpace, Memory > &  tensor)
inlinestatic

Definition at line 449 of file Stokhos_CrsProductTensor.hpp.

template<typename ValueType, class ExecutionSpace, class Memory = void>
template<class DstDevice , class DstMemory >
static void Stokhos::CrsProductTensor< ValueType, ExecutionSpace, Memory >::deep_copy ( const CrsProductTensor< ValueType, DstDevice, DstMemory > &  dst,
const CrsProductTensor< ValueType, ExecutionSpace, Memory > &  src 
)
inlinestatic

Definition at line 469 of file Stokhos_CrsProductTensor.hpp.

Friends And Related Function Documentation

template<typename ValueType, class ExecutionSpace, class Memory = void>
template<class , class , class >
friend class CrsProductTensor
friend

Definition at line 85 of file Stokhos_CrsProductTensor.hpp.

Member Data Documentation

template<typename ValueType, class ExecutionSpace, class Memory = void>
const size_type Stokhos::CrsProductTensor< ValueType, ExecutionSpace, Memory >::host_vectorsize = 2
static

Definition at line 67 of file Stokhos_CrsProductTensor.hpp.

template<typename ValueType, class ExecutionSpace, class Memory = void>
const bool Stokhos::CrsProductTensor< ValueType, ExecutionSpace, Memory >::use_intrinsics = false
static

Definition at line 68 of file Stokhos_CrsProductTensor.hpp.

template<typename ValueType, class ExecutionSpace, class Memory = void>
const size_type Stokhos::CrsProductTensor< ValueType, ExecutionSpace, Memory >::num_entry_align = 1
static

Definition at line 69 of file Stokhos_CrsProductTensor.hpp.

template<typename ValueType, class ExecutionSpace, class Memory = void>
const size_type Stokhos::CrsProductTensor< ValueType, ExecutionSpace, Memory >::cuda_vectorsize = 32
static

Definition at line 71 of file Stokhos_CrsProductTensor.hpp.

template<typename ValueType, class ExecutionSpace, class Memory = void>
const bool Stokhos::CrsProductTensor< ValueType, ExecutionSpace, Memory >::is_cuda
static
Initial value:
=
false

Definition at line 72 of file Stokhos_CrsProductTensor.hpp.

template<typename ValueType, class ExecutionSpace, class Memory = void>
const size_type Stokhos::CrsProductTensor< ValueType, ExecutionSpace, Memory >::vectorsize = is_cuda ? cuda_vectorsize : host_vectorsize
static

Definition at line 78 of file Stokhos_CrsProductTensor.hpp.

template<typename ValueType, class ExecutionSpace, class Memory = void>
const size_type Stokhos::CrsProductTensor< ValueType, ExecutionSpace, Memory >::tensor_align = vectorsize
static

Definition at line 81 of file Stokhos_CrsProductTensor.hpp.

template<typename ValueType, class ExecutionSpace, class Memory = void>
coord_array_type Stokhos::CrsProductTensor< ValueType, ExecutionSpace, Memory >::m_coord
private

Definition at line 94 of file Stokhos_CrsProductTensor.hpp.

template<typename ValueType, class ExecutionSpace, class Memory = void>
coord2_array_type Stokhos::CrsProductTensor< ValueType, ExecutionSpace, Memory >::m_coord2
private

Definition at line 95 of file Stokhos_CrsProductTensor.hpp.

template<typename ValueType, class ExecutionSpace, class Memory = void>
value_array_type Stokhos::CrsProductTensor< ValueType, ExecutionSpace, Memory >::m_value
private

Definition at line 96 of file Stokhos_CrsProductTensor.hpp.

template<typename ValueType, class ExecutionSpace, class Memory = void>
entry_array_type Stokhos::CrsProductTensor< ValueType, ExecutionSpace, Memory >::m_num_entry
private

Definition at line 97 of file Stokhos_CrsProductTensor.hpp.

template<typename ValueType, class ExecutionSpace, class Memory = void>
row_map_array_type Stokhos::CrsProductTensor< ValueType, ExecutionSpace, Memory >::m_row_map
private

Definition at line 98 of file Stokhos_CrsProductTensor.hpp.

template<typename ValueType, class ExecutionSpace, class Memory = void>
size_type Stokhos::CrsProductTensor< ValueType, ExecutionSpace, Memory >::m_dim
private

Definition at line 99 of file Stokhos_CrsProductTensor.hpp.

template<typename ValueType, class ExecutionSpace, class Memory = void>
size_type Stokhos::CrsProductTensor< ValueType, ExecutionSpace, Memory >::m_entry_max
private

Definition at line 100 of file Stokhos_CrsProductTensor.hpp.

template<typename ValueType, class ExecutionSpace, class Memory = void>
size_type Stokhos::CrsProductTensor< ValueType, ExecutionSpace, Memory >::m_nnz
private

Definition at line 101 of file Stokhos_CrsProductTensor.hpp.

template<typename ValueType, class ExecutionSpace, class Memory = void>
size_type Stokhos::CrsProductTensor< ValueType, ExecutionSpace, Memory >::m_flops
private

Definition at line 102 of file Stokhos_CrsProductTensor.hpp.

template<typename ValueType, class ExecutionSpace, class Memory = void>
size_type Stokhos::CrsProductTensor< ValueType, ExecutionSpace, Memory >::m_avg_entries_per_row
private

Definition at line 103 of file Stokhos_CrsProductTensor.hpp.


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