Belos  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
List of all members
Belos::MultiVecTraits< ScalarType, MV > Class Template Reference

Traits class which defines basic operations on multivectors. More...

#include <BelosMultiVecTraits.hpp>

Static Public Member Functions

Creation methods
static Teuchos::RCP< MV > Clone (const MV &mv, const int numvecs)
 Creates a new empty MV containing numvecs columns. More...
 
static Teuchos::RCP< MV > CloneCopy (const MV &mv)
 Creates a new MV and copies contents of mv into the new vector (deep copy). More...
 
static Teuchos::RCP< MV > CloneCopy (const MV &mv, const std::vector< int > &index)
 Creates a new MV and copies the selected contents of mv into the new vector (deep copy). More...
 
static Teuchos::RCP< MV > CloneCopy (const MV &mv, const Teuchos::Range1D &index)
 Deep copy of specified columns of mv. More...
 
static Teuchos::RCP< MV > CloneViewNonConst (MV &mv, const std::vector< int > &index)
 Creates a new MV that shares the selected contents of mv (shallow copy). More...
 
static Teuchos::RCP< MV > CloneViewNonConst (MV &mv, const Teuchos::Range1D &index)
 Non-const view of specified columns of mv. More...
 
static Teuchos::RCP< const MV > CloneView (const MV &mv, const std::vector< int > &index)
 Creates a new const MV that shares the selected contents of mv (shallow copy). More...
 
static Teuchos::RCP< MV > CloneView (MV &mv, const Teuchos::Range1D &index)
 Const view of specified columns of mv. More...
 
Attribute methods
static ptrdiff_t GetGlobalLength (const MV &mv)
 Return the number of rows in the given multivector mv. More...
 
static int GetNumberVecs (const MV &mv)
 Obtain the number of vectors in mv. More...
 
static bool HasConstantStride (const MV &mv)
 Whether the given multivector mv has constant stride. More...
 
Update methods
static void MvTimesMatAddMv (const ScalarType alpha, const MV &A, const Teuchos::SerialDenseMatrix< int, ScalarType > &B, const ScalarType beta, MV &mv)
 Update mv with $ \alpha AB + \beta mv $. More...
 
static void MvAddMv (const ScalarType alpha, const MV &A, const ScalarType beta, const MV &B, MV &mv)
 Replace mv with $\alpha A + \beta B$. More...
 
static void MvScale (MV &mv, const ScalarType alpha)
 Scale each element of the vectors in mv with alpha. More...
 
static void MvScale (MV &mv, const std::vector< ScalarType > &alpha)
 Scale each element of the i-th vector in mv with alpha[i]. More...
 
static void MvTransMv (const ScalarType alpha, const MV &A, const MV &mv, Teuchos::SerialDenseMatrix< int, ScalarType > &B)
 Compute a dense matrix B through the matrix-matrix multiply $ \alpha A^Hmv $. More...
 
static void MvDot (const MV &mv, const MV &A, std::vector< ScalarType > &b)
 Compute a vector b where the components are the individual dot-products of the i-th columns of A and mv, i.e. $b[i] = A[i]^Hmv[i]$. More...
 
Norm method
static void MvNorm (const MV &mv, std::vector< typename Teuchos::ScalarTraits< ScalarType >::magnitudeType > &normvec, NormType type=TwoNorm)
 Compute the norm of each individual vector of mv. Upon return, normvec[i] holds the value of $||mv_i||$, the i-th column of mv. More...
 
Initialization methods
static void SetBlock (const MV &A, const std::vector< int > &index, MV &mv)
 Copy the vectors in A to a set of vectors in mv indicated by the indices given in index. More...
 
static void SetBlock (const MV &A, const Teuchos::Range1D &index, MV &mv)
 Deep copy of A into specified columns of mv. More...
 
static void Assign (const MV &A, MV &mv)
 mv := A More...
 
static void MvRandom (MV &mv)
 Replace the vectors in mv with random vectors. More...
 
static void MvInit (MV &mv, const ScalarType alpha=Teuchos::ScalarTraits< ScalarType >::zero())
 Replace each element of the vectors in mv with alpha. More...
 
Print method
static void MvPrint (const MV &mv, std::ostream &os)
 Print the mv multi-vector to the os output stream. More...
 

Detailed Description

template<class ScalarType, class MV>
class Belos::MultiVecTraits< ScalarType, MV >

Traits class which defines basic operations on multivectors.

Template Parameters
ScalarTypeThe type of the entries in the multivectors.
MVThe type of the multivectors themselves.

