| 
    Anasazi Version of the Day
    
   | 
 
Specialized adapter class for Anasazi::MultiVec that uses Epetra_MultiVector and Epetra_Operator to define the inner-product. More...
#include <AnasaziSpecializedEpetraAdapter.hpp>
  
Public Member Functions | |
Constructors/Destructors  | |
| EpetraOpMultiVec (const Teuchos::RCP< Epetra_Operator > &Op, const Epetra_BlockMap &Map_in, const int numvecs) | |
| Basic EpetraOpMultiVec constructor.   | |
| EpetraOpMultiVec (const Teuchos::RCP< Epetra_Operator > &Op, const Epetra_BlockMap &Map_in, double *array, const int numvecs, const int stride=0) | |
| Create multi-vector with values from two dimensional array.   | |
| EpetraOpMultiVec (const Teuchos::RCP< Epetra_Operator > &Op, Epetra_DataAccess CV, const Epetra_MultiVector &P_vec, const std::vector< int > &index) | |
| Create multi-vector from list of vectors in an existing EpetraOpMultiVec.   | |
| EpetraOpMultiVec (const EpetraOpMultiVec &P_vec) | |
| Copy constructor.   | |
| virtual | ~EpetraOpMultiVec () | 
| Destructor.   | |
Creation methods  | |
| MultiVec< double > * | Clone (const int numvecs) const | 
Creates a new empty EpetraOpMultiVec containing numvecs columns.   | |
| MultiVec< double > * | CloneCopy () const | 
Creates a new EpetraOpMultiVec and copies contents of *this into the new vector (deep copy).   | |
| MultiVec< double > * | CloneCopy (const std::vector< int > &index) const | 
Creates a new EpetraOpMultiVec and copies the selected contents of *this into the new vector (deep copy).  | |
| MultiVec< double > * | CloneViewNonConst (const std::vector< int > &index) | 
Creates a new EpetraOpMultiVec that shares the selected contents of *this.   | |
| const MultiVec< double > * | CloneView (const std::vector< int > &index) const | 
Creates a new EpetraOpMultiVec that shares the selected contents of *this.   | |
Accessor methods  | |
Attribute methods  | |
| ptrdiff_t | GetGlobalLength () const | 
| The number of rows in the multivector.   | |
| int | GetNumberVecs () const | 
| Obtain the vector length of *this.   | |
Update methods  | |
| void | MvTimesMatAddMv (double alpha, const MultiVec< double > &A, const Teuchos::SerialDenseMatrix< int, double > &B, double beta) | 
Update *this with  | |
| void | MvAddMv (double alpha, const MultiVec< double > &A, double beta, const MultiVec< double > &B) | 
Replace *this with  | |
| void | MvTransMv (double alpha, const MultiVec< double > &A, Teuchos::SerialDenseMatrix< int, double > &B) const | 
Compute a dense matrix B through the matrix-matrix multiply  | |
| void | MvDot (const MultiVec< double > &A, std::vector< double > &b) const | 
Compute a vector b where the components are the individual dot-products, i.e. A[i] is the i-th column of A.   | |
| void | MvScale (double alpha) | 
Scale each element of the vectors in *this with alpha.   | |
| void | MvScale (const std::vector< double > &alpha) | 
Scale each element of the i-th vector in *this with alpha[i].   | |
Norm method  | |
| void | MvNorm (std::vector< double > &normvec) const | 
Compute the 2-norm of each individual vector of *this. Upon return, normvec[i] holds the 2-norm of the i-th vector of *this.   | |
Initialization methods  | |
| void | SetBlock (const MultiVec< double > &A, const std::vector< int > &index) | 
Copy the vectors in A to a set of vectors in *this.  | |
| void | MvRandom () | 
Fill the vectors in *this with random numbers.   | |
| void | MvInit (double alpha) | 
Replace each element of the vectors in *this with alpha.   | |
Accessor methods (inherited from EpetraMultiVecAccessor)  | |
| Epetra_MultiVector * | GetEpetraMultiVec () | 
| Return the pointer to the Epetra_MultiVector object.   | |
| const Epetra_MultiVector * | GetEpetraMultiVec () const | 
| Return the pointer to the Epetra_MultiVector object.   | |
  Public Member Functions inherited from Anasazi::MultiVec< double > | |
