42 #ifndef BELOS_MULTI_VEC_HPP 
   43 #define BELOS_MULTI_VEC_HPP 
   90 template <
class ScalarType>
 
  161   virtual void MvScale ( 
const ScalarType alpha ) = 0;
 
  164   virtual void MvScale ( 
const std::vector<ScalarType>& alpha ) = 0;
 
  203   virtual void MvInit ( 
const ScalarType alpha ) = 0;
 
  210   virtual void MvPrint ( std::ostream& os ) 
const = 0;
 
  213 #ifdef HAVE_BELOS_TSQR 
  242       const bool forceNonnegativeDiagonal=
false)
 
  246       "are using does not implement the TSQR-related method factorExplicit().");
 
  289       "are using does not implement the TSQR-related method revealRank().");
 
  293 #endif // HAVE_BELOS_TSQR 
  315 template<
class ScalarType>
 
  330       const bool forceNonnegativeDiagonal=
false)
 
  332     A.factorExplicit (Q, R, forceNonnegativeDiagonal);
 
  341     return Q.revealRank (R, tol);
 
  349     return Teuchos::parameterList ();
 
  363   template<
class ScalarType>
 
  394       std::vector<int> indVec (index.
size ());
 
  395       for (
int k = 0; k < index.
size (); ++k) {
 
  413       std::vector<int> indVec (index.
size ());
 
  414       for (
int k = 0; k < index.
size (); ++k) {
 
  433     { mv.
MvAddMv(alpha, A, beta, B); }
 
  445     { mv.
MvDot( A, b ); }
 
  448     { mv.
MvNorm(normvec,type); }
 
  463         numVecsLhs != numVecsRhs, std::invalid_argument, 
 
  464   "Belos::MultiVecTraits::Assign: Input multivector A has " << numVecsRhs
 
  465   << 
" columns, which differs from the number of columns " << numVecsLhs
 
  466   << 
" in the output multivector mv.");
 
  470       std::vector<int> index (numVecsRhs);
 
  471       for (
int k = 0; k < numVecsRhs; ++k) {
 
  487 #ifdef HAVE_BELOS_TSQR 
  496 #endif // HAVE_BELOS_TSQR 
Collection of types and exceptions used within the Belos solvers. 
 
static void MvTimesMatAddMv(ScalarType alpha, const MultiVec< ScalarType > &A, const Teuchos::SerialDenseMatrix< int, ScalarType > &B, ScalarType beta, MultiVec< ScalarType > &mv)
 
virtual MultiVec< ScalarType > * Clone(const int numvecs) const =0
Create a new MultiVec with numvecs columns. 
 
virtual MultiVec< ScalarType > * CloneCopy() const =0
Create a new MultiVec and copy contents of *this into it (deep copy). 
 
static Teuchos::RCP< const MultiVec< ScalarType > > CloneView(const MultiVec< ScalarType > &mv, const std::vector< int > &index)
 
static void MvTransMv(const ScalarType alpha, const MultiVec< ScalarType > &A, const MultiVec< ScalarType > &mv, Teuchos::SerialDenseMatrix< int, ScalarType > &B)
 
static Teuchos::RCP< MultiVec< ScalarType > > Clone(const MultiVec< ScalarType > &mv, const int numvecs)
Create a new empty MultiVec containing numvecs columns. 
 
static void MvScale(MultiVec< ScalarType > &mv, const ScalarType alpha)
 
virtual int GetNumberVecs() const =0
The number of vectors (i.e., columns) in the multivector. 
 
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 ve...
 
Teuchos::SerialDenseMatrix< ordinal_type, scalar_type > dense_matrix_type
 
MultiVec()
Default constructor. 
 
static void Assign(const MultiVec< ScalarType > &A, MultiVec< ScalarType > &mv)
 
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
 
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). 
 
Declaration of basic traits for the multivector type. 
 
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. 
 
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 ve...
 
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). 
 
virtual void MvScale(const ScalarType alpha)=0
Scale each element of the vectors in *this with alpha. 
 
