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

Interface for multivectors used by Belos' linear solvers. More...

#include <BelosMultiVec.hpp>

Inheritance diagram for Belos::MultiVec< ScalarType >:
Inheritance graph
[legend]

Public Member Functions

Constructor/Destructor
 MultiVec ()
 Default constructor. More...
 
virtual ~MultiVec ()
 Destructor (virtual for memory safety of derived classes). More...
 
Creation methods for new multivectors
virtual MultiVec< ScalarType > * Clone (const int numvecs) const =0
 Create a new MultiVec with numvecs columns. More...
 
virtual MultiVec< ScalarType > * CloneCopy () const =0
 Create a new MultiVec and copy contents of *this into it (deep copy). More...
 
virtual MultiVec< ScalarType > * CloneCopy (const std::vector< int > &index) const =0
 Creates a new Belos::MultiVec and copies the selected contents of *this into the new multivector (deep copy). The copied vectors from *this are indicated by the index.size() indices in index. More...
 
virtual MultiVec< ScalarType > * CloneViewNonConst (const std::vector< int > &index)=0
 Creates a new Belos::MultiVec that shares the selected contents of *this. The index of the numvecs vectors copied from *this are indicated by the indices given in index. More...
 
virtual const MultiVec
< ScalarType > * 
CloneView (const std::vector< int > &index) const =0
 Creates a new Belos::MultiVec that shares the selected contents of *this. The index of the numvecs vectors copied from *this are indicated by the indices given in index. More...
 
Dimension information methods
virtual ptrdiff_t GetGlobalLength () const =0
 The number of rows in the multivector. More...
 
virtual int GetNumberVecs () const =0
 The number of vectors (i.e., columns) in the multivector. More...
 
Update methods
virtual void MvTimesMatAddMv (const ScalarType alpha, const MultiVec< ScalarType > &A, const Teuchos::SerialDenseMatrix< int, ScalarType > &B, const ScalarType beta)=0
 Update *this with alpha * A * B + beta * (*this). More...
 
virtual void MvAddMv (const ScalarType alpha, const MultiVec< ScalarType > &A, const ScalarType beta, const MultiVec< ScalarType > &B)=0
 Replace *this with alpha * A + beta * B. More...
 
virtual void MvScale (const ScalarType alpha)=0
 Scale each element of the vectors in *this with alpha. More...
 
virtual void MvScale (const std::vector< ScalarType > &alpha)=0
 Scale each element of the i-th vector in *this with alpha[i]. More...
 
virtual void MvTransMv (const ScalarType alpha, const MultiVec< ScalarType > &A, Teuchos::SerialDenseMatrix< int, ScalarType > &B) const =0
 Compute a dense matrix B through the matrix-matrix multiply alpha * A^T * (*this). More...
 
virtual void MvDot (const MultiVec< ScalarType > &A, std::vector< ScalarType > &b) const =0
 Compute the dot product of each column of *this with the corresponding column of A. More...
 
Norm method
virtual void MvNorm (std::vector< typename Teuchos::ScalarTraits< ScalarType >::magnitudeType > &normvec, NormType type=TwoNorm) const =0
 Compute the norm of each vector in *this. More...
 
Initialization methods
virtual void SetBlock (const MultiVec< ScalarType > &A, const std::vector< int > &index)=0
 Copy the vectors in A to a set of vectors in *this. More...
 
virtual void MvRandom ()=0
 Fill all the vectors in *this with random numbers. More...
 
virtual void MvInit (const ScalarType alpha)=0
 Replace each element of the vectors in *this with alpha. More...
 
Print method
virtual void MvPrint (std::ostream &os) const =0
 Print *this multivector to the os output stream. More...
 

Detailed Description

template<class ScalarType>
class Belos::MultiVec< ScalarType >

Interface for multivectors used by Belos' linear solvers.

Author
Michael Heroux, Rich Lehoucq, and Heidi Thornquist
Template Parameters
ScalarTypeThe type of entries of the multivector.

