Belos  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Types | Public Member Functions | List of all members
Belos::details::LocalDenseMatrixOps< Scalar > Class Template Reference

Low-level operations on non-distributed dense matrices. More...

#include <BelosProjectedLeastSquaresSolver.hpp>

Public Types

typedef Scalar scalar_type
 The template parameter of this class. More...
 
typedef Teuchos::ScalarTraits
< Scalar >::magnitudeType 
magnitude_type
 The type of the magnitude of a scalar_type value. More...
 
typedef
Teuchos::SerialDenseMatrix
< int, Scalar > 
mat_type
 The type of a dense matrix (or vector) of scalar_type. More...
 

Public Member Functions

void conjugateTranspose (mat_type &A_star, const mat_type &A) const
 A_star := (conjugate) transpose of A. More...
 
void conjugateTransposeOfUpperTriangular (mat_type &L, const mat_type &R) const
 L := (conjugate) transpose of R (upper triangular). More...
 
void zeroOutStrictLowerTriangle (mat_type &A) const
 Zero out everything below the diagonal of A. More...
 
void partition (Teuchos::RCP< mat_type > &A_11, Teuchos::RCP< mat_type > &A_21, Teuchos::RCP< mat_type > &A_12, Teuchos::RCP< mat_type > &A_22, mat_type &A, const int numRows1, const int numRows2, const int numCols1, const int numCols2)
 A -> [A_11, A_21, A_12, A_22]. More...
 
void matScale (mat_type &A, const scalar_type &alpha) const
 A := alpha * A. More...
 
void axpy (mat_type &Y, const scalar_type &alpha, const mat_type &X) const
 Y := Y + alpha * X. More...
 
void matAdd (mat_type &A, const mat_type &B) const
 A := A + B. More...
 
void matSub (mat_type &A, const mat_type &B) const
 A := A - B. More...
 
void rightUpperTriSolve (mat_type &B, const mat_type &R) const
 In Matlab notation: B = B / R, where R is upper triangular. More...
 
void matMatMult (const scalar_type &beta, mat_type &C, const scalar_type &alpha, const mat_type &A, const mat_type &B) const
 C := beta*C + alpha*A*B. More...
 
int infNaNCount (const mat_type &A, const bool upperTriangular=false) const
 Return the number of Inf or NaN entries in the matrix A. More...
 
std::pair< bool, std::pair
< magnitude_type,
magnitude_type > > 
isUpperTriangular (const mat_type &A) const
 Is the matrix A upper triangular / trapezoidal? More...
 
std::pair< bool, std::pair
< magnitude_type,
magnitude_type > > 
isUpperHessenberg (const mat_type &A) const
 Is the matrix A upper Hessenberg? More...
 
void ensureUpperTriangular (const mat_type &A, const char *const matrixName) const
 Throw an exception if A is not upper triangular / trapezoidal. More...
 
void ensureUpperHessenberg (const mat_type &A, const char *const matrixName) const
 Throw an exception if A is not (strictly) upper Hessenberg. More...
 
void ensureUpperHessenberg (const mat_type &A, const char *const matrixName, const magnitude_type relativeTolerance) const
 Throw an exception if A is not "approximately" upper Hessenberg. More...
 
void ensureMinimumDimensions (const mat_type &A, const char *const matrixName, const int minNumRows, const int minNumCols) const
 Ensure that the matrix A is at least minNumRows by minNumCols. More...
 
void ensureEqualDimensions (const mat_type &A, const char *const matrixName, const int numRows, const int numCols) const
 Ensure that the matrix A is exactly numRows by numCols. More...
 

Detailed Description

template<class Scalar>
class Belos::details::LocalDenseMatrixOps< Scalar >

Low-level operations on non-distributed dense matrices.

Author
Mark Hoemmen

This class provides a convenient wrapper around some BLAS operations, operating on non-distributed (hence "local") dense matrices.

Definition at line 349 of file BelosProjectedLeastSquaresSolver.hpp.

Member Typedef Documentation

template<class Scalar>
Belos::details::LocalDenseMatrixOps< Scalar >::scalar_type

The template parameter of this class.

Definition at line 353 of file BelosProjectedLeastSquaresSolver.hpp.

template<class Scalar>
Belos::details::LocalDenseMatrixOps< Scalar >::magnitude_type

The type of the magnitude of a scalar_type value.

