Anasazi  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
List of all members
Anasazi::MultiVecTraits< double, Epetra_MultiVector > Class Template Reference

Template specialization of Anasazi::MultiVecTraits class using the Epetra_MultiVector class. More...

#include <AnasaziEpetraAdapter.hpp>

Static Public Member Functions

Creation methods
static Teuchos::RCP
< Epetra_MultiVector
Clone (const Epetra_MultiVector &mv, const int outNumVecs)
 Creates a new empty Epetra_MultiVector containing numVecs columns. More...
 
static Teuchos::RCP
< Epetra_MultiVector
CloneCopy (const Epetra_MultiVector &mv)
 Creates a new Epetra_MultiVector and copies contents of mv into the new vector (deep copy). More...
 
static Teuchos::RCP
< Epetra_MultiVector
CloneCopy (const Epetra_MultiVector &mv, const std::vector< int > &index)
 Creates a new Epetra_MultiVector and copies the selected contents of mv into the new vector (deep copy). More...
 
static Teuchos::RCP
< Epetra_MultiVector
CloneCopy (const Epetra_MultiVector &mv, const Teuchos::Range1D &index)
 
static Teuchos::RCP
< Epetra_MultiVector
CloneViewNonConst (Epetra_MultiVector &mv, const std::vector< int > &index)
 Creates a new Epetra_MultiVector that shares the selected contents of mv (shallow copy). More...
 
static Teuchos::RCP
< Epetra_MultiVector
CloneViewNonConst (Epetra_MultiVector &mv, const Teuchos::Range1D &index)
 
static Teuchos::RCP< const
Epetra_MultiVector
CloneView (const Epetra_MultiVector &mv, const std::vector< int > &index)
 Creates a new const Epetra_MultiVector that shares the selected contents of mv (shallow copy). More...
 
static Teuchos::RCP
< Epetra_MultiVector
CloneView (const Epetra_MultiVector &mv, const Teuchos::Range1D &index)
 
Attribute methods
static ptrdiff_t GetGlobalLength (const Epetra_MultiVector &mv)
 Obtain the vector length of mv. More...
 
static int GetNumberVecs (const Epetra_MultiVector &mv)
 Obtain the number of vectors in mv. More...
 
static bool HasConstantStride (const Epetra_MultiVector &mv)
 
Update methods
static void MvTimesMatAddMv (double alpha, const Epetra_MultiVector &A, const Teuchos::SerialDenseMatrix< int, double > &B, double beta, Epetra_MultiVector &mv)
 Update mv with $ \alpha AB + \beta mv $. More...
 
static void MvAddMv (double alpha, const Epetra_MultiVector &A, double beta, const Epetra_MultiVector &B, Epetra_MultiVector &mv)
 Replace mv with $\alpha A + \beta B$. More...
 
static void MvTransMv (double alpha, const Epetra_MultiVector &A, const Epetra_MultiVector &mv, Teuchos::SerialDenseMatrix< int, double > &B)
 Compute a dense matrix B through the matrix-matrix multiply $ \alpha A^Tmv $. More...
 
static void MvDot (const Epetra_MultiVector &A, const Epetra_MultiVector &B, std::vector< double > &b)
 Compute a vector b where the components are the individual dot-products of the i-th columns of A and mv, i.e. $ b[i] = A[i]^Tmv[i]$. More...
 
Norm method
static void MvNorm (const Epetra_MultiVector &mv, std::vector< double > &normvec)
 Compute the 2-norm of each individual vector of mv. Upon return, normvec[i] holds the value of $||mv_i||_2$, the i-th column of mv. More...
 
Initialization methods
static void SetBlock (const Epetra_MultiVector &A, const std::vector< int > &index, Epetra_MultiVector &mv)
 Copy the vectors in A to a set of vectors in mv indicated by the indices given in index. More...
 
static void SetBlock (const Epetra_MultiVector &A, const Teuchos::Range1D &index, Epetra_MultiVector &mv)
 
static void Assign (const Epetra_MultiVector &A, Epetra_MultiVector &mv)
 
static void MvScale (Epetra_MultiVector &mv, double alpha)
 Scale each element of the vectors in mv with alpha. More...
 
static void MvScale (Epetra_MultiVector &mv, const std::vector< double > &alpha)
 Scale each element of the i-th vector in mv with alpha[i]. More...
 
