42 #ifndef __Belos_ProjectedLeastSquaresSolver_hpp
43 #define __Belos_ProjectedLeastSquaresSolver_hpp
75 template<
class Scalar>
77 printMatrix (std::ostream& out,
78 const std::string& name,
83 const int numRows =
A.numRows();
84 const int numCols =
A.numCols();
86 out << name <<
" = " << endl <<
'[';
89 for (
int i = 0; i < numRows; ++i) {
96 for (
int i = 0; i < numRows; ++i) {
97 for (
int j = 0; j < numCols; ++j) {
101 }
else if (i < numRows-1) {
112 template<
class Scalar>
114 print (std::ostream& out,
116 const std::string& linePrefix)
120 const int numRows =
A.numRows();
121 const int numCols =
A.numCols();
123 out << linePrefix <<
'[';
126 for (
int i = 0; i < numRows; ++i) {
133 for (
int i = 0; i < numRows; ++i) {
134 for (
int j = 0; j < numCols; ++j) {
136 out << linePrefix <<
" ";
141 }
else if (i < numRows-1) {
147 out << linePrefix <<
']' << endl;
158 template<
class Scalar>
238 H (maxNumIterations+1, maxNumIterations),
239 R (maxNumIterations+1, maxNumIterations),
240 y (maxNumIterations+1, 1),
241 z (maxNumIterations+1, 1),
274 const Scalar initialResidualNorm (beta);
275 z(0,0) = initialResidualNorm;
293 const int maxNumIterations)
299 "ProjectedLeastSquaresProblem::reset: initial "
300 "residual beta = " << beta <<
" < 0.");
302 "ProjectedLeastSquaresProblem::reset: maximum number "
303 "of iterations " << maxNumIterations <<
" <= 0.");
306 const int errcode =
H.
reshape (maxNumIterations+1, maxNumIterations);
308 "Failed to reshape H into a " << (maxNumIterations+1)
309 <<
" x " << maxNumIterations <<
" matrix.");
314 const int errcode =
R.
reshape (maxNumIterations+1, maxNumIterations);
316 "Failed to reshape R into a " << (maxNumIterations+1)
317 <<
" x " << maxNumIterations <<
" matrix.");
322 const int errcode =
y.
reshape (maxNumIterations+1, 1);
324 "Failed to reshape y into a " << (maxNumIterations+1)
325 <<
" x " << 1 <<
" matrix.");
330 const int errcode =
z.
reshape (maxNumIterations+1, 1);
332 "Failed to reshape z into a " << (maxNumIterations+1)
333 <<
" x " << 1 <<
" matrix.");
348 template<
class Scalar>
372 for (
int i = 0; i <
A.numRows(); ++i) {
373 for (
int j = 0; j <
A.numCols(); ++j) {
385 for (
int j = 0; j < N; ++j) {
386 for (
int i = 0; i <= j; ++i) {
396 const int N = std::min (
A.numRows(),
A.numCols());
398 for (
int j = 0; j < N; ++j) {
399 for (
int i = j+1; i <
A.numRows(); ++i) {
430 A_22 =
rcp (
new mat_type (
View,
A, numRows2, numCols2, numRows1, numCols1));
438 const int numRows =
A.numRows();
439 const int numCols =
A.numCols();
441 if (numRows == 0 || numCols == 0) {
444 for (
int j = 0; j < numCols; ++j) {
447 for (
int i = 0; i < numRows; ++i) {
463 const int numRows = Y.
numRows();
464 const int numCols = Y.
numCols();
467 std::invalid_argument,
"Dimensions of X and Y don't "
469 <<
", and Y is " << numRows <<
" x " << numCols <<
".");
470 for (
int j = 0; j < numCols; ++j) {
471 for (
int i = 0; i < numRows; ++i) {
472 Y(i,j) += alpha * X(i,j);
481 const int numRows =
A.numRows();
482 const int numCols =
A.numCols();
485 B.numRows() != numRows ||
B.numCols() != numCols,
486 std::invalid_argument,
487 "matAdd: The input matrices A and B have incompatible dimensions. "
488 "A is " << numRows <<
" x " << numCols <<
", but B is " <<
489 B.numRows () <<
" x " <<
B.numCols () <<
".");
490 if (numRows == 0 || numCols == 0) {
493 for (
int j = 0; j < numCols; ++j) {
497 for (
int i = 0; i < numRows; ++i) {
508 const int numRows =
A.numRows();
509 const int numCols =
A.numCols();
512 B.numRows() != numRows ||
B.numCols() != numCols,
513 std::invalid_argument,
514 "matSub: The input matrices A and B have incompatible dimensions. "
515 "A is " << numRows <<
" x " << numCols <<
", but B is " <<
516 B.numRows () <<
" x " <<
B.numCols () <<
".");
517 if (numRows == 0 || numCols == 0) {
520 for (
int j = 0; j < numCols; ++j) {
524 for (
int i = 0; i < numRows; ++i) {
540 std::invalid_argument,
541 "rightUpperTriSolve: R and B have incompatible "
542 "dimensions. B has " <<
B.numCols() <<
" columns, "
543 "but R has " << R.
numRows() <<
" rows.");
549 B.values(),
B.stride());
567 std::invalid_argument,
568 "matMatMult: The input matrices A and B have "
569 "incompatible dimensions. A is " <<
A.numRows()
570 <<
" x " <<
A.numCols() <<
", but B is "
571 <<
B.numRows() <<
" x " <<
B.numCols() <<
".");
573 std::invalid_argument,
574 "matMatMult: The input matrix A and the output "
575 "matrix C have incompatible dimensions. A has "
576 <<
A.numRows() <<
" rows, but C has " <<
C.numRows()
579 std::invalid_argument,
580 "matMatMult: The input matrix B and the output "
581 "matrix C have incompatible dimensions. B has "
582 <<
B.numCols() <<
" columns, but C has "
583 <<
C.numCols() <<
" columns.");
586 alpha,
A.values(),
A.stride(),
B.values(),
B.stride(),
587 beta,
C.values(),
C.stride());
601 for (
int j = 0; j <
A.numCols(); ++j) {
602 if (upperTriangular) {
603 for (
int i = 0; i <= j && i <
A.numRows(); ++i) {
609 for (
int i = 0; i <
A.numRows(); ++i) {
624 std::pair<bool, std::pair<magnitude_type, magnitude_type> >
631 for (
int j = 0; j <
A.numCols(); ++j) {
635 for (
int i = 0; i <= j && i <
A.numRows(); ++i) {
637 upperTri += A_ij_mag * A_ij_mag;
640 for (
int i = j+1; i <
A.numRows(); ++i) {
642 lowerTri += A_ij_mag * A_ij_mag;
648 return std::make_pair (count == 0, std::make_pair (lowerTri, upperTri));
657 std::pair<bool, std::pair<magnitude_type, magnitude_type> >
664 for (
int j = 0; j <
A.numCols(); ++j) {
668 for (
int i = 0; i <= j+1 && i <
A.numRows(); ++i) {
670 upper += A_ij_mag * A_ij_mag;
673 for (
int i = j+2; i <
A.numRows(); ++i) {
675 lower += A_ij_mag * A_ij_mag;
681 return std::make_pair (count == 0, std::make_pair (lower, upper));
692 const char*
const matrixName)
const
694 std::pair<bool, std::pair<magnitude_type, magnitude_type> > result =
698 "The " <<
A.numRows() <<
" x " <<
A.numCols()
699 <<
" matrix " << matrixName <<
" is not upper "
700 "triangular. ||tril(A)||_F = "
701 << result.second.first <<
" and ||A||_F = "
702 << result.second.second <<
".");
713 const char*
const matrixName)
const
715 std::pair<bool, std::pair<magnitude_type, magnitude_type> > result =
719 "The " <<
A.numRows() <<
" x " <<
A.numCols()
720 <<
" matrix " << matrixName <<
" is not upper "
721 "triangular. ||tril(A(2:end, :))||_F = "
722 << result.second.first <<
" and ||A||_F = "
723 << result.second.second <<
".");
743 const char*
const matrixName,
746 std::pair<bool, std::pair<magnitude_type, magnitude_type> > result =
752 result.second.first :
753 result.second.first / result.second.second);
755 "The " <<
A.numRows() <<
" x " <<
A.numCols()
756 <<
" matrix " << matrixName <<
" is not upper "
757 "triangular. ||tril(A(2:end, :))||_F "
758 << (result.second.second ==
STM::zero() ?
"" :
" / ||A||_F")
759 <<
" = " << err <<
" > " << relativeTolerance <<
".");
776 const char*
const matrixName,
777 const int minNumRows,
778 const int minNumCols)
const
781 std::invalid_argument,
782 "The matrix " << matrixName <<
" is " <<
A.numRows()
783 <<
" x " <<
A.numCols() <<
", and therefore does not "
784 "satisfy the minimum dimensions " << minNumRows
785 <<
" x " << minNumCols <<
".");
801 const char*
const matrixName,
803 const int numCols)
const
806 std::invalid_argument,
807 "The matrix " << matrixName <<
" is supposed to be "
808 << numRows <<
" x " << numCols <<
", but is "
809 <<
A.numRows() <<
" x " <<
A.numCols() <<
" instead.");
843 const char* strings[] = {
"None",
"Some",
"Lots"};
845 std::invalid_argument,
846 "Invalid enum value " << x <<
".");
847 return std::string (strings[x]);
854 const char* strings[] = {
"None",
"Some",
"Lots"};
856 if (x == strings[r]) {
857 return static_cast<ERobustness> (r);
861 "Invalid robustness string " << x <<
".");
877 using Teuchos::stringToIntegralParameterEntryValidator;
884 docs[0] =
"Use the BLAS' triangular solve. This may result in Inf or "
885 "NaN output if the triangular matrix is rank deficient.";
886 docs[1] =
"Robustness somewhere between \"None\" and \"Lots\".";
887 docs[2] =
"Solve the triangular system in a least-squares sense, using "
888 "an SVD-based algorithm. This will always succeed, though the "
889 "solution may not make sense for GMRES.";
894 const std::string pname (
"Robustness of Projected Least-Squares Solve");
896 return stringToIntegralParameterEntryValidator<ERobustness> (strs, docs,
962 template<
class Scalar>
1076 std::pair<int, bool>
1123 std::pair<int, bool>
1131 "The output X and right-hand side B have different "
1132 "numbers of rows. X has " << X.
numRows() <<
" rows"
1133 ", and B has " <<
B.numRows() <<
" rows.");
1139 "The output X has more columns than the "
1140 "right-hand side B. X has " << X.
numCols()
1141 <<
" columns and B has " <<
B.numCols()
1159 std::pair<int, bool>
1174 std::pair<int, bool>
1189 "The input matrix R has fewer columns than rows. "
1190 "R is " << M <<
" x " << N <<
".");
1196 "The input/output matrix X has only "
1197 << X.
numRows() <<
" rows, but needs at least "
1198 << N <<
" rows to match the matrix for a "
1199 "left-side solve R \\ X.");
1202 "The input/output matrix X has only "
1203 << X.
numCols() <<
" columns, but needs at least "
1204 << N <<
" columns to match the matrix for a "
1205 "right-side solve X / R.");
1209 std::invalid_argument,
1210 "Invalid robustness value " << robustness <<
".");
1218 "There " << (count != 1 ?
"are" :
"is")
1219 <<
" " << count <<
" Inf or NaN entr"
1220 << (count != 1 ?
"ies" :
"y")
1221 <<
" in the upper triangle of R.");
1224 "There " << (count != 1 ?
"are" :
"is")
1225 <<
" " << count <<
" Inf or NaN entr"
1226 << (count != 1 ?
"ies" :
"y") <<
" in the "
1227 "right-hand side(s) X.");
1232 bool foundRankDeficiency =
false;
1260 warn_ <<
"Upper triangular solve: Found Infs and/or NaNs in the "
1261 "solution after using the fast algorithm. Retrying using a more "
1262 "robust algorithm." << std::endl;
1303 foundRankDeficiency =
true;
1308 "Should never get here! Invalid robustness value "
1309 << robustness <<
". Please report this bug to the "
1310 "Belos developers.");
1312 return std::make_pair (rank, foundRankDeficiency);
1330 out <<
"Testing Givens rotations:" << endl;
1333 out <<
" x = " << x <<
", y = " << y << endl;
1335 Scalar theCosine, theSine, result;
1338 out <<
"-- After computing rotation:" << endl;
1339 out <<
"---- cos,sin = " << theCosine <<
"," << theSine << endl;
1340 out <<
"---- x = " << x <<
", y = " << y
1341 <<
", result = " << result << endl;
1343 blas.
ROT (1, &x, 1, &y, 1, &theCosine, &theSine);
1344 out <<
"-- After applying rotation:" << endl;
1345 out <<
"---- cos,sin = " << theCosine <<
"," << theSine << endl;
1346 out <<
"---- x = " << x <<
", y = " << y << endl;
1387 const bool testBlockGivens=
false,
1388 const bool extraVerbose=
false)
1394 "numCols = " << numCols <<
" <= 0.");
1395 const int numRows = numCols + 1;
1400 mat_type R_givens (numRows, numCols);
1403 Array<Scalar> theCosines (numCols);
1404 Array<Scalar> theSines (numCols);
1406 mat_type R_blockGivens (numRows, numCols);
1407 mat_type y_blockGivens (numRows, 1);
1408 mat_type z_blockGivens (numRows, 1);
1409 Array<Scalar> blockCosines (numCols);
1410 Array<Scalar> blockSines (numCols);
1411 const int panelWidth = std::min (3, numCols);
1413 mat_type R_lapack (numRows, numCols);
1420 printMatrix<Scalar> (out,
"H", H);
1421 printMatrix<Scalar> (out,
"z", z);
1427 if (testBlockGivens) {
1428 z_blockGivens.
assign (z);
1438 for (
int curCol = 0; curCol < numCols; ++curCol) {
1440 theCosines, theSines, curCol);
1442 solveGivens (y_givens, R_givens, z_givens, numCols-1);
1447 if (testBlockGivens) {
1448 const bool testBlocksAtATime =
true;
1449 if (testBlocksAtATime) {
1451 for (
int startCol = 0; startCol < numCols; startCol += panelWidth) {
1452 int endCol = std::min (startCol + panelWidth - 1, numCols - 1);
1453 residualNormBlockGivens =
1455 blockCosines, blockSines, startCol, endCol);
1461 for (
int startCol = 0; startCol < numCols; ++startCol) {
1462 residualNormBlockGivens =
1464 blockCosines, blockSines, startCol, startCol);
1471 solveGivens (y_blockGivens, R_blockGivens, z_blockGivens, numCols-1);
1476 solveLapack (H, R_lapack, y_lapack, z_lapack, numCols-1);
1500 const magnitude_type blockGivensSolutionError = testBlockGivens ?
1508 mat_type R_factorFromGivens (numCols, numCols);
1509 mat_type R_factorFromBlockGivens (numCols, numCols);
1510 mat_type R_factorFromLapack (numCols, numCols);
1512 for (
int j = 0; j < numCols; ++j) {
1513 for (
int i = 0; i <= j; ++i) {
1514 R_factorFromGivens(i,j) = R_givens(i,j);
1515 if (testBlockGivens) {
1516 R_factorFromBlockGivens(i,j) = R_blockGivens(i,j);
1518 R_factorFromLapack(i,j) = R_lapack(i,j);
1522 printMatrix<Scalar> (out,
"R_givens", R_factorFromGivens);
1523 printMatrix<Scalar> (out,
"y_givens", y_givens_view);
1524 printMatrix<Scalar> (out,
"z_givens", z_givens);
1526 if (testBlockGivens) {
1527 printMatrix<Scalar> (out,
"R_blockGivens", R_factorFromBlockGivens);
1528 printMatrix<Scalar> (out,
"y_blockGivens", y_blockGivens_view);
1529 printMatrix<Scalar> (out,
"z_blockGivens", z_blockGivens);
1532 printMatrix<Scalar> (out,
"R_lapack", R_factorFromLapack);
1533 printMatrix<Scalar> (out,
"y_lapack", y_lapack_view);
1534 printMatrix<Scalar> (out,
"z_lapack", z_lapack);
1540 out <<
"||H||_F = " << H_norm << endl;
1542 out <<
"||H y_givens - z||_2 / ||H||_F = "
1544 if (testBlockGivens) {
1545 out <<
"||H y_blockGivens - z||_2 / ||H||_F = "
1548 out <<
"||H y_lapack - z||_2 / ||H||_F = "
1551 out <<
"||y_givens - y_lapack||_2 / ||y_lapack||_2 = "
1552 << givensSolutionError << endl;
1553 if (testBlockGivens) {
1554 out <<
"||y_blockGivens - y_lapack||_2 / ||y_lapack||_2 = "
1555 << blockGivensSolutionError << endl;
1558 out <<
"Least-squares condition number = "
1559 << leastSquaresCondNum << endl;
1577 wiggleFactor * leastSquaresCondNum;
1579 solutionErrorBoundFactor *
STS::eps();
1580 out <<
"Solution error bound: " << solutionErrorBoundFactor
1581 <<
" * eps = " << solutionErrorBound << endl;
1594 }
else if (givensSolutionError > solutionErrorBound) {
1596 }
else if (testBlockGivens) {
1599 }
else if (blockGivensSolutionError > solutionErrorBound) {
1627 const int testProblemSize,
1629 const bool verbose=
false)
1637 std::ostream& verboseOut = verbose ? out : blackHole;
1639 verboseOut <<
"Testing upper triangular solves" << endl;
1643 verboseOut <<
"-- Generating test matrix" << endl;
1644 const int N = testProblemSize;
1647 for (
int j = 0; j < N; ++j) {
1648 for (
int i = 0; i <= j; ++i) {
1668 verboseOut <<
"-- Solving RX=B" << endl;
1676 verboseOut <<
"---- ||R*X - B||_F = " << Resid.
normFrobenius() << endl;
1677 verboseOut <<
"---- ||R||_F ||X||_F + ||B||_F = "
1692 B_star_copy.
assign (B_star);
1694 verboseOut <<
"-- Solving YR=B^*" << endl;
1699 Resid2.
assign (B_star_copy);
1701 verboseOut <<
"---- ||Y*R - B^*||_F = " << Resid2.
normFrobenius() << endl;
1702 verboseOut <<
"---- ||Y||_F ||R||_F + ||B^*||_F = "
1745 "solveLeastSquaresUsingSVD: The input matrix A "
1746 "contains " << count <<
"Inf and/or NaN entr"
1747 << (count != 1 ?
"ies" :
"y") <<
".");
1750 "solveLeastSquaresUsingSVD: The input matrix X "
1751 "contains " << count <<
"Inf and/or NaN entr"
1752 << (count != 1 ?
"ies" :
"y") <<
".");
1754 const int N = std::min (
A.numRows(),
A.numCols());
1773 Array<magnitude_type> rwork (1);
1775 rwork.resize (std::max (1, 5 * N));
1782 lapack.GELSS (
A.numRows(),
A.numCols(), X.
numCols(),
1785 &lworkScalar, -1, &rwork[0], &info);
1787 "_GELSS workspace query returned INFO = "
1788 << info <<
" != 0.");
1789 const int lwork = static_cast<int> (
STS::real (lworkScalar));
1791 "_GELSS workspace query returned LWORK = "
1792 << lwork <<
" < 0.");
1794 Array<Scalar> work (std::max (1, lwork));
1796 lapack.GELSS (
A.numRows(),
A.numCols(), X.
numCols(),
1799 &work[0], lwork, &rwork[0], &info);
1801 "_GELSS returned INFO = " << info <<
" != 0.");
1817 const int numRows = curCol + 2;
1841 for (
int j = 0; j < H.
numCols(); ++j) {
1842 for (
int i = j+2; i < H.
numRows(); ++i) {
1854 const int numTrials = 1000;
1856 for (
int trial = 0; trial < numTrials && z_init ==
STM::zero(); ++trial) {
1860 "After " << numTrials <<
" trial"
1861 << (numTrials != 1 ?
"s" :
"")
1862 <<
", we were unable to generate a nonzero pseudo"
1863 "random real number. This most likely indicates a "
1864 "broken pseudorandom number generator.");
1885 const bool useLartg =
false;
1890 lapack.
LARTG (x, y, &theCosine, &theSine, &result);
1899 blas.
ROTG (&x_temp, &y_temp, &theCosine, &theSine);
1912 const int numRows =
A.numRows();
1913 const int numCols =
A.numCols();
1915 std::invalid_argument,
1916 "The sigmas array is only of length " << sigmas.
size()
1917 <<
", but must be of length at least "
1918 << std::min (numRows, numCols)
1919 <<
" in order to hold all the singular values of the "
1925 mat_type A_copy (numRows, numCols);
1932 Array<magnitude_type> rwork (std::max (std::min (numRows, numCols) - 1, 1));
1933 lapack.
GESVD (
'N',
'N', numRows, numCols,
1934 A_copy.values(), A_copy.stride(), &sigmas[0],
1935 (Scalar*) NULL, 1, (Scalar*) NULL, 1,
1936 &lworkScalar, -1, &rwork[0], &info);
1939 "LAPACK _GESVD workspace query failed with INFO = "
1941 const int lwork = static_cast<int> (
STS::real (lworkScalar));
1943 "LAPACK _GESVD workspace query returned LWORK = "
1944 << lwork <<
" < 0.");
1950 lapack.
GESVD (
'N',
'N', numRows, numCols,
1951 A_copy.values(), A_copy.stride(), &sigmas[0],
1952 (Scalar*) NULL, 1, (Scalar*) NULL, 1,
1953 &work[0], lwork, &rwork[0], &info);
1955 "LAPACK _GESVD failed with INFO = " << info <<
".");
1965 std::pair<magnitude_type, magnitude_type>
1970 const int numRows =
A.numRows();
1971 const int numCols =
A.numCols();
1973 Array<magnitude_type> sigmas (std::min (numRows, numCols));
1975 return std::make_pair (sigmas[0], sigmas[std::min(numRows, numCols) - 1]);
1994 const std::pair<magnitude_type, magnitude_type> sigmaMaxMin =
2000 "The test matrix is rank deficient; LAPACK's _GESVD "
2001 "routine reports that its smallest singular value is "
2005 const magnitude_type A_cond = sigmaMaxMin.first / sigmaMaxMin.second;
2032 return 2 * A_cond / cosTheta + tanTheta * A_cond * A_cond;
2047 return r.normFrobenius ();
2058 const int numRows = x_exact.
numRows();
2059 const int numCols = x_exact.
numCols();
2061 mat_type x_diff (numRows, numCols);
2062 for (
int j = 0; j < numCols; ++j) {
2063 for (
int i = 0; i < numRows; ++i) {
2064 x_diff(i,j) = x_exact(i,j) - x_approx(i,j);
2132 const int numRows = curCol + 2;
2133 const int LDR = R.
stride();
2134 const bool extraDebug =
false;
2137 cerr <<
"updateColumnGivens: curCol = " << curCol << endl;
2154 printMatrix<Scalar> (cerr,
"R_col before ", R_col);
2160 for (
int j = 0; j < curCol; ++j) {
2163 Scalar theCosine = theCosines[j];
2164 Scalar theSine = theSines[j];
2167 cerr <<
" j = " << j <<
": (cos,sin) = "
2168 << theCosines[j] <<
"," << theSines[j] << endl;
2170 blas.
ROT (1, &R_col(j,0), LDR, &R_col(j+1,0), LDR,
2171 &theCosine, &theSine);
2173 if (extraDebug && curCol > 0) {
2174 printMatrix<Scalar> (cerr,
"R_col after applying previous "
2175 "Givens rotations", R_col);
2180 Scalar theCosine, theSine, result;
2182 theCosine, theSine, result);
2183 theCosines[curCol] = theCosine;
2184 theSines[curCol] = theSine;
2187 cerr <<
" New cos,sin = " << theCosine <<
"," << theSine << endl;
2193 R_col(curCol, 0) = result;
2197 printMatrix<Scalar> (cerr,
"R_col after applying current "
2198 "Givens rotation", R_col);
2206 const int LDZ = z.
stride();
2208 &z(curCol,0), LDZ, &z(curCol+1,0), LDZ,
2209 &theCosine, &theSine);
2215 printMatrix<Scalar> (cerr,
"z_after", z);
2264 const int numRows = curCol + 2;
2265 const int numCols = curCol + 1;
2266 const int LDR = R.
stride();
2283 lapack.
GELS (
'N', numRows, numCols, y_view.numCols(),
2284 NULL, LDR, NULL, y_view.stride(),
2285 &lworkScalar, -1, &info);
2287 "LAPACK _GELS workspace query failed with INFO = "
2288 << info <<
", for a " << numRows <<
" x " << numCols
2289 <<
" matrix with " << y_view.numCols()
2290 <<
" right hand side"
2291 << ((y_view.numCols() != 1) ?
"s" :
"") <<
".");
2294 "LAPACK _GELS workspace query returned an LWORK with "
2295 "negative real part: LWORK = " << lworkScalar
2296 <<
". That should never happen. Please report this "
2297 "to the Belos developers.");
2300 "LAPACK _GELS workspace query returned an LWORK with "
2301 "nonzero imaginary part: LWORK = " << lworkScalar
2302 <<
". That should never happen. Please report this "
2303 "to the Belos developers.");
2309 const int lwork = std::max (1, static_cast<int> (
STS::real (lworkScalar)));
2317 lapack.
GELS (
'N', numRows, numCols, y_view.numCols(),
2319 y_view.values(), y_view.stride(),
2320 (lwork > 0 ? &work[0] : (Scalar*) NULL),
2324 "Solving projected least-squares problem with LAPACK "
2325 "_GELS failed with INFO = " << info <<
", for a "
2326 << numRows <<
" x " << numCols <<
" matrix with "
2327 << y_view.numCols() <<
" right hand side"
2328 << (y_view.numCols() != 1 ?
"s" :
"") <<
".");
2356 "updateColumnGivens: startCol = " << startCol
2357 <<
" > endCol = " << endCol <<
".");
2360 for (
int curCol = startCol; curCol <= endCol; ++curCol) {
2388 const int numRows = endCol + 2;
2389 const int numColsToUpdate = endCol - startCol + 1;
2390 const int LDR = R.
stride();
2405 for (
int j = 0; j < startCol; ++j) {
2406 blas.
ROT (numColsToUpdate,
2407 &R(j, startCol), LDR, &R(j+1, startCol), LDR,
2408 &theCosines[j], &theSines[j]);
2412 for (
int curCol = startCol; curCol < endCol; ++curCol) {
2415 for (
int j = startCol; j < curCol; ++j) {
2416 blas.
ROT (1, &R(j, curCol), LDR, &R(j+1, curCol), LDR,
2417 &theCosines[j], &theSines[j]);
2421 Scalar theCosine, theSine, result;
2423 theCosine, theSine, result);
2424 theCosines[curCol] = theCosine;
2425 theSines[curCol] = theSine;
2430 R(curCol+1, curCol) = result;
2438 const int LDZ = z.
stride();
2440 &z(curCol,0), LDZ, &z(curCol+1,0), LDZ,
2441 &theCosine, &theSine);
2453 #endif // __Belos_ProjectedLeastSquaresSolver_hpp