Anasazi  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
Public Types | Public Member Functions | List of all members
Anasazi::TsqrMatOrthoManager< Scalar, MV, OP > Class Template Reference

MatOrthoManager subclass using TSQR or SVQB. More...

#include <AnasaziTsqrOrthoManager.hpp>

Inheritance diagram for Anasazi::TsqrMatOrthoManager< Scalar, MV, OP >:
Anasazi::MatOrthoManager< Scalar, MV, OP > Anasazi::OutOfPlaceNormalizerMixin< Scalar, MV > Teuchos::ParameterListAcceptorDefaultBase Anasazi::OrthoManager< ScalarType, MV > Teuchos::ParameterListAcceptor

Public Types

typedef MV multivector_type
 Multivector type with which this class was specialized. More...
 
typedef OP operator_type
 Operator type with which this class was specialized. More...
 
- Public Types inherited from Anasazi::OutOfPlaceNormalizerMixin< Scalar, MV >
typedef MV multivector_type
 Multivector type with which this class was specialized. More...
 

Public Member Functions

 TsqrMatOrthoManager (const Teuchos::RCP< Teuchos::ParameterList > &params, const std::string &label="Belos", Teuchos::RCP< const OP > Op=Teuchos::null)
 Constructor (that sets user-specified parameters). More...
 
 TsqrMatOrthoManager (const std::string &label="Belos", Teuchos::RCP< const OP > Op=Teuchos::null)
 Constructor (that sets default parameters). More...
 
virtual ~TsqrMatOrthoManager ()
 Destructor (declared virtual for memory safety of derived classes). More...
 
Teuchos::RCP< const
Teuchos::ParameterList
getValidParameters () const
 Get default parameters for TsqrMatOrthoManager. More...
 
Teuchos::RCP< const
Teuchos::ParameterList
getFastParameters ()
 Get "fast" parameters for TsqrMatOrthoManager. More...
 
virtual void setOp (Teuchos::RCP< const OP > Op)
 Set operator used for inner product. More...
 
Teuchos::RCP< const OP > getOp () const
 Get operator used for inner product. More...
 
void projectMat (MV &X, Teuchos::Array< Teuchos::RCP< const MV > > Q, Teuchos::Array< Teuchos::RCP< mat_type > > C=Teuchos::tuple(Teuchos::RCP< mat_type >(Teuchos::null)), Teuchos::RCP< MV > MX=Teuchos::null, Teuchos::Array< Teuchos::RCP< const MV > > MQ=Teuchos::tuple(Teuchos::null)) const
 Provides matrix-based projection method. More...
 
int projectAndNormalizeMat (MV &X, Teuchos::Array< Teuchos::RCP< const MV > > Q, Teuchos::Array< Teuchos::RCP< mat_type > > C=Teuchos::tuple(Teuchos::RCP< mat_type >(Teuchos::null)), Teuchos::RCP< mat_type > B=Teuchos::null, Teuchos::RCP< MV > MX=Teuchos::null, Teuchos::Array< Teuchos::RCP< const MV > > MQ=Teuchos::tuple(Teuchos::RCP< const MV >(Teuchos::null))) const
 Provides matrix-based projection/orthonormalization method. More...
 
int normalizeOutOfPlace (MV &X, MV &Q, mat_ptr B) const
 Normalize X into Q*B. More...
 
int projectAndNormalizeOutOfPlace (MV &X_in, MV &X_out, Teuchos::Array< mat_ptr > C, mat_ptr B, Teuchos::ArrayView< Teuchos::RCP< const MV > > Q) const
 Project and normalize X_in into X_out. More...
 
magnitude_type orthonormErrorMat (const MV &X, Teuchos::RCP< const MV > MX=Teuchos::null) const
 This method computes the error in orthonormality of a multivector. More...
 
magnitude_type orthogErrorMat (const MV &X, const MV &Y, Teuchos::RCP< const MV > MX=Teuchos::null, Teuchos::RCP< const MV > MY=Teuchos::null) const
 This method computes the error in orthogonality of two multivectors. More...
 