| MultiVec () | |
| Default constructor.   | |
| virtual | ~MultiVec () | 
| Destructor (virtual for memory safety of derived classes).   | |
| virtual void | MvNorm (std::vector< typename Teuchos::ScalarTraits< double >::magnitudeType > &normvec) const=0 | 
Compute the 2-norm of each vector in *this.  | |
  Public Member Functions inherited from Anasazi::EpetraMultiVecAccessor | |
| virtual | ~EpetraMultiVecAccessor () | 
| Destructor.   | |
Print method | |
| void | MvPrint (std::ostream &os) const | 
Print *this EpetraOpMultiVec.   | |
Specialized adapter class for Anasazi::MultiVec that uses Epetra_MultiVector and Epetra_Operator to define the inner-product.
Definition at line 66 of file AnasaziSpecializedEpetraAdapter.hpp.
| Anasazi::EpetraOpMultiVec::EpetraOpMultiVec | ( | const Teuchos::RCP< Epetra_Operator > & | Op, | 
| const Epetra_BlockMap & | Map_in, | ||
| const int | numvecs | ||
| ) | 
Basic EpetraOpMultiVec constructor.
| Op | [in] A reference-counted pointer to an existing fully constructed Epetra_Operator. | 
| Map | [in] An Epetra_LocalMap, Epetra_Map or Epetra_BlockMap. | 
| numvecs | [in] Number of vectors in multi-vector. | 
Definition at line 27 of file AnasaziSpecializedEpetraAdapter.cpp.
| Anasazi::EpetraOpMultiVec::EpetraOpMultiVec | ( | const Teuchos::RCP< Epetra_Operator > & | Op, | 
| const Epetra_BlockMap & | Map_in, | ||
| double * | array, | ||
| const int | numvecs, | ||
| const int | stride = 0  | 
        ||
| ) | 
Create multi-vector with values from two dimensional array.
| Op | [in] A reference-counted pointer to an existing fully constructed Epetra_Operator. | 
| Map | [in] An Epetra_LocalMap, Epetra_Map or Epetra_BlockMap | 
| array | [in] Pointer to an array of double precision numbers. The first vector starts at array, the second at array+stride, and so on. This array is copied.  | 
| numvecs | [in] Number of vectors in the multi-vector. | 
| stride | [in] The stride between vectors in memory of array. | 
Definition at line 34 of file AnasaziSpecializedEpetraAdapter.cpp.
| Anasazi::EpetraOpMultiVec::EpetraOpMultiVec | ( | const Teuchos::RCP< Epetra_Operator > & | Op, | 
| Epetra_DataAccess | CV, | ||
| const Epetra_MultiVector & | P_vec, | ||
| const std::vector< int > & | index | ||
| ) | 
Create multi-vector from list of vectors in an existing EpetraOpMultiVec.
| Op | [in] A reference-counted pointer to an existing fully constructed Epetra_Operator. | 
| P_vec | [in] A reference-counted pointer to an existing fully constructed Epetra_MultiVector. | 
| index | [in] A integer vector containing the indices of the vectors to copy out of P_vec. | 
Definition at line 42 of file AnasaziSpecializedEpetraAdapter.cpp.
| Anasazi::EpetraOpMultiVec::EpetraOpMultiVec | ( | const EpetraOpMultiVec & | P_vec | ) | 
Copy constructor.
Definition at line 50 of file AnasaziSpecializedEpetraAdapter.cpp.
      
  | 
  inlinevirtual | 
Destructor.
Definition at line 105 of file AnasaziSpecializedEpetraAdapter.hpp.
      
  | 
  virtual | 
Creates a new empty EpetraOpMultiVec containing numvecs columns. 
Implements Anasazi::MultiVec< double >.
Definition at line 66 of file AnasaziSpecializedEpetraAdapter.cpp.
      
  | 
  virtual | 
Creates a new EpetraOpMultiVec and copies contents of *this into the new vector (deep copy). 
Implements Anasazi::MultiVec< double >.
Definition at line 77 of file AnasaziSpecializedEpetraAdapter.cpp.
      
  | 
  virtual | 
Creates a new EpetraOpMultiVec and copies the selected contents of *this into the new vector (deep copy). 
 