Definition at line 356 of file BelosProjectedLeastSquaresSolver.hpp.

template<class Scalar>
Belos::details::LocalDenseMatrixOps< Scalar >::mat_type

The type of a dense matrix (or vector) of scalar_type.

Definition at line 359 of file BelosProjectedLeastSquaresSolver.hpp.

Member Function Documentation

template<class Scalar>
void Belos::details::LocalDenseMatrixOps< Scalar >::conjugateTranspose ( mat_type A_star,
const mat_type A 
) const
inline

A_star := (conjugate) transpose of A.

Definition at line 370 of file BelosProjectedLeastSquaresSolver.hpp.

template<class Scalar>
void Belos::details::LocalDenseMatrixOps< Scalar >::conjugateTransposeOfUpperTriangular ( mat_type L,
const mat_type R 
) const
inline

L := (conjugate) transpose of R (upper triangular).

Definition at line 381 of file BelosProjectedLeastSquaresSolver.hpp.

template<class Scalar>
void Belos::details::LocalDenseMatrixOps< Scalar >::zeroOutStrictLowerTriangle ( mat_type A) const
inline

Zero out everything below the diagonal of A.

Definition at line 394 of file BelosProjectedLeastSquaresSolver.hpp.

template<class Scalar>
void Belos::details::LocalDenseMatrixOps< Scalar >::partition ( Teuchos::RCP< mat_type > &  A_11,
Teuchos::RCP< mat_type > &  A_21,
Teuchos::RCP< mat_type > &  A_12,
Teuchos::RCP< mat_type > &  A_22,
mat_type A,
const int  numRows1,
const int  numRows2,
const int  numCols1,
const int  numCols2 
)
inline

A -> [A_11, A_21, A_12, A_22].

The first four arguments are the output arguments. They are views of their respective submatrices of A.

Note
SerialDenseMatrix's operator= and copy constructor both always copy the matrix deeply, so I can't make the output arguments "mat_type&".

Definition at line 414 of file BelosProjectedLeastSquaresSolver.hpp.

template<class Scalar>
void Belos::details::LocalDenseMatrixOps< Scalar >::matScale ( mat_type A,
const scalar_type alpha 
) const
inline

A := alpha * A.

Definition at line 435 of file BelosProjectedLeastSquaresSolver.hpp.

template<class Scalar>
void Belos::details::LocalDenseMatrixOps< Scalar >::axpy ( mat_type Y,
const scalar_type alpha,
const mat_type X 
) const
inline

Y := Y + alpha * X.

"AXPY" stands for "alpha times X plus y," and is the traditional abbreviation for this operation.

Definition at line 459 of file BelosProjectedLeastSquaresSolver.hpp.

template<class Scalar>
void Belos::details::LocalDenseMatrixOps< Scalar >::matAdd ( mat_type A,
const mat_type B 
) const
inline

A := A + B.

Definition at line 479 of file BelosProjectedLeastSquaresSolver.hpp.

template<class Scalar>
void Belos::details::LocalDenseMatrixOps< Scalar >::matSub ( mat_type A,
const mat_type B 
) const
inline

A := A - B.

Definition at line 506 of file BelosProjectedLeastSquaresSolver.hpp.

template<class Scalar>
void Belos::details::LocalDenseMatrixOps< Scalar >::rightUpperTriSolve ( mat_type B,
const mat_type R 
) const
inline

In Matlab notation: B = B / R, where R is upper triangular.

This method only looks at the upper left R.numCols() by R.numCols() part of R.

Definition at line 536 of file BelosProjectedLeastSquaresSolver.hpp.

template<class Scalar>
void Belos::details::LocalDenseMatrixOps< Scalar >::matMatMult ( const scalar_type beta,
mat_type C,
const scalar_type alpha,
const mat_type A,
const mat_type B 
) const
inline

C := beta*C + alpha*A*B.

This method is a thin wrapper around the BLAS' _GEMM routine. The matrix C is NOT allowed to alias the matrices A or B. This method makes no effort to check for aliasing.

Definition at line 558 of file BelosProjectedLeastSquaresSolver.hpp.

template<class Scalar>
int Belos::details::LocalDenseMatrixOps< Scalar >::infNaNCount ( const mat_type A,
const bool  upperTriangular = false 
) const
inline

Return the number of Inf or NaN entries in the matrix A.