- Public Member Functions inherited from Anasazi::MatOrthoManager< Scalar, MV, OP >
 MatOrthoManager (Teuchos::RCP< const OP > Op=Teuchos::null)
 Default constructor. More...
 
virtual ~MatOrthoManager ()
 Destructor. More...
 
int getOpCounter () const
 Retrieve operator counter. More...
 
void resetOpCounter ()
 Reset the operator counter to zero. More...
 
void innerProdMat (const MV &X, const MV &Y, Teuchos::SerialDenseMatrix< int, Scalar > &Z, Teuchos::RCP< const MV > MX=Teuchos::null, Teuchos::RCP< const MV > MY=Teuchos::null) const
 Provides a matrix-based inner product. More...
 
void normMat (const MV &X, std::vector< typename Teuchos::ScalarTraits< Scalar >::magnitudeType > &normvec, Teuchos::RCP< const MV > MX=Teuchos::null) const
 Provides the norm induced by the matrix-based inner product. More...
 
virtual int normalizeMat (MV &X, Teuchos::RCP< Teuchos::SerialDenseMatrix< int, Scalar > > B=Teuchos::null, Teuchos::RCP< MV > MX=Teuchos::null) const =0
 Provides matrix-based orthonormalization method. More...
 
void innerProd (const MV &X, const MV &Y, Teuchos::SerialDenseMatrix< int, Scalar > &Z) const
 Implements the interface OrthoManager::innerProd(). More...
 
void norm (const MV &X, std::vector< typename Teuchos::ScalarTraits< Scalar >::magnitudeType > &normvec) const
 Implements the interface OrthoManager::norm(). More...
 
void project (MV &X, Teuchos::Array< Teuchos::RCP< const MV > > Q, Teuchos::Array< Teuchos::RCP< Teuchos::SerialDenseMatrix< int, Scalar > > > C=Teuchos::tuple(Teuchos::RCP< Teuchos::SerialDenseMatrix< int, Scalar > >(Teuchos::null))) const
 Implements the interface OrthoManager::project(). More...
 
int normalize (MV &X, Teuchos::RCP< Teuchos::SerialDenseMatrix< int, Scalar > > B=Teuchos::null) const
 Implements the interface OrthoManager::normalize(). More...
 
int projectAndNormalize (MV &X, Teuchos::Array< Teuchos::RCP< const MV > > Q, Teuchos::Array< Teuchos::RCP< Teuchos::SerialDenseMatrix< int, Scalar > > > C=Teuchos::tuple(Teuchos::RCP< Teuchos::SerialDenseMatrix< int, Scalar > >(Teuchos::null)), Teuchos::RCP< Teuchos::SerialDenseMatrix< int, Scalar > > B=Teuchos::null) const
 Implements the interface OrthoManager::projectAndNormalize(). More...
 
Teuchos::ScalarTraits< Scalar >
::magnitudeType 
orthonormError (const MV &X) const
 Implements the interface OrthoManager::orthonormError(). More...
 
Teuchos::ScalarTraits< Scalar >
::magnitudeType 
orthogError (const MV &X1, const MV &X2) const
 Implements the interface OrthoManager::orthogError(). More...
 
- Public Member Functions inherited from Anasazi::OrthoManager< ScalarType, MV >
 OrthoManager ()
 Default constructor. More...
 
virtual ~OrthoManager ()
 Destructor. More...
 
virtual void innerProd (const MV &X, const MV &Y, Teuchos::SerialDenseMatrix< int, ScalarType > &Z) const =0
 Provides the inner product defining the orthogonality concepts. More...
 
virtual void norm (const MV &X, std::vector< typename Teuchos::ScalarTraits< ScalarType >::magnitudeType > &normvec) const =0
 Provides the norm induced by innerProd(). More...
 
virtual void project (MV &X, Teuchos::Array< Teuchos::RCP< const MV > > Q, Teuchos::Array< Teuchos::RCP< Teuchos::SerialDenseMatrix< int, ScalarType > > > C=Teuchos::tuple(Teuchos::RCP< Teuchos::SerialDenseMatrix< int, ScalarType > >(Teuchos::null))) const =0
 Given a list of mutually orthogonal and internally orthonormal bases Q, this method projects a multivector X onto the space orthogonal to the individual Q[i], optionally returning the coefficients of X for the individual Q[i]. All of this is done with respect to the inner product innerProd(). More...
 
