Anasazi Version of the Day
Public Member Functions | Protected Member Functions
TSQR::Trilinos::TsqrAdaptor< S, LO, GO, MV > Class Template Reference

Abstract interface between multivector and TSQR. More...

#include <TsqrAdaptor.hpp>

List of all members.

Public Member Functions

virtual factor_output_type factor (multivector_type &A, dense_matrix_type &R, const bool contiguousCacheBlocks=false)
 Compute QR factorization of the multivector A.
virtual void explicitQ (const multivector_type &Q_in, const factor_output_type &factorOutput, multivector_type &Q_out, const bool contiguousCacheBlocks=false)
 Compute the explicit Q factor.
local_ordinal_type revealRank (multivector_type &Q, dense_matrix_type &R, const magnitude_type relativeTolerance, const bool contiguousCacheBlocks=false) const
 Rank-revealing decomposition.
virtual void cacheBlock (const multivector_type &A_in, multivector_type &A_out)
 Cache-block A_in into A_out.
virtual void unCacheBlock (const multivector_type &A_in, multivector_type &A_out)
 Un-cache-block A_in into A_out.
virtual std::pair
< magnitude_type,
magnitude_type > 
verify (const multivector_type &A, const multivector_type &Q, const Teuchos::SerialDenseMatrix< local_ordinal_type, scalar_type > &R)

Protected Member Functions

void init (const multivector_type &mv, const Teuchos::ParameterList &plist)

Detailed Description

template<class S, class LO, class GO, class MV>
class TSQR::Trilinos::TsqrAdaptor< S, LO, GO, MV >

Abstract interface between multivector and TSQR.

Child classes of TsqrAdaptor tell TSQR how to compute a factorization of a specific Trilinos multivector class MV. Currently, Tpetra::MultiVector< S, LO, GO, NodeType > for any NodeType is supported, via TsqrEpetraAdaptor (see AnasaziEpetraAdaptor.hpp) resp. TsqrTpetraAdaptor (include AnasaziTpetraAdaptor.hpp). At the moment, the latter will only be efficient if NodeType is not a GPU node. Support for Epetra_MultiVector and Thyra multivectors may be added on request.

TsqrAdaptor uses the appropriate specialization of TsqrTypeAdaptor to figure out which variant of TSQR to use on the given multivector type. For example, with Tpetra::MultiVector< S, LO, GO, NodeType >, if NodeType is Kokkos::TBBNode, the TBB-parallel intranode variant of TSQR will be used. The caller is responsible for constructing the intranode and internode TSQR objects.

Implementers who want to support TSQR with a new MultiVector (MV) type must create a subclass of that type, using e.g., TsqrTpetraAdaptor as a model. They must then create a new TsqrTypeAdaptor specialization (with the appropriate typedefs), and a new TsqrCommFactory subclass. The TsqrCommFactory subclass gets the underlying communicator object (e.g., Teuchos::Comm<int>) from a "prototype" multivector and wraps it into TSQR::MessengerBase< S > and TSQR::MessengerBase< LO > objects for TSQR.

Implementers who wish to change which TSQR implementation is used for a particular MultiVector type (for which a TsqrAdaptor child class exists) should change the corresponding (possibly partial) specialization of TsqrTypeAdaptor. Certainly the node_tsqr_type (and perhaps also the dist_tsqr_type) typedef(s) in the TsqrTypeAdaptor specialization must be changed. If no corresponding TsqrFactory subclass exists for that combination of node_tsqr_type and dist_tsqr_type, a new TsqrFactory subclass may also have to be created, to tell Anasazi how to instantiate those node_tsqr_type and dist_tsqr_type objects.

Implementers who wish to add a new TSQR factorization must create a new TsqrFactory subclass.

Definition at line 105 of file TsqrAdaptor.hpp.


Member Function Documentation

template<class S, class LO, class GO, class MV>
virtual factor_output_type TSQR::Trilinos::TsqrAdaptor< S, LO, GO, MV >::factor ( multivector_type &  A,
dense_matrix_type R,
const bool  contiguousCacheBlocks = false 
) [inline, virtual]

Compute QR factorization of the multivector A.

Compute the QR factorization in place of the multivector A. The Q factor is represented implicitly; part of that is stored in place in A (overwriting the input), and the other part is returned. The returned object as well as the representation in A are both inputs of explicitQ(). The R factor is copied into R.

Parameters:
A[in/out] On input, the multivector whose QR factorization is to be computed. Overwritten on output with part of the implicit representation of the Q factor.
R[out] On output, the R factor from the QR factorization of A. Represented as a square dense matrix (not in packed form) with the same number of columns as A. The lower triangle of R is overwritten with zeros on output.
contiguousCacheBlocks[in] Whether the data in A has been reorganized so that the elements of each cache block are stored contiguously (i.e., via the output of cacheBlock()). The default is false, which means that each process' row block of A is stored as a matrix in column-major order, with leading dimension >= the number of rows in the row block.
Returns:
Additional information that, together with the A output, encodes the implicitly represented Q factor from the QR factorization of the A input.
Note:
Virtual but implemented, because this default implementation is correct for all multivector_type types, but not necessarily efficient. It should be efficient if fetchNonConstView(A) does not require copying the contents of A (e.g., from GPU memory to CPU memory).