This traits class tells Belos' solvers how to perform multivector operations for the multivector type MV. These operations include creating copies or views, finding the number of rows or columns (i.e., vectors) in a given multivector, and computing inner products, norms, and vector sums. (Belos' solvers use the OperatorTraits traits class to apply operators to multivectors.)

Belos gives users two different ways to tell its solvers how to compute with multivectors of a given type MV. The first and preferred way is for users to specialize MultiVecTraits, this traits class, for their given MV type. Belos provides specializations for MV = Epetra_MultiVector and Tpetra::MultiVector, and Stratimikos provides a specialization for Thyra::MultiVectorBase. The second way is for users to make their multivector type (or a wrapper thereof) inherit from MultiVec. This works because Belos provides a specialization of MultiVecTraits for MultiVec. Specializing MultiVecTraits is more flexible because it does not require a multivector type to inherit from MultiVec; this is possible even if you do not have control over the interface of a class.

If you have a different multivector type MV that you would like to use with Belos, and if that type does not inherit from MultiVec, then you must implement a specialization of MultiVecTraits for MV. Otherwise, this traits class will report a compile-time error (relating to UndefinedMultiVecTraits). Specializing MultiVecTraits for your MV type is not hard. Just look at the examples for Epetra_MultiVector (in belos/epetra/src/BelosEpetraAdapter.hpp) and Tpetra::MultiVector (in belos/tpetra/src/BelosTpetraAdapter.hpp).

Note
You do not need to write a specialization of MultiVecTraits if you are using Epetra, Tpetra, or Thyra multivectors. Belos already provides specializations for the Epetra and Tpetra multivector types, and Stratimikos provides a specialization for the Thyra type. Just relax and enjoy using the solvers!
Examples:
BlockCG/BlockCGEpetraExFile.cpp, BlockCG/BlockPrecCGEpetraExFile.cpp, BlockCG/PseudoBlockCGEpetraExFile.cpp, BlockCG/PseudoBlockPrecCGEpetraExFile.cpp, BlockGmres/BlockFlexGmresEpetraExFile.cpp, BlockGmres/BlockGmresEpetraExFile.cpp, BlockGmres/BlockGmresPolyEpetraExFile.cpp, BlockGmres/BlockPrecGmresEpetraExFile.cpp, BlockGmres/PseudoBlockGmresEpetraExFile.cpp, BlockGmres/PseudoBlockPrecGmresEpetraExFile.cpp, GCRODR/GCRODREpetraExFile.cpp, GCRODR/PrecGCRODREpetraExFile.cpp, PCPG/PCPGEpetraExFile.cpp, TFQMR/PseudoBlockTFQMREpetraExFile.cpp, and TFQMR/TFQMREpetraExFile.cpp.

Definition at line 129 of file BelosMultiVecTraits.hpp.

Member Function Documentation

template<class ScalarType , class MV >
static Teuchos::RCP<MV> Belos::MultiVecTraits< ScalarType, MV >::Clone ( const MV &  mv,
const int  numvecs 
)
inlinestatic

Creates a new empty MV containing numvecs columns.

Returns
Reference-counted pointer to the new multivector of type MV.

Definition at line 138 of file BelosMultiVecTraits.hpp.

template<class ScalarType , class MV >
static Teuchos::RCP<MV> Belos::MultiVecTraits< ScalarType, MV >::CloneCopy ( const MV &  mv)
inlinestatic

Creates a new MV and copies contents of mv into the new vector (deep copy).

Returns
Reference-counted pointer to the new multivector of type MV.

Definition at line 145 of file BelosMultiVecTraits.hpp.

template<class ScalarType , class MV >
static Teuchos::RCP<MV> Belos::MultiVecTraits< ScalarType, MV >::CloneCopy ( const MV &  mv,
const std::vector< int > &  index 
)
inlinestatic

Creates a new MV and copies the selected contents of mv into the new vector (deep copy).

The copied vectors from mv are indicated by the index.size() indices in index.

Returns
Reference-counted pointer to the new multivector of type MV.

Definition at line 153 of file BelosMultiVecTraits.hpp.

template<class ScalarType , class MV >
static Teuchos::RCP<MV> Belos::MultiVecTraits< ScalarType, MV >::CloneCopy ( const MV &  mv,
const Teuchos::Range1D index 
)
inlinestatic

Deep copy of specified columns of mv.

Create a new MV, and copy (deep copy) the columns of mv specified by the given inclusive index range into the new multivector.

Parameters
mv[in] Multivector to copy
index[in] Inclusive index range of columns of mv
Returns
Reference-counted pointer to the new multivector of type MV.

Definition at line 165 of file BelosMultiVecTraits.hpp.

