|
Anasazi
Version of the Day
|
Basic adapter class for Anasazi::MultiVec that uses Epetra_MultiVector. More...
#include <AnasaziEpetraAdapter.hpp>
Public Member Functions | |
| ptrdiff_t | GetGlobalLength () const |
| The number of rows in the multivector. More... | |
Constructors/Destructors | |
| EpetraMultiVec (const Epetra_BlockMap &Map_in, const int numvecs) | |
| Basic EpetraMultiVec constructor. More... | |
| EpetraMultiVec (const Epetra_MultiVector &P_vec) | |
| Copy constructor. More... | |
| EpetraMultiVec (const Epetra_BlockMap &Map_in, double *array, const int numvecs, const int stride=0) | |
| Create multi-vector with values from two dimensional array. More... | |
| EpetraMultiVec (Epetra_DataAccess CV, const Epetra_MultiVector &P_vec, const std::vector< int > &index) | |
| Create multi-vector from list of vectors in an existing EpetraMultiVec. More... | |
| virtual | ~EpetraMultiVec () |
| Destructor. More... | |
Creation methods | |
| MultiVec< double > * | Clone (const int numvecs) const |
Creates a new empty EpetraMultiVec containing numvecs columns. More... | |
| MultiVec< double > * | CloneCopy () const |
Creates a new EpetraMultiVec and copies contents of *this into the new vector (deep copy). More... | |
| MultiVec< double > * | CloneCopy (const std::vector< int > &index) const |
Creates a new EpetraMultiVec and copies the selected contents of *this into the new vector (deep copy). More... | |
| MultiVec< double > * | CloneViewNonConst (const std::vector< int > &index) |
Creates a new EpetraMultiVec that shares the selected contents of *this. More... | |
| const MultiVec< double > * | CloneView (const std::vector< int > &index) const |
Creates a new EpetraMultiVec that shares the selected contents of *this. More... | |
Attribute methods | |
Obtain the vector length of *this. | |
| int | GetNumberVecs () const |
| The number of vectors (i.e., columns) in the multivector. More... | |
Update methods | |
| void | MvTimesMatAddMv (double alpha, const MultiVec< double > &A, const Teuchos::SerialDenseMatrix< int, double > &B, double beta) |
Update *this with . More... | |
| void | MvAddMv (double alpha, const MultiVec< double > &A, double beta, const MultiVec< double > &B) |
Replace *this with . More... | |
| void | MvTransMv (double alpha, const MultiVec< double > &A, Teuchos::SerialDenseMatrix< int, double > &B) const |
Compute a dense matrix B through the matrix-matrix multiply . More... | |
| 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. where A[i] is the i-th column of A. More... | |
| void | MvScale (double alpha) |
Scale each element of the vectors in *this with alpha. More... | |
| void | MvScale (const std::vector< double > &alpha) |
Scale each element of the i-th vector in *this with alpha[i]. More... | |
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. More... | |
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. More... | |
| void | MvRandom () |
Fill the vectors in *this with random numbers. More... | |
| void | MvInit (double alpha) |
Replace each element of the vectors in *this with alpha. More... | |
Accessor methods (inherited from EpetraMultiVecAccessor) | |
| Epetra_MultiVector * | GetEpetraMultiVec () |
| Return the pointer to the Epetra_MultiVector object. More... | |
| const Epetra_MultiVector * | GetEpetraMultiVec () const |
| Return the pointer to the Epetra_MultiVector object. More... | |
Print method | |
| void | MvPrint (std::ostream &os) const |
Print *this EpetraMultiVec. More... | |
Public Member Functions inherited from Anasazi::MultiVec< double > | |
| MultiVec () | |
| Default constructor. More... | |
| virtual | ~MultiVec () |
| Destructor (virtual for memory safety of derived classes). More... | |
| virtual void | MvNorm (std::vector< typename Teuchos::ScalarTraits< double >::magnitudeType > &normvec) const =0 |
Compute the 2-norm of each vector in *this. More... | |
Public Member Functions inherited from Anasazi::EpetraMultiVecAccessor | |
| virtual | ~EpetraMultiVecAccessor () |
| Destructor. More... | |
Basic adapter class for Anasazi::MultiVec that uses Epetra_MultiVector.
Definition at line 96 of file AnasaziEpetraAdapter.hpp.
| Anasazi::EpetraMultiVec::EpetraMultiVec | ( | const Epetra_BlockMap & | Map_in, |
| const int | numvecs | ||
| ) |
Basic EpetraMultiVec constructor.
| Map | [in] An Epetra_LocalMap, Epetra_Map or Epetra_BlockMap. |
| numvecs | [in] Number of vectors in multi-vector. |
Definition at line 33 of file AnasaziEpetraAdapter.cpp.
| Anasazi::EpetraMultiVec::EpetraMultiVec | ( | const Epetra_MultiVector & | P_vec | ) |
Copy constructor.
Definition at line 47 of file AnasaziEpetraAdapter.cpp.
| Anasazi::EpetraMultiVec::EpetraMultiVec | ( | const Epetra_BlockMap & | Map_in, |
| double * | array, | ||
| const int | numvecs, | ||
| const int | stride = 0 |
||
| ) |
Create multi-vector with values from two dimensional array.
| 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 26 of file AnasaziEpetraAdapter.cpp.
| Anasazi::EpetraMultiVec::EpetraMultiVec | ( | Epetra_DataAccess | CV, |
| const Epetra_MultiVector & | P_vec, | ||
| const std::vector< int > & | index | ||
| ) |
Create multi-vector from list of vectors in an existing EpetraMultiVec.
| CV | [in] Enumerated type set to Copy or View. |
| P_vec | [in] 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 39 of file AnasaziEpetraAdapter.cpp.
|
inlinevirtual |
Destructor.
Definition at line 133 of file AnasaziEpetraAdapter.hpp.
|
virtual |
Creates a new empty EpetraMultiVec containing numvecs columns.
Implements Anasazi::MultiVec< double >.
Definition at line 62 of file AnasaziEpetraAdapter.cpp.
|
virtual |
Creates a new EpetraMultiVec and copies contents of *this into the new vector (deep copy).
Implements Anasazi::MultiVec< double >.
Definition at line 73 of file AnasaziEpetraAdapter.cpp.
|
virtual |
Creates a new EpetraMultiVec 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 80 of file AnasaziEpetraAdapter.cpp.
|
virtual |
Creates a new EpetraMultiVec 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 87 of file AnasaziEpetraAdapter.cpp.
|
virtual |
Creates a new EpetraMultiVec 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 93 of file AnasaziEpetraAdapter.cpp.
|
inlinevirtual |
The number of rows in the multivector.
Implements Anasazi::MultiVec< double >.
Definition at line 183 of file AnasaziEpetraAdapter.hpp.
|
inlinevirtual |
The number of vectors (i.e., columns) in the multivector.
Implements Anasazi::MultiVec< double >.
Definition at line 193 of file AnasaziEpetraAdapter.hpp.
|
virtual |
Update *this with
.
Implements Anasazi::MultiVec< double >.
Definition at line 126 of file AnasaziEpetraAdapter.cpp.
|
virtual |
Replace *this with
.
Implements Anasazi::MultiVec< double >.
Definition at line 146 of file AnasaziEpetraAdapter.cpp.
|
virtual |
Compute a dense matrix B through the matrix-matrix multiply
.
Implements Anasazi::MultiVec< double >.
Definition at line 165 of file AnasaziEpetraAdapter.cpp.
|
virtual |
Compute a vector b where the components are the individual dot-products, i.e.
where A[i] is the i-th column of A.
Implements Anasazi::MultiVec< double >.
Definition at line 190 of file AnasaziEpetraAdapter.cpp.
|
inlinevirtual |
Scale each element of the vectors in *this with alpha.
Implements Anasazi::MultiVec< double >.
Definition at line 228 of file AnasaziEpetraAdapter.hpp.
|
virtual |
Scale each element of the i-th vector in *this with alpha[i].
Implements Anasazi::MultiVec< double >.
Definition at line 211 of file AnasaziEpetraAdapter.cpp.
|
inline |
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 244 of file AnasaziEpetraAdapter.hpp.
|
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 100 of file AnasaziEpetraAdapter.cpp.
|
inlinevirtual |
Fill the vectors in *this with random numbers.
Implements Anasazi::MultiVec< double >.
Definition at line 263 of file AnasaziEpetraAdapter.hpp.
|
inlinevirtual |
Replace each element of the vectors in *this with alpha.
Implements Anasazi::MultiVec< double >.
Definition at line 270 of file AnasaziEpetraAdapter.hpp.
|
inlinevirtual |
Return the pointer to the Epetra_MultiVector object.
Reimplemented from Anasazi::EpetraMultiVecAccessor.
Definition at line 279 of file AnasaziEpetraAdapter.hpp.
|
inlinevirtual |
Return the pointer to the Epetra_MultiVector object.
Reimplemented from Anasazi::EpetraMultiVecAccessor.
Definition at line 282 of file AnasaziEpetraAdapter.hpp.
|
inlinevirtual |
Print *this EpetraMultiVec.
Implements Anasazi::MultiVec< double >.
Definition at line 291 of file AnasaziEpetraAdapter.hpp.
1.8.5