Stokhos Package Browser (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Stokhos_SymmetricDiagonalSpec.hpp
Go to the documentation of this file.
1 // @HEADER
2 // *****************************************************************************
3 // Stokhos Package
4 //
5 // Copyright 2009 NTESS and the Stokhos contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 // @HEADER
9 
10 #ifndef STOKHOS_SYMMETRIC_DIAGONAL_SPEC_HPP
11 #define STOKHOS_SYMMETRIC_DIAGONAL_SPEC_HPP
12 
13 #include "Kokkos_StaticCrsGraph.hpp"
14 
15 namespace Stokhos {
16 
30 template< class ExecutionSpace >
32 public:
33 
34  typedef unsigned size_type;
35 
37  KOKKOS_INLINE_FUNCTION
38  unsigned dimension() const { return m_dimension ; }
39 
41  KOKKOS_INLINE_FUNCTION
42  unsigned matrix_offset( const unsigned row , const unsigned column ) const
43  {
44  const int diag_count = 1 + ( m_dimension >> 1 );
45  const int diag = (int) column - (int) row ;
46 
47  unsigned offset = 0 ;
48 
49  if ( ( 0 <= diag && diag < diag_count ) || ( diag <= - diag_count ) ) {
50  offset = row + m_dimension * ( ( m_dimension + diag ) % m_dimension );
51  }
52  else {
53  offset = column + m_dimension * ( ( m_dimension - diag ) % m_dimension );
54  }
55 
56  return offset ;
57  }
58 
60  KOKKOS_INLINE_FUNCTION
61  unsigned matrix_size() const
62  { return ( m_dimension * ( m_dimension + 1 ) ) >> 1 ; }
63 
65  : m_dimension( 0 ) {}
66 
68  : m_dimension( rhs.m_dimension ) {}
69 
70  SymmetricDiagonalSpec & operator =
71  ( const SymmetricDiagonalSpec & rhs )
72  { m_dimension = rhs.m_dimension ; return *this ; }
73 
74  explicit
75  SymmetricDiagonalSpec( const unsigned dim )
76  : m_dimension( dim ) {}
77 
78 private:
79  unsigned m_dimension ;
80 };
81 
82 template < typename Device >
84 public:
85  typedef Device execution_space ;
86  typedef typename execution_space::size_type size_type ;
88 
89  template< typename MatrixValue , typename VectorValue >
90  KOKKOS_INLINE_FUNCTION
91  static void apply( const block_type & block ,
92  const MatrixValue * a ,
93  const VectorValue * const x ,
94  VectorValue * const y )
95  {
96  const size_type dimension = block.dimension();
97  const size_type dim_half = ( dimension + 1 ) >> 1 ;
98 
99  // Multiply the main diagonal (first diagonal)
100  for ( size_type j = 0 ; j < dimension ; ++j ) {
101  y[j] += a[j] * x[j] ; // Contiguous access
102  }
103 
104  // Multiply remaining full diagionals, each diagonal is accessed twice
105  for ( size_type d = 1 ; d < dim_half ; ++d ) {
106  size_type kx = d ;
107  size_type kxr = dimension - d ;
108 
109  a += dimension ; // next diagonal
110 
111  for ( size_type j = 0 ; j < dimension ; ++j ) {
112  y[j] += a[j] * x[kx] + a[kxr] * x[kxr]; // Contiguous access
113  if ( dimension == ++kx ) kx = 0 ;
114  if ( dimension == ++kxr ) kxr = 0 ;
115  }
116  }
117 
118  // If even number of diagonals then the last diagonal is half-length
119  if ( ! ( dimension & 01 ) ) {
120  size_type kx = dim_half ;
121 
122  a += dimension ; // next diagonal
123 
124  for ( size_type j = 0 ; j < dim_half ; ++j , ++kx ) {
125  y[j] += a[j] * x[kx] ; // Contiguous access
126  y[kx] += a[j] * x[j] ; // Contiguous access
127  }
128  }
129  }
130 
131  KOKKOS_INLINE_FUNCTION
132  static size_type matrix_size( const block_type & block )
133  { return block.matrix_size(); }
134 };
135 
136 } // namespace Stokhos
137 
138 #endif /* #ifndef STOKHOS_SYMMETRIC_DIAGONAL_SPEC_HPP */
static KOKKOS_INLINE_FUNCTION size_type matrix_size(const block_type &block)
static KOKKOS_INLINE_FUNCTION void apply(const block_type &block, const MatrixValue *a, const VectorValue *const x, VectorValue *const y)
KOKKOS_INLINE_FUNCTION unsigned matrix_size() const
Storage size for block coefficients.
KOKKOS_INLINE_FUNCTION unsigned matrix_offset(const unsigned row, const unsigned column) const
Storage location for the (row,column) entry.
SymmetricDiagonalSpec(const SymmetricDiagonalSpec &rhs)
Symmetric diagonal storage for a dense matrix.
KOKKOS_INLINE_FUNCTION unsigned dimension() const
Dimension of vector block.