Belos Package Browser (Single Doxygen Collection)
Development
|
Methods for solving GMRES' projected least-squares problem. 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 | |
ProjectedLeastSquaresSolver (std::ostream &warnStream, const ERobustness defaultRobustness=ROBUSTNESS_NONE) | |
Constructor. More... | |
magnitude_type | updateColumn (ProjectedLeastSquaresProblem< Scalar > &problem, const int curCol) |
Update column curCol of the projected least-squares problem. More... | |
magnitude_type | updateColumns (ProjectedLeastSquaresProblem< Scalar > &problem, const int startCol, const int endCol) |
Update columns [startCol,endCol] of the projected least-squares problem. More... | |
void | solve (ProjectedLeastSquaresProblem< Scalar > &problem, const int curCol) |
Solve the projected least-squares problem. More... | |
std::pair< int, bool > | solveUpperTriangularSystem (Teuchos::ESide side, mat_type &X, const mat_type &R, const mat_type &B) |
Solve the given square upper triangular linear system(s). More... | |
std::pair< int, bool > | solveUpperTriangularSystem (Teuchos::ESide side, mat_type &X, const mat_type &R, const mat_type &B, const ERobustness robustness) |
Solve the given square upper triangular linear system(s). More... | |
std::pair< int, bool > | solveUpperTriangularSystemInPlace (Teuchos::ESide side, mat_type &X, const mat_type &R) |
Solve square upper triangular linear system(s) in place. More... | |
std::pair< int, bool > | solveUpperTriangularSystemInPlace (Teuchos::ESide side, mat_type &X, const mat_type &R, const ERobustness robustness) |
Solve square upper triangular linear system(s) in place. More... | |
bool | testGivensRotations (std::ostream &out) |
Test Givens rotations. More... | |
bool | testUpdateColumn (std::ostream &out, const int numCols, const bool testBlockGivens=false, const bool extraVerbose=false) |
Test update and solve using Givens rotations. More... | |
bool | testTriangularSolves (std::ostream &out, const int testProblemSize, const ERobustness robustness, const bool verbose=false) |
Test upper triangular solves. More... | |
Private Types | |
typedef Teuchos::ScalarTraits < scalar_type > | STS |
typedef Teuchos::ScalarTraits < magnitude_type > | STM |
typedef Teuchos::BLAS< int, scalar_type > | blas_type |
typedef Teuchos::LAPACK< int, scalar_type > | lapack_type |
Private Member Functions | |
int | solveLeastSquaresUsingSVD (mat_type &A, mat_type &X) |
Solve using the SVD. More... | |
void | solveGivens (mat_type &y, mat_type &R, const mat_type &z, const int curCol) |
Solve the projected least-squares problem, assuming Givens rotations updates. More... | |
void | makeRandomProblem (mat_type &H, mat_type &z) |
Make a random projected least-squares problem. More... | |
void | computeGivensRotation (const Scalar &x, const Scalar &y, Scalar &theCosine, Scalar &theSine, Scalar &result) |
Compute the Givens rotation corresponding to [x; y]. More... | |
void | singularValues (const mat_type &A, Teuchos::ArrayView< magnitude_type > sigmas) |
Compute the singular values of A. Store them in the sigmas array. More... | |
std::pair< magnitude_type, magnitude_type > | extremeSingularValues (const mat_type &A) |
The (largest, smallest) singular values of the given matrix. More... | |
magnitude_type | leastSquaresConditionNumber (const mat_type &A, const mat_type &b, const magnitude_type &residualNorm) |
Normwise 2-norm condition number of the least-squares problem. More... | |
magnitude_type | leastSquaresResidualNorm (const mat_type &A, const mat_type &x, const mat_type &b) |
(Frobenius norm if b has more than one column). More... | |
magnitude_type | solutionError (const mat_type &x_approx, const mat_type &x_exact) |
||x_approx - x_exact||_2 // ||x_exact||_2. More... | |
magnitude_type | updateColumnGivens (const mat_type &H, mat_type &R, mat_type &y, mat_type &z, Teuchos::ArrayView< scalar_type > theCosines, Teuchos::ArrayView< scalar_type > theSines, const int curCol) |
Update current column using Givens rotations. More... | |
magnitude_type | solveLapack (const mat_type &H, mat_type &R, mat_type &y, mat_type &z, const int curCol) |
Solve the least-squares problem using LAPACK. More... | |
magnitude_type | updateColumnsGivens (const mat_type &H, mat_type &R, mat_type &y, mat_type &z, Teuchos::ArrayView< scalar_type > theCosines, Teuchos::ArrayView< scalar_type > theSines, const int startCol, const int endCol) |
Update columns [startCol,endCol] of the projected least-squares problem. More... | |
magnitude_type | updateColumnsGivensBlock (const mat_type &H, mat_type &R, mat_type &y, mat_type &z, Teuchos::ArrayView< scalar_type > theCosines, Teuchos::ArrayView< scalar_type > theSines, const int startCol, const int endCol) |
Update columns [startCol,endCol] of the projected least-squares problem. More... | |
Private Attributes | |
std::ostream & | warn_ |
Stream to which to output warnings. More... | |
ERobustness | defaultRobustness_ |
Default robustness level, for things like triangular solves. More... | |
Methods for solving GMRES' projected least-squares problem.
Scalar | The type of the matrix and vector entries in the least-squares problem. |
Expected use of this class:
ProjectedLeastSquaresProblem
struct instance to store the projected problem in your GMRES solver.You can defer Step 4 as long as you want. Step 4 must always follow Step 3.
Purposes of this class:
"Robust" here means regularizing the least-squares solve, so that the solution is well-defined even if the problem is ill-conditioned. Many distributed-memory iterative solvers, including those in Belos, currently solve the projected least-squares problem redundantly on different processes. If those processes are heterogeneous or implement the BLAS and LAPACK themselves in parallel (via multithreading, for example), then different BLAS or LAPACK calls on different processes may result in different answers. The answers may be significantly different if the projected problem is singular or ill-conditioned. This is bad because GMRES variants use the projected problem's solution as the coefficients for the solution update. The solution update coefficients must be (almost nearly) the same on all processes. Regularizing the projected problem is one way to ensure that different processes compute (almost) the same solution.
Definition at line 963 of file BelosProjectedLeastSquaresSolver.hpp.
Belos::details::ProjectedLeastSquaresSolver< Scalar >::scalar_type |
The template parameter of this class.
The type of the matrix and vector entries in the projected least-squares problem to solve, and the type of the resulting GMRES solution update coefficients.
Definition at line 971 of file BelosProjectedLeastSquaresSolver.hpp.
Belos::details::ProjectedLeastSquaresSolver< Scalar >::magnitude_type |
The type of the magnitude of a scalar_type
value.
If scalar_type is complex-valued, then magnitude_type is real.
Definition at line 976 of file BelosProjectedLeastSquaresSolver.hpp.
Belos::details::ProjectedLeastSquaresSolver< Scalar >::mat_type |
The type of a dense matrix (or vector) of scalar_type
.
Definition at line 979 of file BelosProjectedLeastSquaresSolver.hpp.
|
private |
Definition at line 982 of file BelosProjectedLeastSquaresSolver.hpp.
|
private |
Definition at line 983 of file BelosProjectedLeastSquaresSolver.hpp.
|
private |
Definition at line 984 of file BelosProjectedLeastSquaresSolver.hpp.
|
private |
Definition at line 985 of file BelosProjectedLeastSquaresSolver.hpp.
|
inline |
Constructor.
warnStream | [out] Stream to which to output warnings. Set to a Teuchos::oblackholestream if you don't want to display warnings. |
defaultRobustness | [in] Default robustness level for operations like triangular solves. For example, at a low robustness level, triangular solves might fail if the triangular matrix is singular or ill-conditioned in a particular way. At a high robustness level, triangular solves will always succeed, but may only be solved in a least-squares sense. |
Definition at line 1001 of file BelosProjectedLeastSquaresSolver.hpp.
|
inline |
Update column curCol of the projected least-squares problem.
The upper Hessenberg matrix H is read but not touched. The R factor, the cosines and sines, and the right-hand side z are updated. This method does not compute the solution of the least-squares problem; call solve()
for that.
problem | [in/out] The projected least-squares problem. |
curCol | [in] Zero-based index of the current column to update. |
Definition at line 1021 of file BelosProjectedLeastSquaresSolver.hpp.
|
inline |
Update columns [startCol,endCol] of the projected least-squares problem.
The upper Hessenberg matrix H is read but not touched. The R factor, the cosines and sines, and the right-hand side z are updated. This method does not compute the solution of the least-squares problem; call solve()
for that.
problem | [in/out] The projected least-squares problem. |
startCol | [in] Zero-based index of the first column to update. |
endCol | [in] Zero-based index of the last column (inclusive) to update. |
Definition at line 1043 of file BelosProjectedLeastSquaresSolver.hpp.
|
inline |
Solve the projected least-squares problem.
Call this method only after calling updateColumn() or updateColumns(). If you call updateColumn(), use the same column index when calling this method as you did when calling updateColumn(). If you call updateColumns(), use its endCol
argument as the column index for calling this method.
problem | [in/out] The projected least-squares problem. |
curCol | [in] Zero-based index of the most recently updated column of the least-squares problem. |
Definition at line 1066 of file BelosProjectedLeastSquaresSolver.hpp.
|
inline |
Solve the given square upper triangular linear system(s).
This method does the same thing as the five-argument overload, but uses the default robustness level.
Definition at line 1077 of file BelosProjectedLeastSquaresSolver.hpp.
|
inline |
Solve the given square upper triangular linear system(s).
side | [in] Whether to solve with the triangular matrix R on the left side or right side of the matrix X of unknowns. side = Teuchos::LEFT_SIDE or Teuchos::RIGHT_SIDE. |
X | [out] The matrix of unknowns. (We consider a matrix with one column a vector.) |
R | [in] The square upper triangular matrix. If R has more rows than columns, this method uses the number of columns of R as the dimension of the square matrix, and ignores the "extra" rows of R. If R has more columns than rows, this method throws an exception. |
B | [in] The matrix of right-hand sides. B must have at least as many columns as X, otherwise this method will throw an exception. If B has more columns than X, those columns will be ignored. |
robustness | [in] The robustness level to use for the solve. ROBUSTNESS_NONE is fastest but does not attempt to detect rank deficiency in R. ROBUSTNESS_LOTS is slowest, but can "solve" the problem (in a least-squares minimum-norm-solution sense) even if R is singular. ROBUSTNESS_SOME may do something in between. |
Depending on the specified robustness level, this method may also attempt to compute the numerical rank and detect rank deficiency in R. The first return value is the computed numerical rank, and the second return value is a Boolean indicating whether any rank deficiency was detected. If the method does not attempt to compute the numerical rank or detect rank deficiency, it will return (N, false), where N is the number of columns in R.
Definition at line 1124 of file BelosProjectedLeastSquaresSolver.hpp.
|
inline |
Solve square upper triangular linear system(s) in place.
This method does the same thing as its four-argument overload, but uses the default robustness level.
Definition at line 1160 of file BelosProjectedLeastSquaresSolver.hpp.
|
inline |
Solve square upper triangular linear system(s) in place.
This is the "in-place" version of solveUpperTriangularSystem()
. The difference between that method and this one is that this method assumes that the right-hand side(s) is/are already stored in X on input. It then overwrites X with the solution.
Definition at line 1175 of file BelosProjectedLeastSquaresSolver.hpp.
|
inline |
Test Givens rotations.
This routine tests both computing Givens rotations (via computeGivensRotation()
) and applying them.
out | [out] Stream to which to write test output. |
Definition at line 1326 of file BelosProjectedLeastSquaresSolver.hpp.
|
inline |
Test update and solve using Givens rotations.
out | [out] Output stream to which to print results. |
numCols | [in] Number of columns in the projected least-squares problem to test. (The number of rows is one plus the number of columns.) |
testBlockGivens | [in] Whether to test the "block" (i.e., panel) version of the Givens rotations update. |
extraVerbose | [in] Whether to print extra-verbose output (e.g., the test problem and results). |
Test updating and solving the least-squares problem using Givens rotations by comparison against LAPACK's least-squares solver. First generate a random least-squares problem that looks like it comes from GMRES. The matrix is upper Hessenberg, and the right-hand side starts out with the first entry being nonzero with nonnegative real part and zero imaginary part, and all the other entries being zero. Then compare the results of updateColumnGivens()
(applied to each column in turn) and solveGivens() (applied at the end) with the results of solveLapack()
(applied at the end). Print out the results to the given output stream.
Definition at line 1385 of file BelosProjectedLeastSquaresSolver.hpp.
|
inline |
Test upper triangular solves.
out | [out] Stream to which this method writes output. |
testProblemSize | [in] Dimension (number of rows and columns) of the upper triangular matrices tested by this method. |
robustness | [in] Robustness level of the upper triangular solves. |
verbose | [in] Whether to print more verbose output. |
Definition at line 1626 of file BelosProjectedLeastSquaresSolver.hpp.
|
inlineprivate |
Solve using the SVD.
X | [in/out] On input: the right-hand side(s) B. On output: the solution to the least-squares problem. |
A | [in/out] On input: the matrix in the least-squares problem. On output, A will be overwritten with factorization information. |
Definition at line 1737 of file BelosProjectedLeastSquaresSolver.hpp.
|
inlineprivate |
Solve the projected least-squares problem, assuming Givens rotations updates.
Call this method after invoking either updateColumnGivens()
with the same curCol, or updateColumnsGivens()
with curCol = endCol.
Definition at line 1812 of file BelosProjectedLeastSquaresSolver.hpp.
|
inlineprivate |
Make a random projected least-squares problem.
Definition at line 1832 of file BelosProjectedLeastSquaresSolver.hpp.
|
inlineprivate |
Compute the Givens rotation corresponding to [x; y].
The result of applying the rotation is [result; 0].
Definition at line 1877 of file BelosProjectedLeastSquaresSolver.hpp.
|
inlineprivate |
Compute the singular values of A. Store them in the sigmas array.
Definition at line 1906 of file BelosProjectedLeastSquaresSolver.hpp.
|
inlineprivate |
The (largest, smallest) singular values of the given matrix.
We use these for computing the 2-norm condition number of the matrix A. We separate out the singular values rather than returning their quotient, so that you can see the value of the largest singular value, even if the smallest singular value is zero.
Definition at line 1966 of file BelosProjectedLeastSquaresSolver.hpp.
|
inlineprivate |
Normwise 2-norm condition number of the least-squares problem.
A | [in] The matrix in the least-squares problem. |
b | [in] The right-hand side in the least-squares problem. |
residualNorm | [in] Residual norm from solving the least-squares problem using a known good method. |
For details on the condition number formula, see Section 3.3 of J. W. Demmel, "Applied Numerical Linear Algebra," SIAM Press.
Definition at line 1989 of file BelosProjectedLeastSquaresSolver.hpp.
|
inlineprivate |
(Frobenius norm if b has more than one column).
Definition at line 2037 of file BelosProjectedLeastSquaresSolver.hpp.
|
inlineprivate |
||x_approx - x_exact||_2 // ||x_exact||_2.
Use the Frobenius norm if more than one column. Don't scale if ||x_exact|| == 0.
Definition at line 2055 of file BelosProjectedLeastSquaresSolver.hpp.
|
inlineprivate |
Update current column using Givens rotations.
Update the current column of the projected least-squares problem, using Givens rotations. This updates the QR factorization of the upper Hessenberg matrix H. The resulting R factor is stored in the matrix R. The Q factor is stored implicitly in the list of cosines and sines, representing the Givens rotations applied to the problem. These Givens rotations are also applied to the right-hand side z.
H | [in] The upper Hessenberg matrix from GMRES. We only view H( 1:curCol+2, 1:curCol+1 ). It's copied and not overwritten, so as not to disturb any backtracking or other features of GMRES. |
R | [in/out] Matrix with the same dimensions as H, used for storing the incrementally computed upper triangular factor from the QR factorization of H. |
y | [out] Vector (one-column matrix) with the same number of rows as H. On output: the solution of the projected least-squares problem. The vector should have one more entry than necessary for the solution, because of the way we solve the least-squares problem. |
z | [in/out] Vector (one-column matrix) with the same number of rows as H. On input: the current right-hand side of the projected least-squares problem. On output: the updated right-hand side. |
theCosines | [in/out] On input: Array of cosines from the previously computed Givens rotations. On output: the same cosines, with the new Givens rotation's cosine appended. |
theSines | [in/out] On input: Array of sines from the previously computed Givens rotations. On output: the same sines, with the new Givens rotation's sine appended. |
curCol | [in] Index of the current column to update. |
Definition at line 2121 of file BelosProjectedLeastSquaresSolver.hpp.
|
inlineprivate |
Solve the least-squares problem using LAPACK.
In particular, this method uses LAPACK's _GELS dense least-squares solver. _GELS assumes that the input matrix has full rank, so this method might fail. This method is inefficient compared to the Givens rotations approach, but useful for testing the latter.
H | [in] The upper Hessenberg matrix from GMRES. We only view H( 1:curCol+2, 1:curCol+1 ). It's copied and not overwritten, so as not to disturb any backtracking or other features of GMRES. |
R | [out] Matrix with the same dimensions as H. On output: the upper triangle of R contains the R factor in the QR factorization of H. |
y | [out] Vector (one-column matrix) with the same number of rows as H. On output: the solution of the projected least-squares problem. The vector should have one more entry than necessary for the solution, because of the way we solve the least-squares problem. |
z | [in] Vector (one-column matrix) with the same number of rows as H. On input: the current right-hand side of the projected least-squares problem. |
curCol | [in] Index of the current column of H. |
Definition at line 2258 of file BelosProjectedLeastSquaresSolver.hpp.
|
inlineprivate |
Update columns [startCol,endCol] of the projected least-squares problem.
This method implements a "left-looking panel QR factorization" of the upper Hessenberg matrix in the projected least-squares problem. It's "left-looking" because we don't updating anything to the right of columns [startCol, endCol], which is the "panel."
Definition at line 2346 of file BelosProjectedLeastSquaresSolver.hpp.
|
inlineprivate |
Update columns [startCol,endCol] of the projected least-squares problem.
This method implements a "left-looking panel QR factorization" of the upper Hessenberg matrix in the projected least-squares problem. It's "left-looking" because we don't updating anything to the right of columns [startCol, endCol], which is the "panel."
Definition at line 2379 of file BelosProjectedLeastSquaresSolver.hpp.
|
private |
Stream to which to output warnings.
Definition at line 1716 of file BelosProjectedLeastSquaresSolver.hpp.
|
private |
Default robustness level, for things like triangular solves.
Definition at line 1722 of file BelosProjectedLeastSquaresSolver.hpp.