The copied vectors from *this are indicated by the index.size() indices in index.
Implements Anasazi::MultiVec< double >.
Definition at line 84 of file AnasaziSpecializedEpetraAdapter.cpp.
      
  | 
  virtual | 
Creates a new EpetraOpMultiVec that shares the selected contents of *this. 
The index of the numvecs vectors shallow copied from *this are indicated by the indices given in index.
Implements Anasazi::MultiVec< double >.
Definition at line 91 of file AnasaziSpecializedEpetraAdapter.cpp.
      
  | 
  virtual | 
Creates a new EpetraOpMultiVec that shares the selected contents of *this. 
The index of the numvecs vectors shallow copied from *this are indicated by the indices given in index.
Implements Anasazi::MultiVec< double >.
Definition at line 97 of file AnasaziSpecializedEpetraAdapter.cpp.
      
  | 
  inlinevirtual | 
The number of rows in the multivector.
Implements Anasazi::MultiVec< double >.
Definition at line 161 of file AnasaziSpecializedEpetraAdapter.hpp.
      
  | 
  inlinevirtual | 
Obtain the vector length of *this.
Implements Anasazi::MultiVec< double >.
Definition at line 170 of file AnasaziSpecializedEpetraAdapter.hpp.
      
  | 
  virtual | 
Update *this with 
Implements Anasazi::MultiVec< double >.
Definition at line 130 of file AnasaziSpecializedEpetraAdapter.cpp.
      
  | 
  virtual | 
Replace *this with 
Implements Anasazi::MultiVec< double >.
Definition at line 150 of file AnasaziSpecializedEpetraAdapter.cpp.
      
  | 
  virtual | 
Compute a dense matrix B through the matrix-matrix multiply 
Implements Anasazi::MultiVec< double >.
Definition at line 169 of file AnasaziSpecializedEpetraAdapter.cpp.
      
  | 
  virtual | 
Compute a vector b where the components are the individual dot-products, i.e. ![$ b[i] = A[i]^H(this[i])$](form_119.png)
A[i] is the i-th column of A. 
Implements Anasazi::MultiVec< double >.
Definition at line 198 of file AnasaziSpecializedEpetraAdapter.cpp.
      
  | 
  inlinevirtual | 
Scale each element of the vectors in *this with alpha. 
Implements Anasazi::MultiVec< double >.
Definition at line 205 of file AnasaziSpecializedEpetraAdapter.hpp.
      
  | 
  virtual | 
Scale each element of the i-th vector in *this with alpha[i]. 
Implements Anasazi::MultiVec< double >.
Definition at line 245 of file AnasaziSpecializedEpetraAdapter.cpp.
| void Anasazi::EpetraOpMultiVec::MvNorm | ( | std::vector< double > & | normvec | ) | const | 
Compute the 2-norm of each individual vector of *this. 
 Upon return, normvec[i] holds the 2-norm of the i-th vector of *this. 
Definition at line 224 of file AnasaziSpecializedEpetraAdapter.cpp.
      
  | 
  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. 
Implements Anasazi::MultiVec< double >.
Definition at line 104 of file AnasaziSpecializedEpetraAdapter.cpp.
      
  | 
  inlinevirtual | 
Fill the vectors in *this with random numbers. 
Implements Anasazi::MultiVec< double >.
Definition at line 236 of file AnasaziSpecializedEpetraAdapter.hpp.
      
  | 
  inlinevirtual | 
Replace each element of the vectors in *this with alpha. 
Implements Anasazi::MultiVec< double >.
Definition at line 243 of file AnasaziSpecializedEpetraAdapter.hpp.
      
  | 
  inlinevirtual | 
Return the pointer to the Epetra_MultiVector object.
Reimplemented from Anasazi::EpetraMultiVecAccessor.
Definition at line 252 of file AnasaziSpecializedEpetraAdapter.hpp.
      
  | 
  inlinevirtual | 
Return the pointer to the Epetra_MultiVector object.
Reimplemented from Anasazi::EpetraMultiVecAccessor.
Definition at line 255 of file AnasaziSpecializedEpetraAdapter.hpp.
      
  | 
  inlinevirtual | 
Print *this EpetraOpMultiVec. 
Implements Anasazi::MultiVec< double >.
Definition at line 264 of file AnasaziSpecializedEpetraAdapter.hpp.