Belos Package Browser (Single Doxygen Collection)  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Types | Public Member Functions | Private Types | Private Member Functions | Private Attributes | List of all members
Belos::TsqrOrthoManagerImpl< Scalar, MV > Class Template Reference

TSQR-based OrthoManager subclass implementation. More...

#include <BelosTsqrOrthoManagerImpl.hpp>

Inheritance diagram for Belos::TsqrOrthoManagerImpl< Scalar, MV >:
Inheritance graph
[legend]

Public Types

typedef Scalar scalar_type
 
typedef Teuchos::ScalarTraits
< Scalar >::magnitudeType 
magnitude_type
 
typedef MV multivector_type
 
typedef
Teuchos::SerialDenseMatrix
< int, Scalar > 
mat_type
 Type of the projection and normalization coefficients. More...
 
typedef Teuchos::RCP< mat_typemat_ptr
 

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 setReorthogonalizationCallback (const Teuchos::RCP< ReorthogonalizationCallback< Scalar > > &callback)
 Set callback to be invoked on reorthogonalization. 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
 Compute the 2-norm of each column j of X. More...
 
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
 Relative tolerance for triggering a block reorthogonalization. More...
 
magnitude_type relativeRankTolerance () const
 Relative tolerance for determining (via the SVD) whether a block is of full numerical rank. More...
 
- Public Member Functions inherited from Teuchos::ParameterListAcceptorDefaultBase
RCP< ParameterListgetNonconstParameterList ()
 
RCP< ParameterListunsetParameterList ()
 
RCP< const ParameterListgetParameterList () const
 

Private Types

typedef Teuchos::ScalarTraits
< Scalar > 
SCT
 
typedef Teuchos::ScalarTraits
< magnitude_type
SCTM
 
typedef MultiVecTraits< Scalar,
MV > 
MVT
 
typedef MVT::tsqr_adaptor_type tsqr_adaptor_type
 

Private Member Functions

void raiseReorthogFault (const std::vector< magnitude_type > &normsAfterFirstPass, const std::vector< magnitude_type > &normsAfterSecondPass, const std::vector< int > &faultIndices)
 Throw an exception indicating a reorthgonalization fault. More...
 
void checkProjectionDims (int &ncols_X, int &num_Q_blocks, int &ncols_Q_total, const MV &X, Teuchos::ArrayView< Teuchos::RCP< const MV > > Q) const
 Return through output arguments some relevant dimension information about X and Q. More...
 
void allocateProjectionCoefficients (Teuchos::Array< mat_ptr > &C, Teuchos::ArrayView< Teuchos::RCP< const MV > > Q, const MV &X, const bool attemptToRecycle=true) const
 Allocate projection coefficients. More...
 
int projectAndNormalizeImpl (MV &X_in, MV &X_out, const bool outOfPlace, Teuchos::Array< mat_ptr > C, mat_ptr B, Teuchos::ArrayView< Teuchos::RCP< const MV > > Q)
 Implementation of projection and normalization. More...
 
void rawProject (MV &X, Teuchos::ArrayView< Teuchos::RCP< const MV > > Q, Teuchos::ArrayView< mat_ptr > C) const
 One projection pass of X against the Q[i] blocks. More...
 
void rawProject (MV &X, const Teuchos::RCP< const MV > &Q, const mat_ptr &C) const
 Overload of rawProject() for one Q block. More...
 
int rawNormalize (MV &X, MV &Q, mat_type &B)
 One out-of-place normalization pass. More...
 
int normalizeOne (MV &X, mat_ptr B) const
 Normalize a "multivector" of only one column. More...
 
int normalizeImpl (MV &X, MV &Q, mat_ptr B, const bool outOfPlace)
 Normalize X into Q*B, with out-of-place option. More...
 

Private Attributes

Teuchos::RCP
< Teuchos::ParameterList
params_
 Configuration parameters. More...
 
Teuchos::RCP< const
Teuchos::ParameterList
defaultParams_
 Default configuration parameters. More...
 
std::string label_
 Label for timers (if timers are used). More...
 
tsqr_adaptor_type tsqrAdaptor_
 Interface to TSQR implementation. More...
 
Teuchos::RCP< MV > Q_
 Scratch space for TSQR. More...
 
magnitude_type eps_
 Machine precision for Scalar. More...
 
bool randomizeNullSpace_
 Whether to fill null space vectors with random data. More...
 
bool reorthogonalizeBlocks_
 Whether to reorthogonalize blocks at all. More...
 
bool throwOnReorthogFault_
 Whether to throw an exception on a orthogonalization fault. More...
 
magnitude_type blockReorthogThreshold_
 Relative reorthogonalization threshold in Block Gram-Schmidt. More...
 
magnitude_type relativeRankTolerance_
 Relative tolerance for measuring the numerical rank of a matrix. More...
 
bool forceNonnegativeDiagonal_
 Force R factor of normalization to have a nonnegative diagonal. More...
 
Teuchos::RCP
< ReorthogonalizationCallback
< Scalar > > 
reorthogCallback_
 Callback invoked if reorthogonalization is necessary. More...
 

Additional Inherited Members

- Protected Member Functions inherited from Teuchos::ParameterListAcceptorDefaultBase
void setMyParamList (const RCP< ParameterList > &paramList)
 
RCP< ParameterListgetMyNonconstParamList ()
 
RCP< const ParameterListgetMyParamList () const
 

Detailed Description

template<class Scalar, class MV>
class Belos::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 192 of file BelosTsqrOrthoManagerImpl.hpp.

Member Typedef Documentation

template<class Scalar , class MV >
typedef Scalar Belos::TsqrOrthoManagerImpl< Scalar, MV >::scalar_type

Definition at line 195 of file BelosTsqrOrthoManagerImpl.hpp.

template<class Scalar , class MV >
typedef Teuchos::ScalarTraits<Scalar>::magnitudeType Belos::TsqrOrthoManagerImpl< Scalar, MV >::magnitude_type

Definition at line 196 of file BelosTsqrOrthoManagerImpl.hpp.

template<class Scalar , class MV >
typedef MV Belos::TsqrOrthoManagerImpl< Scalar, MV >::multivector_type

Definition at line 197 of file BelosTsqrOrthoManagerImpl.hpp.

template<class Scalar , class MV >
typedef Teuchos::SerialDenseMatrix<int, Scalar> Belos::TsqrOrthoManagerImpl< Scalar, MV >::mat_type

Type of the projection and normalization coefficients.

Definition at line 199 of file BelosTsqrOrthoManagerImpl.hpp.

template<class Scalar , class MV >
typedef Teuchos::RCP<mat_type> Belos::TsqrOrthoManagerImpl< Scalar, MV >::mat_ptr

Definition at line 200 of file BelosTsqrOrthoManagerImpl.hpp.

template<class Scalar , class MV >
typedef Teuchos::ScalarTraits<Scalar> Belos::TsqrOrthoManagerImpl< Scalar, MV >::SCT
private

Definition at line 203 of file BelosTsqrOrthoManagerImpl.hpp.

template<class Scalar , class MV >
typedef Teuchos::ScalarTraits<magnitude_type> Belos::TsqrOrthoManagerImpl< Scalar, MV >::SCTM
private

Definition at line 204 of file BelosTsqrOrthoManagerImpl.hpp.

template<class Scalar , class MV >
typedef MultiVecTraits<Scalar, MV> Belos::TsqrOrthoManagerImpl< Scalar, MV >::MVT
private

Definition at line 205 of file BelosTsqrOrthoManagerImpl.hpp.

template<class Scalar , class MV >
typedef MVT::tsqr_adaptor_type Belos::TsqrOrthoManagerImpl< Scalar, MV >::tsqr_adaptor_type
private

Definition at line 206 of file BelosTsqrOrthoManagerImpl.hpp.

Constructor & Destructor Documentation

template<class Scalar , class MV >
Belos::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 800 of file BelosTsqrOrthoManagerImpl.hpp.

template<class Scalar , class MV >
Belos::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 823 of file BelosTsqrOrthoManagerImpl.hpp.

Member Function Documentation

template<class Scalar , class MV >
Teuchos::RCP< const Teuchos::ParameterList > Belos::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 1445 of file BelosTsqrOrthoManagerImpl.hpp.

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

Set parameters from the given parameter list.

Implements Teuchos::ParameterListAcceptor.

Definition at line 741 of file BelosTsqrOrthoManagerImpl.hpp.

template<class Scalar , class MV >
Teuchos::RCP< const Teuchos::ParameterList > Belos::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 1510 of file BelosTsqrOrthoManagerImpl.hpp.

template<class Scalar , class MV >
void Belos::TsqrOrthoManagerImpl< Scalar, MV >::setReorthogonalizationCallback ( const Teuchos::RCP< ReorthogonalizationCallback< Scalar > > &  callback)
inline

Set callback to be invoked on reorthogonalization.

This callback is invoked right after the first projection step, and only if reorthogonalization will be necessary. It is called before actually reorthogonalizing. The first argument is a Teuchos::ArrayView of the norms of the columns of the input multivector before the first projection pass, and the second argument is a Teuchos::ArrayView of their norms after the first projection pass.

The callback is null by default. If the callback is null, no callback will be invoked.

For details and suggested uses, please refer to the documentation of ReorthogonalizationCallback.

Warning
We assume that the input arguments of the callback's operator() are only valid views within the scope of the function. Your callback should not keep the views.

Definition at line 278 of file BelosTsqrOrthoManagerImpl.hpp.

template<class Scalar , class MV >
void Belos::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 290 of file BelosTsqrOrthoManagerImpl.hpp.

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

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

Definition at line 307 of file BelosTsqrOrthoManagerImpl.hpp.

template<class Scalar , class MV >
void Belos::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 318 of file BelosTsqrOrthoManagerImpl.hpp.

template<class Scalar , class MV >
void Belos::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 846 of file BelosTsqrOrthoManagerImpl.hpp.

template<class Scalar , class MV >
void Belos::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 858 of file BelosTsqrOrthoManagerImpl.hpp.

template<class Scalar , class MV >
int Belos::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 937 of file BelosTsqrOrthoManagerImpl.hpp.

template<class Scalar , class MV >
int Belos::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 1055 of file BelosTsqrOrthoManagerImpl.hpp.

template<class Scalar , class MV >
int Belos::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 407 of file BelosTsqrOrthoManagerImpl.hpp.

template<class Scalar , class MV >
int Belos::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. On output: 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
Note
We expose this interface to applications for the same reason that we expose normalizeOutOfPlace().

Definition at line 437 of file BelosTsqrOrthoManagerImpl.hpp.

template<class Scalar , class MV >
magnitude_type Belos::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 453 of file BelosTsqrOrthoManagerImpl.hpp.

template<class Scalar , class MV >
magnitude_type Belos::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 467 of file BelosTsqrOrthoManagerImpl.hpp.

template<class Scalar , class MV >
magnitude_type Belos::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 480 of file BelosTsqrOrthoManagerImpl.hpp.

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

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

Definition at line 484 of file BelosTsqrOrthoManagerImpl.hpp.

template<class Scalar , class MV >
void Belos::TsqrOrthoManagerImpl< Scalar, MV >::raiseReorthogFault ( const std::vector< magnitude_type > &  normsAfterFirstPass,
const std::vector< magnitude_type > &  normsAfterSecondPass,
const std::vector< int > &  faultIndices 
)
private

Throw an exception indicating a reorthgonalization fault.

Definition at line 1424 of file BelosTsqrOrthoManagerImpl.hpp.

template<class Scalar , class MV >
void Belos::TsqrOrthoManagerImpl< Scalar, MV >::checkProjectionDims ( int &  ncols_X,
int &  num_Q_blocks,
int &  ncols_Q_total,
const MV &  X,
Teuchos::ArrayView< Teuchos::RCP< const MV > >  Q 
) const
private

Return through output arguments some relevant dimension information about X and Q.

Parameters
ncols_X[out] Number of columns in X
num_Q_blocks[out] Number of entries in the Q array
ncols_Q_total[out] Total number of columns in all the entries of Q
X[in] Multivector to project against the Q[i]
Q[in] Array of multivectors against which to project X

Definition at line 1913 of file BelosTsqrOrthoManagerImpl.hpp.

template<class Scalar , class MV >
void Belos::TsqrOrthoManagerImpl< Scalar, MV >::allocateProjectionCoefficients ( Teuchos::Array< mat_ptr > &  C,
Teuchos::ArrayView< Teuchos::RCP< const MV > >  Q,
const MV &  X,
const bool  attemptToRecycle = true 
) const
private

Allocate projection coefficients.

Parameters
C[out] Array of projection coefficient matrices
Q[in] Array of MV against which to project
X[in] MV to project against the entries of Q
attemptToRecycle[in] Hint whether to check the existing entries of C to see if they have already been allocated and have the right dimensions. This function will do the right thing regardless, but the hint might improve performance by avoiding unnecessary allocations or checks.

Definition at line 1015 of file BelosTsqrOrthoManagerImpl.hpp.

template<class Scalar , class MV >
int Belos::TsqrOrthoManagerImpl< Scalar, MV >::projectAndNormalizeImpl ( MV &  X_in,
MV &  X_out,
const bool  outOfPlace,
Teuchos::Array< mat_ptr C,
mat_ptr  B,
Teuchos::ArrayView< Teuchos::RCP< const MV > >  Q 
)
private

Implementation of projection and normalization.

Implementation of projectAndNormalize() (in which case X_out is not read or written, so it may alias X_in, and outOfPlace==false) and projectAndNormalizeOutOfPlace() (in which case X_out is written, and outOfPlace==true).

Returns
Rank of X_in after projection

Definition at line 1088 of file BelosTsqrOrthoManagerImpl.hpp.

template<class Scalar , class MV >
void Belos::TsqrOrthoManagerImpl< Scalar, MV >::rawProject ( MV &  X,
Teuchos::ArrayView< Teuchos::RCP< const MV > >  Q,
Teuchos::ArrayView< mat_ptr C 
) const
private

One projection pass of X against the Q[i] blocks.

Perform one projection pass (Modified Block Gram-Schmidt) of X against the Q[i] blocks. Does not allocate C[i] coefficients, and does not reorthogonalize.

template<class Scalar , class MV >
void Belos::TsqrOrthoManagerImpl< Scalar, MV >::rawProject ( MV &  X,
const Teuchos::RCP< const MV > &  Q,
const mat_ptr C 
) const
private

Overload of rawProject() for one Q block.

template<class Scalar , class MV >
int Belos::TsqrOrthoManagerImpl< Scalar, MV >::rawNormalize ( MV &  X,
MV &  Q,
mat_type B 
)
private

One out-of-place normalization pass.

Compute one normalization pass of X into Q*B. Overwrite X with invalid values.

Parameters
X[in/out] On input: multivector whose columns are to be orthogonalized ("normalized"). On output: overwritten with invalid values.
Q[out] The orthogonalized ("normalized") columns of X. If X on input had (numerical) rank r, the first r columns are a column space basis for X, and the remaining columns are a null space basis for X.
B[out] Normalization coefficients: X = Q*B. If X on input had full (numerical) rank, B is upper triangular. Otherwise, B may not be upper triangular, but the factorization X = Q*B is still valid.
Returns
The rank of the input X
Warning
Q must have exactly as many columns as X.
B must have been allocated and must have the right dimensions (square, with number of rows/columns equal to the number of columns in X).

Definition at line 1538 of file BelosTsqrOrthoManagerImpl.hpp.

template<class Scalar , class MV >
int Belos::TsqrOrthoManagerImpl< Scalar, MV >::normalizeOne ( MV &  X,
mat_ptr  B 
) const
private

Normalize a "multivector" of only one column.

Special case of normalize() when X has only one column. The operation is done in place in X. We assume that B(0,0) makes sense, since we're going to assign to it.

Parameters
X[in/out] On input: multivector of one column. On output: if that column has nonzero 2-norm, the column scaled by its 2-norm; otherwise, if the column has zero 2-norm, it is not modified.
B[out] Matrix of dimension 1 x 1. On output, the 2-norm of X.
Returns
One if X has nonzero 2-norm, else zero.
Warning
We do no checking of the dimensions of X or B, and we do not resize B if it has dimensions other than 1 x 1.

Definition at line 1560 of file BelosTsqrOrthoManagerImpl.hpp.

template<class Scalar , class MV >
int Belos::TsqrOrthoManagerImpl< Scalar, MV >::normalizeImpl ( MV &  X,
MV &  Q,
mat_ptr  B,
const bool  outOfPlace 
)
private

Normalize X into Q*B, with out-of-place option.

If outOfPlace is true, write the normalized vectors to Q, leaving the contents of X invalid. Otherwise, write the normalized vectors to X, leaving the contents of Q invalid. Regardless, if X on input had (numerical) rank r, the first r normalized vectors are a column space basis for X, and the remaining vectors are a null space basis for X.

Parameters
X[in/out] On input: multivector whose columns are to be orthogonalized ("normalized"). On output: if outOfPlace, overwritten with invalid values; else, the normalized vector(s).
Q[out] If outOfPlace, overwritten with invalid values; else, the normalized vector(s).
B[out] Normalization coefficients: X = Q*B. If X on input had full (numerical) rank, B is upper triangular. Otherwise, B may not be upper triangular, but the factorization X = Q*B is still valid.
Returns
The rank of X
Note
Q must have at least as many columns as X. It may have more columns than X. This routine doesn't try to allocate space for Q if it is too small.