virtual int normalize (MV &X, Teuchos::RCP< Teuchos::SerialDenseMatrix< int, ScalarType > > B=Teuchos::null) const =0
 This method takes a multivector X and attempts to compute a basis for $ colspan(X)$. This basis is orthonormal with respect to innerProd(). More...
 
virtual int projectAndNormalize (MV &X, Teuchos::Array< Teuchos::RCP< const MV > > Q, Teuchos::Array< Teuchos::RCP< Teuchos::SerialDenseMatrix< int, ScalarType > > > C=Teuchos::tuple(Teuchos::RCP< Teuchos::SerialDenseMatrix< int, ScalarType > >(Teuchos::null)), Teuchos::RCP< Teuchos::SerialDenseMatrix< int, ScalarType > > B=Teuchos::null) const =0
 Given a set of bases Q[i] and a multivector X, this method computes an orthonormal basis for $ colspan(X) - \sum_i colspan(Q[i])$. More...
 
- Public Member Functions inherited from Anasazi::OutOfPlaceNormalizerMixin< Scalar, MV >
virtual ~OutOfPlaceNormalizerMixin ()
 Trivial virtual destructor, to silence compiler warnings. More...
 

Detailed Description

template<class Scalar, class MV, class OP>
class Anasazi::TsqrMatOrthoManager< Scalar, MV, OP >

MatOrthoManager subclass using TSQR or SVQB.

When the inner product matrix has not been set, this class uses TSQR + Block Gram-Schmidt (via TsqrOrthoManagerImpl). If the inner product matrix has been set, then this class uses the SVQB algorithm (Stathopoulos and Wu 2002: CholeskyQR + SVD) for orthogonalization.

TSQR uses multivector scratch space. However, scratch space initialization is "lazy," so scratch space will not be allocated if TSQR is not used.

Definition at line 332 of file AnasaziTsqrOrthoManager.hpp.

Member Typedef Documentation

template<class Scalar , class MV , class OP >
typedef MV Anasazi::TsqrMatOrthoManager< Scalar, MV, OP >::multivector_type

Multivector type with which this class was specialized.

Definition at line 341 of file AnasaziTsqrOrthoManager.hpp.

template<class Scalar , class MV , class OP >
typedef OP Anasazi::TsqrMatOrthoManager< Scalar, MV, OP >::operator_type

Operator type with which this class was specialized.

Definition at line 343 of file AnasaziTsqrOrthoManager.hpp.

Constructor & Destructor Documentation

template<class Scalar , class MV , class OP >
Anasazi::TsqrMatOrthoManager< Scalar, MV, OP >::TsqrMatOrthoManager ( const Teuchos::RCP< Teuchos::ParameterList > &  params,
const std::string &  label = "Belos",
Teuchos::RCP< const OP >  Op = Teuchos::null 
)
inline

Constructor (that sets user-specified parameters).

Parameters
params[in/out] Configuration parameters, both for this orthogonalization manager, and for TSQR itself (as the "TSQR implementation" sublist). This can be null, in which case default parameters will be set for now; you can always call setParameterList() later to change these.
label[in] Label for timers. This only matters if the compile-time option for enabling timers is set.
Op[in] Inner product with respect to which to orthogonalize vectors. If Teuchos::null, use the Euclidean inner product.

Call getValidParameters() for default parameters and their documentation, including TSQR implementation parameters. Call getFastParameters() to get documented parameters for faster computation, possibly at the expense of accuracy and robustness.

Definition at line 392 of file AnasaziTsqrOrthoManager.hpp.

template<class Scalar , class MV , class OP >
Anasazi::TsqrMatOrthoManager< Scalar, MV, OP >::TsqrMatOrthoManager ( const std::string &  label = "Belos",
Teuchos::RCP< const OP >  Op = Teuchos::null 
)
inline

Constructor (that sets default parameters).

Parameters
Op[in] Inner product with respect to which to orthogonalize vectors. If Teuchos::null, use the Euclidean inner product.
label[in] Label for timers. This only matters if the compile-time option for enabling timers is set.