template<class ScalarType , class MV >
static Teuchos::RCP<MV> Belos::MultiVecTraits< ScalarType, MV >::CloneViewNonConst ( MV &  mv,
const std::vector< int > &  index 
)
inlinestatic

Creates a new MV that shares the selected contents of mv (shallow copy).

The index of the numvecs vectors shallow copied from mv are indicated by the indices given in index.

Returns
Reference-counted pointer to the new multivector of type MV.

Definition at line 173 of file BelosMultiVecTraits.hpp.

template<class ScalarType , class MV >
static Teuchos::RCP<MV> Belos::MultiVecTraits< ScalarType, MV >::CloneViewNonConst ( MV &  mv,
const Teuchos::Range1D index 
)
inlinestatic

Non-const view of specified columns of mv.

Return a non-const view of the columns of mv specified by the given inclusive index range.

Parameters
mv[in] Multivector to view (shallow non-const copy)
index[in] Inclusive index range of columns of mv
Returns
Reference-counted pointer to the non-const view of specified columns of mv

Definition at line 184 of file BelosMultiVecTraits.hpp.

template<class ScalarType , class MV >
static Teuchos::RCP<const MV> Belos::MultiVecTraits< ScalarType, MV >::CloneView ( const MV &  mv,
const std::vector< int > &  index 
)
inlinestatic

Creates a new const MV that shares the selected contents of mv (shallow copy).

The index of the numvecs vectors shallow copied from mv are indicated by the indices given in index.

Returns
Reference-counted pointer to the new const multivector of type MV.

Definition at line 192 of file BelosMultiVecTraits.hpp.

template<class ScalarType , class MV >
static Teuchos::RCP<MV> Belos::MultiVecTraits< ScalarType, MV >::CloneView ( MV &  mv,
const Teuchos::Range1D index 
)
inlinestatic

Const view of specified columns of mv.

Return a const view of the columns of mv specified by the given inclusive index range.

Parameters
mv[in] Multivector to view (shallow const copy)
index[in] Inclusive index range of columns of mv
Returns
Reference-counted pointer to the const view of specified columns of mv

Definition at line 203 of file BelosMultiVecTraits.hpp.

template<class ScalarType , class MV >
static ptrdiff_t Belos::MultiVecTraits< ScalarType, MV >::GetGlobalLength ( const MV &  mv)
inlinestatic

Return the number of rows in the given multivector mv.

Obtain the vector length of mv.

Definition at line 214 of file BelosMultiVecTraits.hpp.

template<class ScalarType , class MV >
static int Belos::MultiVecTraits< ScalarType, MV >::GetNumberVecs ( const MV &  mv)
inlinestatic

Obtain the number of vectors in mv.

Definition at line 218 of file BelosMultiVecTraits.hpp.

template<class ScalarType , class MV >
static bool Belos::MultiVecTraits< ScalarType, MV >::HasConstantStride ( const MV &  mv)
inlinestatic

Whether the given multivector mv has constant stride.

Parameters
mv[in] Multivector to check

Knowing whether mv has constant stride is useful for certain orthogonalization methods, for example.

Note
(mfh 13 Jan 2011) This is really a hack for TSQR, which currently can only process multivectors with constant stride. Fixing this can be done in a few different ways:
  • Copy the entire multivector into a constant-stride multivector (performance penalty for copying in and out, storage penalty of one multivector, but easy to implement)
  • Copying "cache blocks" in and out of constant-stride storage (some performance penalty, storage penalty << one multivector, requires a new variant of TSQR::SequentialTsqr and perhaps a restructuring of the TSQR has-a hierarchy)
  • Writing an implementation / back-end of TSQR::Combine that can handle cache blocks not of constant stride. Like the above, a bit more work, but without the copy in/out performance penalty for cache blocks, and doesn't require an extra cache block's worth of storage.

Definition at line 243 of file BelosMultiVecTraits.hpp.

template<class ScalarType , class MV >
static void Belos::MultiVecTraits< ScalarType, MV >::MvTimesMatAddMv ( const ScalarType  alpha,
const MV &  A,
const Teuchos::SerialDenseMatrix< int, ScalarType > &  B,
const ScalarType  beta,
MV &  mv 
)
inlinestatic

Update mv with $ \alpha AB + \beta mv $.

Definition at line 253 of file BelosMultiVecTraits.hpp.

template<class ScalarType , class MV >
static void Belos::MultiVecTraits< ScalarType, MV >::MvAddMv ( const ScalarType  alpha,
const MV &  A,
const ScalarType  beta,
const MV &  B,
MV &  mv 
)
inlinestatic

Replace mv with $\alpha A + \beta B$.

Definition at line 260 of file BelosMultiVecTraits.hpp.