Definition at line 1677 of file BelosTsqrOrthoManagerImpl.hpp.

Member Data Documentation

template<class Scalar , class MV >
Teuchos::RCP<Teuchos::ParameterList> Belos::TsqrOrthoManagerImpl< Scalar, MV >::params_
private

Configuration parameters.

Definition at line 488 of file BelosTsqrOrthoManagerImpl.hpp.

template<class Scalar , class MV >
Teuchos::RCP<const Teuchos::ParameterList> Belos::TsqrOrthoManagerImpl< Scalar, MV >::defaultParams_
mutableprivate

Default configuration parameters.

Definition at line 491 of file BelosTsqrOrthoManagerImpl.hpp.

template<class Scalar , class MV >
std::string Belos::TsqrOrthoManagerImpl< Scalar, MV >::label_
private

Label for timers (if timers are used).

Definition at line 494 of file BelosTsqrOrthoManagerImpl.hpp.

template<class Scalar , class MV >
tsqr_adaptor_type Belos::TsqrOrthoManagerImpl< Scalar, MV >::tsqrAdaptor_
private

Interface to TSQR implementation.

Definition at line 497 of file BelosTsqrOrthoManagerImpl.hpp.

template<class Scalar , class MV >
Teuchos::RCP<MV> Belos::TsqrOrthoManagerImpl< Scalar, MV >::Q_
private

Scratch space for TSQR.

This multivector scratch space is allocated lazily, only if normalize() is called with a multivector input having more than one column. We do our best to avoid reallocation and recycle this space whenever possible. The normalizeOutOfPlace() method does not allocate Q_, which you can use to your advantage if you already have scratch space allocated.

Definition at line 508 of file BelosTsqrOrthoManagerImpl.hpp.

template<class Scalar , class MV >
magnitude_type Belos::TsqrOrthoManagerImpl< Scalar, MV >::eps_
private

Machine precision for Scalar.

Definition at line 511 of file BelosTsqrOrthoManagerImpl.hpp.

template<class Scalar , class MV >
bool Belos::TsqrOrthoManagerImpl< Scalar, MV >::randomizeNullSpace_
private

Whether to fill null space vectors with random data.

If so, this happens after normalization.

Definition at line 516 of file BelosTsqrOrthoManagerImpl.hpp.

template<class Scalar , class MV >
bool Belos::TsqrOrthoManagerImpl< Scalar, MV >::reorthogonalizeBlocks_
private

Whether to reorthogonalize blocks at all.

Reorthogonalization is conditional, based on the block reorthogonalization threshold. Tests for reorthogonalization only happen if this Boolean is set.

Definition at line 523 of file BelosTsqrOrthoManagerImpl.hpp.

template<class Scalar , class MV >
bool Belos::TsqrOrthoManagerImpl< Scalar, MV >::throwOnReorthogFault_
private

Whether to throw an exception on a orthogonalization fault.

Recovery is possible, but expensive.

Definition at line 528 of file BelosTsqrOrthoManagerImpl.hpp.

template<class Scalar , class MV >
magnitude_type Belos::TsqrOrthoManagerImpl< Scalar, MV >::blockReorthogThreshold_
private

Relative reorthogonalization threshold in Block Gram-Schmidt.

Definition at line 531 of file BelosTsqrOrthoManagerImpl.hpp.

template<class Scalar , class MV >
magnitude_type Belos::TsqrOrthoManagerImpl< Scalar, MV >::relativeRankTolerance_
private

Relative tolerance for measuring the numerical rank of a matrix.

Definition at line 534 of file BelosTsqrOrthoManagerImpl.hpp.

template<class Scalar , class MV >
bool Belos::TsqrOrthoManagerImpl< Scalar, MV >::forceNonnegativeDiagonal_
private

Force R factor of normalization to have a nonnegative diagonal.

If true, then (if necessary) do extra work (modifying both the Q and R factors) in the normalization step in order to force the R factor of the current block to have a nonnegative diagonal.

Definition at line 542 of file BelosTsqrOrthoManagerImpl.hpp.

template<class Scalar , class MV >
Teuchos::RCP<ReorthogonalizationCallback<Scalar> > Belos::TsqrOrthoManagerImpl< Scalar, MV >::reorthogCallback_
private

Callback invoked if reorthogonalization is necessary.

Definition at line 556 of file BelosTsqrOrthoManagerImpl.hpp.


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