static void MvRandom (Epetra_MultiVector &mv)
 Replace the vectors in mv with random vectors. More...
 
static void MvInit (Epetra_MultiVector &mv, double alpha=Teuchos::ScalarTraits< double >::zero())
 Replace each element of the vectors in mv with alpha. More...
 
Print method
static void MvPrint (const Epetra_MultiVector &mv, std::ostream &os)
 Print the mv multi-vector to the os output stream. More...
 

Detailed Description

template<>
class Anasazi::MultiVecTraits< double, Epetra_MultiVector >

Template specialization of Anasazi::MultiVecTraits class using the Epetra_MultiVector class.

This interface will ensure that any Epetra_MultiVector will be accepted by the Anasazi templated solvers.

Note
The Epetra package performs double-precision arithmetic, so the use of Epetra with Anasazi will only provide a double-precision eigensolver.
Examples:
BlockDavidson/BlockDavidsonEpetraEx.cpp, BlockDavidson/BlockDavidsonEpetraExGen.cpp, BlockDavidson/BlockDavidsonEpetraExGenPrecIfpack.cpp, LOBPCGEpetra.cpp, LOBPCGEpetraEx.cpp, LOBPCGEpetraExGen.cpp, LOBPCGEpetraExGenPrecIfpack.cpp, and LOBPCGEpetraExSimple.cpp.

Definition at line 716 of file AnasaziEpetraAdapter.hpp.

Member Function Documentation

static Teuchos::RCP<Epetra_MultiVector> Anasazi::MultiVecTraits< double, Epetra_MultiVector >::Clone ( const Epetra_MultiVector mv,
const int  outNumVecs 
)
inlinestatic

Creates a new empty Epetra_MultiVector containing numVecs columns.

Returns
Reference-counted pointer to the new Epetra_MultiVector.

Definition at line 728 of file AnasaziEpetraAdapter.hpp.

static Teuchos::RCP<Epetra_MultiVector> Anasazi::MultiVecTraits< double, Epetra_MultiVector >::CloneCopy ( const Epetra_MultiVector mv)
inlinestatic

Creates a new Epetra_MultiVector and copies contents of mv into the new vector (deep copy).

Returns
Reference-counted pointer to the new Epetra_MultiVector.

Definition at line 746 of file AnasaziEpetraAdapter.hpp.

static Teuchos::RCP<Epetra_MultiVector> Anasazi::MultiVecTraits< double, Epetra_MultiVector >::CloneCopy ( const Epetra_MultiVector mv,
const std::vector< int > &  index 
)
inlinestatic

Creates a new Epetra_MultiVector and copies the selected contents of mv into the new vector (deep copy).

The copied vectors from mv are indicated by the indeX.size() indices in index.

Returns
Reference-counted pointer to the new Epetra_MultiVector.

Definition at line 757 of file AnasaziEpetraAdapter.hpp.

static Teuchos::RCP<Epetra_MultiVector> Anasazi::MultiVecTraits< double, Epetra_MultiVector >::CloneViewNonConst ( Epetra_MultiVector mv,
const std::vector< int > &  index 
)
inlinestatic

Creates a new Epetra_MultiVector that shares the selected contents of mv (shallow copy).

The index of the numvecs vectors shallow copied from mv are indicated by the indices given in index.

Returns
Reference-counted pointer to the new Epetra_MultiVector.

Definition at line 847 of file AnasaziEpetraAdapter.hpp.

static Teuchos::RCP<const Epetra_MultiVector> Anasazi::MultiVecTraits< double, Epetra_MultiVector >::CloneView ( const Epetra_MultiVector mv,
const std::vector< int > &  index 
)
inlinestatic

Creates a new const Epetra_MultiVector that shares the selected contents of mv (shallow copy).

The index of the numvecs vectors shallow copied from mv are indicated by the indices given in index.

Returns
Reference-counted pointer to the new const Epetra_MultiVector.

Definition at line 938 of file AnasaziEpetraAdapter.hpp.

static ptrdiff_t Anasazi::MultiVecTraits< double, Epetra_MultiVector >::GetGlobalLength ( const Epetra_MultiVector mv)
inlinestatic

Obtain the vector length of mv.