virtual ptrdiff_t GetGlobalLength() const =0
The number of rows in the multivector. 
 
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). 
 
static Teuchos::RCP< const MultiVec< ScalarType > > CloneView(const MultiVec< ScalarType > &mv, const Teuchos::Range1D &index)
 
static void MvDot(const MultiVec< ScalarType > &mv, const MultiVec< ScalarType > &A, std::vector< ScalarType > &b)
 
MultiVec< ScalarType > MV
 
static int GetNumberVecs(const MV &mv)
Obtain the number of vectors in mv. 
 
virtual void MvNorm(std::vector< typename Teuchos::ScalarTraits< ScalarType >::magnitudeType > &normvec, NormType type=TwoNorm) const =0
Compute the norm of each vector in *this. 
 
int revealRank(MV &Q, dense_matrix_type &R, const magnitude_type &tol)
Compute rank-revealing decomposition using results of factorExplicit(). 
 
Traits class which defines basic operations on multivectors. 
 
static void MvRandom(MultiVec< ScalarType > &mv)
 
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. 
 
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...
 
static ptrdiff_t GetGlobalLength(const MultiVec< ScalarType > &mv)
 
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). 
 
static Teuchos::RCP< MV > Clone(const MV &mv, const int numvecs)
Creates a new empty MV containing numvecs columns. 
 
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
 
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const 
 
static void MvScale(MultiVec< ScalarType > &mv, const std::vector< ScalarType > &alpha)
 
static Teuchos::RCP< MV > CloneCopy(const MV &mv)
Creates a new MV and copies contents of mv into the new vector (deep copy). 
 
TSQR adapter for MultiVec. 
 
static Teuchos::RCP< MultiVec< ScalarType > > CloneViewNonConst(MultiVec< ScalarType > &mv, const std::vector< int > &index)
 
static void MvPrint(const MultiVec< ScalarType > &mv, std::ostream &os)
 
static void MvNorm(const MultiVec< ScalarType > &mv, std::vector< typename Teuchos::ScalarTraits< ScalarType >::magnitudeType > &normvec, NormType type=TwoNorm)
 
static void MvAddMv(ScalarType alpha, const MultiVec< ScalarType > &A, ScalarType beta, const MultiVec< ScalarType > &B, MultiVec< ScalarType > &mv)
 
NormType
The type of vector norm to compute. 
 
void factorExplicit(MV &A, MV &Q, dense_matrix_type &R, const bool forceNonnegativeDiagonal=false)
Compute QR factorization A = QR, using TSQR. 
 
virtual void MvInit(const ScalarType alpha)=0
Replace each element of the vectors in *this with alpha. 
 
static int GetNumberVecs(const MultiVec< ScalarType > &mv)
 
Teuchos::ScalarTraits< scalar_type >::magnitudeType magnitude_type
 
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. 
 
static void MvInit(MultiVec< ScalarType > &mv, ScalarType alpha=Teuchos::ScalarTraits< ScalarType >::zero())
 
Interface for multivectors used by Belos' linear solvers. 
 
static Teuchos::RCP< MultiVec< ScalarType > > CloneCopy(const MultiVec< ScalarType > &mv)
 
virtual ~MultiVec()
Destructor (virtual for memory safety of derived classes). 
 
virtual void MvPrint(std::ostream &os) const =0
Print *this multivector to the os output stream. 
 
static Teuchos::RCP< MultiVec< ScalarType > > CloneViewNonConst(MultiVec< ScalarType > &mv, const Teuchos::Range1D &index)
 
static Teuchos::RCP< MultiVec< ScalarType > > CloneCopy(const MultiVec< ScalarType > &mv, const std::vector< int > &index)
 
virtual void MvRandom()=0
Fill all the vectors in *this with random numbers. 
 
Belos header file which uses auto-configuration information to include necessary C++ headers...
 
static void SetBlock(const MultiVec< ScalarType > &A, const std::vector< int > &index, MultiVec< ScalarType > &mv)
 
void setParameterList(const Teuchos::RCP< Teuchos::ParameterList > ¶ms)