Anasazi  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
Public Member Functions | List of all members
Anasazi::EpetraMultiVec Class Reference

Basic adapter class for Anasazi::MultiVec that uses Epetra_MultiVector. More...

#include <AnasaziEpetraAdapter.hpp>

Inheritance diagram for Anasazi::EpetraMultiVec:
Anasazi::MultiVec< double > Epetra_MultiVector Anasazi::EpetraMultiVecAccessor Epetra_DistObject Epetra_CompObject Epetra_BLAS Epetra_Object Epetra_SrcDistObject

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 $\alpha AB + \beta (*this)$. More...
 
void MvAddMv (double alpha, const MultiVec< double > &A, double beta, const MultiVec< double > &B)
 Replace *this with $\alpha A + \beta B$. 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 $\alpha A^T(*this)$. 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. $ b[i] = A[i]^H(this[i])$ 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_MultiVectorGetEpetraMultiVec ()
 Return the pointer to the Epetra_MultiVector object. More...
 
const Epetra_MultiVectorGetEpetraMultiVec () 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...
 

Detailed Description

Basic adapter class for Anasazi::MultiVec that uses Epetra_MultiVector.

Note
The Epetra package performs double-precision arithmetic, so the use of Epetra with Anasazi will only provide a double-precision eigensolver.
Examples:
BlockKrylovSchur/BlockKrylovSchurEpetraExSVD.cpp, and MVOPTester/MVOPTesterEx.cpp.

Definition at line 128 of file AnasaziEpetraAdapter.hpp.

Constructor & Destructor Documentation

Anasazi::EpetraMultiVec::EpetraMultiVec ( const Epetra_BlockMap Map_in,
const int  numvecs 
)

Basic EpetraMultiVec constructor.

Parameters
Map[in] An Epetra_LocalMap, Epetra_Map or Epetra_BlockMap.
numvecs[in] Number of vectors in multi-vector.
Returns
Pointer to an EpetraMultiVec

Definition at line 65 of file AnasaziEpetraAdapter.cpp.

Anasazi::EpetraMultiVec::EpetraMultiVec ( const Epetra_MultiVector P_vec)

Copy constructor.

Definition at line 79 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.

Parameters
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.
Returns
Pointer to an EpetraMultiVec

Definition at line 58 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.

Parameters
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.
Returns
Pointer to an EpetraMultiVec

Definition at line 71 of file AnasaziEpetraAdapter.cpp.

virtual Anasazi::EpetraMultiVec::~EpetraMultiVec ( )
inlinevirtual

Destructor.

Definition at line 165 of file AnasaziEpetraAdapter.hpp.

Member Function Documentation

MultiVec< double > * Anasazi::EpetraMultiVec::Clone ( const int  numvecs) const
virtual

Creates a new empty EpetraMultiVec containing numvecs columns.

Returns
Pointer to an EpetraMultiVec

Implements Anasazi::MultiVec< double >.

Definition at line 94 of file AnasaziEpetraAdapter.cpp.

MultiVec< double > * Anasazi::EpetraMultiVec::CloneCopy ( ) const
virtual

Creates a new EpetraMultiVec and copies contents of *this into the new vector (deep copy).

Returns
Pointer to an EpetraMultiVec

Implements Anasazi::MultiVec< double >.

Definition at line 105 of file AnasaziEpetraAdapter.cpp.

MultiVec< double > * Anasazi::EpetraMultiVec::CloneCopy ( const std::vector< int > &  index) const
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.

Returns
Pointer to an EpetraMultiVec

Implements Anasazi::MultiVec< double >.

Definition at line 112 of file AnasaziEpetraAdapter.cpp.

MultiVec< double > * Anasazi::EpetraMultiVec::CloneViewNonConst ( const std::vector< int > &  index)
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.

Returns
Pointer to an EpetraMultiVec

Implements Anasazi::MultiVec< double >.

Definition at line 119 of file AnasaziEpetraAdapter.cpp.

const MultiVec< double > * Anasazi::EpetraMultiVec::CloneView ( const std::vector< int > &  index) const
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.

Returns
Pointer to an EpetraMultiVec

Implements Anasazi::MultiVec< double >.

Definition at line 125 of file AnasaziEpetraAdapter.cpp.

ptrdiff_t Anasazi::EpetraMultiVec::GetGlobalLength ( ) const
inlinevirtual