Definition at line 173 of file TsqrAdaptor.hpp.

template<class S, class LO, class GO, class MV>
virtual void TSQR::Trilinos::TsqrAdaptor< S, LO, GO, MV >::explicitQ ( const multivector_type &  Q_in,
const factor_output_type &  factorOutput,
multivector_type &  Q_out,
const bool  contiguousCacheBlocks = false 
) [inline, virtual]

Compute the explicit Q factor.

Compute the explicit (multivector) "thin" (same number of columns as the input) representation of the Q factor computed by factor(), using the implicit representation returned by factor().

Parameters:
Q_in[in] Same as the "A" input of factor()
factorOutput[in] Return value of factor() corresponding to Q_in
Q_out[out] Explicit "thin" representation of the Q factor. "Explicit" means as a regular matrix (in the same multivector storage format as the "A" input of factor()). "Thin" (terminology used by Golub and Van Loan) means that the dimensions of Q_out are the same as the dimensions of the "A" input of factor().
contiguousCacheBlocks[in] See the epinonymous argument of factor(). In this case, it applies to both Q_in and Q_out, which must have the same data layout.
Note:
Virtual but implemented, because this default implementation is correct for all multivector_type types, but not necessarily efficient. It should be efficient if fetchNonConstView(Q_out) and fetchConstView(Q_in) do not require copying (e.g., from GPU memory to CPU memory) the contents of their respective multivector inputs.

Definition at line 222 of file TsqrAdaptor.hpp.

template<class S, class LO, class GO, class MV>
local_ordinal_type TSQR::Trilinos::TsqrAdaptor< S, LO, GO, MV >::revealRank ( multivector_type &  Q,
dense_matrix_type R,
const magnitude_type  relativeTolerance,
const bool  contiguousCacheBlocks = false 
) const [inline]

Rank-revealing decomposition.

Using the R factor from factor() and the explicit Q factor from explicitQ(), compute the SVD of R ( $R = U \Sigma V^*$). R. If R is full rank (with respect to the given relative tolerance), 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: the explicit Q factor computed by explicitQ(). On output: unchanged if R has full (numerical) rank, else $Q := Q \cdot U$, where $U$ is the ncols by ncols matrix of R's left singular vectors.
R[in/out] On input: ncols by ncols upper triangular matrix stored in column-major order. On output: if input has full (numerical) rank, R is unchanged on output. Otherwise, if $R = U \Sigma V^*$ is the SVD of R, on output R is overwritten with $\Sigma \cdot V^*$. This is also an ncols by ncols matrix, but may not necessarily be upper triangular.
Returns:
Rank $r$ of R: $ 0 \leq r \leq ncols$.

Definition at line 276 of file TsqrAdaptor.hpp.

template<class S, class LO, class GO, class MV>
virtual void TSQR::Trilinos::TsqrAdaptor< S, LO, GO, MV >::cacheBlock ( const multivector_type &  A_in,
multivector_type &  A_out 
) [inline, virtual]

Cache-block A_in into A_out.

Copy A_in into A_out, in a reorganized way that improves locality of cache blocks.

Warning:
This may invalidate some invariants of A_out, such as the mapping from index pair (i,j) to element of A_out. Another way to say this is that the multivector object may not be aware that its data has been reorganized underneath it.

Definition at line 305 of file TsqrAdaptor.hpp.

template<class S, class LO, class GO, class MV>
virtual void TSQR::Trilinos::TsqrAdaptor< S, LO, GO, MV >::unCacheBlock ( const multivector_type &  A_in,
multivector_type &  A_out 
) [inline, virtual]

Un-cache-block A_in into A_out.

Undo the transformation performed by cacheBlock(), by copying the contiguously cache blocked data in A_in into the conventionally stored A_out.

Definition at line 343 of file TsqrAdaptor.hpp.

template<class S, class LO, class GO, class MV>
virtual std::pair< magnitude_type, magnitude_type > TSQR::Trilinos::TsqrAdaptor< S, LO, GO, MV >::verify ( const multivector_type &  A,
const multivector_type &  Q,
const Teuchos::SerialDenseMatrix< local_ordinal_type, scalar_type > &  R 
) [inline, virtual]
Returns:
Two magnitudes: first $\| A - QR \|_F / \|A\|_F$, second $\|I - Q^* Q\|_F / \|A\|_F$.

Definition at line 378 of file TsqrAdaptor.hpp.

template<class S, class LO, class GO, class MV>
void TSQR::Trilinos::TsqrAdaptor< S, LO, GO, MV >::init ( const multivector_type &  mv,
const Teuchos::ParameterList plist 
) [inline, protected]

Like the constructor, except you're not supposed to call the constructor of a pure virtual class.

Definition at line 409 of file TsqrAdaptor.hpp.


The documentation for this class was generated from the following file:
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends