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::TsqrOrthoManagerImpl< Scalar, MV > Class Template Reference

TSQR-based OrthoManager subclass implementation. More...

#include <AnasaziTsqrOrthoManagerImpl.hpp>

Inheritance diagram for Anasazi::TsqrOrthoManagerImpl< Scalar, MV >:
Teuchos::ParameterListAcceptorDefaultBase Teuchos::ParameterListAcceptor

Public Types

typedef
Teuchos::SerialDenseMatrix
< int, Scalar > 
mat_type
 

Public Member Functions

Teuchos::RCP< const
Teuchos::ParameterList
getValidParameters () const
 Default valid parameter list. More...
 
void setParameterList (const Teuchos::RCP< Teuchos::ParameterList > &params)
 Set parameters from the given parameter list. More...
 
Teuchos::RCP< const
Teuchos::ParameterList
getFastParameters ()
 Get "fast" parameters for TsqrOrthoManagerImpl. More...
 
 TsqrOrthoManagerImpl (const Teuchos::RCP< Teuchos::ParameterList > &params, const std::string &label)
 Constructor (that sets user-specified parameters). More...
 
 TsqrOrthoManagerImpl (const std::string &label)
 Constructor (that sets default parameters). More...
 
void setLabel (const std::string &label)
 Set the label for timers. More...
 
const std::string & getLabel () const
 Get the label for timers (if timers are enabled). More...
 
void innerProd (const MV &X, const MV &Y, mat_type &Z) const
 Euclidean inner product. More...
 
void norm (const MV &X, std::vector< magnitude_type > &normvec) const
 
void project (MV &X, Teuchos::Array< mat_ptr > C, Teuchos::ArrayView< Teuchos::RCP< const MV > > Q)
 Compute $C := Q^* X$ and $X := X - Q C$. More...
 
int normalize (MV &X, mat_ptr B)
 Orthogonalize the columns of X in place. More...
 
int normalizeOutOfPlace (MV &X, MV &Q, mat_ptr B)
 Normalize X into Q*B, overwriting X. More...
 
int projectAndNormalize (MV &X, Teuchos::Array< mat_ptr > C, mat_ptr B, Teuchos::ArrayView< Teuchos::RCP< const MV > > Q)
 Project X against Q and normalize X. More...
 
int projectAndNormalizeOutOfPlace (MV &X_in, MV &X_out, Teuchos::Array< mat_ptr > C, mat_ptr B, Teuchos::ArrayView< Teuchos::RCP< const MV > > Q)
 Project and normalize X_in into X_out; overwrite X_in. More...
 
magnitude_type orthonormError (const MV &X) const
 Return $ \| I - X^* \cdot X \|_F $. More...
 
magnitude_type orthogError (const MV &X1, const MV &X2) const
 Return the Frobenius norm of the inner product of X1 with itself. More...
 
magnitude_type blockReorthogThreshold () const
 
magnitude_type relativeRankTolerance () const
 

Detailed Description

template<class Scalar, class MV>
class Anasazi::TsqrOrthoManagerImpl< Scalar, MV >

TSQR-based OrthoManager subclass implementation.

Author
Mark Hoemmen

TsqrOrthoManagerImpl implements the interface defined by OrthoManager, as well as the interface defined by OutOfPlaceNormalizerMixin. We use TsqrOrthoManagerImpl to implement TsqrOrthoManager and TsqrMatOrthoManager.

Template Parameters
ScalarThe type of matrix and (multi)vector entries.
MVThe type of (multi)vector inputs and outputs.