The number of rows in the multivector.

Implements Anasazi::MultiVec< double >.

Definition at line 215 of file AnasaziEpetraAdapter.hpp.

int Anasazi::EpetraMultiVec::GetNumberVecs ( ) const
inlinevirtual

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

Implements Anasazi::MultiVec< double >.

Definition at line 225 of file AnasaziEpetraAdapter.hpp.

void Anasazi::EpetraMultiVec::MvTimesMatAddMv ( double  alpha,
const MultiVec< double > &  A,
const Teuchos::SerialDenseMatrix< int, double > &  B,
double  beta 
)
virtual

Update *this with $\alpha AB + \beta (*this)$.

Implements Anasazi::MultiVec< double >.

Examples:
BlockKrylovSchur/BlockKrylovSchurEpetraExSVD.cpp.

Definition at line 158 of file AnasaziEpetraAdapter.cpp.

void Anasazi::EpetraMultiVec::MvAddMv ( double  alpha,
const MultiVec< double > &  A,
double  beta,
const MultiVec< double > &  B 
)
virtual

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

Implements Anasazi::MultiVec< double >.

Definition at line 178 of file AnasaziEpetraAdapter.cpp.

void Anasazi::EpetraMultiVec::MvTransMv ( double  alpha,
const MultiVec< double > &  A,
Teuchos::SerialDenseMatrix< int, double > &  B 
) const
virtual

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

Implements Anasazi::MultiVec< double >.

Definition at line 197 of file AnasaziEpetraAdapter.cpp.

void Anasazi::EpetraMultiVec::MvDot ( const MultiVec< double > &  A,
std::vector< double > &  b 
) const
virtual

Compute a vector b where the components are the individual dot-products, i.e. $ b[i] = A[i]^H(this[i])$ where A[i] is the i-th column of A.

Implements Anasazi::MultiVec< double >.

Definition at line 222 of file AnasaziEpetraAdapter.cpp.

void Anasazi::EpetraMultiVec::MvScale ( double  alpha)
inlinevirtual

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

Implements Anasazi::MultiVec< double >.

Definition at line 260 of file AnasaziEpetraAdapter.hpp.

void Anasazi::EpetraMultiVec::MvScale ( const std::vector< double > &  alpha)
virtual

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

Implements Anasazi::MultiVec< double >.

Definition at line 243 of file AnasaziEpetraAdapter.cpp.

void Anasazi::EpetraMultiVec::MvNorm ( std::vector< double > &  normvec) const
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 276 of file AnasaziEpetraAdapter.hpp.

void Anasazi::EpetraMultiVec::SetBlock ( const MultiVec< double > &  A,
const std::vector< int > &  index 
)
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 132 of file AnasaziEpetraAdapter.cpp.

void Anasazi::EpetraMultiVec::MvRandom ( )
inlinevirtual

Fill the vectors in *this with random numbers.

Implements Anasazi::MultiVec< double >.

Examples:
BlockKrylovSchur/BlockKrylovSchurEpetraExSVD.cpp.

Definition at line 295 of file AnasaziEpetraAdapter.hpp.

void Anasazi::EpetraMultiVec::MvInit ( double  alpha)
inlinevirtual

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

Implements Anasazi::MultiVec< double >.

Definition at line 302 of file AnasaziEpetraAdapter.hpp.

Epetra_MultiVector* Anasazi::EpetraMultiVec::GetEpetraMultiVec ( )
inlinevirtual

Return the pointer to the Epetra_MultiVector object.

Reimplemented from Anasazi::EpetraMultiVecAccessor.

Definition at line 311 of file AnasaziEpetraAdapter.hpp.

const Epetra_MultiVector* Anasazi::EpetraMultiVec::GetEpetraMultiVec ( ) const
inlinevirtual

Return the pointer to the Epetra_MultiVector object.

Reimplemented from Anasazi::EpetraMultiVecAccessor.

Definition at line 314 of file AnasaziEpetraAdapter.hpp.

void Anasazi::EpetraMultiVec::MvPrint ( std::ostream &  os) const
inlinevirtual

Print *this EpetraMultiVec.

Implements Anasazi::MultiVec< double >.

Definition at line 323 of file AnasaziEpetraAdapter.hpp.


The documentation for this class was generated from the following files: