Tpetra parallel linear algebra  Version of the Day
Public Member Functions | List of all members
Tpetra::Experimental::LittleBlock< Scalar, LO > Class Template Reference

Nonowning view of a square dense block in a block matrix. More...

#include <Tpetra_Experimental_BlockView.hpp>

Public Member Functions

 LittleBlock (Scalar *const A, const LO blockSize, const LO strideX, const LO strideY)
 Constructor. More...
 
template<class T >
 LittleBlock (T *const A, const LO blockSize, const LO strideX, const LO strideY, typename std::enable_if< !std::is_same< Scalar, T >::value &&std::is_convertible< Scalar, T >::value &&sizeof(Scalar)==sizeof(T), int * >::type ignoreMe=NULL)
 Constructor that takes an impl_scalar_type pointer. More...
 
LO getBlockSize () const
 The block size (number of rows, and number of columns). More...
 
Scalar * getRawPtr () const
 Pointer to the block's entries, as Scalar*. More...
 
impl_scalar_type & operator() (const LO i, const LO j) const
 Reference to entry (i,j) of the block. More...
 
template<class LittleBlockType >
void update (const Scalar &alpha, const LittleBlockType &X) const
 *this := *this + alpha * X. More...
 
template<class LittleBlockType >
void assign (const LittleBlockType &X) const
 *this := X. More...
 
void scale (const Scalar &alpha) const
 (*this)(i,j) := alpha * (*this)(i,j) for all (i,j). More...
 
void fill (const Scalar &alpha) const
 (*this)(i,j) := alpha for all (i,j). More...
 
template<class LittleBlockType >
void absmax (const LittleBlockType &X) const
 (*this)(i,j) := max(abs((*this)(i,j)), abs(X(i,j))) for all (i,j). More...
 

Detailed Description

template<class Scalar, class LO>
class Tpetra::Experimental::LittleBlock< Scalar, LO >

Nonowning view of a square dense block in a block matrix.

Template Parameters
ScalarThe type of entries in the block.
LOThe type of local indices. See the documentation of the first template parameter of Map for requirements.

"Little" means local (not distributed over multiple MPI processes; stored to maximize locality) and small (think 3x3, not 1000x1000).

The Scalar template parameter may be const or nonconst. This is one reason why instance methods below that take a LittleBlock accept it as a template parameter: that lets you add a const LittleBlock (e.g., LittleBlock<const double, int>) to a nonconst LittleBlock (e.g., LittleBlock<double, int>).

Definition at line 131 of file Tpetra_Experimental_BlockView.hpp.

Constructor & Destructor Documentation

template<class Scalar , class LO >
Tpetra::Experimental::LittleBlock< Scalar, LO >::LittleBlock ( Scalar *const  A,
const LO  blockSize,
const LO  strideX,
const LO  strideY 
)
inline

Constructor.

Parameters
A[in] Pointer to the block's entries
blockSize[in] Dimension of the block (all blocks are square)
strideX[in] Stride between consecutive entries in a column
strideY[in] Stride between consecutive entries in a row

Definition at line 145 of file Tpetra_Experimental_BlockView.hpp.

template<class Scalar , class LO >
template<class T >
Tpetra::Experimental::LittleBlock< Scalar, LO >::LittleBlock ( T *const  A,
const LO  blockSize,
const LO  strideX,
const LO  strideY,
typename std::enable_if< !std::is_same< Scalar, T >::value &&std::is_convertible< Scalar, T >::value &&sizeof(Scalar)==sizeof(T), int * >::type  ignoreMe = NULL 
)
inline

Constructor that takes an impl_scalar_type pointer.

Parameters
A[in] Pointer to the block's entries, as impl_scalar_type* rather than Scalar*
blockSize[in] Dimension of the block (all blocks are square)
strideX[in] Stride between consecutive entries in a column
strideY[in] Stride between consecutive entries in a row

While this constructor is templated on a type T, the intent is that T == impl_scalar_type. (We must template on T rather than using impl_scalar_type directly, because of how std::enable_if works.) The long, complicated std::enable_if expression ensures that this constructor only exists if Scalar differs from impl_scalar_type, but the two types are mutually compatible and have the same size. (They must be bitwise compatible, so that reinterpret_cast makes sense between them.)

Definition at line 173 of file Tpetra_Experimental_BlockView.hpp.

Member Function Documentation

template<class Scalar , class LO >
LO Tpetra::Experimental::LittleBlock< Scalar, LO >::getBlockSize ( ) const
inline

The block size (number of rows, and number of columns).

Definition at line 189 of file Tpetra_Experimental_BlockView.hpp.

template<class Scalar , class LO >
Scalar* Tpetra::Experimental::LittleBlock< Scalar, LO >::getRawPtr ( ) const
inline

Pointer to the block's entries, as Scalar*.

Definition at line 194 of file Tpetra_Experimental_BlockView.hpp.

template<class Scalar , class LO >
impl_scalar_type& Tpetra::Experimental::LittleBlock< Scalar, LO >::operator() ( const LO  i,
const LO  j 
) const
inline

Reference to entry (i,j) of the block.

Note
To Tpetra developers: This is returned as impl_scalar_type and not as Scalar, in order to avoid a lot of reinterpret_cast calls in the inner loop of the sparse matrix-vector multiply kernel of Tpetra::Experimental::BlockCrsMatrix. Any pair of types impl_scalar_type, Scalar used here should always be convertible in either direction, so the return type should not pose any issues in practice.

Definition at line 208 of file Tpetra_Experimental_BlockView.hpp.

template<class Scalar , class LO >
template<class LittleBlockType >
void Tpetra::Experimental::LittleBlock< Scalar, LO >::update ( const Scalar &  alpha,
const LittleBlockType &  X 
) const
inline

*this := *this + alpha * X.

Definition at line 214 of file Tpetra_Experimental_BlockView.hpp.

template<class Scalar , class LO >
template<class LittleBlockType >
void Tpetra::Experimental::LittleBlock< Scalar, LO >::assign ( const LittleBlockType &  X) const
inline

*this := X.

Definition at line 225 of file Tpetra_Experimental_BlockView.hpp.

template<class Scalar , class LO >
void Tpetra::Experimental::LittleBlock< Scalar, LO >::scale ( const Scalar &  alpha) const
inline

(*this)(i,j) := alpha * (*this)(i,j) for all (i,j).

Definition at line 234 of file Tpetra_Experimental_BlockView.hpp.

template<class Scalar , class LO >
void Tpetra::Experimental::LittleBlock< Scalar, LO >::fill ( const Scalar &  alpha) const
inline

(*this)(i,j) := alpha for all (i,j).

Definition at line 244 of file Tpetra_Experimental_BlockView.hpp.

template<class Scalar , class LO >
template<class LittleBlockType >
void Tpetra::Experimental::LittleBlock< Scalar, LO >::absmax ( const LittleBlockType &  X) const
inline

(*this)(i,j) := max(abs((*this)(i,j)), abs(X(i,j))) for all (i,j).

Tpetra uses this operation to implement the ABSMAX CombineMode.

Definition at line 258 of file Tpetra_Experimental_BlockView.hpp.


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