Parameters
A[in] The matrix to check.
upperTriangular[in] If true, only check the upper triangle / trapezoid of A. Otherwise, check all entries of A.

Definition at line 598 of file BelosProjectedLeastSquaresSolver.hpp.

template<class Scalar>
std::pair<bool, std::pair<magnitude_type, magnitude_type> > Belos::details::LocalDenseMatrixOps< Scalar >::isUpperTriangular ( const mat_type A) const
inline

Is the matrix A upper triangular / trapezoidal?

Returns
(is upper triangular?, (squared Frobenius norm of strict lower triangle, squared Frobenius norm of the whole matrix))

Definition at line 625 of file BelosProjectedLeastSquaresSolver.hpp.

template<class Scalar>
std::pair<bool, std::pair<magnitude_type, magnitude_type> > Belos::details::LocalDenseMatrixOps< Scalar >::isUpperHessenberg ( const mat_type A) const
inline

Is the matrix A upper Hessenberg?

Returns
(is upper Hessenberg?, (squared Frobenius norm of the part of A that should be zero if A is upper Hessenberg, squared Frobenius norm of the whole matrix))

Definition at line 658 of file BelosProjectedLeastSquaresSolver.hpp.

template<class Scalar>
void Belos::details::LocalDenseMatrixOps< Scalar >::ensureUpperTriangular ( const mat_type A,
const char *const  matrixName 
) const
inline

Throw an exception if A is not upper triangular / trapezoidal.

Parameters
A[in] The matrix to test.
matrixName[in] Name of the matrix. Used only to make the exception message more informative.

Definition at line 691 of file BelosProjectedLeastSquaresSolver.hpp.

template<class Scalar>
void Belos::details::LocalDenseMatrixOps< Scalar >::ensureUpperHessenberg ( const mat_type A,
const char *const  matrixName 
) const
inline

Throw an exception if A is not (strictly) upper Hessenberg.

Parameters
A[in] The matrix to test.
matrixName[in] Name of the matrix. Used only to make the exception message more informative.

Definition at line 712 of file BelosProjectedLeastSquaresSolver.hpp.

template<class Scalar>
void Belos::details::LocalDenseMatrixOps< Scalar >::ensureUpperHessenberg ( const mat_type A,
const char *const  matrixName,
const magnitude_type  relativeTolerance 
) const
inline

Throw an exception if A is not "approximately" upper Hessenberg.

"Approximately" in this case means that the Frobenius norm of the part of A that should be zero, divided by the Frobenius norm of all of A, is less than or equal to the given relative tolerance.

Parameters
A[in] The matrix to test.
matrixName[in] Name of the matrix. Used only to make the exception message more informative.
relativeTolerance[in] Amount by which we allow the norm of the part of A that should be zero to deviate from zero, relative to the norm of A.

Definition at line 742 of file BelosProjectedLeastSquaresSolver.hpp.

template<class Scalar>
void Belos::details::LocalDenseMatrixOps< Scalar >::ensureMinimumDimensions ( const mat_type A,
const char *const  matrixName,
const int  minNumRows,
const int  minNumCols 
) const
inline

Ensure that the matrix A is at least minNumRows by minNumCols.

If A has fewer rows than minNumRows, or fewer columns than minNumCols, this method throws an informative exception.

Parameters
A[in] The matrix whose dimensions to check.
matrixName[in] Name of the matrix; used to make the exception message more informative.
minNumRows[in] Minimum number of rows allowed in A.
minNumCols[in] Minimum number of columns allowed in A.

Definition at line 775 of file BelosProjectedLeastSquaresSolver.hpp.

template<class Scalar>
void Belos::details::LocalDenseMatrixOps< Scalar >::ensureEqualDimensions ( const mat_type A,
const char *const  matrixName,
const int  numRows,
const int  numCols 
) const
inline

Ensure that the matrix A is exactly numRows by numCols.

If A has a different number of rows than numRows, or a different number of columns than numCols, this method throws an informative exception.

Parameters
A[in] The matrix whose dimensions to check.
matrixName[in] Name of the matrix; used to make the exception message more informative.
numRows[in] Number of rows that A must have.
numCols[in] Number of columns that A must have.

Definition at line 800 of file BelosProjectedLeastSquaresSolver.hpp.


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

Generated on Thu Apr 25 2024 09:27:25 for Belos by doxygen 1.8.5