|
Belos Version of the Day
|
Low-level operations on non-distributed dense matrices. More...
#include <BelosProjectedLeastSquaresSolver.hpp>
Public Types | |
| typedef Scalar | scalar_type |
| The template parameter of this class. | |
| typedef Teuchos::ScalarTraits < Scalar >::magnitudeType | magnitude_type |
The type of the magnitude of a scalar_type value. | |
| typedef Teuchos::SerialDenseMatrix < int, Scalar > | mat_type |
The type of a dense matrix (or vector) of scalar_type. | |
Public Member Functions | |
| void | conjugateTranspose (mat_type &A_star, const mat_type &A) const |
| A_star := (conjugate) transpose of A. | |
| void | conjugateTransposeOfUpperTriangular (mat_type &L, const mat_type &R) const |
| L := (conjugate) transpose of R (upper triangular). | |
| void | zeroOutStrictLowerTriangle (mat_type &A) const |
| Zero out everything below the diagonal of A. | |
| 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]. | |
| void | matScale (mat_type &A, const scalar_type &alpha) const |
| A := alpha * A. | |
| void | axpy (mat_type &Y, const scalar_type &alpha, const mat_type &X) const |
| Y := Y + alpha * X. | |
| void | matAdd (mat_type &A, const mat_type &B) const |
| A := A + B. | |
| void | matSub (mat_type &A, const mat_type &B) const |
| A := A - B. | |
| void | rightUpperTriSolve (mat_type &B, const mat_type &R) const |
| In Matlab notation: B = B / R, where R is upper triangular. | |
| 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. | |
| int | infNaNCount (const mat_type &A, const bool upperTriangular=false) const |
| Return the number of Inf or NaN entries in the matrix A. | |
| std::pair< bool, std::pair < magnitude_type, magnitude_type > > | isUpperTriangular (const mat_type &A) const |
| Is the matrix A upper triangular / trapezoidal? | |
| std::pair< bool, std::pair < magnitude_type, magnitude_type > > | isUpperHessenberg (const mat_type &A) const |
| Is the matrix A upper Hessenberg? | |
| void | ensureUpperTriangular (const mat_type &A, const char *const matrixName) const |
| Throw an exception if A is not upper triangular / trapezoidal. | |
| void | ensureUpperHessenberg (const mat_type &A, const char *const matrixName) const |
| Throw an exception if A is not (strictly) upper Hessenberg. | |
| 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. | |
| 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. | |
| 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. | |
Low-level operations on non-distributed dense matrices.
This class provides a convenient wrapper around some BLAS operations, operating on non-distributed (hence "local") dense matrices.
Definition at line 348 of file BelosProjectedLeastSquaresSolver.hpp.
| Belos::details::LocalDenseMatrixOps< Scalar >::scalar_type |
The template parameter of this class.
Definition at line 352 of file BelosProjectedLeastSquaresSolver.hpp.
| Belos::details::LocalDenseMatrixOps< Scalar >::magnitude_type |
The type of the magnitude of a scalar_type value.
Definition at line 355 of file BelosProjectedLeastSquaresSolver.hpp.
| Belos::details::LocalDenseMatrixOps< Scalar >::mat_type |
The type of a dense matrix (or vector) of scalar_type.
Definition at line 358 of file BelosProjectedLeastSquaresSolver.hpp.
| 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 369 of file BelosProjectedLeastSquaresSolver.hpp.
| 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 380 of file BelosProjectedLeastSquaresSolver.hpp.
| void Belos::details::LocalDenseMatrixOps< Scalar >::zeroOutStrictLowerTriangle | ( | mat_type & | A | ) | const [inline] |
Zero out everything below the diagonal of A.
Definition at line 393 of file BelosProjectedLeastSquaresSolver.hpp.
| 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.
Definition at line 413 of file BelosProjectedLeastSquaresSolver.hpp.
| void Belos::details::LocalDenseMatrixOps< Scalar >::matScale | ( | mat_type & | A, |
| const scalar_type & | alpha | ||
| ) | const [inline] |
A := alpha * A.
Definition at line 434 of file BelosProjectedLeastSquaresSolver.hpp.
| 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 458 of file BelosProjectedLeastSquaresSolver.hpp.
| void Belos::details::LocalDenseMatrixOps< Scalar >::matAdd | ( | mat_type & | A, |
| const mat_type & | B | ||
| ) | const [inline] |
A := A + B.
Definition at line 478 of file BelosProjectedLeastSquaresSolver.hpp.
| void Belos::details::LocalDenseMatrixOps< Scalar >::matSub | ( | mat_type & | A, |
| const mat_type & | B | ||
| ) | const [inline] |
A := A - B.
Definition at line 507 of file BelosProjectedLeastSquaresSolver.hpp.
| 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 539 of file BelosProjectedLeastSquaresSolver.hpp.
| 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 561 of file BelosProjectedLeastSquaresSolver.hpp.
| 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.
| 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 601 of file BelosProjectedLeastSquaresSolver.hpp.
| 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?
Definition at line 628 of file BelosProjectedLeastSquaresSolver.hpp.
| 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?
Definition at line 661 of file BelosProjectedLeastSquaresSolver.hpp.
| 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.
| A | [in] The matrix to test. |
| matrixName | [in] Name of the matrix. Used only to make the exception message more informative. |
Definition at line 694 of file BelosProjectedLeastSquaresSolver.hpp.
| 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.
| A | [in] The matrix to test. |
| matrixName | [in] Name of the matrix. Used only to make the exception message more informative. |
Definition at line 715 of file BelosProjectedLeastSquaresSolver.hpp.
| 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.
| 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 745 of file BelosProjectedLeastSquaresSolver.hpp.
| 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.
| 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 778 of file BelosProjectedLeastSquaresSolver.hpp.
| 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.
| 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 803 of file BelosProjectedLeastSquaresSolver.hpp.
1.7.4