Definition at line 408 of file AnasaziTsqrOrthoManager.hpp.

template<class Scalar , class MV , class OP >
virtual Anasazi::TsqrMatOrthoManager< Scalar, MV, OP >::~TsqrMatOrthoManager ( )
inlinevirtual

Destructor (declared virtual for memory safety of derived classes).

Definition at line 416 of file AnasaziTsqrOrthoManager.hpp.

Member Function Documentation

template<class Scalar , class MV , class OP >
Teuchos::RCP<const Teuchos::ParameterList> Anasazi::TsqrMatOrthoManager< Scalar, MV, OP >::getValidParameters ( ) const
inlinevirtual

Get default parameters for TsqrMatOrthoManager.

Get a (pointer to a) default list of parameters for configuring a TsqrMatOrthoManager instance.

Note
TSQR implementation configuration options are stored under "TSQR implementation" as a sublist.

Reimplemented from Teuchos::ParameterListAcceptor.

Definition at line 425 of file AnasaziTsqrOrthoManager.hpp.

template<class Scalar , class MV , class OP >
Teuchos::RCP<const Teuchos::ParameterList> Anasazi::TsqrMatOrthoManager< Scalar, MV, OP >::getFastParameters ( )
inline

Get "fast" parameters for TsqrMatOrthoManager.

Get a (pointer to a) list of parameters for configuring a TsqrMatOrthoManager instance for maximum speed, at the cost of accuracy (no block reorthogonalization) and robustness to rank deficiency (no randomization of the null space basis).

Note
TSQR implementation configuration options are stored under "TSQR implementation" as a sublist.

Definition at line 438 of file AnasaziTsqrOrthoManager.hpp.

template<class Scalar , class MV , class OP >
virtual void Anasazi::TsqrMatOrthoManager< Scalar, MV, OP >::setOp ( Teuchos::RCP< const OP >  Op)
inlinevirtual

Set operator used for inner product.

Reimplemented from Anasazi::MatOrthoManager< Scalar, MV, OP >.

Definition at line 447 of file AnasaziTsqrOrthoManager.hpp.

template<class Scalar , class MV , class OP >
Teuchos::RCP<const OP> Anasazi::TsqrMatOrthoManager< Scalar, MV, OP >::getOp ( ) const
inlinevirtual

Get operator used for inner product.

Reimplemented from Anasazi::MatOrthoManager< Scalar, MV, OP >.

Definition at line 462 of file AnasaziTsqrOrthoManager.hpp.

template<class Scalar , class MV , class OP >
void Anasazi::TsqrMatOrthoManager< Scalar, MV, OP >::projectMat ( MV &  X,
Teuchos::Array< Teuchos::RCP< const MV > >  Q,
Teuchos::Array< Teuchos::RCP< mat_type > >  C = Teuchos::tuple (Teuchos::RCP<mat_type> (Teuchos::null)),
Teuchos::RCP< MV >  MX = Teuchos::null,
Teuchos::Array< Teuchos::RCP< const MV > >  MQ = Teuchos::tuple (Teuchos::null) 
) const
inlinevirtual

Provides matrix-based projection method.

This method optionally allows the provision of $ M X$ and/or the $ M Q[i]$. See OrthoManager::project() for more details.

Parameters
X,Q,C[in/out] As in OrthoManager::project()
MX[in/out] If specified by the user, on input MX is required to be the image of X under the operator getOp(). On output, MX will be updated to reflect the changes in X.
MQ[in] If specified by the user, on MQ[i] is required to be the image of Q[i] under the operator getOp().

Implements Anasazi::MatOrthoManager< Scalar, MV, OP >.

Definition at line 469 of file AnasaziTsqrOrthoManager.hpp.

template<class Scalar , class MV , class OP >
int Anasazi::TsqrMatOrthoManager< Scalar, MV, OP >::projectAndNormalizeMat ( MV &  X,
Teuchos::Array< Teuchos::RCP< const MV > >  Q,
Teuchos::Array< Teuchos::RCP< mat_type > >  C = Teuchos::tuple (Teuchos::RCP<mat_type> (Teuchos::null)),
Teuchos::RCP< mat_type B = Teuchos::null,
Teuchos::RCP< MV >  MX = Teuchos::null,
Teuchos::Array< Teuchos::RCP< const MV > >  MQ = Teuchos::tuple (Teuchos::RCP<const MV> (Teuchos::null)) 
) const
inlinevirtual

