45 #ifndef __AnasaziTsqrOrthoManagerImpl_hpp 46 #define __AnasaziTsqrOrthoManagerImpl_hpp 52 #include "Teuchos_as.hpp" 53 #include "Teuchos_LAPACK.hpp" 54 #include "Teuchos_ParameterList.hpp" 55 #include "Teuchos_ParameterListAcceptorDefaultBase.hpp" 127 template<
class Scalar,
class MV>
129 public Teuchos::ParameterListAcceptorDefaultBase {
131 typedef Scalar scalar_type;
132 typedef typename Teuchos::ScalarTraits<Scalar>::magnitudeType magnitude_type;
133 typedef MV multivector_type;
136 typedef Teuchos::SerialDenseMatrix<int, Scalar>
mat_type;
137 typedef Teuchos::RCP<mat_type> mat_ptr;
140 typedef Teuchos::ScalarTraits<Scalar> SCT;
141 typedef Teuchos::ScalarTraits<magnitude_type> SCTM;
143 typedef typename MVT::tsqr_adaptor_type tsqr_adaptor_type;
153 Teuchos::RCP<const Teuchos::ParameterList> getValidParameters ()
const;
156 void setParameterList (
const Teuchos::RCP<Teuchos::ParameterList>& params);
168 Teuchos::RCP<const Teuchos::ParameterList> getFastParameters ();
187 const std::string& label);
203 if (label != label_) {
209 const std::string&
getLabel ()
const {
return label_; }
222 MVT::MvTransMv (SCT::one(), X, Y, Z);
243 norm (
const MV& X, std::vector<magnitude_type>& normvec)
const;
256 Teuchos::Array<mat_ptr> C,
257 Teuchos::ArrayView<Teuchos::RCP<const MV> > Q);
272 int normalize (MV& X, mat_ptr B);
293 normalizeOutOfPlace (MV& X, MV& Q, mat_ptr B);
310 Teuchos::Array<mat_ptr> C,
312 Teuchos::ArrayView<Teuchos::RCP<const MV> > Q)
316 return projectAndNormalizeImpl (X, X,
false, C, B, Q);
341 Teuchos::Array<mat_ptr> C,
343 Teuchos::ArrayView<Teuchos::RCP<const MV> > Q)
347 return projectAndNormalizeImpl (X_in, X_out,
true, C, B, Q);
357 const Scalar ONE = SCT::one();
358 const int ncols = MVT::GetNumberVecs(X);
359 mat_type XTX (ncols, ncols);
360 innerProd (X, X, XTX);
361 for (
int k = 0; k < ncols; ++k) {
364 return XTX.normFrobenius();
372 const int ncols_X1 = MVT::GetNumberVecs (X1);
373 const int ncols_X2 = MVT::GetNumberVecs (X2);
374 mat_type X1_T_X2 (ncols_X1, ncols_X2);
375 innerProd (X1, X2, X1_T_X2);
376 return X1_T_X2.normFrobenius();
390 Teuchos::RCP<Teuchos::ParameterList> params_;
393 mutable Teuchos::RCP<const Teuchos::ParameterList> defaultParams_;
399 tsqr_adaptor_type tsqrAdaptor_;
419 bool randomizeNullSpace_;
426 bool reorthogonalizeBlocks_;
431 bool throwOnReorthogFault_;
434 magnitude_type blockReorthogThreshold_;
437 magnitude_type relativeRankTolerance_;
445 bool forceNonnegativeDiagonal_;
449 raiseReorthogFault (
const std::vector<magnitude_type>& normsAfterFirstPass,
450 const std::vector<magnitude_type>& normsAfterSecondPass,
451 const std::vector<int>& faultIndices);
463 checkProjectionDims (
int& ncols_X,
467 Teuchos::ArrayView<Teuchos::RCP<const MV> > Q)
const;
480 allocateProjectionCoefficients (Teuchos::Array<mat_ptr>& C,
481 Teuchos::ArrayView<Teuchos::RCP<const MV> > Q,
483 const bool attemptToRecycle =
true)
const;
494 projectAndNormalizeImpl (MV& X_in,
496 const bool outOfPlace,
497 Teuchos::Array<mat_ptr> C,
499 Teuchos::ArrayView<Teuchos::RCP<const MV> > Q);
508 Teuchos::ArrayView<Teuchos::RCP<const MV> > Q,
509 Teuchos::ArrayView<mat_ptr> C)
const;
514 const Teuchos::RCP<const MV>& Q,
515 const mat_ptr& C)
const;
543 int rawNormalize (MV& X, MV& Q, mat_type& B);
562 int normalizeOne (MV& X, mat_ptr B)
const;
591 int normalizeImpl (MV& X, MV& Q, mat_ptr B,
const bool outOfPlace);
594 template<
class Scalar,
class MV>
599 using Teuchos::ParameterList;
600 using Teuchos::parameterList;
602 using Teuchos::sublist;
603 typedef magnitude_type M;
605 RCP<const ParameterList> defaultParams = getValidParameters ();
607 RCP<ParameterList> tsqrParams;
609 RCP<ParameterList> theParams;
610 if (params.is_null()) {
611 theParams = parameterList (*defaultParams);
619 randomizeNullSpace_ =
620 theParams->get<
bool> (
"randomizeNullSpace",
621 defaultParams->get<
bool> (
"randomizeNullSpace"));
622 reorthogonalizeBlocks_ =
623 theParams->get<
bool> (
"reorthogonalizeBlocks",
624 defaultParams->get<
bool> (
"reorthogonalizeBlocks"));
625 throwOnReorthogFault_ =
626 theParams->get<
bool> (
"throwOnReorthogFault",
627 defaultParams->get<
bool> (
"throwOnReorthogFault"));
628 blockReorthogThreshold_ =
629 theParams->get<M> (
"blockReorthogThreshold",
630 defaultParams->get<M> (
"blockReorthogThreshold"));
631 relativeRankTolerance_ =
632 theParams->get<M> (
"relativeRankTolerance",
633 defaultParams->get<M> (
"relativeRankTolerance"));
634 forceNonnegativeDiagonal_ =
635 theParams->get<
bool> (
"forceNonnegativeDiagonal",
636 defaultParams->get<
bool> (
"forceNonnegativeDiagonal"));
640 if (! theParams->isSublist (
"TSQR implementation")) {
641 theParams->set (
"TSQR implementation",
642 defaultParams->sublist (
"TSQR implementation"));
644 tsqrParams = sublist (theParams,
"TSQR implementation",
true);
648 tsqrAdaptor_.setParameterList (tsqrParams);
651 setMyParamList (theParams);
654 template<
class Scalar,
class MV>
657 const std::string& label) :
661 randomizeNullSpace_ (true),
662 reorthogonalizeBlocks_ (true),
663 throwOnReorthogFault_ (false),
664 blockReorthogThreshold_ (0),
665 relativeRankTolerance_ (0),
666 forceNonnegativeDiagonal_ (false)
671 template<
class Scalar,
class MV>
677 randomizeNullSpace_ (true),
678 reorthogonalizeBlocks_ (true),
679 throwOnReorthogFault_ (false),
680 blockReorthogThreshold_ (0),
681 relativeRankTolerance_ (0),
682 forceNonnegativeDiagonal_ (false)
687 template<
class Scalar,
class MV>
690 norm (
const MV& X, std::vector<magnitude_type>& normVec)
const 695 if (normVec.size() <
static_cast<size_t>(numCols))
696 normVec.resize (numCols);
700 template<
class Scalar,
class MV>
703 Teuchos::Array<mat_ptr> C,
704 Teuchos::ArrayView<Teuchos::RCP<const MV> > Q)
706 int ncols_X, num_Q_blocks, ncols_Q_total;
707 checkProjectionDims (ncols_X, num_Q_blocks, ncols_Q_total, X, Q);
711 if (ncols_X == 0 || num_Q_blocks == 0 || ncols_Q_total == 0)
715 allocateProjectionCoefficients (C, Q, X,
true);
719 std::vector<magnitude_type> columnNormsBefore (ncols_X, magnitude_type(0));
720 if (reorthogonalizeBlocks_)
725 rawProject (X, Q, C);
729 if (reorthogonalizeBlocks_)
731 std::vector<magnitude_type> columnNormsAfter (ncols_X, magnitude_type(0));
741 bool reorthogonalize =
false;
742 for (
int j = 0; j < ncols_X; ++j)
743 if (columnNormsAfter[j] < relThres * columnNormsBefore[j])
745 reorthogonalize =
true;
751 Teuchos::Array<mat_ptr> C2;
752 allocateProjectionCoefficients (C2, Q, X,
false);
756 rawProject (X, Q, C2);
758 for (
int k = 0; k < num_Q_blocks; ++k)
766 template<
class Scalar,
class MV>
770 using Teuchos::Range1D;
785 return normalizeOne (X, B);
827 return normalizeImpl (X, *Q_, B,
false);
830 return normalizeImpl (X, *Q_view, B,
false);
834 template<
class Scalar,
class MV>
838 Teuchos::ArrayView<Teuchos::RCP<const MV> > Q,
840 const bool attemptToRecycle)
const 842 const int num_Q_blocks = Q.size();
844 C.resize (num_Q_blocks);
845 if (attemptToRecycle)
847 for (
int i = 0; i < num_Q_blocks; ++i)
853 C[i] = rcp (
new mat_type (ncols_Qi, ncols_X));
857 if (Ci.numRows() != ncols_Qi || Ci.numCols() != ncols_X)
858 Ci.shape (ncols_Qi, ncols_X);
860 Ci.putScalar (SCT::zero());
866 for (
int i = 0; i < num_Q_blocks; ++i)
869 C[i] = rcp (
new mat_type (ncols_Qi, ncols_X));
874 template<
class Scalar,
class MV>
882 }
else if (numVecs == 1) {
884 using Teuchos::Range1D;
889 const int rank = normalizeOne (X, B);
898 return normalizeImpl (X, Q, B,
true);
902 template<
class Scalar,
class MV>
907 const bool outOfPlace,
908 Teuchos::Array<mat_ptr> C,
910 Teuchos::ArrayView<Teuchos::RCP<const MV> > Q)
912 using Teuchos::Range1D;
919 std::invalid_argument,
920 "Belos::TsqrOrthoManagerImpl::" 921 "projectAndNormalizeOutOfPlace(...):" 923 <<
" columns, but X_in has " 928 int ncols_X, num_Q_blocks, ncols_Q_total;
929 checkProjectionDims (ncols_X, num_Q_blocks, ncols_Q_total, X_in, Q);
936 if (num_Q_blocks == 0 || ncols_Q_total == 0) {
947 allocateProjectionCoefficients (C, Q, X_in,
true);
952 std::vector<magnitude_type> normsBeforeFirstPass (ncols_X, SCTM::zero());
953 if (reorthogonalizeBlocks_) {
958 rawProject (X_in, Q, C);
971 B_out = rcp (
new mat_type (ncols_X, ncols_X));
974 TEUCHOS_TEST_FOR_EXCEPTION(B->numRows() < ncols_X || B->numCols() < ncols_X,
975 std::invalid_argument,
976 "normalizeOne: Input matrix B must be at " 977 "least " << ncols_X <<
" x " << ncols_X
978 <<
", but is instead " << B->numRows()
979 <<
" x " << B->numCols() <<
".");
983 B_out = rcp (
new mat_type (Teuchos::View, *B, ncols_X, ncols_X));
989 const int firstPassRank = outOfPlace ?
1000 int rank = firstPassRank;
1016 if (firstPassRank < ncols_X && randomizeNullSpace_) {
1017 const int numNullSpaceCols = ncols_X - firstPassRank;
1018 const Range1D nullSpaceIndices (firstPassRank, ncols_X - 1);
1021 Teuchos::Array<mat_ptr> C_null (num_Q_blocks);
1022 for (
int k = 0; k < num_Q_blocks; ++k) {
1024 C_null[k] = rcp (
new mat_type (numColsQk, numNullSpaceCols));
1027 RCP<mat_type> B_null (
new mat_type (numNullSpaceCols, numNullSpaceCols));
1029 int randomVectorsRank;
1042 rawProject (*X_in_null, Q, C_null);
1050 rawProject (*X_null, Q, C_null);
1051 randomVectorsRank =
normalize (*X_null, B_null);
1057 TEUCHOS_TEST_FOR_EXCEPTION(randomVectorsRank != numNullSpaceCols,
1059 "Belos::TsqrOrthoManagerImpl::projectAndNormalize" 1060 "OutOfPlace(): After projecting and normalizing the " 1061 "random vectors (used to replace the null space " 1062 "basis vectors from normalizing X), they have rank " 1063 << randomVectorsRank <<
", but should have full " 1064 "rank " << numNullSpaceCols <<
".");
1069 if (reorthogonalizeBlocks_) {
1072 std::vector<magnitude_type>
1073 normsAfterFirstPass (firstPassRank, SCTM::zero());
1074 std::vector<magnitude_type>
1075 normsAfterSecondPass (firstPassRank, SCTM::zero());
1089 Teuchos::BLAS<int, Scalar> blas;
1090 for (
int j = 0; j < firstPassRank; ++j) {
1091 const Scalar*
const B_j = &(*B_out)(0,j);
1094 normsAfterFirstPass[j] = blas.NRM2 (firstPassRank, B_j, 1);
1098 bool reorthogonalize =
false;
1099 for (
int j = 0; j < firstPassRank; ++j) {
1106 const magnitude_type curThreshold =
1108 if (normsAfterFirstPass[j] < curThreshold) {
1109 reorthogonalize =
true;
1124 bool reorthogFault =
false;
1126 std::vector<int> faultIndices;
1127 if (reorthogonalize) {
1128 using Teuchos::Copy;
1129 using Teuchos::NO_TRANS;
1140 Teuchos::Array<mat_ptr> C2;
1141 allocateProjectionCoefficients (C2, Q, X_in,
false);
1146 rawProject (X_in, Q, C2);
1149 RCP<mat_type> B2 (
new mat_type (ncols_X, ncols_X));
1152 const int secondPassRank = outOfPlace ?
1155 rank = secondPassRank;
1160 mat_type B_copy (Copy, *B_out, B_out->numRows(), B_out->numCols());
1162 const int err = B_out->multiply (NO_TRANS, NO_TRANS, SCT::one(),
1163 *B2, B_copy, SCT::zero());
1164 TEUCHOS_TEST_FOR_EXCEPTION(err != 0, std::logic_error,
1165 "Teuchos::SerialDenseMatrix::multiply " 1166 "returned err = " << err <<
" != 0");
1170 for (
int k = 0; k < num_Q_blocks; ++k) {
1171 mat_type C_k_copy (Copy, *C[k], C[k]->numRows(), C[k]->numCols());
1174 const int err2 = C[k]->multiply (NO_TRANS, NO_TRANS, SCT::one(),
1175 *C2[k], B_copy, SCT::one());
1176 TEUCHOS_TEST_FOR_EXCEPTION(err2 != 0, std::logic_error,
1177 "Teuchos::SerialDenseMatrix::multiply " 1178 "returned err = " << err <<
" != 0");
1183 for (
int j = 0; j < rank; ++j) {
1184 const Scalar*
const B2_j = &(*B2)(0,j);
1185 normsAfterSecondPass[j] = blas.NRM2 (rank, B2_j, 1);
1190 reorthogFault =
false;
1191 for (
int j = 0; j < rank; ++j) {
1192 const magnitude_type relativeLowerBound =
1194 if (normsAfterSecondPass[j] < relativeLowerBound) {
1195 reorthogFault =
true;
1196 faultIndices.push_back (j);
1201 if (reorthogFault) {
1202 if (throwOnReorthogFault_) {
1203 raiseReorthogFault (normsAfterFirstPass,
1204 normsAfterSecondPass,
1212 TEUCHOS_TEST_FOR_EXCEPTION(
true, std::logic_error,
1213 "TsqrOrthoManagerImpl has not yet implemented" 1214 " recovery from an orthogonalization fault.");
1222 template<
class Scalar,
class MV>
1226 const std::vector<magnitude_type>& normsAfterSecondPass,
1227 const std::vector<int>& faultIndices)
1230 typedef std::vector<int>::size_type size_type;
1231 std::ostringstream os;
1233 os <<
"Orthogonalization fault at the following column(s) of X:" << endl;
1234 os <<
"Column\tNorm decrease factor" << endl;
1235 for (size_type k = 0; k < faultIndices.size(); ++k)
1237 const int index = faultIndices[k];
1238 const magnitude_type decreaseFactor =
1239 normsAfterSecondPass[index] / normsAfterFirstPass[index];
1240 os << index <<
"\t" << decreaseFactor << endl;
1245 template<
class Scalar,
class MV>
1246 Teuchos::RCP<const Teuchos::ParameterList>
1249 using Teuchos::ParameterList;
1250 using Teuchos::parameterList;
1253 if (defaultParams_.is_null()) {
1254 RCP<ParameterList> params = parameterList (
"TsqrOrthoManagerImpl");
1258 params->set (
"TSQR implementation", *(tsqrAdaptor_.getValidParameters()),
1259 "TSQR implementation parameters.");
1263 const bool defaultRandomizeNullSpace =
true;
1264 params->set (
"randomizeNullSpace", defaultRandomizeNullSpace,
1265 "Whether to fill in null space vectors with random data.");
1267 const bool defaultReorthogonalizeBlocks =
true;
1268 params->set (
"reorthogonalizeBlocks", defaultReorthogonalizeBlocks,
1269 "Whether to do block reorthogonalization as necessary.");
1273 const magnitude_type defaultBlockReorthogThreshold =
1274 magnitude_type(10) * SCTM::squareroot (SCTM::eps());
1275 params->set (
"blockReorthogThreshold", defaultBlockReorthogThreshold,
1276 "If reorthogonalizeBlocks==true, and if the norm of " 1277 "any column within a block decreases by this much or " 1278 "more after orthogonalization, we reorthogonalize.");
1282 const magnitude_type defaultRelativeRankTolerance =
1283 Teuchos::as<magnitude_type>(10) * SCTM::eps();
1288 params->set (
"relativeRankTolerance", defaultRelativeRankTolerance,
1289 "Relative tolerance to determine the numerical rank of a " 1290 "block when normalizing.");
1294 const bool defaultThrowOnReorthogFault =
true;
1295 params->set (
"throwOnReorthogFault", defaultThrowOnReorthogFault,
1296 "Whether to throw an exception if an orthogonalization " 1297 "fault occurs. This only matters if reorthogonalization " 1298 "is enabled (reorthogonalizeBlocks==true).");
1300 const bool defaultForceNonnegativeDiagonal =
false;
1301 params->set (
"forceNonnegativeDiagonal", defaultForceNonnegativeDiagonal,
1302 "Whether to force the R factor produced by the normalization " 1303 "step to have a nonnegative diagonal.");
1305 defaultParams_ = params;
1307 return defaultParams_;
1310 template<
class Scalar,
class MV>
1311 Teuchos::RCP<const Teuchos::ParameterList>
1314 using Teuchos::ParameterList;
1320 RCP<ParameterList> params = rcp (
new ParameterList (*defaultParams));
1329 const bool randomizeNullSpace =
false;
1330 params->set (
"randomizeNullSpace", randomizeNullSpace);
1331 const bool reorthogonalizeBlocks =
false;
1332 params->set (
"reorthogonalizeBlocks", reorthogonalizeBlocks);
1337 template<
class Scalar,
class MV>
1342 Teuchos::SerialDenseMatrix<int, Scalar>& B)
1349 tsqrAdaptor_.factorExplicit (X, Q, B, forceNonnegativeDiagonal_);
1352 rank = tsqrAdaptor_.revealRank (Q, B, relativeRankTolerance_);
1353 }
catch (std::exception& e) {
1359 template<
class Scalar,
class MV>
1363 Teuchos::RCP<Teuchos::SerialDenseMatrix<int, Scalar> > B)
const 1372 const int numRows = B->numRows();
1373 const int numCols = B->numCols();
1374 TEUCHOS_TEST_FOR_EXCEPTION(numRows < 1 || numCols < 1,
1375 std::invalid_argument,
1376 "normalizeOne: Input matrix B must be at " 1377 "least 1 x 1, but is instead " << numRows
1378 <<
" x " << numCols <<
".");
1380 B_out = rcp (
new mat_type (Teuchos::View, *B, 1, 1));
1384 std::vector<magnitude_type> theNorm (1, SCTM::zero());
1386 (*B_out)(0,0) = theNorm[0];
1400 if (theNorm[0] == SCTM::zero()) {
1403 if (randomizeNullSpace_) {
1406 if (theNorm[0] == SCTM::zero()) {
1412 "vector has norm zero!");
1417 const Scalar alpha = SCT::one() / theNorm[0];
1425 const Scalar alpha = SCT::one() / theNorm[0];
1432 template<
class Scalar,
class MV>
1436 Teuchos::ArrayView<Teuchos::RCP<const MV> > Q,
1437 Teuchos::ArrayView<Teuchos::RCP<Teuchos::SerialDenseMatrix<int, Scalar> > > C)
const 1440 const int num_Q_blocks = Q.size();
1441 for (
int i = 0; i < num_Q_blocks; ++i)
1450 const MV& Qi = *Q[i];
1457 template<
class Scalar,
class MV>
1461 const Teuchos::RCP<const MV>& Q,
1462 const Teuchos::RCP<Teuchos::SerialDenseMatrix<int, Scalar> >& C)
const 1469 template<
class Scalar,
class MV>
1474 Teuchos::RCP<Teuchos::SerialDenseMatrix<int, Scalar> > B,
1475 const bool outOfPlace)
1477 using Teuchos::Range1D;
1480 using Teuchos::ScalarTraits;
1481 using Teuchos::tuple;
1486 const bool extraDebug =
false;
1496 std::invalid_argument,
1497 "TsqrOrthoManagerImpl::normalizeImpl(X,Q,B): Q has " 1499 "few, since X has " << numCols <<
" columns.");
1511 B_out = rcp (
new mat_type (numCols, numCols));
1514 TEUCHOS_TEST_FOR_EXCEPTION(B->numRows() < numCols || B->numCols() < numCols,
1515 std::invalid_argument,
1516 "normalizeOne: Input matrix B must be at " 1517 "least " << numCols <<
" x " << numCols
1518 <<
", but is instead " << B->numRows()
1519 <<
" x " << B->numCols() <<
".");
1522 B_out = rcp (
new mat_type (Teuchos::View, *B, numCols, numCols));
1526 std::vector<magnitude_type> norms (numCols);
1528 cerr <<
"Column norms of X before orthogonalization: ";
1529 typedef typename std::vector<magnitude_type>::const_iterator iter_type;
1530 for (iter_type iter = norms.begin(); iter != norms.end(); ++iter) {
1532 if (iter+1 != norms.end())
1546 const int rank = rawNormalize (X, *Q_view, *B_out);
1557 std::vector<magnitude_type> norms (numCols);
1559 cerr <<
"Column norms of Q_view after orthogonalization: ";
1560 for (
typename std::vector<magnitude_type>::const_iterator iter = norms.begin();
1561 iter != norms.end(); ++iter) {
1563 if (iter+1 != norms.end())
1567 TEUCHOS_TEST_FOR_EXCEPTION(rank < 0 || rank > numCols, std::logic_error,
1568 "Belos::TsqrOrthoManagerImpl::normalizeImpl: " 1569 "rawNormalize() returned rank = " << rank <<
" for a " 1570 "matrix X with " << numCols <<
" columns. Please report" 1571 " this bug to the Belos developers.");
1572 if (extraDebug && rank == 0) {
1576 std::vector<magnitude_type> norms (B_ref.numCols());
1577 for (
typename mat_type::ordinalType j = 0; j < B_ref.numCols(); ++j) {
1578 typedef typename mat_type::scalarType mat_scalar_type;
1579 mat_scalar_type sumOfSquares = ScalarTraits<mat_scalar_type>::zero();
1580 for (
typename mat_type::ordinalType i = 0; i <= j; ++i) {
1581 const mat_scalar_type B_ij =
1582 ScalarTraits<mat_scalar_type>::magnitude (B_ref(i,j));
1583 sumOfSquares += B_ij*B_ij;
1585 norms[j] = ScalarTraits<mat_scalar_type>::squareroot (sumOfSquares);
1589 cerr <<
"Norms of columns of B after orthogonalization: ";
1590 for (
typename mat_type::ordinalType j = 0; j < B_ref.numCols(); ++j) {
1592 if (j != B_ref.numCols() - 1)
1600 if (rank == numCols || ! randomizeNullSpace_) {
1609 if (randomizeNullSpace_ && rank < numCols) {
1616 const int nullSpaceNumCols = numCols - rank;
1619 Range1D nullSpaceIndices (rank, numCols-1);
1637 bool anyZero =
false;
1638 typedef typename std::vector<magnitude_type>::const_iterator iter_type;
1639 for (iter_type it = norms.begin(); it != norms.end(); ++it) {
1640 if (*it == SCTM::zero()) {
1645 std::ostringstream os;
1646 os <<
"TsqrOrthoManagerImpl::normalizeImpl: " 1647 "We are being asked to randomize the null space, for a matrix " 1648 "with " << numCols <<
" columns and reported column rank " 1649 << rank <<
". The inclusive range of columns to fill with " 1650 "random data is [" << nullSpaceIndices.lbound() <<
"," 1651 << nullSpaceIndices.ubound() <<
"]. After filling the null " 1652 "space vectors with random numbers, at least one of the vectors" 1653 " has norm zero. Here are the norms of all the null space " 1655 for (iter_type it = norms.begin(); it != norms.end(); ++it) {
1657 if (it+1 != norms.end())
1660 os <<
"].) There is a tiny probability that this could happen " 1661 "randomly, but it is likely a bug. Please report it to the " 1662 "Belos developers, especially if you are able to reproduce the " 1675 RCP<const MV> Q_col =
MVT::CloneView (*Q_view, Range1D(0, rank-1));
1680 mat_ptr C_null (
new mat_type (rank, nullSpaceNumCols));
1681 rawProject (*Q_null, Q_col, C_null);
1693 mat_type B_null (nullSpaceNumCols, nullSpaceNumCols);
1695 const int nullSpaceBasisRank = rawNormalize (*Q_null, *X_null, B_null);
1709 if (nullSpaceBasisRank < nullSpaceNumCols) {
1712 std::ostringstream os;
1713 os <<
"TsqrOrthoManagerImpl::normalizeImpl: " 1714 <<
"We are being asked to randomize the null space, " 1715 <<
"for a matrix with " << numCols <<
" columns and " 1716 <<
"column rank " << rank <<
". After projecting and " 1717 <<
"normalizing the generated random vectors, they " 1718 <<
"only have rank " << nullSpaceBasisRank <<
". They" 1719 <<
" should have full rank " << nullSpaceNumCols
1720 <<
". (The inclusive range of columns to fill with " 1721 <<
"random data is [" << nullSpaceIndices.lbound()
1722 <<
"," << nullSpaceIndices.ubound() <<
"]. The " 1723 <<
"column norms of the resulting Q factor are: [";
1724 for (
typename std::vector<magnitude_type>::size_type k = 0;
1725 k < norms.size(); ++k) {
1727 if (k != norms.size()-1) {
1731 os <<
"].) There is a tiny probability that this could " 1732 <<
"happen randomly, but it is likely a bug. Please " 1733 <<
"report it to the Belos developers, especially if " 1734 <<
"you are able to reproduce the behavior.";
1736 TEUCHOS_TEST_FOR_EXCEPTION(nullSpaceBasisRank < nullSpaceNumCols,
1748 }
else if (rank > 0) {
1750 RCP<const MV> Q_col =
MVT::CloneView (*Q_view, Range1D(0, rank-1));
1759 template<
class Scalar,
class MV>
1766 Teuchos::ArrayView<Teuchos::RCP<const MV> > Q)
const 1774 int the_ncols_X, the_num_Q_blocks, the_ncols_Q_total;
1775 the_num_Q_blocks = Q.size();
1779 the_ncols_Q_total = 0;
1782 using Teuchos::ArrayView;
1784 typedef typename ArrayView<RCP<const MV> >::const_iterator iter_type;
1785 for (iter_type it = Q.begin(); it != Q.end(); ++it) {
1786 const MV& Qi = **it;
1791 ncols_X = the_ncols_X;
1792 num_Q_blocks = the_num_Q_blocks;
1793 ncols_Q_total = the_ncols_Q_total;
1798 #endif // __AnasaziTsqrOrthoManagerImpl_hpp static ptrdiff_t GetGlobalLength(const MV &mv)
Return the number of rows in the given multivector mv.
static void MvRandom(MV &mv)
Replace the vectors in mv with random vectors.
static void MvTimesMatAddMv(const ScalarType alpha, const MV &A, const Teuchos::SerialDenseMatrix< int, ScalarType > &B, const ScalarType beta, MV &mv)
Update mv with .
void norm(const MV &X, std::vector< magnitude_type > &normvec) const
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
Default valid parameter list.
magnitude_type relativeRankTolerance() const
TSQR-based OrthoManager subclass implementation.
magnitude_type orthonormError(const MV &X) const
Return .
magnitude_type blockReorthogThreshold() const
magnitude_type orthogError(const MV &X1, const MV &X2) const
Return the Frobenius norm of the inner product of X1 with itself.
Declaration of basic traits for the multivector type.
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.
Teuchos::SerialDenseMatrix< int, Scalar > mat_type
static void Assign(const MV &A, MV &mv)
mv := A
void setParameterList(const Teuchos::RCP< Teuchos::ParameterList > ¶ms)
Set parameters from the given parameter list.
Namespace Anasazi contains the classes, structs, enums and utilities used by the Anasazi package...
TsqrOrthoManagerImpl(const Teuchos::RCP< Teuchos::ParameterList > ¶ms, const std::string &label)
Constructor (that sets user-specified parameters).
TsqrOrthoManager(Impl) error.
const std::string & getLabel() const
Get the label for timers (if timers are enabled).
void setLabel(const std::string &label)
Set the label for timers.
Templated virtual class for providing orthogonalization/orthonormalization methods.
Traits class which defines basic operations on multivectors.
static void MvScale(MV &mv, const ScalarType alpha)
Scale each element of the vectors in mv with alpha.
Teuchos::RCP< const Teuchos::ParameterList > getFastParameters()
Get "fast" parameters for TsqrOrthoManagerImpl.
void project(MV &X, Teuchos::Array< mat_ptr > C, Teuchos::ArrayView< Teuchos::RCP< const MV > > Q)
Compute and .
Anasazi header file which uses auto-configuration information to include necessary C++ headers...
static Teuchos::RCP< MV > CloneViewNonConst(MV &mv, const std::vector< int > &index)
Creates a new MV that shares the selected contents of mv (shallow copy).
Exception thrown to signal error in an orthogonalization manager method.
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.
static void MvNorm(const MV &mv, std::vector< typename Teuchos::ScalarTraits< ScalarType >::magnitudeType > &normvec)
Compute the 2-norm of each individual vector of mv. Upon return, normvec[i] holds the value of ...
int normalizeOutOfPlace(MV &X, MV &Q, mat_ptr B)
Normalize X into Q*B, overwriting X.
static Teuchos::RCP< const MV > CloneView(const MV &mv, const std::vector< int > &index)
Creates a new const MV that shares the selected contents of mv (shallow copy).
void innerProd(const MV &X, const MV &Y, mat_type &Z) const
Euclidean inner product.
static int GetNumberVecs(const MV &mv)
Obtain the number of vectors in mv.
int normalize(MV &X, mat_ptr B)
Orthogonalize the columns of X in place.
static Teuchos::RCP< MV > Clone(const MV &mv, const int numvecs)
Creates a new empty MV containing numvecs columns.