Tpetra Matrix/Vector Services  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Defines
Public Member Functions
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>

List of all members.

Public Member Functions

 LittleBlock (Scalar *const A, const LO blockSize, const LO strideX, const LO strideY)
 Constructor.
LO getBlockSize () const
 The block size (number of rows, and number of columns).
Scalar * getRawPtr () const
 Pointer to the block's entries.
Scalar & operator() (const LO i, const LO j) const
 Reference to entry (i,j) of the block.
template<class LittleBlockType >
void update (const Scalar &alpha, const LittleBlockType &X) const
 *this := *this + alpha * X.
template<class LittleBlockType >
void assign (const LittleBlockType &X) const
 *this := X.
void scale (const Scalar &alpha) const
 (*this)(i,j) := alpha * (*this)(i,j) for all (i,j).
void fill (const Scalar &alpha) const
 (*this)(i,j) := alpha for all (i,j).
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).

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 85 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 95 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 104 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.

Definition at line 109 of file Tpetra_Experimental_BlockView.hpp.

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

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

Definition at line 114 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 120 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 130 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 139 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 148 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 161 of file Tpetra_Experimental_BlockView.hpp.


The documentation for this class was generated from the following file:
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Defines