Provides matrix-based projection/orthonormalization method.

This method optionally allows the provision of $ M X$ and/or the $ M Q[i]$. See orthoManager::projectAndNormalize() for more details.

Parameters
X,Q,C,B[in/out] As in OrthoManager::projectAndNormalize()
MX[in/out] If specified by the user, on input MX is required to be the image of X under the operator getOp(). On output, MX will be updated to reflect the changes in X.
MQ[in] If specified by the user, on MQ[i] is required to be the image of Q[i] under the operator getOp().
Returns
Rank of the basis computed by this method, less than or equal to the number of columns in X. This specifies how many columns in the returned X and MX and rows in the returned B are valid.

Implements Anasazi::MatOrthoManager< Scalar, MV, OP >.

Definition at line 522 of file AnasaziTsqrOrthoManager.hpp.

template<class Scalar , class MV , class OP >
int Anasazi::TsqrMatOrthoManager< Scalar, MV, OP >::normalizeOutOfPlace ( MV &  X,
MV &  Q,
mat_ptr  B 
) const
inlinevirtual

Normalize X into Q*B.

Parameters
X[in/out] On input: Multivector to normalize. On output: Possibly overwritten with invalid values.
Q[out] On output: Normalized multivector.
B[out] On output: Normalization coefficients.
Returns
Rank of the input multivector X.

Implements Anasazi::OutOfPlaceNormalizerMixin< Scalar, MV >.

Definition at line 553 of file AnasaziTsqrOrthoManager.hpp.

template<class Scalar , class MV , class OP >
int Anasazi::TsqrMatOrthoManager< Scalar, MV, OP >::projectAndNormalizeOutOfPlace ( MV &  X_in,
MV &  X_out,
Teuchos::Array< mat_ptr C,
mat_ptr  B,
Teuchos::ArrayView< Teuchos::RCP< const MV > >  Q 
) const
inlinevirtual

Project and normalize X_in into X_out.

Project X_in against Q, storing projection coefficients in C, and normalize X_in into X_out, storing normalization coefficients in B. On output, X_out has the resulting orthogonal vectors. X_in may be overwritten with invalid values.

Parameters
X_in[in/out] On input: The vectors to project against Q and normalize. On output: possibly overwritten with invalid values.
X_out[out] The normalized input vectors after projection against Q.
C[out] Projection coefficients
B[out] Normalization coefficients
Q[in] The orthogonal basis against which to project
Returns
Rank of X_in after projection

Implements Anasazi::OutOfPlaceNormalizerMixin< Scalar, MV >.

Definition at line 567 of file AnasaziTsqrOrthoManager.hpp.

template<class Scalar , class MV , class OP >
magnitude_type Anasazi::TsqrMatOrthoManager< Scalar, MV, OP >::orthonormErrorMat ( const MV &  X,
Teuchos::RCP< const MV >  MX = Teuchos::null 
) const
inlinevirtual

This method computes the error in orthonormality of a multivector.

This method optionally allows optionally exploits a caller-provided MX.

Implements Anasazi::MatOrthoManager< Scalar, MV, OP >.

Definition at line 591 of file AnasaziTsqrOrthoManager.hpp.

template<class Scalar , class MV , class OP >
magnitude_type Anasazi::TsqrMatOrthoManager< Scalar, MV, OP >::orthogErrorMat ( const MV &  X,
const MV &  Y,
Teuchos::RCP< const MV >  MX = Teuchos::null,
Teuchos::RCP< const MV >  MY = Teuchos::null 
) const
inlinevirtual

This method computes the error in orthogonality of two multivectors.

This method optionally allows optionally exploits a caller-provided MX and/or MY.

Implements Anasazi::MatOrthoManager< Scalar, MV, OP >.

Definition at line 603 of file AnasaziTsqrOrthoManager.hpp.


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