Belos accesses multivectors through a traits interface called MultiVecTraits. If you want to use Belos with your own multivector class MV, you may either specialize MultiVecTraits for MV, or you may wrap MV in your own class that implements MultiVec. Specializing MultiVecTraits works via compile-time polymorphism, whereas implementing the MultiVec interface works via run-time polymorphism. You may pick whichever option you like. However, specializing MultiVecTraits is the preferred method. This is because Belos' linear solvers always use a specialization of MultiVecTraits to access multivector operations. They only use MultiVec through a specialization of the MultiVecTraits traits class, which is implemented below in this header file.

If you want your multivector class (or a wrapper thereof) to implement the MultiVec interface, you should inherit from MultiVec<ScalarType>, where ScalarType is the type of entries in the multivector. For example, a multivector with entries of type double would inherit from MultiVec<double>.

Definition at line 91 of file BelosMultiVec.hpp.

Constructor & Destructor Documentation

template<class ScalarType>
Belos::MultiVec< ScalarType >::MultiVec ( )
inline

Default constructor.

Definition at line 96 of file BelosMultiVec.hpp.

template<class ScalarType>
virtual Belos::MultiVec< ScalarType >::~MultiVec ( )
inlinevirtual

Destructor (virtual for memory safety of derived classes).

Definition at line 99 of file BelosMultiVec.hpp.

Member Function Documentation

template<class ScalarType>
virtual MultiVec<ScalarType>* Belos::MultiVec< ScalarType >::Clone ( const int  numvecs) const
pure virtual

Create a new MultiVec with numvecs columns.

Returns
Pointer to the new multivector with uninitialized values.

Implemented in Belos::GmresPolyMv< ScalarType, MV >.

template<class ScalarType>
virtual MultiVec<ScalarType>* Belos::MultiVec< ScalarType >::CloneCopy ( ) const
pure virtual

Create a new MultiVec and copy contents of *this into it (deep copy).

Returns
Pointer to the new multivector

Implemented in Belos::GmresPolyMv< ScalarType, MV >.

template<class ScalarType>
virtual MultiVec<ScalarType>* Belos::MultiVec< ScalarType >::CloneCopy ( const std::vector< int > &  index) const
pure virtual

Creates a new Belos::MultiVec and copies the selected contents of *this into the new multivector (deep copy). The copied vectors from *this are indicated by the index.size() indices in index.

Returns
Pointer to the new multivector

Implemented in Belos::GmresPolyMv< ScalarType, MV >.

template<class ScalarType>
virtual MultiVec<ScalarType>* Belos::MultiVec< ScalarType >::CloneViewNonConst ( const std::vector< int > &  index)
pure virtual

Creates a new Belos::MultiVec that shares the selected contents of *this. The index of the numvecs vectors copied from *this are indicated by the indices given in index.

Returns
Pointer to the new multivector

Implemented in Belos::GmresPolyMv< ScalarType, MV >.

template<class ScalarType>
virtual const MultiVec<ScalarType>* Belos::MultiVec< ScalarType >::CloneView ( const std::vector< int > &  index) const
pure virtual

Creates a new Belos::MultiVec that shares the selected contents of *this. The index of the numvecs vectors copied from *this are indicated by the indices given in index.

Returns
Pointer to the new multivector

Implemented in Belos::GmresPolyMv< ScalarType, MV >.

template<class ScalarType>
virtual ptrdiff_t Belos::MultiVec< ScalarType >::GetGlobalLength ( ) const
pure virtual

The number of rows in the multivector.

Implemented in Belos::GmresPolyMv< ScalarType, MV >.

template<class ScalarType>
virtual int Belos::MultiVec< ScalarType >::GetNumberVecs ( ) const
pure virtual

The number of vectors (i.e., columns) in the multivector.

Implemented in Belos::GmresPolyMv< ScalarType, MV >.