This class uses a combination of Tall Skinny QR (TSQR) and Block Gram-Schmidt (BGS) to orthogonalize multivectors. The Block Gram-Schmidt procedure used here is inspired by that of G. W. Stewart ("Block Gram-Schmidt Orthogonalization", SISC vol 31 #1 pp. 761–775, 2008). The difference is that we use TSQR+SVD instead of Stewart's careful Gram-Schmidt with reorthogonalization to handle the current block. "Orthogonalization faults" (as defined by Stewart) may still happen, but we do not handle them by default. Rather, we make one BGS pass, do TSQR+SVD, check the resulting column norms, and make a second BGS pass (+ TSQR+SVD) if necessary. If we then detect an orthogonalization fault, we throw TsqrOrthoFault.

Note
Despite the "Impl" part of the name of this class, we don't actually use it for the "pImpl" C++ idiom. We just separate out the TSQR implementation to make it easier to implement the OrthoManager and MatOrthoManager interfaces for the case where the inner product operator is not the identity matrix.

Definition at line 128 of file AnasaziTsqrOrthoManagerImpl.hpp.

Member Typedef Documentation

template<class Scalar , class MV >
Anasazi::TsqrOrthoManagerImpl< Scalar, MV >::mat_type

Type of the projection and normalization coefficients

Definition at line 136 of file AnasaziTsqrOrthoManagerImpl.hpp.

Constructor & Destructor Documentation

template<class Scalar , class MV >
Anasazi::TsqrOrthoManagerImpl< Scalar, MV >::TsqrOrthoManagerImpl ( const Teuchos::RCP< Teuchos::ParameterList > &  params,
const std::string &  label 
)

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.

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 656 of file AnasaziTsqrOrthoManagerImpl.hpp.

template<class Scalar , class MV >
Anasazi::TsqrOrthoManagerImpl< Scalar, MV >::TsqrOrthoManagerImpl ( const std::string &  label)

Constructor (that sets default parameters).

Parameters
label[in] Label for timers. This only matters if the compile-time option for enabling timers is set.

Definition at line 673 of file AnasaziTsqrOrthoManagerImpl.hpp.

Member Function Documentation

template<class Scalar , class MV >
Teuchos::RCP< const Teuchos::ParameterList > Anasazi::TsqrOrthoManagerImpl< Scalar, MV >::getValidParameters ( ) const
virtual

Default valid parameter list.

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

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

Reimplemented from Teuchos::ParameterListAcceptor.

Definition at line 1247 of file AnasaziTsqrOrthoManagerImpl.hpp.

template<class Scalar , class MV >
void Anasazi::TsqrOrthoManagerImpl< Scalar, MV >::setParameterList ( const Teuchos::RCP< Teuchos::ParameterList > &  params)
virtual

Set parameters from the given parameter list.

Implements Teuchos::ParameterListAcceptor.

Definition at line 597 of file AnasaziTsqrOrthoManagerImpl.hpp.

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

Get "fast" parameters for TsqrOrthoManagerImpl.

Get a (pointer to a) list of parameters for configuring a TsqrOrthoManager or 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 1312 of file AnasaziTsqrOrthoManagerImpl.hpp.

template<class Scalar , class MV >
void Anasazi::TsqrOrthoManagerImpl< Scalar, MV >::setLabel ( const std::string &  label)
inline

Set the label for timers.

This only matters if timers are enabled. If timers are enabled and the label changes, this method will clear the old timers and replace them with new ones. The old timers will not appear in the list of timers shown by Teuchos::TimeMonitor::summarize().

Definition at line 202 of file AnasaziTsqrOrthoManagerImpl.hpp.

template<class Scalar , class MV >
const std::string& Anasazi::TsqrOrthoManagerImpl< Scalar, MV >::getLabel ( ) const
inline

Get the label for timers (if timers are enabled).

Definition at line 209 of file AnasaziTsqrOrthoManagerImpl.hpp.

template<class Scalar , class MV >
void Anasazi::TsqrOrthoManagerImpl< Scalar, MV >::innerProd ( const MV &  X,
const MV &  Y,
mat_type Z 
) const
inline

Euclidean inner product.

Compute the Euclidean block inner product X^* Y, and store the result in Z.

Parameters
X[in]
Y[in]
Z[out] On output, $X^* Y$

Definition at line 220 of file AnasaziTsqrOrthoManagerImpl.hpp.

template<class Scalar , class MV >
void Anasazi::TsqrOrthoManagerImpl< Scalar, MV >::norm ( const MV &  X,
std::vector< magnitude_type > &  normvec 
) const

Compute the 2-norm of each column j of X.

Parameters
X[in] Multivector for which to compute column norms.
normVec[out] On output: normvec[j] is the 2-norm of column j of X. normVec is resized if necessary so that it has at least as many entries as there are columns of X.
Note
Performance of this method depends on how MultiVecTraits implements column norm computation for the given multivector type MV. It may or may not be the case that a reduction is performed for every column of X. Furthermore, whether or not the columns of X are contiguous (as opposed to a view of noncontiguous columns) may also affect performance. The computed results should be the same regardless, except perhaps for small rounding differences due to a different order of operations.

Definition at line 690 of file AnasaziTsqrOrthoManagerImpl.hpp.

template<class Scalar , class MV >
void Anasazi::TsqrOrthoManagerImpl< Scalar, MV >::project ( MV &  X,
Teuchos::Array< mat_ptr C,
Teuchos::ArrayView< Teuchos::RCP< const MV > >  Q 
)

Compute $C := Q^* X$ and $X := X - Q C$.

Project X against the span of the (Euclidean) orthogonal vectors Q, and store the resulting coefficients in C.

Parameters
X[in/out] On input: the vectors to project. On output: $X := X - Q C$ where $C := Q^* X$.
C[out] The projection coefficients $C := Q^* X$
Q[in] The orthogonal basis against which to project

Definition at line 702 of file AnasaziTsqrOrthoManagerImpl.hpp.

template<class Scalar , class MV >
int Anasazi::TsqrOrthoManagerImpl< Scalar, MV >::normalize ( MV &  X,
mat_ptr  B 
)

Orthogonalize the columns of X in place.

Orthogonalize the columns of X in place, storing the resulting coefficients in B. Return the rank of X. If X is full rank, then X*B on output is a QR factorization of X on input. If X is not full rank, then the first rank columns of X on output form a basis for the column space of X (on input). Additional options control randomization of the null space basis.

Parameters
X[in/out]
B[out]
Returns
Rank of X

Definition at line 768 of file AnasaziTsqrOrthoManagerImpl.hpp.

template<class Scalar , class MV >
int Anasazi::TsqrOrthoManagerImpl< Scalar, MV >::normalizeOutOfPlace ( MV &  X,
MV &  Q,
mat_ptr  B 
)

Normalize X into Q*B, overwriting X.

Normalize X into Q*B, overwriting X with invalid values.

Parameters
X[in/out] Vector(s) to normalize
Q[out] Normalized vector(s)
B[out] Normalization coefficients
Returns
Rank of X
Note
Q must have at least as many columns as X. It may have more columns than X; those columns are ignored.
We expose this interface to applications because TSQR is not able to compute an orthogonal basis in place; it needs scratch space. Applications can exploit this interface to avoid excessive copying of vectors when using TSQR for orthogonalization.

Definition at line 877 of file AnasaziTsqrOrthoManagerImpl.hpp.

template<class Scalar , class MV >
int Anasazi::TsqrOrthoManagerImpl< Scalar, MV >::projectAndNormalize ( MV &  X,
Teuchos::Array< mat_ptr C,
mat_ptr  B,
Teuchos::ArrayView< Teuchos::RCP< const MV > >  Q 
)
inline

Project X against Q and normalize X.

This method is equivalent (in exact arithmetic) to project(X,C,Q) followed by normalize(X,B). However, the interface allows this method to implement reorthogonalization more efficiently and accurately.

Parameters
X[in/out] The vectors to project against Q and normalize
C[out] The projection coefficients
B[out] The normalization coefficients
Q[in] The orthogonal basis against which to project
Returns
Rank of X after projection

Definition at line 309 of file AnasaziTsqrOrthoManagerImpl.hpp.

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

Project and normalize X_in into X_out; overwrite X_in.

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 and X_in is overwritten with invalid values.

Parameters
X_in[in/out] On input: The vectors to project against Q and normalize. Overwritten with invalid values on output.
X_out[out] On output: the normalized input vectors after projection against Q.
C[out] The projection coefficients
B[out] The normalization coefficients
Q[in] The orthogonal basis against which to project
Returns
Rank of X_in after projection
Note
We expose this interface to applications for the same reason that we expose normalizeOutOfPlace().

Definition at line 339 of file AnasaziTsqrOrthoManagerImpl.hpp.

template<class Scalar , class MV >
magnitude_type Anasazi::TsqrOrthoManagerImpl< Scalar, MV >::orthonormError ( const MV &  X) const
inline

Return $ \| I - X^* \cdot X \|_F $.

Return the Frobenius norm of I - X^* X, which is an absolute measure of the orthogonality of the columns of X.

Definition at line 355 of file AnasaziTsqrOrthoManagerImpl.hpp.

template<class Scalar , class MV >
magnitude_type Anasazi::TsqrOrthoManagerImpl< Scalar, MV >::orthogError ( const MV &  X1,
const MV &  X2 
) const
inline

Return the Frobenius norm of the inner product of X1 with itself.

Definition at line 369 of file AnasaziTsqrOrthoManagerImpl.hpp.

template<class Scalar , class MV >
magnitude_type Anasazi::TsqrOrthoManagerImpl< Scalar, MV >::blockReorthogThreshold ( ) const
inline

Relative tolerance for triggering a block reorthogonalization. If any column norm in a block decreases by this amount, then we reorthogonalize.

Definition at line 382 of file AnasaziTsqrOrthoManagerImpl.hpp.

template<class Scalar , class MV >
magnitude_type Anasazi::TsqrOrthoManagerImpl< Scalar, MV >::relativeRankTolerance ( ) const
inline

Relative tolerance for determining (via the SVD) whether a block is of full numerical rank.

Definition at line 386 of file AnasaziTsqrOrthoManagerImpl.hpp.


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