Belos  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
BelosMultiVec.hpp
Go to the documentation of this file.
1 // @HEADER
2 // *****************************************************************************
3 // Belos: Block Linear Solvers Package
4 //
5 // Copyright 2004-2016 NTESS and the Belos contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 // @HEADER
9 
10 #ifndef BELOS_MULTI_VEC_HPP
11 #define BELOS_MULTI_VEC_HPP
12 
27 
28 #include "BelosMultiVecTraits.hpp"
29 #include "BelosTypes.hpp"
30 #include "BelosConfigDefs.hpp"
31 
32 namespace Belos {
33 
58 template <class ScalarType>
59 class MultiVec {
60 public:
62 
63  MultiVec() {};
65 
67  virtual ~MultiVec () {};
68 
70 
72 
75  virtual MultiVec<ScalarType> * Clone ( const int numvecs ) const = 0;
76 
79  virtual MultiVec<ScalarType> * CloneCopy () const = 0;
80 
87  virtual MultiVec<ScalarType> * CloneCopy ( const std::vector<int>& index ) const = 0;
88 
95  virtual MultiVec<ScalarType> * CloneViewNonConst ( const std::vector<int>& index ) = 0;
96 
103  virtual const MultiVec<ScalarType> * CloneView ( const std::vector<int>& index ) const = 0;
104 
106 
108 
110  virtual ptrdiff_t GetGlobalLength () const = 0;
111 
113  virtual int GetNumberVecs () const = 0;
114 
116 
118 
120  virtual void
121  MvTimesMatAddMv (const ScalarType alpha,
122  const MultiVec<ScalarType>& A,
123  const Teuchos::SerialDenseMatrix<int,ScalarType>& B, const ScalarType beta) = 0;
124 
126  virtual void MvAddMv ( const ScalarType alpha, const MultiVec<ScalarType>& A, const ScalarType beta, const MultiVec<ScalarType>& B ) = 0;
127 
129  virtual void MvScale ( const ScalarType alpha ) = 0;
130 
132  virtual void MvScale ( const std::vector<ScalarType>& alpha ) = 0;
133 
137  virtual void MvTransMv ( const ScalarType alpha, const MultiVec<ScalarType>& A, Teuchos::SerialDenseMatrix<int,ScalarType>& B) const = 0;
138 
144  virtual void MvDot ( const MultiVec<ScalarType>& A, std::vector<ScalarType>& b ) const = 0;
145 
147 
149 
155  virtual void MvNorm ( std::vector<typename Teuchos::ScalarTraits<ScalarType>::magnitudeType>& normvec, NormType type = TwoNorm ) const = 0;
156 
158 
160 
165  virtual void SetBlock ( const MultiVec<ScalarType>& A, const std::vector<int>& index ) = 0;
166 
168  virtual void MvRandom () = 0;
169 
171  virtual void MvInit ( const ScalarType alpha ) = 0;
172 
174 
176 
178  virtual void MvPrint ( std::ostream& os ) const = 0;
180 
181 #ifdef HAVE_BELOS_TSQR
182 
184 
207  virtual void
208  factorExplicit (MultiVec<ScalarType>& Q,
210  const bool forceNonnegativeDiagonal=false)
211  {
212  TEUCHOS_TEST_FOR_EXCEPTION(true, std::logic_error, "The Belos::MultiVec<"
213  << Teuchos::TypeNameTraits<ScalarType>::name() << "> subclass which you "
214  "are using does not implement the TSQR-related method factorExplicit().");
215  }
216 
251  virtual int
254  {
255  TEUCHOS_TEST_FOR_EXCEPTION(true, std::logic_error, "The Belos::MultiVec<"
256  << Teuchos::TypeNameTraits<ScalarType>::name() << "> subclass which you "
257  "are using does not implement the TSQR-related method revealRank().");
258  }
259 
261 #endif // HAVE_BELOS_TSQR
262 };
263 
264 
265 namespace details {
283 template<class ScalarType>
285 public:
287  typedef ScalarType scalar_type;
288  typedef int ordinal_type; // This doesn't matter either
289  typedef int node_type; // Nor does this
292 
294  void
296  MV& Q,
298  const bool forceNonnegativeDiagonal=false)
299  {
300  A.factorExplicit (Q, R, forceNonnegativeDiagonal);
301  }
302 
304  int
307  const magnitude_type& tol)
308  {
309  return Q.revealRank (R, tol);
310  }
311 
313  (void) params;
314  }
315 
317  return Teuchos::parameterList ();
318  }
319 };
320 } // namespace details
321 
331  template<class ScalarType>
332  class MultiVecTraits<ScalarType,MultiVec<ScalarType> > {
333  public:
335 
336 
340  Clone (const MultiVec<ScalarType>& mv, const int numvecs) {
341  return Teuchos::rcp (const_cast<MultiVec<ScalarType>&> (mv).Clone (numvecs));
342  }
345  { return Teuchos::rcp( const_cast<MultiVec<ScalarType>&>(mv).CloneCopy() ); }
347  static Teuchos::RCP<MultiVec<ScalarType> > CloneCopy( const MultiVec<ScalarType>& mv, const std::vector<int>& index )
348  { return Teuchos::rcp( const_cast<MultiVec<ScalarType>&>(mv).CloneCopy(index) ); }
351  CloneViewNonConst (MultiVec<ScalarType>& mv, const std::vector<int>& index)
352  {
353  return Teuchos::rcp( mv.CloneViewNonConst(index) );
354  }
355 
358  {
359  // mfh 02 Mar 2013: For now, we'll just use the above index
360  // vector version of CloneViewNonConst to implement this, since
361  // that doesn't require adding to the MultiVec interface.
362  std::vector<int> indVec (index.size ());
363  for (int k = 0; k < index.size (); ++k) {
364  indVec[k] = k;
365  }
366  return CloneViewNonConst (mv, indVec);
367  }
368 
371  CloneView (const MultiVec<ScalarType>& mv, const std::vector<int>& index) {
372  return Teuchos::rcp( const_cast<MultiVec<ScalarType>&>(mv).CloneView(index) );
373  }
374 
377  {
378  // mfh 02 Mar 2013: For now, we'll just use the above index
379  // vector version of CloneView to implement this, since that
380  // doesn't require adding to the MultiVec interface.
381  std::vector<int> indVec (index.size ());
382  for (int k = 0; k < index.size (); ++k) {
383  indVec[k] = k;
384  }
385  return CloneView (mv, indVec);
386  }
387 
389  static ptrdiff_t GetGlobalLength( const MultiVec<ScalarType>& mv )
390  { return mv.GetGlobalLength(); }
392  static int GetNumberVecs( const MultiVec<ScalarType>& mv )
393  { return mv.GetNumberVecs(); }
395  static void MvTimesMatAddMv( ScalarType alpha, const MultiVec<ScalarType>& A,
397  ScalarType beta, MultiVec<ScalarType>& mv )
398  { mv.MvTimesMatAddMv(alpha, A, B, beta); }
400  static void MvAddMv( ScalarType alpha, const MultiVec<ScalarType>& A, ScalarType beta, const MultiVec<ScalarType>& B, MultiVec<ScalarType>& mv )
401  { mv.MvAddMv(alpha, A, beta, B); }
403  static void MvScale ( MultiVec<ScalarType>& mv, const ScalarType alpha )
404  { mv.MvScale( alpha ); }
405 
406  static void MvScale ( MultiVec<ScalarType>& mv, const std::vector<ScalarType>& alpha )
407  { mv.MvScale(alpha); }
409  static void MvTransMv( const ScalarType alpha, const MultiVec<ScalarType>& A, const MultiVec<ScalarType>& mv, Teuchos::SerialDenseMatrix<int,ScalarType>& B )
410  { mv.MvTransMv(alpha, A, B); }
412  static void MvDot( const MultiVec<ScalarType>& mv, const MultiVec<ScalarType>& A, std::vector<ScalarType>& b )
413  { mv.MvDot( A, b ); }
415  static void MvNorm( const MultiVec<ScalarType>& mv, std::vector<typename Teuchos::ScalarTraits<ScalarType>::magnitudeType>& normvec, NormType type = TwoNorm )
416  { mv.MvNorm(normvec,type); }
418  static void SetBlock( const MultiVec<ScalarType>& A, const std::vector<int>& index, MultiVec<ScalarType>& mv )
419  { mv.SetBlock(A, index); }
420 
421  static void
424  {
425  // mfh 02 Mar 2013: For now, we'll just use SetBlock to implement this,
426  // since that doesn't require adding to the MultiVec interface.
427  const int numVecsRhs = GetNumberVecs (A);
428  const int numVecsLhs = GetNumberVecs (mv);
429 
431  numVecsLhs != numVecsRhs, std::invalid_argument,
432  "Belos::MultiVecTraits::Assign: Input multivector A has " << numVecsRhs
433  << " columns, which differs from the number of columns " << numVecsLhs
434  << " in the output multivector mv.");
435 
436  // mfh 02 Mar 2013: It's pretty silly to build this each time.
437  // However, at least that makes the code correct.
438  std::vector<int> index (numVecsRhs);
439  for (int k = 0; k < numVecsRhs; ++k) {
440  index[k] = k;
441  }
442  SetBlock (A, index, mv);
443  }
444 
446  static void MvRandom( MultiVec<ScalarType>& mv )
447  { mv.MvRandom(); }
449  static void MvInit( MultiVec<ScalarType>& mv, ScalarType alpha = Teuchos::ScalarTraits<ScalarType>::zero() )
450  { mv.MvInit(alpha); }
452  static void MvPrint( const MultiVec<ScalarType>& mv, std::ostream& os )
453  { mv.MvPrint(os); }
454 
455 #ifdef HAVE_BELOS_TSQR
456  typedef details::MultiVecTsqrAdapter<ScalarType> tsqr_adaptor_type;
464 #endif // HAVE_BELOS_TSQR
465  };
466 
467 
468 } // namespace Belos
469 
470 #endif
471 
472 // end of file BelosMultiVec.hpp
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)
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.
Definition: BelosTypes.hpp:65
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
Ordinal size() const
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&#39; 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 > &params)

Generated on Wed Oct 23 2024 09:24:26 for Belos by doxygen 1.8.5