template<class ScalarType>
virtual void Belos::MultiVec< ScalarType >::MvTimesMatAddMv ( const ScalarType  alpha,
const MultiVec< ScalarType > &  A,
const Teuchos::SerialDenseMatrix< int, ScalarType > &  B,
const ScalarType  beta 
)
pure virtual

Update *this with alpha * A * B + beta * (*this).

Implemented in Belos::GmresPolyMv< ScalarType, MV >.

template<class ScalarType>
virtual void Belos::MultiVec< ScalarType >::MvAddMv ( const ScalarType  alpha,
const MultiVec< ScalarType > &  A,
const ScalarType  beta,
const MultiVec< ScalarType > &  B 
)
pure virtual

Replace *this with alpha * A + beta * B.

Implemented in Belos::GmresPolyMv< ScalarType, MV >.

template<class ScalarType>
virtual void Belos::MultiVec< ScalarType >::MvScale ( const ScalarType  alpha)
pure virtual

Scale each element of the vectors in *this with alpha.

Implemented in Belos::GmresPolyMv< ScalarType, MV >.

template<class ScalarType>
virtual void Belos::MultiVec< ScalarType >::MvScale ( const std::vector< ScalarType > &  alpha)
pure virtual

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

Implemented in Belos::GmresPolyMv< ScalarType, MV >.

template<class ScalarType>
virtual void Belos::MultiVec< ScalarType >::MvTransMv ( const ScalarType  alpha,
const MultiVec< ScalarType > &  A,
Teuchos::SerialDenseMatrix< int, ScalarType > &  B 
) const
pure virtual

Compute a dense matrix B through the matrix-matrix multiply alpha * A^T * (*this).

Implemented in Belos::GmresPolyMv< ScalarType, MV >.

template<class ScalarType>
virtual void Belos::MultiVec< ScalarType >::MvDot ( const MultiVec< ScalarType > &  A,
std::vector< ScalarType > &  b 
) const
pure virtual

Compute the dot product of each column of *this with the corresponding column of A.

Compute a vector b whose entries are the individual dot-products. That is, b[i] = A[i]^H * (*this)[i] where A[i] is the i-th column of A.

Implemented in Belos::GmresPolyMv< ScalarType, MV >.

template<class ScalarType>
virtual void Belos::MultiVec< ScalarType >::MvNorm ( std::vector< typename Teuchos::ScalarTraits< ScalarType >::magnitudeType > &  normvec,
NormType  type = TwoNorm 
) const
pure virtual

Compute the norm of each vector in *this.

Parameters
normvec[out] On output, normvec[i] holds the norm of the i-th vector of *this.
type[in] The type of norm to compute. The 2-norm is the default.

Implemented in Belos::GmresPolyMv< ScalarType, MV >.

template<class ScalarType>
virtual void Belos::MultiVec< ScalarType >::SetBlock ( const MultiVec< ScalarType > &  A,
const std::vector< int > &  index 
)
pure virtual

Copy the vectors in A to a set of vectors in *this.

The numvecs vectors in A are copied to a subset of vectors in *this indicated by the indices given in index.

Implemented in Belos::GmresPolyMv< ScalarType, MV >.

template<class ScalarType>
virtual void Belos::MultiVec< ScalarType >::MvRandom ( )
pure virtual

Fill all the vectors in *this with random numbers.

Implemented in Belos::GmresPolyMv< ScalarType, MV >.

template<class ScalarType>
virtual void Belos::MultiVec< ScalarType >::MvInit ( const ScalarType  alpha)
pure virtual

Replace each element of the vectors in *this with alpha.

Implemented in Belos::GmresPolyMv< ScalarType, MV >.

template<class ScalarType>
virtual void Belos::MultiVec< ScalarType >::MvPrint ( std::ostream &  os) const
pure virtual

Print *this multivector to the os output stream.

Implemented in Belos::GmresPolyMv< ScalarType, MV >.


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