Stratimikos Package Browser (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Types | Public Member Functions | Private Member Functions | Static Private Member Functions | List of all members
Thyra::TsqrAdaptor< Scalar > Class Template Reference

Stub adaptor from Thyra::MultiVectorBase to TSQR. More...

#include <Thyra_TsqrAdaptor.hpp>

Inheritance diagram for Thyra::TsqrAdaptor< Scalar >:
Inheritance graph
[legend]

Public Types

typedef Thyra::MultiVectorBase
< Scalar > 
MV
 
typedef Scalar scalar_type
 
typedef int ordinal_type
 
typedef
Teuchos::SerialDenseMatrix
< ordinal_type, scalar_type
dense_matrix_type
 
typedef Teuchos::ScalarTraits
< scalar_type >::magnitudeType 
magnitude_type
 

Public Member Functions

 TsqrAdaptor (const Teuchos::RCP< Teuchos::ParameterList > &)
 Constructor (that accepts a parameter list). More...
 
 TsqrAdaptor ()
 Constructor (that uses default parameters). More...
 
Teuchos::RCP< const
Teuchos::ParameterList
getValidParameters () const
 
void setParameterList (const Teuchos::RCP< Teuchos::ParameterList > &)
 
void factorExplicit (MV &, MV &, dense_matrix_type &, const bool=false)
 Compute QR factorization [Q,R] = qr(A,0). More...
 
int revealRank (MV &, dense_matrix_type &, const magnitude_type &)
 Rank-revealing decomposition. More...
 
- Public Member Functions inherited from Teuchos::ParameterListAcceptorDefaultBase
RCP< ParameterListgetNonconstParameterList ()
 
RCP< ParameterListunsetParameterList ()
 
RCP< const ParameterListgetParameterList () const
 

Private Member Functions

void prepareDistTsqr (const MV &)
 Finish interprocess TSQR initialization. More...
 
void prepareTsqr (const MV &)
 Finish TSQR initialization. More...
 

Static Private Member Functions

static Teuchos::RCP< const
Teuchos::Comm< int > > 
getComm (const MV &X)
 Attempt to get a communicator out of the given multivector. 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 Thyra::TsqrAdaptor< Scalar >

Stub adaptor from Thyra::MultiVectorBase to TSQR.

TSQR (Tall Skinny QR factorization) is an orthogonalization kernel that is as accurate as Householder QR, yet requires only $2 \log P$ messages between $P$ MPI processes, independently of the number of columns in the multivector.

TSQR works independently of the particular multivector implementation, and interfaces to the latter via an adaptor class. This class is the adaptor class for MultiVectorBase. It templates on the MultiVector (MV) type so that it can pick up that class' typedefs. In particular, TSQR chooses its intranode implementation based on the Kokkos Node type of the multivector.

Warning
This is a stub adaptor that just placates the compiler and does nothing. It's not hard to implement a Thyra adaptor, but in order for the adaptor to be efficient, it requires special cases for extracting the actual multivector implementation (e.g., Epetra_MultiVector or Tpetra::MultiVector) out of the Thyra wrapper.

Definition at line 86 of file Thyra_TsqrAdaptor.hpp.

Member Typedef Documentation

template<class Scalar >
typedef Thyra::MultiVectorBase<Scalar> Thyra::TsqrAdaptor< Scalar >::MV

Definition at line 88 of file Thyra_TsqrAdaptor.hpp.

template<class Scalar >
typedef Scalar Thyra::TsqrAdaptor< Scalar >::scalar_type

Definition at line 89 of file Thyra_TsqrAdaptor.hpp.

template<class Scalar >
typedef int Thyra::TsqrAdaptor< Scalar >::ordinal_type

Definition at line 90 of file Thyra_TsqrAdaptor.hpp.

Definition at line 91 of file Thyra_TsqrAdaptor.hpp.

template<class Scalar >
typedef Teuchos::ScalarTraits<scalar_type>::magnitudeType Thyra::TsqrAdaptor< Scalar >::magnitude_type

Definition at line 92 of file Thyra_TsqrAdaptor.hpp.

Constructor & Destructor Documentation

template<class Scalar >
Thyra::TsqrAdaptor< Scalar >::TsqrAdaptor ( const Teuchos::RCP< Teuchos::ParameterList > &  )
inline

Constructor (that accepts a parameter list).

Parameters
plist[in] List of parameters for configuring TSQR. The specific parameter keys that are read depend on the TSQR implementation. For details, call getValidParameters() and examine the documentation embedded therein.

Definition at line 100 of file Thyra_TsqrAdaptor.hpp.

template<class Scalar >
Thyra::TsqrAdaptor< Scalar >::TsqrAdaptor ( )
inline

Constructor (that uses default parameters).

Definition at line 106 of file Thyra_TsqrAdaptor.hpp.

Member Function Documentation

template<class Scalar >
Teuchos::RCP<const Teuchos::ParameterList> Thyra::TsqrAdaptor< Scalar >::getValidParameters ( ) const
inlinevirtual

Reimplemented from Teuchos::ParameterListAcceptor.

Definition at line 112 of file Thyra_TsqrAdaptor.hpp.

template<class Scalar >
void Thyra::TsqrAdaptor< Scalar >::setParameterList ( const Teuchos::RCP< Teuchos::ParameterList > &  )
inlinevirtual

Implements Teuchos::ParameterListAcceptor.

Definition at line 118 of file Thyra_TsqrAdaptor.hpp.

template<class Scalar >
void Thyra::TsqrAdaptor< Scalar >::factorExplicit ( MV ,
MV ,
dense_matrix_type ,
const bool  = false 
)
inline

Compute QR factorization [Q,R] = qr(A,0).

Parameters
A[in/out] On input: the multivector to factor. Overwritten with garbage on output.
Q[out] On output: the (explicitly stored) Q factor in the QR factorization of the (input) multivector A.
R[out] On output: the R factor in the QR factorization of the (input) multivector A.
forceNonnegativeDiagonal[in] If true, then (if necessary) do extra work (modifying both the Q and R factors) in order to force the R factor to have a nonnegative diagonal.
Warning
Currently, this method only works if A and Q have the same communicator and row distribution ("map," in Petra terms) as those of the multivector given to this TsqrAdaptor instance's constructor. Otherwise, the result of this method is undefined.

Definition at line 145 of file Thyra_TsqrAdaptor.hpp.

template<class Scalar >
int Thyra::TsqrAdaptor< Scalar >::revealRank ( MV ,
dense_matrix_type ,
const magnitude_type  
)
inline

Rank-revealing decomposition.

Using the R factor and explicit Q factor from factorExplicit(), compute the singular value decomposition (SVD) of R ( $R = U \Sigma V^*$). If R is full rank (with respect to the given relative tolerance tol), don't change Q or R. Otherwise, compute $Q := Q \cdot U$ and $R := \Sigma V^*$ in place (the latter may be no longer upper triangular).

Parameters
Q[in/out] On input: explicit Q factor computed by factorExplicit(). (Must be an orthogonal resp. unitary matrix.) On output: If R is of full numerical rank with respect to the tolerance tol, Q is unmodified. Otherwise, Q is updated so that the first rank columns of Q are a basis for the column space of A (the original matrix whose QR factorization was computed by factorExplicit()). The remaining columns of Q are a basis for the null space of A.
R[in/out] On input: ncols by ncols upper triangular matrix with leading dimension ldr >= ncols. On output: if input is full rank, R is unchanged on output. Otherwise, if $R = U \Sigma V^*$ is the SVD of R, on output R is overwritten with $ V^*$. This is also an ncols by ncols matrix, but may not necessarily be upper triangular.
tol[in] Relative tolerance for computing the numerical rank of the matrix R.
Returns
Rank $r$ of R: $ 0 \leq r \leq ncols$.

Definition at line 184 of file Thyra_TsqrAdaptor.hpp.

template<class Scalar >
static Teuchos::RCP<const Teuchos::Comm<int> > Thyra::TsqrAdaptor< Scalar >::getComm ( const MV X)
inlinestaticprivate

Attempt to get a communicator out of the given multivector.

This only works if the multivector's range (VectorSpaceBase) is actually an SpmdVectorSpaceBase object, and if that object's Comm is either an MpiComm (in an MPI build) or a SerialComm (in either an MPI build or a no-MPI build).

If the attempt does not succeed, this method throws std::runtime_error. If it does succeed, it returns the (suitably wrapped) communicator.

Definition at line 203 of file Thyra_TsqrAdaptor.hpp.

template<class Scalar >
void Thyra::TsqrAdaptor< Scalar >::prepareDistTsqr ( const MV )
inlineprivate

Finish interprocess TSQR initialization.

Input X is a valid Thyra::MultiVectorBase instance whose communicator wrapper we will use to prepare TSQR. It is not modified.

Note
It's OK to call this method more than once; it is idempotent.

This method may fail if MV is not the right kind of multivector, that is, if it does not have a communicator or if we don't know how to extract a communicator from it. If it fails in this way, it will throw std::runtime_error.

Definition at line 281 of file Thyra_TsqrAdaptor.hpp.

template<class Scalar >
void Thyra::TsqrAdaptor< Scalar >::prepareTsqr ( const MV )
inlineprivate

Finish TSQR initialization.

The intranode and internode TSQR implementations both have a two-stage initialization procedure: first, setting parameters (which may happen at construction), and second, getting information they need from the multivector input in order to finish initialization. For intranode TSQR, this may include the Kokkos Node instance; for internode TSQR, this includes the communicator. The second stage of initialization happens in this class' computational routines; all of those routines accept one or more multivector inputs, which this method can use for finishing initialization. Thus, users of this class never need to see the two-stage initialization.

Parameters
X[in] Multivector object, used only to access the underlying communicator object (in this case, the Teuchos::Comm<int>) and (possibly) the Kokkos Node instance. All multivector objects used with this adapter must have the same communicator and Kokkos Node instance (if applicable).

Definition at line 303 of file Thyra_TsqrAdaptor.hpp.


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