Belos Package Browser (Single Doxygen Collection)  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
BelosTsqrOrthoManagerImpl.hpp
Go to the documentation of this file.
1 // @HEADER
2 // *****************************************************************************
3 // Belos: Block Linear Solvers Package
4 //
5 // Copyright 2004-2016 NTESS and the Belos contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 // @HEADER
9 
13 #ifndef __BelosTsqrOrthoManagerImpl_hpp
14 #define __BelosTsqrOrthoManagerImpl_hpp
15 
16 #include "BelosConfigDefs.hpp" // HAVE_BELOS_TSQR
17 #include "BelosMultiVecTraits.hpp"
18 #include "BelosOrthoManager.hpp" // OrthoError, etc.
19 
20 #include "Teuchos_as.hpp"
23 #ifdef BELOS_TEUCHOS_TIME_MONITOR
24 # include "Teuchos_TimeMonitor.hpp"
25 #endif // BELOS_TEUCHOS_TIME_MONITOR
26 #include <algorithm>
27 #include <functional>
28 
29 namespace Belos {
30 
34  class TsqrOrthoError : public OrthoError {
35  public:
36  TsqrOrthoError (const std::string& what_arg) :
37  OrthoError (what_arg) {}
38  };
39 
59  class TsqrOrthoFault : public OrthoError {
60  public:
61  TsqrOrthoFault (const std::string& what_arg) :
62  OrthoError (what_arg) {}
63  };
64 
97  template<class Scalar>
99  {
100  public:
102  typedef Scalar scalar_type;
108 
110  virtual ~ReorthogonalizationCallback ();
111 
116  virtual void
117  operator() (Teuchos::ArrayView<magnitude_type> normsBeforeFirstPass,
118  Teuchos::ArrayView<magnitude_type> normsAfterFirstPass) = 0;
119  };
120 
121  template<class Scalar>
123 
155  template<class Scalar, class MV>
158  public:
159  typedef Scalar scalar_type;
161  typedef MV multivector_type;
165 
166  private:
170  typedef typename MVT::tsqr_adaptor_type tsqr_adaptor_type;
171 
172  public:
181 
184 
196 
214  const std::string& label);
215 
220  TsqrOrthoManagerImpl (const std::string& label);
221 
241  void
243  {
244  reorthogCallback_ = callback;
245  }
246 
254  void setLabel (const std::string& label) {
255  if (label != label_) {
256  label_ = label;
257 
258 #ifdef BELOS_TEUCHOS_TIME_MONITOR
259  clearTimer (label, "All orthogonalization");
260  clearTimer (label, "Projection");
261  clearTimer (label, "Normalization");
262 
263  timerOrtho_ = makeTimer (label, "All orthogonalization");
264  timerProject_ = makeTimer (label, "Projection");
265  timerNormalize_ = makeTimer (label, "Normalization");
266 #endif // BELOS_TEUCHOS_TIME_MONITOR
267  }
268  }
269 
271  const std::string& getLabel () const { return label_; }
272 
281  void
282  innerProd (const MV& X, const MV& Y, mat_type& Z) const
283  {
284  MVT::MvTransMv (SCT::one(), X, Y, Z);
285  }
286 
304  void
305  norm (const MV& X, std::vector<magnitude_type>& normVec) const;
306 
316  void
317  project (MV& X,
320 
334  int normalize (MV& X, mat_ptr B);
335 
354  int
355  normalizeOutOfPlace (MV& X, MV& Q, mat_ptr B);
356 
370  int
373  mat_ptr B,
375  {
376  // "false" means we work on X in place. The second argument is
377  // not read or written in that case.
378  return projectAndNormalizeImpl (X, X, false, C, B, Q);
379  }
380 
400  int
402  MV& X_out,
404  mat_ptr B,
406  {
407  // "true" means we work on X_in out of place, writing the
408  // results into X_out.
409  return projectAndNormalizeImpl (X_in, X_out, true, C, B, Q);
410  }
411 
417  orthonormError (const MV &X) const
418  {
419  const Scalar ONE = SCT::one();
420  const int ncols = MVT::GetNumberVecs(X);
421  mat_type XTX (ncols, ncols);
422  innerProd (X, X, XTX);
423  for (int k = 0; k < ncols; ++k) {
424  XTX(k,k) -= ONE;
425  }
426  return XTX.normFrobenius();
427  }
428 
431  orthogError (const MV &X1,
432  const MV &X2) const
433  {
434  const int ncols_X1 = MVT::GetNumberVecs (X1);
435  const int ncols_X2 = MVT::GetNumberVecs (X2);
436  mat_type X1_T_X2 (ncols_X1, ncols_X2);
437  innerProd (X1, X2, X1_T_X2);
438  return X1_T_X2.normFrobenius();
439  }
440 
445 
449 
450  private:
453 
456 
458  std::string label_;
459 
462 
473 
476 
481 
488 
493 
496 
499 
507 
508 #ifdef BELOS_TEUCHOS_TIME_MONITOR
509  Teuchos::RCP<Teuchos::Time> timerOrtho_;
511 
513  Teuchos::RCP<Teuchos::Time> timerProject_;
514 
516  Teuchos::RCP<Teuchos::Time> timerNormalize_;
517 #endif // BELOS_TEUCHOS_TIME_MONITOR
518 
521 
522 #ifdef BELOS_TEUCHOS_TIME_MONITOR
532  makeTimer (const std::string& prefix,
533  const std::string& timerName)
534  {
535  const std::string timerLabel =
536  prefix.empty() ? timerName : (prefix + ": " + timerName);
537  return Teuchos::TimeMonitor::getNewCounter (timerLabel);
538  }
539 
545  void
546  clearTimer (const std::string& prefix,
547  const std::string& timerName)
548  {
549  const std::string timerLabel =
550  prefix.empty() ? timerName : (prefix + ": " + timerName);
552  }
553 #endif // BELOS_TEUCHOS_TIME_MONITOR
554 
556  void
557  raiseReorthogFault (const std::vector<magnitude_type>& normsAfterFirstPass,
558  const std::vector<magnitude_type>& normsAfterSecondPass,
559  const std::vector<int>& faultIndices);
560 
570  void
571  checkProjectionDims (int& ncols_X,
572  int& num_Q_blocks,
573  int& ncols_Q_total,
574  const MV& X,
576 
587  void
590  const MV& X,
591  const bool attemptToRecycle = true) const;
592 
601  int
602  projectAndNormalizeImpl (MV& X_in,
603  MV& X_out,
604  const bool outOfPlace,
606  mat_ptr B,
608 
614  void
615  rawProject (MV& X,
618 
620  void
621  rawProject (MV& X,
622  const Teuchos::RCP<const MV>& Q,
623  const mat_ptr& C) const;
624 
651  int rawNormalize (MV& X, MV& Q, mat_type& B);
652 
670  int normalizeOne (MV& X, mat_ptr B) const;
671 
699  int normalizeImpl (MV& X, MV& Q, mat_ptr B, const bool outOfPlace);
700  };
701 
702  template<class Scalar, class MV>
703  void
706  {
708  using Teuchos::parameterList;
709  using Teuchos::RCP;
710  using Teuchos::sublist;
711  typedef magnitude_type M; // abbreviation.
712 
713  RCP<const ParameterList> defaultParams = getValidParameters ();
714  // Sublist of TSQR implementation parameters; to get below.
715  RCP<ParameterList> tsqrParams;
716 
717  RCP<ParameterList> theParams;
718  if (params.is_null()) {
719  theParams = parameterList (*defaultParams);
720  } else {
721  theParams = params;
722 
723  // Don't call validateParametersAndSetDefaults(); we prefer to
724  // ignore parameters that we don't recognize, at least for now.
725  // However, we do fill in missing parameters with defaults.
726 
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"));
745 
746  // Get the sublist of TSQR implementation parameters. Use the
747  // default sublist if one isn't provided.
748  if (! theParams->isSublist ("TSQR implementation")) {
749  theParams->set ("TSQR implementation",
750  defaultParams->sublist ("TSQR implementation"));
751  }
752  tsqrParams = sublist (theParams, "TSQR implementation", true);
753  }
754 
755  // Send the TSQR implementation parameters to the TSQR adaptor.
756  tsqrAdaptor_.setParameterList (tsqrParams);
757 
758  // Save the input parameter list.
759  setMyParamList (theParams);
760  }
761 
762  template<class Scalar, class MV>
765  const std::string& label) :
766  label_ (label),
767  Q_ (Teuchos::null), // Initialized on demand
768  eps_ (SCTM::eps()), // Machine precision
769  randomizeNullSpace_ (true),
770  reorthogonalizeBlocks_ (true),
771  throwOnReorthogFault_ (false),
772  blockReorthogThreshold_ (0),
773  relativeRankTolerance_ (0),
774  forceNonnegativeDiagonal_ (false)
775  {
776  setParameterList (params); // This also sets tsqrAdaptor_'s parameters.
777 
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
783  }
784 
785  template<class Scalar, class MV>
787  TsqrOrthoManagerImpl (const std::string& label) :
788  label_ (label),
789  Q_ (Teuchos::null), // Initialized on demand
790  eps_ (SCTM::eps()), // Machine precision
791  randomizeNullSpace_ (true),
792  reorthogonalizeBlocks_ (true),
793  throwOnReorthogFault_ (false),
794  blockReorthogThreshold_ (0),
795  relativeRankTolerance_ (0),
796  forceNonnegativeDiagonal_ (false)
797  {
798  setParameterList (Teuchos::null); // Set default parameters.
799 
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
805  }
806 
807  template<class Scalar, class MV>
808  void
810  norm (const MV& X, std::vector<magnitude_type>& normVec) const
811  {
812  const int numCols = MVT::GetNumberVecs (X);
813  // std::vector<T>::size_type is unsigned; int is signed. Mixed
814  // unsigned/signed comparisons trigger compiler warnings.
815  if (normVec.size() < static_cast<size_t>(numCols))
816  normVec.resize (numCols); // Resize normvec if necessary.
817  MVT::MvNorm (X, normVec);
818  }
819 
820  template<class Scalar, class MV>
821  void
825  {
826 #ifdef BELOS_TEUCHOS_TIME_MONITOR
827  // "Projection" only happens in rawProject(), so we only time
828  // projection inside rawProject(). However, we count the time
829  // spend in project() as part of the whole orthogonalization.
830  //
831  // If project() is called from projectAndNormalize(), the
832  // TimeMonitor won't start timerOrtho_, because it is already
833  // running in projectAndNormalize().
834  Teuchos::TimeMonitor timerMonitorOrtho(*timerOrtho_);
835 #endif // BELOS_TEUCHOS_TIME_MONITOR
836 
837  int ncols_X, num_Q_blocks, ncols_Q_total;
838  checkProjectionDims (ncols_X, num_Q_blocks, ncols_Q_total, X, Q);
839  // Test for quick exit: any dimension of X is zero, or there are
840  // zero Q blocks, or the total number of columns of the Q blocks
841  // is zero.
842  if (ncols_X == 0 || num_Q_blocks == 0 || ncols_Q_total == 0)
843  return;
844 
845  // Make space for first-pass projection coefficients
846  allocateProjectionCoefficients (C, Q, X, true);
847 
848  // We only use columnNormsBefore and compute pre-projection column
849  // norms if doing block reorthogonalization.
850  std::vector<magnitude_type> columnNormsBefore (ncols_X, magnitude_type(0));
851  if (reorthogonalizeBlocks_)
852  MVT::MvNorm (X, columnNormsBefore);
853 
854  // Project (first block orthogonalization step):
855  // C := Q^* X, X := X - Q C.
856  rawProject (X, Q, C);
857 
858  // If we are doing block reorthogonalization, reorthogonalize X if
859  // necessary.
860  if (reorthogonalizeBlocks_) {
861  std::vector<magnitude_type> columnNormsAfter (ncols_X, magnitude_type(0));
862  MVT::MvNorm (X, columnNormsAfter);
863 
864  // Relative block reorthogonalization threshold.
865  const magnitude_type relThres = blockReorthogThreshold();
866  // Reorthogonalize X if any of its column norms decreased by a
867  // factor more than the block reorthogonalization threshold.
868  // Don't bother trying to subset the columns; that will make the
869  // columns noncontiguous and thus hinder BLAS 3 optimizations.
870  bool reorthogonalize = false;
871  for (int j = 0; j < ncols_X; ++j) {
872  if (columnNormsAfter[j] < relThres * columnNormsBefore[j]) {
873  reorthogonalize = true;
874  break;
875  }
876  }
877  if (reorthogonalize) {
878  // Notify the caller via callback about the need for
879  // reorthogonalization.
880  if (! reorthogCallback_.is_null()) {
881  reorthogCallback_->operator() (Teuchos::arrayViewFromVector (columnNormsBefore),
882  Teuchos::arrayViewFromVector (columnNormsAfter));
883  }
884  // Second-pass projection coefficients
886  allocateProjectionCoefficients (C2, Q, X, false);
887 
888  // Perform the second projection pass:
889  // C2 = Q' X, X = X - Q*C2
890  rawProject (X, Q, C2);
891  // Update the projection coefficients
892  for (int k = 0; k < num_Q_blocks; ++k)
893  *C[k] += *C2[k];
894  }
895  }
896  }
897 
898 
899  template<class Scalar, class MV>
900  int
902  {
903  using Teuchos::Range1D;
904  using Teuchos::RCP;
905 
906 #ifdef BELOS_TEUCHOS_TIME_MONITOR
907  Teuchos::TimeMonitor timerMonitorNormalize(*timerNormalize_);
908  // If normalize() is called internally -- i.e., called from
909  // projectAndNormalize() -- the timer will not be started or
910  // stopped, because it is already running. TimeMonitor handles
911  // recursive invocation by doing nothing.
912  Teuchos::TimeMonitor timerMonitorOrtho(*timerOrtho_);
913 #endif // BELOS_TEUCHOS_TIME_MONITOR
914 
915  // MVT returns int for this, even though the "local ordinal
916  // type" of the MV may be some other type (for example,
917  // Tpetra::MultiVector<double, int32_t, int64_t, ...>).
918  const int numCols = MVT::GetNumberVecs (X);
919 
920  // This special case (for X having only one column) makes
921  // TsqrOrthoManagerImpl equivalent to Modified Gram-Schmidt
922  // orthogonalization with conditional full reorthogonalization,
923  // if all multivector inputs have only one column. It also
924  // avoids allocating Q_ scratch space and copying data when we
925  // don't need to invoke TSQR at all.
926  if (numCols == 1) {
927  return normalizeOne (X, B);
928  }
929 
930  // We use Q_ as scratch space for the normalization, since TSQR
931  // requires a scratch multivector (it can't factor in place). Q_
932  // should come from a vector space compatible with X's vector
933  // space, and Q_ should have at least as many columns as X.
934  // Otherwise, we have to reallocate. We also have to allocate
935  // (not "re-") Q_ if we haven't allocated it before. (We can't
936  // allocate Q_ until we have some X, so we need a multivector as
937  // the "prototype.")
938  //
939  // NOTE (mfh 11 Jan 2011) We only increase the number of columsn
940  // in Q_, never decrease. This is OK for typical uses of TSQR,
941  // but you might prefer different behavior in some cases.
942  //
943  // NOTE (mfh 10 Mar 2011) We should only reuse the scratch space
944  // Q_ if X and Q_ have compatible data distributions. However,
945  // Belos' current MultiVecTraits interface does not let us check
946  // for this. Thus, we can only check whether X and Q_ have the
947  // same number of rows. This will behave correctly for the common
948  // case in Belos that all multivectors with the same number of
949  // rows have the same data distribution.
950  //
951  // The specific MV implementation may do more checks than this on
952  // its own and throw an exception if X and Q_ are not compatible,
953  // but it may not. If you find that recycling the Q_ space causes
954  // troubles, you may consider modifying the code below to
955  // reallocate Q_ for every X that comes in.
956  if (Q_.is_null() ||
957  MVT::GetGlobalLength(*Q_) != MVT::GetGlobalLength(X) ||
958  numCols > MVT::GetNumberVecs (*Q_)) {
959  Q_ = MVT::Clone (X, numCols);
960  }
961 
962  // normalizeImpl() wants the second MV argument to have the same
963  // number of columns as X. To ensure this, we pass it a view of
964  // Q_ if Q_ has more columns than X. (This is possible if we
965  // previously called normalize() with a different multivector,
966  // since we never reallocate Q_ if it has more columns than
967  // necessary.)
968  if (MVT::GetNumberVecs(*Q_) == numCols) {
969  return normalizeImpl (X, *Q_, B, false);
970  } else {
971  RCP<MV> Q_view = MVT::CloneViewNonConst (*Q_, Range1D(0, numCols-1));
972  return normalizeImpl (X, *Q_view, B, false);
973  }
974  }
975 
976  template<class Scalar, class MV>
977  void
981  const MV& X,
982  const bool attemptToRecycle) const
983  {
984  const int num_Q_blocks = Q.size();
985  const int ncols_X = MVT::GetNumberVecs (X);
986  C.resize (num_Q_blocks);
987  if (attemptToRecycle)
988  {
989  for (int i = 0; i < num_Q_blocks; ++i)
990  {
991  const int ncols_Qi = MVT::GetNumberVecs (*Q[i]);
992  // Create a new C[i] if necessary, otherwise resize if
993  // necessary, otherwise fill with zeros.
994  if (C[i].is_null())
995  C[i] = Teuchos::rcp (new mat_type (ncols_Qi, ncols_X));
996  else
997  {
998  mat_type& Ci = *C[i];
999  if (Ci.numRows() != ncols_Qi || Ci.numCols() != ncols_X)
1000  Ci.shape (ncols_Qi, ncols_X);
1001  else
1002  Ci.putScalar (SCT::zero());
1003  }
1004  }
1005  }
1006  else
1007  {
1008  for (int i = 0; i < num_Q_blocks; ++i)
1009  {
1010  const int ncols_Qi = MVT::GetNumberVecs (*Q[i]);
1011  C[i] = Teuchos::rcp (new mat_type (ncols_Qi, ncols_X));
1012  }
1013  }
1014  }
1015 
1016  template<class Scalar, class MV>
1017  int
1020  {
1021 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1022  Teuchos::TimeMonitor timerMonitorOrtho(*timerOrtho_);
1023  Teuchos::TimeMonitor timerMonitorNormalize(*timerNormalize_);
1024 #endif // BELOS_TEUCHOS_TIME_MONITOR
1025 
1026  const int numVecs = MVT::GetNumberVecs(X);
1027  if (numVecs == 0) {
1028  return 0; // Nothing to do.
1029  } else if (numVecs == 1) {
1030  // Special case for a single column; scale and copy over.
1031  using Teuchos::Range1D;
1032  using Teuchos::RCP;
1033  using Teuchos::rcp;
1034 
1035  // Normalize X in place (faster than TSQR for one column).
1036  const int rank = normalizeOne (X, B);
1037  // Copy results to first column of Q.
1038  RCP<MV> Q_0 = MVT::CloneViewNonConst (Q, Range1D(0,0));
1039  MVT::Assign (X, *Q_0);
1040  return rank;
1041  } else {
1042  // The "true" argument to normalizeImpl() means the output
1043  // vectors go into Q, and the contents of X are overwritten with
1044  // invalid values.
1045  return normalizeImpl (X, Q, B, true);
1046  }
1047  }
1048 
1049  template<class Scalar, class MV>
1050  int
1053  MV& X_out, // Only written if outOfPlace==false.
1054  const bool outOfPlace,
1056  mat_ptr B,
1058  {
1059  using Teuchos::Range1D;
1060  using Teuchos::RCP;
1061  using Teuchos::rcp;
1062 
1063 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1064  // Projection is only timed in rawProject(), and normalization is
1065  // only timed in normalize() and normalizeOutOfPlace().
1066  Teuchos::TimeMonitor timerMonitorOrtho(*timerOrtho_);
1067 #endif // BELOS_TEUCHOS_TIME_MONITOR
1068 
1069  if (outOfPlace) {
1070  // Make sure that X_out has at least as many columns as X_in.
1071  TEUCHOS_TEST_FOR_EXCEPTION(MVT::GetNumberVecs(X_out) < MVT::GetNumberVecs(X_in),
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.");
1078  }
1079  // Fetch dimensions of X_in and Q, and allocate space for first-
1080  // and second-pass projection coefficients (C resp. C2).
1081  int ncols_X, num_Q_blocks, ncols_Q_total;
1082  checkProjectionDims (ncols_X, num_Q_blocks, ncols_Q_total, X_in, Q);
1083 
1084  // Test for quick exit: if any dimension of X is zero.
1085  if (ncols_X == 0) {
1086  return 0;
1087  }
1088  // If there are zero Q blocks or zero Q columns, just normalize!
1089  if (num_Q_blocks == 0 || ncols_Q_total == 0) {
1090  if (outOfPlace) {
1091  return normalizeOutOfPlace (X_in, X_out, B);
1092  } else {
1093  return normalize (X_in, B);
1094  }
1095  }
1096 
1097  // The typical case is that the entries of C have been allocated
1098  // before, so we attempt to recycle the allocations. The call
1099  // below will reallocate if it cannot recycle.
1100  allocateProjectionCoefficients (C, Q, X_in, true);
1101 
1102  // If we are doing block reorthogonalization, then compute the
1103  // column norms of X before projecting for the first time. This
1104  // will help us decide whether we need to reorthogonalize X.
1105  std::vector<magnitude_type> normsBeforeFirstPass (ncols_X, SCTM::zero());
1106  if (reorthogonalizeBlocks_) {
1107  MVT::MvNorm (X_in, normsBeforeFirstPass);
1108  }
1109 
1110  // First (Modified) Block Gram-Schmidt pass, in place in X_in.
1111  rawProject (X_in, Q, C);
1112 
1113  // Make space for the normalization coefficients. This will
1114  // either be a freshly allocated matrix (if B is null), or a view
1115  // of the appropriately sized upper left submatrix of *B (if B is
1116  // not null).
1117  //
1118  // Note that if we let the normalize() routine allocate (in the
1119  // case that B is null), that storage will go away at the end of
1120  // normalize(). (This is because it passes the RCP by value, not
1121  // by reference.)
1122  mat_ptr B_out;
1123  if (B.is_null()) {
1124  B_out = rcp (new mat_type (ncols_X, ncols_X));
1125  } else {
1126  // Make sure that B is no smaller than numCols x numCols.
1127  TEUCHOS_TEST_FOR_EXCEPTION(B->numRows() < ncols_X || B->numCols() < ncols_X,
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() << ".");
1133  // Create a view of the ncols_X by ncols_X upper left
1134  // submatrix of *B. TSQR will write the normalization
1135  // coefficients there.
1136  B_out = rcp (new mat_type (Teuchos::View, *B, ncols_X, ncols_X));
1137  }
1138 
1139  // Rank of X(_in) after first projection pass. If outOfPlace,
1140  // this overwrites X_in with invalid values, and the results go in
1141  // X_out. Otherwise, it's done in place in X_in.
1142  const int firstPassRank = outOfPlace ?
1143  normalizeOutOfPlace (X_in, X_out, B_out) :
1144  normalize (X_in, B_out);
1145  if (B.is_null()) {
1146  // The input matrix B is null, so assign B_out to it. If B was
1147  // not null on input, then B_out is a view of *B, so we don't
1148  // have to do anything here. Note that SerialDenseMatrix uses
1149  // raw pointers to store data and represent views, so we have to
1150  // be careful about scope.
1151  B = B_out;
1152  }
1153  int rank = firstPassRank; // Current rank of X.
1154 
1155  // If X was not full rank after projection and randomizeNullSpace_
1156  // is true, then normalize(OutOfPlace)() replaced the null space
1157  // basis of X with random vectors, and orthogonalized them against
1158  // the column space basis of X. However, we still need to
1159  // orthogonalize the random vectors against the Q[i], after which
1160  // we will need to renormalize them.
1161  //
1162  // If outOfPlace, then we need to work in X_out (where
1163  // normalizeOutOfPlace() wrote the normalized vectors).
1164  // Otherwise, we need to work in X_in.
1165  //
1166  // Note: We don't need to keep the new projection coefficients,
1167  // since they are multiplied by the "small" part of B
1168  // corresponding to the null space of the original X.
1169  if (firstPassRank < ncols_X && randomizeNullSpace_) {
1170  const int numNullSpaceCols = ncols_X - firstPassRank;
1171  const Range1D nullSpaceIndices (firstPassRank, ncols_X - 1);
1172 
1173  // Space for projection coefficients (will be thrown away)
1174  Teuchos::Array<mat_ptr> C_null (num_Q_blocks);
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));
1178  }
1179  // Space for normalization coefficients (will be thrown away).
1180  RCP<mat_type> B_null (new mat_type (numNullSpaceCols, numNullSpaceCols));
1181 
1182  int randomVectorsRank;
1183  if (outOfPlace) {
1184  // View of the null space basis columns of X.
1185  // normalizeOutOfPlace() wrote them into X_out.
1186  RCP<MV> X_out_null = MVT::CloneViewNonConst (X_out, nullSpaceIndices);
1187  // Use X_in for scratch space. Copy X_out_null into the
1188  // last few columns of X_in (X_in_null) and do projections
1189  // in there. (This saves us a copy wen we renormalize
1190  // (out of place) back into X_out.)
1191  RCP<MV> X_in_null = MVT::CloneViewNonConst (X_in, nullSpaceIndices);
1192  MVT::Assign (*X_out_null, *X_in_null);
1193  // Project the new random vectors against the Q blocks, and
1194  // renormalize the result into X_out_null.
1195  rawProject (*X_in_null, Q, C_null);
1196  randomVectorsRank = normalizeOutOfPlace (*X_in_null, *X_out_null, B_null);
1197  } else {
1198  // View of the null space columns of X.
1199  // They live in X_in.
1200  RCP<MV> X_null = MVT::CloneViewNonConst (X_in, nullSpaceIndices);
1201  // Project the new random vectors against the Q blocks,
1202  // and renormalize the result (in place).
1203  rawProject (*X_null, Q, C_null);
1204  randomVectorsRank = normalize (*X_null, B_null);
1205  }
1206  // While unusual, it is still possible for the random data not
1207  // to be full rank after projection and normalization. In that
1208  // case, we could try another set of random data and recurse as
1209  // necessary, but instead for now we just raise an exception.
1210  TEUCHOS_TEST_FOR_EXCEPTION(randomVectorsRank != numNullSpaceCols,
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 << ".");
1218  }
1219 
1220  // Whether or not X_in was full rank after projection, we still
1221  // might want to reorthogonalize against Q.
1222  if (reorthogonalizeBlocks_) {
1223  // We are only interested in the column space basis of X
1224  // resp. X_out.
1225  std::vector<magnitude_type>
1226  normsAfterFirstPass (firstPassRank, SCTM::zero());
1227  std::vector<magnitude_type>
1228  normsAfterSecondPass (firstPassRank, SCTM::zero());
1229 
1230  // Compute post-first-pass (pre-normalization) norms. We could
1231  // have done that using MVT::MvNorm() on X_in after projecting,
1232  // but before the first normalization. However, that operation
1233  // may be expensive. It is also unnecessary: after calling
1234  // normalize(), the 2-norm of B(:,j) is the 2-norm of X_in(:,j)
1235  // before normalization, in exact arithmetic.
1236  //
1237  // NOTE (mfh 06 Nov 2010) This is one way that combining
1238  // projection and normalization into a single kernel --
1239  // projectAndNormalize() -- pays off. In project(), we have to
1240  // compute column norms of X before and after projection. Here,
1241  // we get them for free from the normalization coefficients.
1243  for (int j = 0; j < firstPassRank; ++j) {
1244  const Scalar* const B_j = &(*B_out)(0,j);
1245  // Teuchos::BLAS::NRM2 returns a magnitude_type result on
1246  // Scalar inputs.
1247  normsAfterFirstPass[j] = blas.NRM2 (firstPassRank, B_j, 1);
1248  }
1249  // Test whether any of the norms dropped below the
1250  // reorthogonalization threshold.
1251  bool reorthogonalize = false;
1252  for (int j = 0; j < firstPassRank; ++j) {
1253  // If any column's norm decreased too much, mark this block
1254  // for reorthogonalization. Note that this test will _not_
1255  // activate reorthogonalization if a column's norm before the
1256  // first project-and-normalize step was zero. It _will_
1257  // activate reorthogonalization if the column's norm before
1258  // was not zero, but is zero now.
1259  const magnitude_type curThreshold =
1260  blockReorthogThreshold() * normsBeforeFirstPass[j];
1261  if (normsAfterFirstPass[j] < curThreshold) {
1262  reorthogonalize = true;
1263  break;
1264  }
1265  }
1266 
1267  // Notify the caller via callback about the need for
1268  // reorthogonalization.
1269  if (! reorthogCallback_.is_null()) {
1270  using Teuchos::arrayViewFromVector;
1271  (*reorthogCallback_) (arrayViewFromVector (normsBeforeFirstPass),
1272  arrayViewFromVector (normsAfterFirstPass));
1273  }
1274 
1275  // Perform another Block Gram-Schmidt pass if necessary. "Twice
1276  // is enough" (Kahan's theorem) for a Krylov method, unless
1277  // (using Stewart's term) there is an "orthogonalization fault"
1278  // (indicated by reorthogFault).
1279  //
1280  // NOTE (mfh 07 Nov 2010) For now, we include the entire block
1281  // of X, including possible random data (that was already
1282  // projected and normalized above). It might make more sense
1283  // just to process the first firstPassRank columns of X.
1284  // However, the resulting reorthogonalization should still be
1285  // correct regardless.
1286  bool reorthogFault = false;
1287  // Indices of X at which there was an orthogonalization fault.
1288  std::vector<int> faultIndices;
1289  if (reorthogonalize) {
1290  using Teuchos::Copy;
1291  using Teuchos::NO_TRANS;
1292 
1293  // If we're using out-of-place normalization, copy X_out
1294  // (results of first project and normalize pass) back into
1295  // X_in, for the second project and normalize pass.
1296  if (outOfPlace) {
1297  MVT::Assign (X_out, X_in);
1298  }
1299 
1300  // C2 is only used internally, so we know that we are
1301  // allocating fresh and not recycling allocations. Stating
1302  // this lets us save time checking dimensions.
1304  allocateProjectionCoefficients (C2, Q, X_in, false);
1305 
1306  // Block Gram-Schmidt (again). Delay updating the block
1307  // coefficients until we have the new normalization
1308  // coefficients, which we need in order to do the update.
1309  rawProject (X_in, Q, C2);
1310 
1311  // Coefficients for (re)normalization of X_in.
1312  RCP<mat_type> B2 (new mat_type (ncols_X, ncols_X));
1313 
1314  // Normalize X_in (into X_out, if working out of place).
1315  const int secondPassRank = outOfPlace ?
1316  normalizeOutOfPlace (X_in, X_out, B2) :
1317  normalize (X_in, B2);
1318  rank = secondPassRank; // Current rank of X
1319 
1320  // Update normalization coefficients. We begin with copying
1321  // B_out, since the BLAS' _GEMM routine doesn't let us alias
1322  // its input and output arguments.
1323  mat_type B_copy (Copy, *B_out, B_out->numRows(), B_out->numCols());
1324  // B_out := B2 * B_out (where input B_out is in B_copy).
1325  const int err = B_out->multiply (NO_TRANS, NO_TRANS, SCT::one(),
1326  *B2, B_copy, SCT::zero());
1327  TEUCHOS_TEST_FOR_EXCEPTION(err != 0, std::logic_error,
1328  "Teuchos::SerialDenseMatrix::multiply "
1329  "returned err = " << err << " != 0");
1330  // Update the block coefficients from the projection step. We
1331  // use B_copy for this (a copy of B_out, the first-pass
1332  // normalization coefficients).
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());
1335 
1336  // C[k] := C2[k]*B_copy + C[k].
1337  const int err1 = C[k]->multiply (NO_TRANS, NO_TRANS, SCT::one(),
1338  *C2[k], B_copy, SCT::one());
1339  TEUCHOS_TEST_FOR_EXCEPTION(err1 != 0, std::logic_error,
1340  "Teuchos::SerialDenseMatrix::multiply "
1341  "returned err = " << err1 << " != 0");
1342  }
1343  // Compute post-second-pass (pre-normalization) norms, using
1344  // B2 (the coefficients from the second normalization step) in
1345  // the same way as with B_out before.
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);
1349  }
1350  // Test whether any of the norms dropped below the
1351  // reorthogonalization threshold. If so, it's an
1352  // orthogonalization fault, which requires expensive recovery.
1353  reorthogFault = false;
1354  for (int j = 0; j < rank; ++j) {
1355  const magnitude_type relativeLowerBound =
1356  blockReorthogThreshold() * normsAfterFirstPass[j];
1357  if (normsAfterSecondPass[j] < relativeLowerBound) {
1358  reorthogFault = true;
1359  faultIndices.push_back (j);
1360  }
1361  }
1362  } // if (reorthogonalize) // reorthogonalization pass
1363 
1364  if (reorthogFault) {
1365  if (throwOnReorthogFault_) {
1366  raiseReorthogFault (normsAfterFirstPass,
1367  normsAfterSecondPass,
1368  faultIndices);
1369  } else {
1370  // NOTE (mfh 19 Jan 2011) We could handle the fault here by
1371  // slowly reorthogonalizing, one vector at a time, the
1372  // offending vectors of X. However, we choose not to
1373  // implement this for now. If it becomes a problem, let us
1374  // know and we will prioritize implementing this.
1375  TEUCHOS_TEST_FOR_EXCEPTION(true, std::logic_error,
1376  "TsqrOrthoManagerImpl has not yet implemented"
1377  " recovery from an orthogonalization fault.");
1378  }
1379  }
1380  } // if (reorthogonalizeBlocks_)
1381  return rank;
1382  }
1383 
1384 
1385  template<class Scalar, class MV>
1386  void
1388  raiseReorthogFault (const std::vector<magnitude_type>& normsAfterFirstPass,
1389  const std::vector<magnitude_type>& normsAfterSecondPass,
1390  const std::vector<int>& faultIndices)
1391  {
1392  using std::endl;
1393  typedef std::vector<int>::size_type size_type;
1394  std::ostringstream os;
1395 
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];
1400  const magnitude_type decreaseFactor =
1401  normsAfterSecondPass[index] / normsAfterFirstPass[index];
1402  os << index << "\t" << decreaseFactor << endl;
1403  }
1404  throw TsqrOrthoFault (os.str());
1405  }
1406 
1407  template<class Scalar, class MV>
1410  {
1411  using Teuchos::ParameterList;
1412  using Teuchos::parameterList;
1413  using Teuchos::RCP;
1414 
1415  if (defaultParams_.is_null()) {
1416  RCP<ParameterList> params = parameterList ("TsqrOrthoManagerImpl");
1417  //
1418  // TSQR parameters (set as a sublist).
1419  //
1420  params->set ("TSQR implementation", *(tsqrAdaptor_.getValidParameters()),
1421  "TSQR implementation parameters.");
1422  //
1423  // Orthogonalization parameters
1424  //
1425  const bool defaultRandomizeNullSpace = true;
1426  params->set ("randomizeNullSpace", defaultRandomizeNullSpace,
1427  "Whether to fill in null space vectors with random data.");
1428 
1429  const bool defaultReorthogonalizeBlocks = true;
1430  params->set ("reorthogonalizeBlocks", defaultReorthogonalizeBlocks,
1431  "Whether to do block reorthogonalization as necessary.");
1432 
1433  // This parameter corresponds to the "blk_tol_" parameter in
1434  // Belos' DGKSOrthoManager. We choose the same default value.
1435  const magnitude_type defaultBlockReorthogThreshold =
1436  magnitude_type(10) * SCTM::squareroot (SCTM::eps());
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.");
1441 
1442  // This parameter corresponds to the "sing_tol_" parameter in
1443  // Belos' DGKSOrthoManager. We choose the same default value.
1444  const magnitude_type defaultRelativeRankTolerance =
1445  Teuchos::as<magnitude_type>(10) * SCTM::eps();
1446 
1447  // If the relative rank tolerance is zero, then we will always
1448  // declare blocks to be numerically full rank, as long as no
1449  // singular values are zero.
1450  params->set ("relativeRankTolerance", defaultRelativeRankTolerance,
1451  "Relative tolerance to determine the numerical rank of a "
1452  "block when normalizing.");
1453 
1454  // See Stewart's 2008 paper on block Gram-Schmidt for a definition
1455  // of "orthogonalization fault."
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).");
1461 
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.");
1466 
1467  defaultParams_ = params;
1468  }
1469  return defaultParams_;
1470  }
1471 
1472  template<class Scalar, class MV>
1475  {
1476  using Teuchos::ParameterList;
1477  using Teuchos::RCP;
1478  using Teuchos::rcp;
1479 
1480  RCP<const ParameterList> defaultParams = getValidParameters();
1481  // Start with a clone of the default parameters.
1482  RCP<ParameterList> params = rcp (new ParameterList (*defaultParams));
1483 
1484  // Disable reorthogonalization and randomization of the null
1485  // space basis. Reorthogonalization tolerances don't matter,
1486  // since we aren't reorthogonalizing blocks in the fast
1487  // settings. We can leave the default values. Also,
1488  // (re)orthogonalization faults may only occur with
1489  // reorthogonalization, so we don't have to worry about the
1490  // "throwOnReorthogFault" setting.
1491  const bool randomizeNullSpace = false;
1492  params->set ("randomizeNullSpace", randomizeNullSpace);
1493  const bool reorthogonalizeBlocks = false;
1494  params->set ("reorthogonalizeBlocks", reorthogonalizeBlocks);
1495 
1496  return params;
1497  }
1498 
1499  template<class Scalar, class MV>
1500  int
1503  MV& Q,
1505  {
1506  int rank;
1507  try {
1508  // This call only computes the QR factorization X = Q B.
1509  // It doesn't compute the rank of X. That comes from
1510  // revealRank() below.
1511  tsqrAdaptor_.factorExplicit (X, Q, B, forceNonnegativeDiagonal_);
1512  // This call will only modify *B if *B on input is not of full
1513  // numerical rank.
1514  rank = tsqrAdaptor_.revealRank (Q, B, relativeRankTolerance_);
1515  } catch (std::exception& e) {
1516  throw TsqrOrthoError (e.what()); // Toss the exception up the chain.
1517  }
1518  return rank;
1519  }
1520 
1521  template<class Scalar, class MV>
1522  int
1526  {
1527  // Make space for the normalization coefficient. This will either
1528  // be a freshly allocated matrix (if B is null), or a view of the
1529  // 1x1 upper left submatrix of *B (if B is not null).
1530  mat_ptr B_out;
1531  if (B.is_null()) {
1532  B_out = Teuchos::rcp (new mat_type (1, 1));
1533  } else {
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 << ".");
1540  // Create a view of the 1x1 upper left submatrix of *B.
1541  B_out = Teuchos::rcp (new mat_type (Teuchos::View, *B, 1, 1));
1542  }
1543 
1544  // Compute the norm of X, and write the result to B_out.
1545  std::vector<magnitude_type> theNorm (1, SCTM::zero());
1546  MVT::MvNorm (X, theNorm);
1547  (*B_out)(0,0) = theNorm[0];
1548 
1549  if (B.is_null()) {
1550  // The input matrix B is null, so assign B_out to it. If B was
1551  // not null on input, then B_out is a view of *B, so we don't
1552  // have to do anything here. Note that SerialDenseMatrix uses
1553  // raw pointers to store data and represent views, so we have to
1554  // be careful about scope.
1555  B = B_out;
1556  }
1557 
1558  // Scale X by its norm, if its norm is zero. Otherwise, do the
1559  // right thing based on whether the user wants us to fill the null
1560  // space with random vectors.
1561  if (theNorm[0] == SCTM::zero()) {
1562  // Make a view of the first column of Q, fill it with random
1563  // data, and normalize it. Throw away the resulting norm.
1564  if (randomizeNullSpace_) {
1565  MVT::MvRandom(X);
1566  MVT::MvNorm (X, theNorm);
1567  if (theNorm[0] == SCTM::zero()) {
1568  // It is possible that a random vector could have all zero
1569  // entries, but unlikely. We could try again, but it's also
1570  // possible that multiple tries could result in zero
1571  // vectors. We choose instead to give up.
1572  throw TsqrOrthoError("normalizeOne: a supposedly random "
1573  "vector has norm zero!");
1574  } else {
1575  // NOTE (mfh 09 Nov 2010) I'm assuming that dividing a
1576  // Scalar by a magnitude_type is defined and that it results
1577  // in a Scalar.
1578  const Scalar alpha = SCT::one() / theNorm[0];
1579  MVT::MvScale (X, alpha);
1580  }
1581  }
1582  return 0; // The rank of the matrix (actually just one vector) X.
1583  } else {
1584  // NOTE (mfh 09 Nov 2010) I'm assuming that dividing a Scalar by
1585  // a magnitude_type is defined and that it results in a Scalar.
1586  const Scalar alpha = SCT::one() / theNorm[0];
1587  MVT::MvScale (X, alpha);
1588  return 1; // The rank of the matrix (actually just one vector) X.
1589  }
1590  }
1591 
1592 
1593  template<class Scalar, class MV>
1594  void
1596  rawProject (MV& X,
1599  {
1600 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1601  Teuchos::TimeMonitor timerMonitorNormalize(*timerProject_);
1602 #endif // BELOS_TEUCHOS_TIME_MONITOR
1603 
1604  // "Modified Gram-Schmidt" version of Block Gram-Schmidt.
1605  const int num_Q_blocks = Q.size();
1606  for (int i = 0; i < num_Q_blocks; ++i)
1607  {
1608  // TEUCHOS_TEST_FOR_EXCEPTION(C[i].is_null(), std::logic_error,
1609  // "TsqrOrthoManagerImpl::rawProject(): C["
1610  // << i << "] is null");
1611  // TEUCHOS_TEST_FOR_EXCEPTION(Q[i].is_null(), std::logic_error,
1612  // "TsqrOrthoManagerImpl::rawProject(): Q["
1613  // << i << "] is null");
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);
1618  }
1619  }
1620 
1621 
1622  template<class Scalar, class MV>
1623  void
1625  rawProject (MV& X,
1626  const Teuchos::RCP<const MV>& Q,
1628  {
1629 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1630  Teuchos::TimeMonitor timerMonitorNormalize(*timerProject_);
1631 #endif // BELOS_TEUCHOS_TIME_MONITOR
1632 
1633  // Block Gram-Schmidt
1634  innerProd (*Q, X, *C);
1635  MVT::MvTimesMatAddMv (-SCT::one(), *Q, *C, SCT::one(), X);
1636  }
1637 
1638  template<class Scalar, class MV>
1639  int
1642  MV& Q,
1644  const bool outOfPlace)
1645  {
1646  using Teuchos::Range1D;
1647  using Teuchos::RCP;
1648  using Teuchos::rcp;
1649  using Teuchos::ScalarTraits;
1650  using Teuchos::tuple;
1651 
1652  const int numCols = MVT::GetNumberVecs (X);
1653  if (numCols == 0) {
1654  return 0; // Fast exit for an empty input matrix.
1655  }
1656 
1657  // We allow Q to have more columns than X. In that case, we only
1658  // touch the first numCols columns of Q.
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.");
1664  // TSQR wants a Q with the same number of columns as X, so have it
1665  // work on a nonconstant view of Q with the same number of columns
1666  // as X.
1667  RCP<MV> Q_view = MVT::CloneViewNonConst (Q, Range1D (0, numCols-1));
1668 
1669  // Make space for the normalization coefficients. This will
1670  // either be a freshly allocated matrix (if B is null), or a view
1671  // of the appropriately sized upper left submatrix of *B (if B is
1672  // not null).
1673  mat_ptr B_out;
1674  if (B.is_null ()) {
1675  B_out = rcp (new mat_type (numCols, numCols));
1676  } else {
1677  // Make sure that B is no smaller than numCols x numCols.
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() << ".");
1683  // Create a view of the numCols x numCols upper left submatrix
1684  // of *B. TSQR will write the normalization coefficients there.
1685  B_out = rcp (new mat_type (Teuchos::View, *B, numCols, numCols));
1686  }
1687 
1688  // Compute rank-revealing decomposition (in this case, TSQR of X
1689  // followed by SVD of the R factor and appropriate updating of the
1690  // resulting Q factor) of X. X is modified in place and filled
1691  // with garbage, and Q_view contains the resulting explicit Q
1692  // factor. Later, we will copy this back into X.
1693  //
1694  // The matrix *B_out will only be upper triangular if X is of full
1695  // numerical rank. Otherwise, the entries below the diagonal may
1696  // be filled in as well.
1697  const int rank = rawNormalize (X, *Q_view, *B_out);
1698  if (B.is_null ()) {
1699  // The input matrix B is null, so assign B_out to it. If B was
1700  // not null on input, then B_out is a view of *B, so we don't
1701  // have to do anything here. Note that SerialDenseMatrix uses
1702  // raw pointers to store data and represent views, so we have to
1703  // be careful about scope.
1704  B = B_out;
1705  }
1706 
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.");
1712 
1713  // If X is full rank or we don't want to replace its null space
1714  // basis with random vectors, then we're done.
1715  if (rank == numCols || ! randomizeNullSpace_) {
1716  // If we're supposed to be working in place in X, copy the
1717  // results back from Q_view into X.
1718  if (! outOfPlace) {
1719  MVT::Assign (*Q_view, X);
1720  }
1721  return rank;
1722  }
1723 
1724  if (randomizeNullSpace_ && rank < numCols) {
1725  // X wasn't full rank. Replace the null space basis of X (in
1726  // the last numCols-rank columns of Q_view) with random data,
1727  // project it against the first rank columns of Q_view, and
1728  // normalize.
1729  //
1730  // Number of columns to fill with random data.
1731  const int nullSpaceNumCols = numCols - rank;
1732  // Inclusive range of indices of columns of X to fill with
1733  // random data.
1734  Range1D nullSpaceIndices (rank, numCols-1);
1735 
1736  // rawNormalize wrote the null space basis vectors into Q_view.
1737  // We continue to work in place in Q_view by writing the random
1738  // data there and (if there is a nontrival column space)
1739  // projecting in place against the column space basis vectors
1740  // (also in Q_view).
1741  RCP<MV> Q_null = MVT::CloneViewNonConst (*Q_view, nullSpaceIndices);
1742  // Replace the null space basis with random data.
1743  MVT::MvRandom (*Q_null);
1744 
1745  // Make sure that the "random" data isn't all zeros. This is
1746  // statistically nearly impossible, but we test for debugging
1747  // purposes.
1748  {
1749  std::vector<magnitude_type> norms (MVT::GetNumberVecs (*Q_null));
1750  MVT::MvNorm (*Q_null, norms);
1751 
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()) {
1756  anyZero = true;
1757  }
1758  }
1759  if (anyZero) {
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 "
1769  "vectors: [";
1770  for (iter_type it = norms.begin(); it != norms.end(); ++it) {
1771  os << *it;
1772  if (it+1 != norms.end())
1773  os << ", ";
1774  }
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 "
1778  "behavior.";
1779  TEUCHOS_TEST_FOR_EXCEPTION(anyZero, TsqrOrthoError, os.str());
1780  }
1781  }
1782 
1783  if (rank > 0) {
1784  // Project the random data against the column space basis of
1785  // X, using a simple block projection ("Block Classical
1786  // Gram-Schmidt"). This is accurate because we've already
1787  // orthogonalized the column space basis of X nearly to
1788  // machine precision via a QR factorization (TSQR) with
1789  // accuracy comparable to Householder QR.
1790  RCP<const MV> Q_col = MVT::CloneView (*Q_view, Range1D (0, rank-1));
1791 
1792  // Temporary storage for projection coefficients. We don't
1793  // need to keep them, since they represent the null space
1794  // basis (for which the coefficients are logically zero).
1795  mat_ptr C_null (new mat_type (rank, nullSpaceNumCols));
1796  rawProject (*Q_null, Q_col, C_null);
1797  }
1798  // Normalize the projected random vectors, so that they are
1799  // mutually orthogonal (as well as orthogonal to the column
1800  // space basis of X). We use X for the output of the
1801  // normalization: for out-of-place normalization (outOfPlace ==
1802  // true), X is overwritten with "invalid values" anyway, and for
1803  // in-place normalization (outOfPlace == false), we want the
1804  // result to be in X anyway.
1805  RCP<MV> X_null = MVT::CloneViewNonConst (X, nullSpaceIndices);
1806  // Normalization coefficients for projected random vectors.
1807  // Will be thrown away.
1808  mat_type B_null (nullSpaceNumCols, nullSpaceNumCols);
1809  // Write the normalized vectors to X_null (in X).
1810  const int nullSpaceBasisRank = rawNormalize (*Q_null, *X_null, B_null);
1811 
1812  // It's possible, but unlikely, that X_null doesn't have full
1813  // rank (after the projection step). We could recursively fill
1814  // in more random vectors until we finally get a full rank
1815  // matrix, but instead we just throw an exception.
1816  //
1817  // NOTE (mfh 08 Nov 2010) Perhaps we should deal with this case
1818  // more elegantly. Recursion might be one way to solve it, but
1819  // be sure to check that the recursion will terminate. We could
1820  // do this by supplying an additional argument to rawNormalize,
1821  // which is the null space basis rank from the previous
1822  // iteration. The rank has to decrease each time, or the
1823  // recursion may go on forever.
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) {
1841  os << norms[k];
1842  if (k != norms.size()-1) {
1843  os << ", ";
1844  }
1845  }
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.";
1850 
1851  TEUCHOS_TEST_FOR_EXCEPTION(nullSpaceBasisRank < nullSpaceNumCols,
1852  TsqrOrthoError, os.str ());
1853  }
1854  // If we're normalizing out of place, copy the X_null
1855  // vectors back into Q_null; the Q_col vectors are already
1856  // where they are supposed to be in that case.
1857  //
1858  // If we're normalizing in place, leave X_null alone (it's
1859  // already where it needs to be, in X), but copy Q_col back
1860  // into the first rank columns of X.
1861  if (outOfPlace) {
1862  MVT::Assign (*X_null, *Q_null);
1863  } else if (rank > 0) {
1864  // MVT::Assign() doesn't accept empty ranges of columns.
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);
1868  }
1869  }
1870  return rank;
1871  }
1872 
1873 
1874  template<class Scalar, class MV>
1875  void
1877  checkProjectionDims (int& ncols_X,
1878  int& num_Q_blocks,
1879  int& ncols_Q_total,
1880  const MV& X,
1882  {
1883  // First assign to temporary values, so the function won't
1884  // commit any externally visible side effects unless it will
1885  // return normally (without throwing an exception). (I'm being
1886  // cautious; MVT::GetNumberVecs() probably shouldn't have any
1887  // externally visible side effects, unless it is logging to a
1888  // file or something.)
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);
1892 
1893  // Compute the total number of columns of all the Q[i] blocks.
1894  the_ncols_Q_total = 0;
1895  // You should be angry if your compiler doesn't support type
1896  // inference ("auto"). That's why I need this awful typedef.
1897  using Teuchos::ArrayView;
1898  using Teuchos::RCP;
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);
1903  }
1904 
1905  // Commit temporary values to the output arguments.
1906  ncols_X = the_ncols_X;
1907  num_Q_blocks = the_num_Q_blocks;
1908  ncols_Q_total = the_ncols_Q_total;
1909  }
1910 
1911 } // namespace Belos
1912 
1913 #endif // __BelosTsqrOrthoManagerImpl_hpp
1914 
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 > &params, 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.
static void clearCounter(const std::string &name)
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)
static RCP< Time > getNewCounter(const std::string &name)
void setParameterList(const Teuchos::RCP< Teuchos::ParameterList > &params)
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 &quot;fast&quot; 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.
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 &quot;multivector&quot; 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...
Orthogonalization fault.
bool randomizeNullSpace_
Whether to fill null space vectors with random data.
OrdinalType numRows() const
MultiVecTraits< Scalar, MV > MVT