45 #ifndef __BelosTsqrOrthoManagerImpl_hpp
46 #define __BelosTsqrOrthoManagerImpl_hpp
56 #ifdef BELOS_TEUCHOS_TIME_MONITOR
58 #endif // BELOS_TEUCHOS_TIME_MONITOR
130 template<
class Scalar>
132 public std::binary_function<Teuchos::ArrayView<typename Teuchos::ScalarTraits<Scalar>::magnitudeType>,
133 Teuchos::ArrayView<typename Teuchos::ScalarTraits<Scalar>::magnitudeType>,
157 template<
class Scalar>
191 template<
class Scalar,
class MV>
250 const std::string& label);
294 #ifdef BELOS_TEUCHOS_TIME_MONITOR
295 clearTimer (label,
"All orthogonalization");
296 clearTimer (label,
"Projection");
297 clearTimer (label,
"Normalization");
299 timerOrtho_ = makeTimer (label,
"All orthogonalization");
300 timerProject_ = makeTimer (label,
"Projection");
301 timerNormalize_ = makeTimer (label,
"Normalization");
302 #endif // BELOS_TEUCHOS_TIME_MONITOR
341 norm (
const MV& X, std::vector<magnitude_type>& normVec)
const;
459 for (
int k = 0; k < ncols; ++k) {
472 mat_type X1_T_X2 (ncols_X1, ncols_X2);
544 #ifdef BELOS_TEUCHOS_TIME_MONITOR
553 #endif // BELOS_TEUCHOS_TIME_MONITOR
558 #ifdef BELOS_TEUCHOS_TIME_MONITOR
568 makeTimer (
const std::string& prefix,
569 const std::string& timerName)
571 const std::string timerLabel =
572 prefix.empty() ? timerName : (prefix +
": " + timerName);
582 clearTimer (
const std::string& prefix,
583 const std::string& timerName)
585 const std::string timerLabel =
586 prefix.empty() ? timerName : (prefix +
": " + timerName);
589 #endif // BELOS_TEUCHOS_TIME_MONITOR
594 const std::vector<magnitude_type>& normsAfterSecondPass,
595 const std::vector<int>& faultIndices);
627 const bool attemptToRecycle =
true)
const;
640 const bool outOfPlace,
738 template<
class Scalar,
class MV>
744 using Teuchos::parameterList;
746 using Teuchos::sublist;
749 RCP<const ParameterList> defaultParams = getValidParameters ();
751 RCP<ParameterList> tsqrParams;
753 RCP<ParameterList> theParams;
755 theParams = parameterList (*defaultParams);
763 randomizeNullSpace_ =
764 theParams->
get<
bool> (
"randomizeNullSpace",
765 defaultParams->get<
bool> (
"randomizeNullSpace"));
766 reorthogonalizeBlocks_ =
767 theParams->get<
bool> (
"reorthogonalizeBlocks",
768 defaultParams->get<
bool> (
"reorthogonalizeBlocks"));
769 throwOnReorthogFault_ =
770 theParams->get<
bool> (
"throwOnReorthogFault",
771 defaultParams->get<
bool> (
"throwOnReorthogFault"));
772 blockReorthogThreshold_ =
773 theParams->get<M> (
"blockReorthogThreshold",
774 defaultParams->get<M> (
"blockReorthogThreshold"));
775 relativeRankTolerance_ =
776 theParams->get<M> (
"relativeRankTolerance",
777 defaultParams->get<M> (
"relativeRankTolerance"));
778 forceNonnegativeDiagonal_ =
779 theParams->get<
bool> (
"forceNonnegativeDiagonal",
780 defaultParams->get<
bool> (
"forceNonnegativeDiagonal"));
784 if (! theParams->isSublist (
"TSQR implementation")) {
785 theParams->set (
"TSQR implementation",
786 defaultParams->sublist (
"TSQR implementation"));
788 tsqrParams = sublist (theParams,
"TSQR implementation",
true);
792 tsqrAdaptor_.setParameterList (tsqrParams);
795 setMyParamList (theParams);
798 template<
class Scalar,
class MV>
801 const std::string& label) :
805 randomizeNullSpace_ (true),
806 reorthogonalizeBlocks_ (true),
807 throwOnReorthogFault_ (false),
808 blockReorthogThreshold_ (0),
809 relativeRankTolerance_ (0),
810 forceNonnegativeDiagonal_ (false)
814 #ifdef BELOS_TEUCHOS_TIME_MONITOR
815 timerOrtho_ = makeTimer (label,
"All orthogonalization");
816 timerProject_ = makeTimer (label,
"Projection");
817 timerNormalize_ = makeTimer (label,
"Normalization");
818 #endif // BELOS_TEUCHOS_TIME_MONITOR
821 template<
class Scalar,
class MV>
827 randomizeNullSpace_ (true),
828 reorthogonalizeBlocks_ (true),
829 throwOnReorthogFault_ (false),
830 blockReorthogThreshold_ (0),
831 relativeRankTolerance_ (0),
832 forceNonnegativeDiagonal_ (false)
836 #ifdef BELOS_TEUCHOS_TIME_MONITOR
837 timerOrtho_ = makeTimer (label,
"All orthogonalization");
838 timerProject_ = makeTimer (label,
"Projection");
839 timerNormalize_ = makeTimer (label,
"Normalization");
840 #endif // BELOS_TEUCHOS_TIME_MONITOR
843 template<
class Scalar,
class MV>
846 norm (
const MV& X, std::vector<magnitude_type>& normVec)
const
848 const int numCols = MVT::GetNumberVecs (X);
851 if (normVec.size() <
static_cast<size_t>(numCols))
852 normVec.resize (numCols);
853 MVT::MvNorm (X, normVec);
856 template<
class Scalar,
class MV>
862 #ifdef BELOS_TEUCHOS_TIME_MONITOR
871 #endif // BELOS_TEUCHOS_TIME_MONITOR
873 int ncols_X, num_Q_blocks, ncols_Q_total;
874 checkProjectionDims (ncols_X, num_Q_blocks, ncols_Q_total, X, Q);
878 if (ncols_X == 0 || num_Q_blocks == 0 || ncols_Q_total == 0)
882 allocateProjectionCoefficients (C, Q, X,
true);
886 std::vector<magnitude_type> columnNormsBefore (ncols_X,
magnitude_type(0));
887 if (reorthogonalizeBlocks_)
888 MVT::MvNorm (X, columnNormsBefore);
892 rawProject (X, Q, C);
896 if (reorthogonalizeBlocks_) {
897 std::vector<magnitude_type> columnNormsAfter (ncols_X,
magnitude_type(0));
898 MVT::MvNorm (X, columnNormsAfter);
906 bool reorthogonalize =
false;
907 for (
int j = 0; j < ncols_X; ++j) {
908 if (columnNormsAfter[j] < relThres * columnNormsBefore[j]) {
909 reorthogonalize =
true;
913 if (reorthogonalize) {
916 if (! reorthogCallback_.is_null()) {
917 reorthogCallback_->operator() (Teuchos::arrayViewFromVector (columnNormsBefore),
918 Teuchos::arrayViewFromVector (columnNormsAfter));
922 allocateProjectionCoefficients (C2, Q, X,
false);
926 rawProject (X, Q, C2);
928 for (
int k = 0; k < num_Q_blocks; ++k)
935 template<
class Scalar,
class MV>
942 #ifdef BELOS_TEUCHOS_TIME_MONITOR
949 #endif // BELOS_TEUCHOS_TIME_MONITOR
954 const int numCols = MVT::GetNumberVecs (X);
963 return normalizeOne (X, B);
993 MVT::GetGlobalLength(*Q_) != MVT::GetGlobalLength(X) ||
994 numCols > MVT::GetNumberVecs (*Q_)) {
995 Q_ = MVT::Clone (X, numCols);
1004 if (MVT::GetNumberVecs(*Q_) == numCols) {
1005 return normalizeImpl (X, *Q_, B,
false);
1007 RCP<MV> Q_view = MVT::CloneViewNonConst (*Q_, Range1D(0, numCols-1));
1008 return normalizeImpl (X, *Q_view, B,
false);
1012 template<
class Scalar,
class MV>
1018 const bool attemptToRecycle)
const
1020 const int num_Q_blocks = Q.size();
1021 const int ncols_X = MVT::GetNumberVecs (X);
1023 if (attemptToRecycle)
1025 for (
int i = 0; i < num_Q_blocks; ++i)
1027 const int ncols_Qi = MVT::GetNumberVecs (*Q[i]);
1036 Ci.
shape (ncols_Qi, ncols_X);
1044 for (
int i = 0; i < num_Q_blocks; ++i)
1046 const int ncols_Qi = MVT::GetNumberVecs (*Q[i]);
1052 template<
class Scalar,
class MV>
1057 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1060 #endif // BELOS_TEUCHOS_TIME_MONITOR
1062 const int numVecs = MVT::GetNumberVecs(X);
1065 }
else if (numVecs == 1) {
1072 const int rank = normalizeOne (X, B);
1074 RCP<MV> Q_0 = MVT::CloneViewNonConst (Q, Range1D(0,0));
1075 MVT::Assign (X, *Q_0);
1081 return normalizeImpl (X, Q, B,
true);
1085 template<
class Scalar,
class MV>
1090 const bool outOfPlace,
1099 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1103 #endif // BELOS_TEUCHOS_TIME_MONITOR
1108 std::invalid_argument,
1109 "Belos::TsqrOrthoManagerImpl::"
1110 "projectAndNormalizeImpl(..., outOfPlace=true, ...):"
1111 "X_out has " << MVT::GetNumberVecs(X_out)
1112 <<
" columns, but X_in has "
1113 << MVT::GetNumberVecs(X_in) <<
" columns.");
1117 int ncols_X, num_Q_blocks, ncols_Q_total;
1118 checkProjectionDims (ncols_X, num_Q_blocks, ncols_Q_total, X_in, Q);
1125 if (num_Q_blocks == 0 || ncols_Q_total == 0) {
1127 return normalizeOutOfPlace (X_in, X_out, B);
1129 return normalize (X_in, B);
1136 allocateProjectionCoefficients (C, Q, X_in,
true);
1141 std::vector<magnitude_type> normsBeforeFirstPass (ncols_X, SCTM::zero());
1142 if (reorthogonalizeBlocks_) {
1143 MVT::MvNorm (X_in, normsBeforeFirstPass);
1147 rawProject (X_in, Q, C);
1164 std::invalid_argument,
1165 "normalizeOne: Input matrix B must be at "
1166 "least " << ncols_X <<
" x " << ncols_X
1167 <<
", but is instead " << B->numRows()
1168 <<
" x " << B->numCols() <<
".");
1178 const int firstPassRank = outOfPlace ?
1179 normalizeOutOfPlace (X_in, X_out, B_out) :
1180 normalize (X_in, B_out);
1189 int rank = firstPassRank;
1205 if (firstPassRank < ncols_X && randomizeNullSpace_) {
1206 const int numNullSpaceCols = ncols_X - firstPassRank;
1207 const Range1D nullSpaceIndices (firstPassRank, ncols_X - 1);
1211 for (
int k = 0; k < num_Q_blocks; ++k) {
1212 const int numColsQk = MVT::GetNumberVecs(*Q[k]);
1213 C_null[k] =
rcp (
new mat_type (numColsQk, numNullSpaceCols));
1216 RCP<mat_type> B_null (
new mat_type (numNullSpaceCols, numNullSpaceCols));
1218 int randomVectorsRank;
1222 RCP<MV> X_out_null = MVT::CloneViewNonConst (X_out, nullSpaceIndices);
1227 RCP<MV> X_in_null = MVT::CloneViewNonConst (X_in, nullSpaceIndices);
1228 MVT::Assign (*X_out_null, *X_in_null);
1231 rawProject (*X_in_null, Q, C_null);
1232 randomVectorsRank = normalizeOutOfPlace (*X_in_null, *X_out_null, B_null);
1236 RCP<MV> X_null = MVT::CloneViewNonConst (X_in, nullSpaceIndices);
1239 rawProject (*X_null, Q, C_null);
1240 randomVectorsRank = normalize (*X_null, B_null);
1248 "Belos::TsqrOrthoManagerImpl::projectAndNormalize"
1249 "OutOfPlace(): After projecting and normalizing the "
1250 "random vectors (used to replace the null space "
1251 "basis vectors from normalizing X), they have rank "
1252 << randomVectorsRank <<
", but should have full "
1253 "rank " << numNullSpaceCols <<
".");
1258 if (reorthogonalizeBlocks_) {
1261 std::vector<magnitude_type>
1262 normsAfterFirstPass (firstPassRank, SCTM::zero());
1263 std::vector<magnitude_type>
1264 normsAfterSecondPass (firstPassRank, SCTM::zero());
1279 for (
int j = 0; j < firstPassRank; ++j) {
1280 const Scalar*
const B_j = &(*B_out)(0,j);
1283 normsAfterFirstPass[j] = blas.
NRM2 (firstPassRank, B_j, 1);
1287 bool reorthogonalize =
false;
1288 for (
int j = 0; j < firstPassRank; ++j) {
1296 blockReorthogThreshold() * normsBeforeFirstPass[j];
1297 if (normsAfterFirstPass[j] < curThreshold) {
1298 reorthogonalize =
true;
1305 if (! reorthogCallback_.is_null()) {
1306 using Teuchos::arrayViewFromVector;
1307 (*reorthogCallback_) (arrayViewFromVector (normsBeforeFirstPass),
1308 arrayViewFromVector (normsAfterFirstPass));
1322 bool reorthogFault =
false;
1324 std::vector<int> faultIndices;
1325 if (reorthogonalize) {
1333 MVT::Assign (X_out, X_in);
1340 allocateProjectionCoefficients (C2, Q, X_in,
false);
1345 rawProject (X_in, Q, C2);
1348 RCP<mat_type>
B2 (
new mat_type (ncols_X, ncols_X));
1351 const int secondPassRank = outOfPlace ?
1352 normalizeOutOfPlace (X_in, X_out, B2) :
1353 normalize (X_in, B2);
1354 rank = secondPassRank;
1359 mat_type B_copy (
Copy, *B_out, B_out->numRows(), B_out->numCols());
1362 *B2, B_copy, SCT::zero());
1364 "Teuchos::SerialDenseMatrix::multiply "
1365 "returned err = " << err <<
" != 0");
1369 for (
int k = 0; k < num_Q_blocks; ++k) {
1370 mat_type C_k_copy (
Copy, *C[k], C[k]->numRows(), C[k]->numCols());
1374 *C2[k], B_copy, SCT::one());
1376 "Teuchos::SerialDenseMatrix::multiply "
1377 "returned err = " << err1 <<
" != 0");
1382 for (
int j = 0; j < rank; ++j) {
1383 const Scalar*
const B2_j = &(*B2)(0,j);
1384 normsAfterSecondPass[j] = blas.
NRM2 (rank, B2_j, 1);
1389 reorthogFault =
false;
1390 for (
int j = 0; j < rank; ++j) {
1392 blockReorthogThreshold() * normsAfterFirstPass[j];
1393 if (normsAfterSecondPass[j] < relativeLowerBound) {
1394 reorthogFault =
true;
1395 faultIndices.push_back (j);
1400 if (reorthogFault) {
1401 if (throwOnReorthogFault_) {
1402 raiseReorthogFault (normsAfterFirstPass,
1403 normsAfterSecondPass,
1412 "TsqrOrthoManagerImpl has not yet implemented"
1413 " recovery from an orthogonalization fault.");
1421 template<
class Scalar,
class MV>
1425 const std::vector<magnitude_type>& normsAfterSecondPass,
1426 const std::vector<int>& faultIndices)
1429 typedef std::vector<int>::size_type size_type;
1430 std::ostringstream os;
1432 os <<
"Orthogonalization fault at the following column(s) of X:" << endl;
1433 os <<
"Column\tNorm decrease factor" << endl;
1434 for (size_type k = 0; k < faultIndices.size(); ++k) {
1435 const int index = faultIndices[k];
1437 normsAfterSecondPass[index] / normsAfterFirstPass[index];
1438 os << index <<
"\t" << decreaseFactor << endl;
1443 template<
class Scalar,
class MV>
1448 using Teuchos::parameterList;
1451 if (defaultParams_.is_null()) {
1452 RCP<ParameterList> params = parameterList (
"TsqrOrthoManagerImpl");
1456 params->set (
"TSQR implementation", *(tsqrAdaptor_.getValidParameters()),
1457 "TSQR implementation parameters.");
1461 const bool defaultRandomizeNullSpace =
true;
1462 params->set (
"randomizeNullSpace", defaultRandomizeNullSpace,
1463 "Whether to fill in null space vectors with random data.");
1465 const bool defaultReorthogonalizeBlocks =
true;
1466 params->set (
"reorthogonalizeBlocks", defaultReorthogonalizeBlocks,
1467 "Whether to do block reorthogonalization as necessary.");
1473 params->set (
"blockReorthogThreshold", defaultBlockReorthogThreshold,
1474 "If reorthogonalizeBlocks==true, and if the norm of "
1475 "any column within a block decreases by this much or "
1476 "more after orthogonalization, we reorthogonalize.");
1481 Teuchos::as<magnitude_type>(10) * SCTM::eps();
1486 params->set (
"relativeRankTolerance", defaultRelativeRankTolerance,
1487 "Relative tolerance to determine the numerical rank of a "
1488 "block when normalizing.");
1492 const bool defaultThrowOnReorthogFault =
true;
1493 params->set (
"throwOnReorthogFault", defaultThrowOnReorthogFault,
1494 "Whether to throw an exception if an orthogonalization "
1495 "fault occurs. This only matters if reorthogonalization "
1496 "is enabled (reorthogonalizeBlocks==true).");
1498 const bool defaultForceNonnegativeDiagonal =
false;
1499 params->set (
"forceNonnegativeDiagonal", defaultForceNonnegativeDiagonal,
1500 "Whether to force the R factor produced by the normalization "
1501 "step to have a nonnegative diagonal.");
1503 defaultParams_ = params;
1505 return defaultParams_;
1508 template<
class Scalar,
class MV>
1516 RCP<const ParameterList> defaultParams = getValidParameters();
1518 RCP<ParameterList> params =
rcp (
new ParameterList (*defaultParams));
1527 const bool randomizeNullSpace =
false;
1528 params->set (
"randomizeNullSpace", randomizeNullSpace);
1529 const bool reorthogonalizeBlocks =
false;
1530 params->set (
"reorthogonalizeBlocks", reorthogonalizeBlocks);
1535 template<
class Scalar,
class MV>
1547 tsqrAdaptor_.factorExplicit (X, Q, B, forceNonnegativeDiagonal_);
1550 rank = tsqrAdaptor_.revealRank (Q, B, relativeRankTolerance_);
1551 }
catch (std::exception& e) {
1557 template<
class Scalar,
class MV>
1570 const int theNumRows =
B->numRows ();
1571 const int theNumCols =
B->numCols ();
1573 theNumRows < 1 || theNumCols < 1, std::invalid_argument,
1574 "normalizeOne: Input matrix B must be at least 1 x 1, but "
1575 "is instead " << theNumRows <<
" x " << theNumCols <<
".");
1581 std::vector<magnitude_type> theNorm (1, SCTM::zero());
1582 MVT::MvNorm (X, theNorm);
1583 (*B_out)(0,0) = theNorm[0];
1597 if (theNorm[0] == SCTM::zero()) {
1600 if (randomizeNullSpace_) {
1602 MVT::MvNorm (X, theNorm);
1603 if (theNorm[0] == SCTM::zero()) {
1609 "vector has norm zero!");
1614 const Scalar alpha = SCT::one() / theNorm[0];
1615 MVT::MvScale (X, alpha);
1622 const Scalar alpha = SCT::one() / theNorm[0];
1623 MVT::MvScale (X, alpha);
1629 template<
class Scalar,
class MV>
1636 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1638 #endif // BELOS_TEUCHOS_TIME_MONITOR
1641 const int num_Q_blocks = Q.size();
1642 for (
int i = 0; i < num_Q_blocks; ++i)
1650 mat_type& Ci = *
C[i];
1651 const MV& Qi = *Q[i];
1652 innerProd (Qi, X, Ci);
1653 MVT::MvTimesMatAddMv (-SCT::one(), Qi, Ci, SCT::one(), X);
1658 template<
class Scalar,
class MV>
1665 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1667 #endif // BELOS_TEUCHOS_TIME_MONITOR
1670 innerProd (*Q, X, *
C);
1671 MVT::MvTimesMatAddMv (-SCT::one(), *Q, *
C, SCT::one(), X);
1674 template<
class Scalar,
class MV>
1680 const bool outOfPlace)
1686 using Teuchos::tuple;
1688 const int numCols = MVT::GetNumberVecs (X);
1696 MVT::GetNumberVecs (Q) < numCols, std::invalid_argument,
1697 "TsqrOrthoManagerImpl::normalizeImpl: Q has "
1698 << MVT::GetNumberVecs(Q) <<
" columns. This is too "
1699 "few, since X has " << numCols <<
" columns.");
1703 RCP<MV> Q_view = MVT::CloneViewNonConst (Q, Range1D (0, numCols-1));
1715 B->numRows() < numCols ||
B->numCols() < numCols, std::invalid_argument,
1716 "TsqrOrthoManagerImpl::normalizeImpl: Input matrix B must be at least "
1717 << numCols <<
" x " << numCols <<
", but is instead " <<
B->numRows ()
1718 <<
" x " <<
B->numCols() <<
".");
1733 const int rank = rawNormalize (X, *Q_view, *B_out);
1744 rank < 0 || rank > numCols, std::logic_error,
1745 "Belos::TsqrOrthoManagerImpl::normalizeImpl: rawNormalize returned rank "
1746 " = " << rank <<
" for a matrix X with " << numCols <<
" columns. "
1747 "Please report this bug to the Belos developers.");
1751 if (rank == numCols || ! randomizeNullSpace_) {
1755 MVT::Assign (*Q_view, X);
1760 if (randomizeNullSpace_ && rank < numCols) {
1767 const int nullSpaceNumCols = numCols - rank;
1770 Range1D nullSpaceIndices (rank, numCols-1);
1777 RCP<MV> Q_null = MVT::CloneViewNonConst (*Q_view, nullSpaceIndices);
1779 MVT::MvRandom (*Q_null);
1785 std::vector<magnitude_type> norms (MVT::GetNumberVecs (*Q_null));
1786 MVT::MvNorm (*Q_null, norms);
1788 bool anyZero =
false;
1789 typedef typename std::vector<magnitude_type>::const_iterator iter_type;
1790 for (iter_type it = norms.begin(); it != norms.end(); ++it) {
1791 if (*it == SCTM::zero()) {
1796 std::ostringstream os;
1797 os <<
"TsqrOrthoManagerImpl::normalizeImpl: "
1798 "We are being asked to randomize the null space, for a matrix "
1799 "with " << numCols <<
" columns and reported column rank "
1800 << rank <<
". The inclusive range of columns to fill with "
1801 "random data is [" << nullSpaceIndices.lbound() <<
","
1802 << nullSpaceIndices.ubound() <<
"]. After filling the null "
1803 "space vectors with random numbers, at least one of the vectors"
1804 " has norm zero. Here are the norms of all the null space "
1806 for (iter_type it = norms.begin(); it != norms.end(); ++it) {
1808 if (it+1 != norms.end())
1811 os <<
"].) There is a tiny probability that this could happen "
1812 "randomly, but it is likely a bug. Please report it to the "
1813 "Belos developers, especially if you are able to reproduce the "
1826 RCP<const MV> Q_col = MVT::CloneView (*Q_view, Range1D (0, rank-1));
1832 rawProject (*Q_null, Q_col, C_null);
1841 RCP<MV> X_null = MVT::CloneViewNonConst (X, nullSpaceIndices);
1844 mat_type B_null (nullSpaceNumCols, nullSpaceNumCols);
1846 const int nullSpaceBasisRank = rawNormalize (*Q_null, *X_null, B_null);
1860 if (nullSpaceBasisRank < nullSpaceNumCols) {
1861 std::vector<magnitude_type> norms (MVT::GetNumberVecs(*X_null));
1862 MVT::MvNorm (*X_null, norms);
1863 std::ostringstream os;
1864 os <<
"TsqrOrthoManagerImpl::normalizeImpl: "
1865 <<
"We are being asked to randomize the null space, "
1866 <<
"for a matrix with " << numCols <<
" columns and "
1867 <<
"column rank " << rank <<
". After projecting and "
1868 <<
"normalizing the generated random vectors, they "
1869 <<
"only have rank " << nullSpaceBasisRank <<
". They"
1870 <<
" should have full rank " << nullSpaceNumCols
1871 <<
". (The inclusive range of columns to fill with "
1872 <<
"random data is [" << nullSpaceIndices.lbound()
1873 <<
"," << nullSpaceIndices.ubound() <<
"]. The "
1874 <<
"column norms of the resulting Q factor are: [";
1875 for (
typename std::vector<magnitude_type>::size_type k = 0;
1876 k < norms.size(); ++k) {
1878 if (k != norms.size()-1) {
1882 os <<
"].) There is a tiny probability that this could "
1883 <<
"happen randomly, but it is likely a bug. Please "
1884 <<
"report it to the Belos developers, especially if "
1885 <<
"you are able to reproduce the behavior.";
1898 MVT::Assign (*X_null, *Q_null);
1899 }
else if (rank > 0) {
1901 RCP<const MV> Q_col = MVT::CloneView (*Q_view, Range1D (0, rank-1));
1902 RCP<MV> X_col = MVT::CloneViewNonConst (X, Range1D (0, rank-1));
1903 MVT::Assign (*Q_col, *X_col);
1910 template<
class Scalar,
class MV>
1925 int the_ncols_X, the_num_Q_blocks, the_ncols_Q_total;
1926 the_num_Q_blocks = Q.size();
1927 the_ncols_X = MVT::GetNumberVecs (X);
1930 the_ncols_Q_total = 0;
1935 typedef typename ArrayView<RCP<const MV> >::const_iterator iter_type;
1936 for (iter_type it = Q.begin(); it != Q.end(); ++it) {
1937 const MV& Qi = **it;
1938 the_ncols_Q_total += MVT::GetNumberVecs (Qi);
1942 ncols_X = the_ncols_X;
1943 num_Q_blocks = the_num_Q_blocks;
1944 ncols_Q_total = the_ncols_Q_total;
1949 #endif // __BelosTsqrOrthoManagerImpl_hpp
MVT::tsqr_adaptor_type tsqr_adaptor_type
Teuchos::ScalarTraits< magnitude_type > SCTM
magnitude_type relativeRankTolerance_
Relative tolerance for measuring the numerical rank of a matrix.
int projectAndNormalize(MV &X, Teuchos::Array< mat_ptr > C, mat_ptr B, Teuchos::ArrayView< Teuchos::RCP< const MV > > Q)
Project X against Q and normalize X.
Interface of callback invoked by TsqrOrthoManager on reorthogonalization.
bool is_null(const boost::shared_ptr< T > &p)
magnitude_type blockReorthogThreshold_
Relative reorthogonalization threshold in Block Gram-Schmidt.
void project(MV &X, Teuchos::Array< mat_ptr > C, Teuchos::ArrayView< Teuchos::RCP< const MV > > Q)
Compute and .
void norm(const MV &X, std::vector< magnitude_type > &normVec) const
Compute the 2-norm of each column j of X.
TsqrOrthoManagerImpl(const Teuchos::RCP< Teuchos::ParameterList > ¶ms, const std::string &label)
Constructor (that sets user-specified parameters).
void checkProjectionDims(int &ncols_X, int &num_Q_blocks, int &ncols_Q_total, const MV &X, Teuchos::ArrayView< Teuchos::RCP< const MV > > Q) const
Return through output arguments some relevant dimension information about X and Q.
magnitude_type orthogError(const MV &X1, const MV &X2) const
Return the Frobenius norm of the inner product of X1 with itself.
Teuchos::SerialDenseMatrix< int, Scalar > mat_type
Type of the projection and normalization coefficients.
Teuchos::RCP< MV > Q_
Scratch space for TSQR.
T & get(ParameterList &l, const std::string &name)
void setParameterList(const Teuchos::RCP< Teuchos::ParameterList > ¶ms)
Set parameters from the given parameter list.
int projectAndNormalizeImpl(MV &X_in, MV &X_out, const bool outOfPlace, Teuchos::Array< mat_ptr > C, mat_ptr B, Teuchos::ArrayView< Teuchos::RCP< const MV > > Q)
Implementation of projection and normalization.
magnitude_type eps_
Machine precision for Scalar.
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
Teuchos::RCP< const Teuchos::ParameterList > getFastParameters()
Get "fast" parameters for TsqrOrthoManagerImpl.
void rawProject(MV &X, Teuchos::ArrayView< Teuchos::RCP< const MV > > Q, Teuchos::ArrayView< mat_ptr > C) const
One projection pass of X against the Q[i] blocks.
TSQR-based OrthoManager subclass implementation.
Declaration of basic traits for the multivector type.
void setReorthogonalizationCallback(const Teuchos::RCP< ReorthogonalizationCallback< Scalar > > &callback)
Set callback to be invoked on reorthogonalization.
Teuchos::RCP< ReorthogonalizationCallback< Scalar > > reorthogCallback_
Callback invoked if reorthogonalization is necessary.
Error in TsqrOrthoManager or TsqrOrthoManagerImpl.
ScalarTraits< ScalarType >::magnitudeType NRM2(const OrdinalType &n, const ScalarType *x, const OrdinalType &incx) const
bool forceNonnegativeDiagonal_
Force R factor of normalization to have a nonnegative diagonal.
int normalize(MV &X, mat_ptr B)
Orthogonalize the columns of X in place.
bool throwOnReorthogFault_
Whether to throw an exception on a orthogonalization fault.
void setLabel(const std::string &label)
Set the label for timers.
Teuchos::RCP< mat_type > mat_ptr
static void MvTransMv(const ScalarType alpha, const MV &A, const MV &mv, Teuchos::SerialDenseMatrix< int, ScalarType > &B)
Compute a dense matrix B through the matrix-matrix multiply .
ScalarTraits< ScalarType >::magnitudeType normFrobenius() const
static int GetNumberVecs(const MV &mv)
Obtain the number of vectors in mv.
Teuchos::ScalarTraits< Scalar >::magnitudeType magnitude_type
The type of a norm result.
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
Default valid parameter list.
tsqr_adaptor_type tsqrAdaptor_
Interface to TSQR implementation.
Traits class which defines basic operations on multivectors.
Teuchos::RCP< const Teuchos::ParameterList > defaultParams_
Default configuration parameters.
Scalar scalar_type
The template parameter of this class; the type of an inner product result.
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Teuchos::ScalarTraits< Scalar > SCT
bool is_null(const RCP< T > &p)
virtual void operator()(Teuchos::ArrayView< magnitude_type > normsBeforeFirstPass, Teuchos::ArrayView< magnitude_type > normsAfterFirstPass)=0
Callback invoked by TsqrOrthoManager on reorthogonalization.
int putScalar(const ScalarType value=Teuchos::ScalarTraits< ScalarType >::zero())
const std::string & getLabel() const
Get the label for timers (if timers are enabled).
magnitude_type orthonormError(const MV &X) const
Return .
Exception thrown to signal error in an orthogonalization manager method.
Teuchos::ScalarTraits< Scalar >::magnitudeType magnitude_type
magnitude_type relativeRankTolerance() const
Relative tolerance for determining (via the SVD) whether a block is of full numerical rank...
void resize(size_type new_size, const value_type &x=value_type())
TsqrOrthoError(const std::string &what_arg)
int normalizeOne(MV &X, mat_ptr B) const
Normalize a "multivector" of only one column.
OrdinalType numCols() const
int normalizeImpl(MV &X, MV &Q, mat_ptr B, const bool outOfPlace)
Normalize X into Q*B, with out-of-place option.
Teuchos::RCP< Teuchos::ParameterList > params_
Configuration parameters.
virtual ~ReorthogonalizationCallback()
Destructor (virtual for memory safety of derived classes)
Templated virtual class for providing orthogonalization/orthonormalization methods.
int rawNormalize(MV &X, MV &Q, mat_type &B)
One out-of-place normalization pass.
bool reorthogonalizeBlocks_
Whether to reorthogonalize blocks at all.
TsqrOrthoFault(const std::string &what_arg)
void raiseReorthogFault(const std::vector< magnitude_type > &normsAfterFirstPass, const std::vector< magnitude_type > &normsAfterSecondPass, const std::vector< int > &faultIndices)
Throw an exception indicating a reorthgonalization fault.
std::string label_
Label for timers (if timers are used).
magnitude_type blockReorthogThreshold() const
Relative tolerance for triggering a block reorthogonalization.
void allocateProjectionCoefficients(Teuchos::Array< mat_ptr > &C, Teuchos::ArrayView< Teuchos::RCP< const MV > > Q, const MV &X, const bool attemptToRecycle=true) const
Allocate projection coefficients.
int shape(OrdinalType numRows, OrdinalType numCols)
int projectAndNormalizeOutOfPlace(MV &X_in, MV &X_out, Teuchos::Array< mat_ptr > C, mat_ptr B, Teuchos::ArrayView< Teuchos::RCP< const MV > > Q)
Project and normalize X_in into X_out; overwrite X_in.
int normalizeOutOfPlace(MV &X, MV &Q, mat_ptr B)
Normalize X into Q*B, overwriting X.
void innerProd(const MV &X, const MV &Y, mat_type &Z) const
Euclidean inner product.
Belos header file which uses auto-configuration information to include necessary C++ headers...
bool randomizeNullSpace_
Whether to fill null space vectors with random data.
OrdinalType numRows() const
MultiVecTraits< Scalar, MV > MVT