Definition at line 1029 of file AnasaziEpetraAdapter.hpp.

static int Anasazi::MultiVecTraits< double, Epetra_MultiVector >::GetNumberVecs ( const Epetra_MultiVector mv)
inlinestatic

Obtain the number of vectors in mv.

Definition at line 1038 of file AnasaziEpetraAdapter.hpp.

static void Anasazi::MultiVecTraits< double, Epetra_MultiVector >::MvTimesMatAddMv ( double  alpha,
const Epetra_MultiVector A,
const Teuchos::SerialDenseMatrix< int, double > &  B,
double  beta,
Epetra_MultiVector mv 
)
inlinestatic

Update mv with $ \alpha AB + \beta mv $.

Definition at line 1050 of file AnasaziEpetraAdapter.hpp.

static void Anasazi::MultiVecTraits< double, Epetra_MultiVector >::MvAddMv ( double  alpha,
const Epetra_MultiVector A,
double  beta,
const Epetra_MultiVector B,
Epetra_MultiVector mv 
)
inlinestatic

Replace mv with $\alpha A + \beta B$.

Definition at line 1063 of file AnasaziEpetraAdapter.hpp.

static void Anasazi::MultiVecTraits< double, Epetra_MultiVector >::MvTransMv ( double  alpha,
const Epetra_MultiVector A,
const Epetra_MultiVector mv,
Teuchos::SerialDenseMatrix< int, double > &  B 
)
inlinestatic

Compute a dense matrix B through the matrix-matrix multiply $ \alpha A^Tmv $.

Definition at line 1120 of file AnasaziEpetraAdapter.hpp.

static void Anasazi::MultiVecTraits< double, Epetra_MultiVector >::MvDot ( const Epetra_MultiVector A,
const Epetra_MultiVector B,
std::vector< double > &  b 
)
inlinestatic

Compute a vector b where the components are the individual dot-products of the i-th columns of A and mv, i.e. $ b[i] = A[i]^Tmv[i]$.

Definition at line 1135 of file AnasaziEpetraAdapter.hpp.

static void Anasazi::MultiVecTraits< double, Epetra_MultiVector >::MvNorm ( const Epetra_MultiVector mv,
std::vector< double > &  normvec 
)
inlinestatic

Compute the 2-norm of each individual vector of mv. Upon return, normvec[i] holds the value of $||mv_i||_2$, the i-th column of mv.

Definition at line 1158 of file AnasaziEpetraAdapter.hpp.

static void Anasazi::MultiVecTraits< double, Epetra_MultiVector >::SetBlock ( const Epetra_MultiVector A,
const std::vector< int > &  index,
Epetra_MultiVector mv 
)
inlinestatic

Copy the vectors in A to a set of vectors in mv indicated by the indices given in index.

Definition at line 1175 of file AnasaziEpetraAdapter.hpp.

static void Anasazi::MultiVecTraits< double, Epetra_MultiVector >::MvScale ( Epetra_MultiVector mv,
double  alpha 
)
inlinestatic

Scale each element of the vectors in mv with alpha.

Definition at line 1321 of file AnasaziEpetraAdapter.hpp.

static void Anasazi::MultiVecTraits< double, Epetra_MultiVector >::MvScale ( Epetra_MultiVector mv,
const std::vector< double > &  alpha 
)
inlinestatic

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

Definition at line 1329 of file AnasaziEpetraAdapter.hpp.

static void Anasazi::MultiVecTraits< double, Epetra_MultiVector >::MvRandom ( Epetra_MultiVector mv)
inlinestatic

Replace the vectors in mv with random vectors.

Definition at line 1345 of file AnasaziEpetraAdapter.hpp.

static void Anasazi::MultiVecTraits< double, Epetra_MultiVector >::MvInit ( Epetra_MultiVector mv,
double  alpha = Teuchos::ScalarTraits<double>::zero() 
)
inlinestatic

Replace each element of the vectors in mv with alpha.

Definition at line 1353 of file AnasaziEpetraAdapter.hpp.

static void Anasazi::MultiVecTraits< double, Epetra_MultiVector >::MvPrint ( const Epetra_MultiVector mv,
std::ostream &  os 
)
inlinestatic

Print the mv multi-vector to the os output stream.

Definition at line 1366 of file AnasaziEpetraAdapter.hpp.


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