template<class ScalarType , class MV >
static void Belos::MultiVecTraits< ScalarType, MV >::MvScale ( MV &  mv,
const ScalarType  alpha 
)
inlinestatic

Scale each element of the vectors in mv with alpha.

Definition at line 265 of file BelosMultiVecTraits.hpp.

template<class ScalarType , class MV >
static void Belos::MultiVecTraits< ScalarType, MV >::MvScale ( MV &  mv,
const std::vector< ScalarType > &  alpha 
)
inlinestatic

Scale each element of the i-th vector in mv with alpha[i].

Definition at line 270 of file BelosMultiVecTraits.hpp.

template<class ScalarType , class MV >
static void Belos::MultiVecTraits< ScalarType, MV >::MvTransMv ( const ScalarType  alpha,
const MV &  A,
const MV &  mv,
Teuchos::SerialDenseMatrix< int, ScalarType > &  B 
)
inlinestatic

Compute a dense matrix B through the matrix-matrix multiply $ \alpha A^Hmv $.

Definition at line 275 of file BelosMultiVecTraits.hpp.

template<class ScalarType , class MV >
static void Belos::MultiVecTraits< ScalarType, MV >::MvDot ( const MV &  mv,
const MV &  A,
std::vector< ScalarType > &  b 
)
inlinestatic

Compute a vector b where the components are the individual dot-products of the i-th columns of A and mv, i.e. $b[i] = A[i]^Hmv[i]$.

Definition at line 280 of file BelosMultiVecTraits.hpp.

template<class ScalarType , class MV >
static void Belos::MultiVecTraits< ScalarType, MV >::MvNorm ( const MV &  mv,
std::vector< typename Teuchos::ScalarTraits< ScalarType >::magnitudeType > &  normvec,
NormType  type = TwoNorm 
)
inlinestatic

Compute the norm of each individual vector of mv. Upon return, normvec[i] holds the value of $||mv_i||$, the i-th column of mv.

Parameters
mv,:multi-vector
normvec,:std::vector containing the norm for each vector in mv
NormType,:norm type (default: two-norm)

Definition at line 294 of file BelosMultiVecTraits.hpp.

template<class ScalarType , class MV >
static void Belos::MultiVecTraits< ScalarType, MV >::SetBlock ( const MV &  A,
const std::vector< int > &  index,
MV &  mv 
)
inlinestatic

Copy the vectors in A to a set of vectors in mv indicated by the indices given in index.

The numvecs vectors in A are copied to a subset of vectors in mv indicated by the indices given in index, i.e. mv[index[i]] = A[i].

Definition at line 306 of file BelosMultiVecTraits.hpp.

template<class ScalarType , class MV >
static void Belos::MultiVecTraits< ScalarType, MV >::SetBlock ( const MV &  A,
const Teuchos::Range1D index,
MV &  mv 
)
inlinestatic

Deep copy of A into specified columns of mv.

(Deeply) copy the first index.size() columns of A into the columns of mv specified by the given index range.

Postcondition: mv[i] = A[i - index.lbound()] for all i in [index.lbound(), index.ubound()]

Parameters
A[in] Source multivector
index[in] Inclusive index range of columns of mv; index set of the target
mv[out] Target multivector

Definition at line 321 of file BelosMultiVecTraits.hpp.

template<class ScalarType , class MV >
static void Belos::MultiVecTraits< ScalarType, MV >::Assign ( const MV &  A,
MV &  mv 
)
inlinestatic

mv := A

Assign (deep copy) A into mv.

Definition at line 327 of file BelosMultiVecTraits.hpp.

template<class ScalarType , class MV >
static void Belos::MultiVecTraits< ScalarType, MV >::MvRandom ( MV &  mv)
inlinestatic

Replace the vectors in mv with random vectors.

Definition at line 332 of file BelosMultiVecTraits.hpp.

template<class ScalarType , class MV >
static void Belos::MultiVecTraits< ScalarType, MV >::MvInit ( MV &  mv,
const ScalarType  alpha = Teuchos::ScalarTraits<ScalarType>::zero() 
)
inlinestatic

Replace each element of the vectors in mv with alpha.

Definition at line 337 of file BelosMultiVecTraits.hpp.

template<class ScalarType , class MV >
static void Belos::MultiVecTraits< ScalarType, MV >::MvPrint ( const MV &  mv,
std::ostream &  os 
)
inlinestatic

Print the mv multi-vector to the os output stream.

Definition at line 347 of file BelosMultiVecTraits.hpp.


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

Generated on Thu Apr 18 2024 09:24:47 for Belos by doxygen 1.8.5