13 #ifndef __BelosTsqrOrthoManagerImpl_hpp
14 #define __BelosTsqrOrthoManagerImpl_hpp
23 #ifdef BELOS_TEUCHOS_TIME_MONITOR
25 #endif // BELOS_TEUCHOS_TIME_MONITOR
97 template<
class Scalar>
121 template<
class Scalar>
155 template<
class Scalar,
class MV>
214 const std::string& label);
258 #ifdef BELOS_TEUCHOS_TIME_MONITOR
259 clearTimer (label,
"All orthogonalization");
260 clearTimer (label,
"Projection");
261 clearTimer (label,
"Normalization");
263 timerOrtho_ = makeTimer (label,
"All orthogonalization");
264 timerProject_ = makeTimer (label,
"Projection");
265 timerNormalize_ = makeTimer (label,
"Normalization");
266 #endif // BELOS_TEUCHOS_TIME_MONITOR
305 norm (
const MV& X, std::vector<magnitude_type>& normVec)
const;
423 for (
int k = 0; k < ncols; ++k) {
436 mat_type X1_T_X2 (ncols_X1, ncols_X2);
508 #ifdef BELOS_TEUCHOS_TIME_MONITOR
517 #endif // BELOS_TEUCHOS_TIME_MONITOR
522 #ifdef BELOS_TEUCHOS_TIME_MONITOR
532 makeTimer (
const std::string& prefix,
533 const std::string& timerName)
535 const std::string timerLabel =
536 prefix.empty() ? timerName : (prefix +
": " + timerName);
546 clearTimer (
const std::string& prefix,
547 const std::string& timerName)
549 const std::string timerLabel =
550 prefix.empty() ? timerName : (prefix +
": " + timerName);
553 #endif // BELOS_TEUCHOS_TIME_MONITOR
558 const std::vector<magnitude_type>& normsAfterSecondPass,
559 const std::vector<int>& faultIndices);
591 const bool attemptToRecycle =
true)
const;
604 const bool outOfPlace,
702 template<
class Scalar,
class MV>
708 using Teuchos::parameterList;
710 using Teuchos::sublist;
713 RCP<const ParameterList> defaultParams = getValidParameters ();
715 RCP<ParameterList> tsqrParams;
717 RCP<ParameterList> theParams;
719 theParams = parameterList (*defaultParams);
727 randomizeNullSpace_ =
728 theParams->
get<
bool> (
"randomizeNullSpace",
729 defaultParams->get<
bool> (
"randomizeNullSpace"));
730 reorthogonalizeBlocks_ =
731 theParams->get<
bool> (
"reorthogonalizeBlocks",
732 defaultParams->get<
bool> (
"reorthogonalizeBlocks"));
733 throwOnReorthogFault_ =
734 theParams->get<
bool> (
"throwOnReorthogFault",
735 defaultParams->get<
bool> (
"throwOnReorthogFault"));
736 blockReorthogThreshold_ =
737 theParams->get<M> (
"blockReorthogThreshold",
738 defaultParams->get<M> (
"blockReorthogThreshold"));
739 relativeRankTolerance_ =
740 theParams->get<M> (
"relativeRankTolerance",
741 defaultParams->get<M> (
"relativeRankTolerance"));
742 forceNonnegativeDiagonal_ =
743 theParams->get<
bool> (
"forceNonnegativeDiagonal",
744 defaultParams->get<
bool> (
"forceNonnegativeDiagonal"));
748 if (! theParams->isSublist (
"TSQR implementation")) {
749 theParams->set (
"TSQR implementation",
750 defaultParams->sublist (
"TSQR implementation"));
752 tsqrParams = sublist (theParams,
"TSQR implementation",
true);
756 tsqrAdaptor_.setParameterList (tsqrParams);
759 setMyParamList (theParams);
762 template<
class Scalar,
class MV>
765 const std::string& label) :
769 randomizeNullSpace_ (true),
770 reorthogonalizeBlocks_ (true),
771 throwOnReorthogFault_ (false),
772 blockReorthogThreshold_ (0),
773 relativeRankTolerance_ (0),
774 forceNonnegativeDiagonal_ (false)
778 #ifdef BELOS_TEUCHOS_TIME_MONITOR
779 timerOrtho_ = makeTimer (label,
"All orthogonalization");
780 timerProject_ = makeTimer (label,
"Projection");
781 timerNormalize_ = makeTimer (label,
"Normalization");
782 #endif // BELOS_TEUCHOS_TIME_MONITOR
785 template<
class Scalar,
class MV>
791 randomizeNullSpace_ (true),
792 reorthogonalizeBlocks_ (true),
793 throwOnReorthogFault_ (false),
794 blockReorthogThreshold_ (0),
795 relativeRankTolerance_ (0),
796 forceNonnegativeDiagonal_ (false)
800 #ifdef BELOS_TEUCHOS_TIME_MONITOR
801 timerOrtho_ = makeTimer (label,
"All orthogonalization");
802 timerProject_ = makeTimer (label,
"Projection");
803 timerNormalize_ = makeTimer (label,
"Normalization");
804 #endif // BELOS_TEUCHOS_TIME_MONITOR
807 template<
class Scalar,
class MV>
810 norm (
const MV& X, std::vector<magnitude_type>& normVec)
const
812 const int numCols = MVT::GetNumberVecs (X);
815 if (normVec.size() <
static_cast<size_t>(numCols))
816 normVec.resize (numCols);
817 MVT::MvNorm (X, normVec);
820 template<
class Scalar,
class MV>
826 #ifdef BELOS_TEUCHOS_TIME_MONITOR
835 #endif // BELOS_TEUCHOS_TIME_MONITOR
837 int ncols_X, num_Q_blocks, ncols_Q_total;
838 checkProjectionDims (ncols_X, num_Q_blocks, ncols_Q_total, X, Q);
842 if (ncols_X == 0 || num_Q_blocks == 0 || ncols_Q_total == 0)
846 allocateProjectionCoefficients (C, Q, X,
true);
850 std::vector<magnitude_type> columnNormsBefore (ncols_X,
magnitude_type(0));
851 if (reorthogonalizeBlocks_)
852 MVT::MvNorm (X, columnNormsBefore);
856 rawProject (X, Q, C);
860 if (reorthogonalizeBlocks_) {
861 std::vector<magnitude_type> columnNormsAfter (ncols_X,
magnitude_type(0));
862 MVT::MvNorm (X, columnNormsAfter);
870 bool reorthogonalize =
false;
871 for (
int j = 0; j < ncols_X; ++j) {
872 if (columnNormsAfter[j] < relThres * columnNormsBefore[j]) {
873 reorthogonalize =
true;
877 if (reorthogonalize) {
880 if (! reorthogCallback_.is_null()) {
881 reorthogCallback_->operator() (Teuchos::arrayViewFromVector (columnNormsBefore),
882 Teuchos::arrayViewFromVector (columnNormsAfter));
886 allocateProjectionCoefficients (C2, Q, X,
false);
890 rawProject (X, Q, C2);
892 for (
int k = 0; k < num_Q_blocks; ++k)
899 template<
class Scalar,
class MV>
906 #ifdef BELOS_TEUCHOS_TIME_MONITOR
913 #endif // BELOS_TEUCHOS_TIME_MONITOR
918 const int numCols = MVT::GetNumberVecs (X);
927 return normalizeOne (X, B);
957 MVT::GetGlobalLength(*Q_) != MVT::GetGlobalLength(X) ||
958 numCols > MVT::GetNumberVecs (*Q_)) {
959 Q_ = MVT::Clone (X, numCols);
968 if (MVT::GetNumberVecs(*Q_) == numCols) {
969 return normalizeImpl (X, *Q_, B,
false);
971 RCP<MV> Q_view = MVT::CloneViewNonConst (*Q_, Range1D(0, numCols-1));
972 return normalizeImpl (X, *Q_view, B,
false);
976 template<
class Scalar,
class MV>
982 const bool attemptToRecycle)
const
984 const int num_Q_blocks = Q.size();
985 const int ncols_X = MVT::GetNumberVecs (X);
987 if (attemptToRecycle)
989 for (
int i = 0; i < num_Q_blocks; ++i)
991 const int ncols_Qi = MVT::GetNumberVecs (*Q[i]);
1000 Ci.
shape (ncols_Qi, ncols_X);
1008 for (
int i = 0; i < num_Q_blocks; ++i)
1010 const int ncols_Qi = MVT::GetNumberVecs (*Q[i]);
1016 template<
class Scalar,
class MV>
1021 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1024 #endif // BELOS_TEUCHOS_TIME_MONITOR
1026 const int numVecs = MVT::GetNumberVecs(X);
1029 }
else if (numVecs == 1) {
1036 const int rank = normalizeOne (X, B);
1038 RCP<MV> Q_0 = MVT::CloneViewNonConst (Q, Range1D(0,0));
1039 MVT::Assign (X, *Q_0);
1045 return normalizeImpl (X, Q, B,
true);
1049 template<
class Scalar,
class MV>
1054 const bool outOfPlace,
1063 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1067 #endif // BELOS_TEUCHOS_TIME_MONITOR
1072 std::invalid_argument,
1073 "Belos::TsqrOrthoManagerImpl::"
1074 "projectAndNormalizeImpl(..., outOfPlace=true, ...):"
1075 "X_out has " << MVT::GetNumberVecs(X_out)
1076 <<
" columns, but X_in has "
1077 << MVT::GetNumberVecs(X_in) <<
" columns.");
1081 int ncols_X, num_Q_blocks, ncols_Q_total;
1082 checkProjectionDims (ncols_X, num_Q_blocks, ncols_Q_total, X_in, Q);
1089 if (num_Q_blocks == 0 || ncols_Q_total == 0) {
1091 return normalizeOutOfPlace (X_in, X_out, B);
1093 return normalize (X_in, B);
1100 allocateProjectionCoefficients (C, Q, X_in,
true);
1105 std::vector<magnitude_type> normsBeforeFirstPass (ncols_X, SCTM::zero());
1106 if (reorthogonalizeBlocks_) {
1107 MVT::MvNorm (X_in, normsBeforeFirstPass);
1111 rawProject (X_in, Q, C);
1128 std::invalid_argument,
1129 "normalizeOne: Input matrix B must be at "
1130 "least " << ncols_X <<
" x " << ncols_X
1131 <<
", but is instead " << B->numRows()
1132 <<
" x " << B->numCols() <<
".");
1142 const int firstPassRank = outOfPlace ?
1143 normalizeOutOfPlace (X_in, X_out, B_out) :
1144 normalize (X_in, B_out);
1153 int rank = firstPassRank;
1169 if (firstPassRank < ncols_X && randomizeNullSpace_) {
1170 const int numNullSpaceCols = ncols_X - firstPassRank;
1171 const Range1D nullSpaceIndices (firstPassRank, ncols_X - 1);
1175 for (
int k = 0; k < num_Q_blocks; ++k) {
1176 const int numColsQk = MVT::GetNumberVecs(*Q[k]);
1177 C_null[k] =
rcp (
new mat_type (numColsQk, numNullSpaceCols));
1180 RCP<mat_type> B_null (
new mat_type (numNullSpaceCols, numNullSpaceCols));
1182 int randomVectorsRank;
1186 RCP<MV> X_out_null = MVT::CloneViewNonConst (X_out, nullSpaceIndices);
1191 RCP<MV> X_in_null = MVT::CloneViewNonConst (X_in, nullSpaceIndices);
1192 MVT::Assign (*X_out_null, *X_in_null);
1195 rawProject (*X_in_null, Q, C_null);
1196 randomVectorsRank = normalizeOutOfPlace (*X_in_null, *X_out_null, B_null);
1200 RCP<MV> X_null = MVT::CloneViewNonConst (X_in, nullSpaceIndices);
1203 rawProject (*X_null, Q, C_null);
1204 randomVectorsRank = normalize (*X_null, B_null);
1212 "Belos::TsqrOrthoManagerImpl::projectAndNormalize"
1213 "OutOfPlace(): After projecting and normalizing the "
1214 "random vectors (used to replace the null space "
1215 "basis vectors from normalizing X), they have rank "
1216 << randomVectorsRank <<
", but should have full "
1217 "rank " << numNullSpaceCols <<
".");
1222 if (reorthogonalizeBlocks_) {
1225 std::vector<magnitude_type>
1226 normsAfterFirstPass (firstPassRank, SCTM::zero());
1227 std::vector<magnitude_type>
1228 normsAfterSecondPass (firstPassRank, SCTM::zero());
1243 for (
int j = 0; j < firstPassRank; ++j) {
1244 const Scalar*
const B_j = &(*B_out)(0,j);
1247 normsAfterFirstPass[j] = blas.
NRM2 (firstPassRank, B_j, 1);
1251 bool reorthogonalize =
false;
1252 for (
int j = 0; j < firstPassRank; ++j) {
1260 blockReorthogThreshold() * normsBeforeFirstPass[j];
1261 if (normsAfterFirstPass[j] < curThreshold) {
1262 reorthogonalize =
true;
1269 if (! reorthogCallback_.is_null()) {
1270 using Teuchos::arrayViewFromVector;
1271 (*reorthogCallback_) (arrayViewFromVector (normsBeforeFirstPass),
1272 arrayViewFromVector (normsAfterFirstPass));
1286 bool reorthogFault =
false;
1288 std::vector<int> faultIndices;
1289 if (reorthogonalize) {
1297 MVT::Assign (X_out, X_in);
1304 allocateProjectionCoefficients (C2, Q, X_in,
false);
1309 rawProject (X_in, Q, C2);
1312 RCP<mat_type>
B2 (
new mat_type (ncols_X, ncols_X));
1315 const int secondPassRank = outOfPlace ?
1316 normalizeOutOfPlace (X_in, X_out, B2) :
1317 normalize (X_in, B2);
1318 rank = secondPassRank;
1323 mat_type B_copy (
Copy, *B_out, B_out->numRows(), B_out->numCols());
1326 *B2, B_copy, SCT::zero());
1328 "Teuchos::SerialDenseMatrix::multiply "
1329 "returned err = " << err <<
" != 0");
1333 for (
int k = 0; k < num_Q_blocks; ++k) {
1334 mat_type C_k_copy (
Copy, *C[k], C[k]->numRows(), C[k]->numCols());
1338 *C2[k], B_copy, SCT::one());
1340 "Teuchos::SerialDenseMatrix::multiply "
1341 "returned err = " << err1 <<
" != 0");
1346 for (
int j = 0; j < rank; ++j) {
1347 const Scalar*
const B2_j = &(*B2)(0,j);
1348 normsAfterSecondPass[j] = blas.
NRM2 (rank, B2_j, 1);
1353 reorthogFault =
false;
1354 for (
int j = 0; j < rank; ++j) {
1356 blockReorthogThreshold() * normsAfterFirstPass[j];
1357 if (normsAfterSecondPass[j] < relativeLowerBound) {
1358 reorthogFault =
true;
1359 faultIndices.push_back (j);
1364 if (reorthogFault) {
1365 if (throwOnReorthogFault_) {
1366 raiseReorthogFault (normsAfterFirstPass,
1367 normsAfterSecondPass,
1376 "TsqrOrthoManagerImpl has not yet implemented"
1377 " recovery from an orthogonalization fault.");
1385 template<
class Scalar,
class MV>
1389 const std::vector<magnitude_type>& normsAfterSecondPass,
1390 const std::vector<int>& faultIndices)
1393 typedef std::vector<int>::size_type size_type;
1394 std::ostringstream os;
1396 os <<
"Orthogonalization fault at the following column(s) of X:" << endl;
1397 os <<
"Column\tNorm decrease factor" << endl;
1398 for (size_type k = 0; k < faultIndices.size(); ++k) {
1399 const int index = faultIndices[k];
1401 normsAfterSecondPass[index] / normsAfterFirstPass[index];
1402 os << index <<
"\t" << decreaseFactor << endl;
1407 template<
class Scalar,
class MV>
1412 using Teuchos::parameterList;
1415 if (defaultParams_.is_null()) {
1416 RCP<ParameterList> params = parameterList (
"TsqrOrthoManagerImpl");
1420 params->set (
"TSQR implementation", *(tsqrAdaptor_.getValidParameters()),
1421 "TSQR implementation parameters.");
1425 const bool defaultRandomizeNullSpace =
true;
1426 params->set (
"randomizeNullSpace", defaultRandomizeNullSpace,
1427 "Whether to fill in null space vectors with random data.");
1429 const bool defaultReorthogonalizeBlocks =
true;
1430 params->set (
"reorthogonalizeBlocks", defaultReorthogonalizeBlocks,
1431 "Whether to do block reorthogonalization as necessary.");
1437 params->set (
"blockReorthogThreshold", defaultBlockReorthogThreshold,
1438 "If reorthogonalizeBlocks==true, and if the norm of "
1439 "any column within a block decreases by this much or "
1440 "more after orthogonalization, we reorthogonalize.");
1445 Teuchos::as<magnitude_type>(10) * SCTM::eps();
1450 params->set (
"relativeRankTolerance", defaultRelativeRankTolerance,
1451 "Relative tolerance to determine the numerical rank of a "
1452 "block when normalizing.");
1456 const bool defaultThrowOnReorthogFault =
true;
1457 params->set (
"throwOnReorthogFault", defaultThrowOnReorthogFault,
1458 "Whether to throw an exception if an orthogonalization "
1459 "fault occurs. This only matters if reorthogonalization "
1460 "is enabled (reorthogonalizeBlocks==true).");
1462 const bool defaultForceNonnegativeDiagonal =
false;
1463 params->set (
"forceNonnegativeDiagonal", defaultForceNonnegativeDiagonal,
1464 "Whether to force the R factor produced by the normalization "
1465 "step to have a nonnegative diagonal.");
1467 defaultParams_ = params;
1469 return defaultParams_;
1472 template<
class Scalar,
class MV>
1480 RCP<const ParameterList> defaultParams = getValidParameters();
1482 RCP<ParameterList> params =
rcp (
new ParameterList (*defaultParams));
1491 const bool randomizeNullSpace =
false;
1492 params->set (
"randomizeNullSpace", randomizeNullSpace);
1493 const bool reorthogonalizeBlocks =
false;
1494 params->set (
"reorthogonalizeBlocks", reorthogonalizeBlocks);
1499 template<
class Scalar,
class MV>
1511 tsqrAdaptor_.factorExplicit (X, Q, B, forceNonnegativeDiagonal_);
1514 rank = tsqrAdaptor_.revealRank (Q, B, relativeRankTolerance_);
1515 }
catch (std::exception& e) {
1521 template<
class Scalar,
class MV>
1534 const int theNumRows =
B->numRows ();
1535 const int theNumCols =
B->numCols ();
1537 theNumRows < 1 || theNumCols < 1, std::invalid_argument,
1538 "normalizeOne: Input matrix B must be at least 1 x 1, but "
1539 "is instead " << theNumRows <<
" x " << theNumCols <<
".");
1545 std::vector<magnitude_type> theNorm (1, SCTM::zero());
1546 MVT::MvNorm (X, theNorm);
1547 (*B_out)(0,0) = theNorm[0];
1561 if (theNorm[0] == SCTM::zero()) {
1564 if (randomizeNullSpace_) {
1566 MVT::MvNorm (X, theNorm);
1567 if (theNorm[0] == SCTM::zero()) {
1573 "vector has norm zero!");
1578 const Scalar alpha = SCT::one() / theNorm[0];
1579 MVT::MvScale (X, alpha);
1586 const Scalar alpha = SCT::one() / theNorm[0];
1587 MVT::MvScale (X, alpha);
1593 template<
class Scalar,
class MV>
1600 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1602 #endif // BELOS_TEUCHOS_TIME_MONITOR
1605 const int num_Q_blocks = Q.size();
1606 for (
int i = 0; i < num_Q_blocks; ++i)
1614 mat_type& Ci = *
C[i];
1615 const MV& Qi = *Q[i];
1616 innerProd (Qi, X, Ci);
1617 MVT::MvTimesMatAddMv (-SCT::one(), Qi, Ci, SCT::one(), X);
1622 template<
class Scalar,
class MV>
1629 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1631 #endif // BELOS_TEUCHOS_TIME_MONITOR
1634 innerProd (*Q, X, *
C);
1635 MVT::MvTimesMatAddMv (-SCT::one(), *Q, *
C, SCT::one(), X);
1638 template<
class Scalar,
class MV>
1644 const bool outOfPlace)
1650 using Teuchos::tuple;
1652 const int numCols = MVT::GetNumberVecs (X);
1660 MVT::GetNumberVecs (Q) < numCols, std::invalid_argument,
1661 "TsqrOrthoManagerImpl::normalizeImpl: Q has "
1662 << MVT::GetNumberVecs(Q) <<
" columns. This is too "
1663 "few, since X has " << numCols <<
" columns.");
1667 RCP<MV> Q_view = MVT::CloneViewNonConst (Q, Range1D (0, numCols-1));
1679 B->numRows() < numCols ||
B->numCols() < numCols, std::invalid_argument,
1680 "TsqrOrthoManagerImpl::normalizeImpl: Input matrix B must be at least "
1681 << numCols <<
" x " << numCols <<
", but is instead " <<
B->numRows ()
1682 <<
" x " <<
B->numCols() <<
".");
1697 const int rank = rawNormalize (X, *Q_view, *B_out);
1708 rank < 0 || rank > numCols, std::logic_error,
1709 "Belos::TsqrOrthoManagerImpl::normalizeImpl: rawNormalize returned rank "
1710 " = " << rank <<
" for a matrix X with " << numCols <<
" columns. "
1711 "Please report this bug to the Belos developers.");
1715 if (rank == numCols || ! randomizeNullSpace_) {
1719 MVT::Assign (*Q_view, X);
1724 if (randomizeNullSpace_ && rank < numCols) {
1731 const int nullSpaceNumCols = numCols - rank;
1734 Range1D nullSpaceIndices (rank, numCols-1);
1741 RCP<MV> Q_null = MVT::CloneViewNonConst (*Q_view, nullSpaceIndices);
1743 MVT::MvRandom (*Q_null);
1749 std::vector<magnitude_type> norms (MVT::GetNumberVecs (*Q_null));
1750 MVT::MvNorm (*Q_null, norms);
1752 bool anyZero =
false;
1753 typedef typename std::vector<magnitude_type>::const_iterator iter_type;
1754 for (iter_type it = norms.begin(); it != norms.end(); ++it) {
1755 if (*it == SCTM::zero()) {
1760 std::ostringstream os;
1761 os <<
"TsqrOrthoManagerImpl::normalizeImpl: "
1762 "We are being asked to randomize the null space, for a matrix "
1763 "with " << numCols <<
" columns and reported column rank "
1764 << rank <<
". The inclusive range of columns to fill with "
1765 "random data is [" << nullSpaceIndices.lbound() <<
","
1766 << nullSpaceIndices.ubound() <<
"]. After filling the null "
1767 "space vectors with random numbers, at least one of the vectors"
1768 " has norm zero. Here are the norms of all the null space "
1770 for (iter_type it = norms.begin(); it != norms.end(); ++it) {
1772 if (it+1 != norms.end())
1775 os <<
"].) There is a tiny probability that this could happen "
1776 "randomly, but it is likely a bug. Please report it to the "
1777 "Belos developers, especially if you are able to reproduce the "
1790 RCP<const MV> Q_col = MVT::CloneView (*Q_view, Range1D (0, rank-1));
1796 rawProject (*Q_null, Q_col, C_null);
1805 RCP<MV> X_null = MVT::CloneViewNonConst (X, nullSpaceIndices);
1808 mat_type B_null (nullSpaceNumCols, nullSpaceNumCols);
1810 const int nullSpaceBasisRank = rawNormalize (*Q_null, *X_null, B_null);
1824 if (nullSpaceBasisRank < nullSpaceNumCols) {
1825 std::vector<magnitude_type> norms (MVT::GetNumberVecs(*X_null));
1826 MVT::MvNorm (*X_null, norms);
1827 std::ostringstream os;
1828 os <<
"TsqrOrthoManagerImpl::normalizeImpl: "
1829 <<
"We are being asked to randomize the null space, "
1830 <<
"for a matrix with " << numCols <<
" columns and "
1831 <<
"column rank " << rank <<
". After projecting and "
1832 <<
"normalizing the generated random vectors, they "
1833 <<
"only have rank " << nullSpaceBasisRank <<
". They"
1834 <<
" should have full rank " << nullSpaceNumCols
1835 <<
". (The inclusive range of columns to fill with "
1836 <<
"random data is [" << nullSpaceIndices.lbound()
1837 <<
"," << nullSpaceIndices.ubound() <<
"]. The "
1838 <<
"column norms of the resulting Q factor are: [";
1839 for (
typename std::vector<magnitude_type>::size_type k = 0;
1840 k < norms.size(); ++k) {
1842 if (k != norms.size()-1) {
1846 os <<
"].) There is a tiny probability that this could "
1847 <<
"happen randomly, but it is likely a bug. Please "
1848 <<
"report it to the Belos developers, especially if "
1849 <<
"you are able to reproduce the behavior.";
1862 MVT::Assign (*X_null, *Q_null);
1863 }
else if (rank > 0) {
1865 RCP<const MV> Q_col = MVT::CloneView (*Q_view, Range1D (0, rank-1));
1866 RCP<MV> X_col = MVT::CloneViewNonConst (X, Range1D (0, rank-1));
1867 MVT::Assign (*Q_col, *X_col);
1874 template<
class Scalar,
class MV>
1889 int the_ncols_X, the_num_Q_blocks, the_ncols_Q_total;
1890 the_num_Q_blocks = Q.size();
1891 the_ncols_X = MVT::GetNumberVecs (X);
1894 the_ncols_Q_total = 0;
1899 typedef typename ArrayView<RCP<const MV> >::const_iterator iter_type;
1900 for (iter_type it = Q.begin(); it != Q.end(); ++it) {
1901 const MV& Qi = **it;
1902 the_ncols_Q_total += MVT::GetNumberVecs (Qi);
1906 ncols_X = the_ncols_X;
1907 num_Q_blocks = the_num_Q_blocks;
1908 ncols_Q_total = the_ncols_Q_total;
1913 #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