$darkmode
Eigen  5.0.1-dev
SparseQR.h
1 // This file is part of Eigen, a lightweight C++ template library
2 // for linear algebra.
3 //
4 // Copyright (C) 2012-2013 Desire Nuentsa <desire.nuentsa_wakam@inria.fr>
5 // Copyright (C) 2012-2014 Gael Guennebaud <gael.guennebaud@inria.fr>
6 //
7 // This Source Code Form is subject to the terms of the Mozilla
8 // Public License v. 2.0. If a copy of the MPL was not distributed
9 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
10 
11 #ifndef EIGEN_SPARSE_QR_H
12 #define EIGEN_SPARSE_QR_H
13 
14 // IWYU pragma: private
15 #include "./InternalHeaderCheck.h"
16 
17 namespace Eigen {
18 
19 template <typename MatrixType, typename OrderingType>
20 class SparseQR;
21 template <typename SparseQRType>
22 struct SparseQRMatrixQReturnType;
23 template <typename SparseQRType>
24 struct SparseQRMatrixQTransposeReturnType;
25 template <typename SparseQRType, typename Derived>
26 struct SparseQR_QProduct;
27 namespace internal {
28 template <typename SparseQRType>
29 struct traits<SparseQRMatrixQReturnType<SparseQRType> > {
30  typedef typename SparseQRType::MatrixType ReturnType;
31  typedef typename ReturnType::StorageIndex StorageIndex;
32  typedef typename ReturnType::StorageKind StorageKind;
33  enum { RowsAtCompileTime = Dynamic, ColsAtCompileTime = Dynamic };
34 };
35 template <typename SparseQRType>
36 struct traits<SparseQRMatrixQTransposeReturnType<SparseQRType> > {
37  typedef typename SparseQRType::MatrixType ReturnType;
38 };
39 template <typename SparseQRType, typename Derived>
40 struct traits<SparseQR_QProduct<SparseQRType, Derived> > {
41  typedef typename Derived::PlainObject ReturnType;
42 };
43 } // End namespace internal
44 
87 template <typename MatrixType_, typename OrderingType_>
88 class SparseQR : public SparseSolverBase<SparseQR<MatrixType_, OrderingType_> > {
89  protected:
90  typedef SparseSolverBase<SparseQR<MatrixType_, OrderingType_> > Base;
91  using Base::m_isInitialized;
92 
93  public:
94  using Base::_solve_impl;
95  typedef MatrixType_ MatrixType;
96  typedef OrderingType_ OrderingType;
97  typedef typename MatrixType::Scalar Scalar;
98  typedef typename MatrixType::RealScalar RealScalar;
99  typedef typename MatrixType::StorageIndex StorageIndex;
100  typedef SparseMatrix<Scalar, ColMajor, StorageIndex> QRMatrixType;
101  typedef Matrix<StorageIndex, Dynamic, 1> IndexVector;
102  typedef Matrix<Scalar, Dynamic, 1> ScalarVector;
103  typedef PermutationMatrix<Dynamic, Dynamic, StorageIndex> PermutationType;
104 
105  enum { ColsAtCompileTime = MatrixType::ColsAtCompileTime, MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime };
106 
107  public:
108  SparseQR()
109  : m_analysisIsok(false), m_lastError(""), m_useDefaultThreshold(true), m_isQSorted(false), m_isEtreeOk(false) {}
110 
117  explicit SparseQR(const MatrixType& mat)
118  : m_analysisIsok(false), m_lastError(""), m_useDefaultThreshold(true), m_isQSorted(false), m_isEtreeOk(false) {
119  compute(mat);
120  }
121 
128  void compute(const MatrixType& mat) {
129  analyzePattern(mat);
130  factorize(mat);
131  }
132  void analyzePattern(const MatrixType& mat);
133  void factorize(const MatrixType& mat);
134 
137  inline Index rows() const { return m_pmat.rows(); }
138 
141  inline Index cols() const { return m_pmat.cols(); }
142 
156  const QRMatrixType& matrixR() const { return m_R; }
157 
162  Index rank() const {
163  eigen_assert(m_isInitialized && "The factorization should be called first, use compute()");
164  return m_nonzeropivots;
165  }
166 
185  SparseQRMatrixQReturnType<SparseQR> matrixQ() const { return SparseQRMatrixQReturnType<SparseQR>(*this); }
186 
191  eigen_assert(m_isInitialized && "Decomposition is not initialized.");
192  return m_outputPerm_c;
193  }
194 
198  std::string lastErrorMessage() const { return m_lastError; }
199 
201  template <typename Rhs, typename Dest>
202  bool _solve_impl(const MatrixBase<Rhs>& B, MatrixBase<Dest>& dest) const {
203  eigen_assert(m_isInitialized && "The factorization should be called first, use compute()");
204  eigen_assert(this->rows() == B.rows() &&
205  "SparseQR::solve() : invalid number of rows in the right hand side matrix");
206 
207  Index rank = this->rank();
208 
209  // Compute Q^* * b;
210  typename Dest::PlainObject y, b;
211  y = this->matrixQ().adjoint() * B;
212  b = y;
213 
214  // Solve with the triangular matrix R
215  y.resize((std::max<Index>)(cols(), y.rows()), y.cols());
216  y.topRows(rank) = this->matrixR().topLeftCorner(rank, rank).template triangularView<Upper>().solve(b.topRows(rank));
217  y.bottomRows(y.rows() - rank).setZero();
218 
219  // Apply the column permutation
220  if (m_perm_c.size())
221  dest = colsPermutation() * y.topRows(cols());
222  else
223  dest = y.topRows(cols());
224 
225  m_info = Success;
226  return true;
227  }
228 
234  void setPivotThreshold(const RealScalar& threshold) {
235  m_useDefaultThreshold = false;
236  m_threshold = threshold;
237  }
238 
243  template <typename Rhs>
244  inline const Solve<SparseQR, Rhs> solve(const MatrixBase<Rhs>& B) const {
245  eigen_assert(m_isInitialized && "The factorization should be called first, use compute()");
246  eigen_assert(this->rows() == B.rows() &&
247  "SparseQR::solve() : invalid number of rows in the right hand side matrix");
248  return Solve<SparseQR, Rhs>(*this, B.derived());
249  }
250  template <typename Rhs>
251  inline const Solve<SparseQR, Rhs> solve(const SparseMatrixBase<Rhs>& B) const {
252  eigen_assert(m_isInitialized && "The factorization should be called first, use compute()");
253  eigen_assert(this->rows() == B.rows() &&
254  "SparseQR::solve() : invalid number of rows in the right hand side matrix");
255  return Solve<SparseQR, Rhs>(*this, B.derived());
256  }
257 
267  eigen_assert(m_isInitialized && "Decomposition is not initialized.");
268  return m_info;
269  }
270 
272  inline void _sort_matrix_Q() {
273  if (this->m_isQSorted) return;
274  // The matrix Q is sorted during the transposition
276  this->m_Q = mQrm;
277  this->m_isQSorted = true;
278  }
279 
280  protected:
281  bool m_analysisIsok;
282  bool m_factorizationIsok;
283  mutable ComputationInfo m_info;
284  std::string m_lastError;
285  QRMatrixType m_pmat; // Temporary matrix
286  QRMatrixType m_R; // The triangular factor matrix
287  QRMatrixType m_Q; // The orthogonal reflectors
288  ScalarVector m_hcoeffs; // The Householder coefficients
289  PermutationType m_perm_c; // Fill-reducing Column permutation
290  PermutationType m_pivotperm; // The permutation for rank revealing
291  PermutationType m_outputPerm_c; // The final column permutation
292  RealScalar m_threshold; // Threshold to determine null Householder reflections
293  bool m_useDefaultThreshold; // Use default threshold
294  Index m_nonzeropivots; // Number of non zero pivots found
295  IndexVector m_etree; // Column elimination tree
296  IndexVector m_firstRowElt; // First element in each row
297  bool m_isQSorted; // whether Q is sorted or not
298  bool m_isEtreeOk; // whether the elimination tree match the initial input matrix
299 
300  template <typename, typename>
301  friend struct SparseQR_QProduct;
302 };
303 
313 template <typename MatrixType, typename OrderingType>
315  eigen_assert(
316  mat.isCompressed() &&
317  "SparseQR requires a sparse matrix in compressed mode. Call .makeCompressed() before passing it to SparseQR");
318  // Copy to a column major matrix if the input is rowmajor
319  std::conditional_t<MatrixType::IsRowMajor, QRMatrixType, const MatrixType&> matCpy(mat);
320  // Compute the column fill reducing ordering
321  OrderingType ord;
322  ord(matCpy, m_perm_c);
323  Index n = mat.cols();
324  Index m = mat.rows();
325  Index diagSize = (std::min)(m, n);
326 
327  if (!m_perm_c.size()) {
328  m_perm_c.resize(n);
329  m_perm_c.indices().setLinSpaced(n, 0, StorageIndex(n - 1));
330  }
331 
332  // Compute the column elimination tree of the permuted matrix
333  m_outputPerm_c = m_perm_c.inverse();
334  internal::coletree(matCpy, m_etree, m_firstRowElt, m_outputPerm_c.indices().data());
335  m_isEtreeOk = true;
336 
337  m_R.resize(m, n);
338  m_Q.resize(m, diagSize);
339 
340  // Allocate space for nonzero elements: rough estimation
341  m_R.reserve(2 * mat.nonZeros()); // FIXME Get a more accurate estimation through symbolic factorization with the
342  // etree
343  m_Q.reserve(2 * mat.nonZeros());
344  m_hcoeffs.resize(diagSize);
345  m_analysisIsok = true;
346 }
347 
355 template <typename MatrixType, typename OrderingType>
357  using std::abs;
358 
359  eigen_assert(m_analysisIsok && "analyzePattern() should be called before this step");
360  StorageIndex m = StorageIndex(mat.rows());
361  StorageIndex n = StorageIndex(mat.cols());
362  StorageIndex diagSize = (std::min)(m, n);
363  IndexVector mark((std::max)(m, n));
364  mark.setConstant(-1); // Record the visited nodes
365  IndexVector Ridx(n), Qidx(m); // Store temporarily the row indexes for the current column of R and Q
366  Index nzcolR, nzcolQ; // Number of nonzero for the current column of R and Q
367  ScalarVector tval(m); // The dense vector used to compute the current column
368 
369  m_R.setZero();
370  m_Q.setZero();
371  m_pmat = mat;
372  if (!m_isEtreeOk) {
373  m_outputPerm_c = m_perm_c.inverse();
374  internal::coletree(m_pmat, m_etree, m_firstRowElt, m_outputPerm_c.indices().data());
375  m_isEtreeOk = true;
376  }
377 
378  m_pmat.uncompress(); // To have the innerNonZeroPtr allocated
379 
380  // Apply the fill-in reducing permutation lazily:
381  {
382  // If the input is row major, copy the original column indices,
383  // otherwise directly use the input matrix
384  //
385  IndexVector originalOuterIndicesCpy;
386  const StorageIndex* originalOuterIndices = mat.outerIndexPtr();
387  if (MatrixType::IsRowMajor) {
388  originalOuterIndicesCpy = IndexVector::Map(m_pmat.outerIndexPtr(), n + 1);
389  originalOuterIndices = originalOuterIndicesCpy.data();
390  }
391 
392  for (int i = 0; i < n; i++) {
393  Index p = m_perm_c.size() ? m_perm_c.indices()(i) : i;
394  m_pmat.outerIndexPtr()[p] = originalOuterIndices[i];
395  m_pmat.innerNonZeroPtr()[p] = originalOuterIndices[i + 1] - originalOuterIndices[i];
396  }
397  }
398 
399  /* Compute the default threshold as in MatLab, see:
400  * Tim Davis, "Algorithm 915, SuiteSparseQR: Multifrontal Multithreaded Rank-Revealing
401  * Sparse QR Factorization, ACM Trans. on Math. Soft. 38(1), 2011, Page 8:3
402  */
403  RealScalar pivotThreshold;
404  if (m_useDefaultThreshold) {
405  RealScalar max2Norm = 0.0;
406  for (int j = 0; j < n; j++) max2Norm = numext::maxi(max2Norm, m_pmat.col(j).norm());
407  if (max2Norm == RealScalar(0)) max2Norm = RealScalar(1);
408  pivotThreshold = 20 * (m + n) * max2Norm * NumTraits<RealScalar>::epsilon();
409  } else {
410  pivotThreshold = m_threshold;
411  }
412 
413  // Initialize the numerical permutation
414  m_pivotperm.setIdentity(n);
415 
416  StorageIndex nonzeroCol = 0; // Record the number of valid pivots
417  m_Q.startVec(0);
418 
419  // Left looking rank-revealing QR factorization: compute a column of R and Q at a time
420  for (StorageIndex col = 0; col < n; ++col) {
421  mark.setConstant(-1);
422  m_R.startVec(col);
423  mark(nonzeroCol) = col;
424  Qidx(0) = nonzeroCol;
425  nzcolR = 0;
426  nzcolQ = 1;
427  bool found_diag = nonzeroCol >= m;
428  tval.setZero();
429 
430  // Symbolic factorization: find the nonzero locations of the column k of the factors R and Q, i.e.,
431  // all the nodes (with indexes lower than rank) reachable through the column elimination tree (etree) rooted at node
432  // k. Note: if the diagonal entry does not exist, then its contribution must be explicitly added, thus the trick
433  // with found_diag that permits to do one more iteration on the diagonal element if this one has not been found.
434  for (typename QRMatrixType::InnerIterator itp(m_pmat, col); itp || !found_diag; ++itp) {
435  StorageIndex curIdx = nonzeroCol;
436  if (itp) curIdx = StorageIndex(itp.row());
437  if (curIdx == nonzeroCol) found_diag = true;
438 
439  // Get the nonzeros indexes of the current column of R
440  StorageIndex st = m_firstRowElt(curIdx); // The traversal of the etree starts here
441  if (st < 0) {
442  m_lastError = "Empty row found during numerical factorization";
443  m_info = InvalidInput;
444  return;
445  }
446 
447  // Traverse the etree
448  Index bi = nzcolR;
449  for (; mark(st) != col; st = m_etree(st)) {
450  Ridx(nzcolR) = st; // Add this row to the list,
451  mark(st) = col; // and mark this row as visited
452  nzcolR++;
453  }
454 
455  // Reverse the list to get the topological ordering
456  Index nt = nzcolR - bi;
457  for (Index i = 0; i < nt / 2; i++) std::swap(Ridx(bi + i), Ridx(nzcolR - i - 1));
458 
459  // Copy the current (curIdx,pcol) value of the input matrix
460  if (itp)
461  tval(curIdx) = itp.value();
462  else
463  tval(curIdx) = Scalar(0);
464 
465  // Compute the pattern of Q(:,k)
466  if (curIdx > nonzeroCol && mark(curIdx) != col) {
467  Qidx(nzcolQ) = curIdx; // Add this row to the pattern of Q,
468  mark(curIdx) = col; // and mark it as visited
469  nzcolQ++;
470  }
471  }
472 
473  // Browse all the indexes of R(:,col) in reverse order
474  for (Index i = nzcolR - 1; i >= 0; i--) {
475  Index curIdx = Ridx(i);
476 
477  // Apply the curIdx-th householder vector to the current column (temporarily stored into tval)
478  Scalar tdot(0);
479 
480  // First compute q' * tval
481  tdot = m_Q.col(curIdx).dot(tval);
482 
483  tdot *= m_hcoeffs(curIdx);
484 
485  // Then update tval = tval - q * tau
486  tval -= tdot * m_Q.col(curIdx);
487 
488  // Detect fill-in for the current column of Q
489  if (m_etree(Ridx(i)) == nonzeroCol) {
490  for (typename QRMatrixType::InnerIterator itq(m_Q, curIdx); itq; ++itq) {
491  StorageIndex iQ = StorageIndex(itq.row());
492  if (mark(iQ) != col) {
493  Qidx(nzcolQ++) = iQ; // Add this row to the pattern of Q,
494  mark(iQ) = col; // and mark it as visited
495  }
496  }
497  }
498  } // End update current column
499 
500  Scalar tau = RealScalar(0);
501  RealScalar beta = 0;
502 
503  if (nonzeroCol < diagSize) {
504  // Compute the Householder reflection that eliminate the current column
505  // FIXME this step should call the Householder module.
506  Scalar c0 = nzcolQ ? tval(Qidx(0)) : Scalar(0);
507 
508  // First, the squared norm of Q((col+1):m, col)
509  RealScalar sqrNorm = 0.;
510  for (Index itq = 1; itq < nzcolQ; ++itq) sqrNorm += numext::abs2(tval(Qidx(itq)));
511  if (sqrNorm == RealScalar(0) && numext::imag(c0) == RealScalar(0)) {
512  beta = numext::real(c0);
513  tval(Qidx(0)) = 1;
514  } else {
515  using std::sqrt;
516  beta = sqrt(numext::abs2(c0) + sqrNorm);
517  if (numext::real(c0) >= RealScalar(0)) beta = -beta;
518  tval(Qidx(0)) = 1;
519  for (Index itq = 1; itq < nzcolQ; ++itq) tval(Qidx(itq)) /= (c0 - beta);
520  tau = numext::conj((beta - c0) / beta);
521  }
522  }
523 
524  // Insert values in R
525  for (Index i = nzcolR - 1; i >= 0; i--) {
526  Index curIdx = Ridx(i);
527  if (curIdx < nonzeroCol) {
528  m_R.insertBackByOuterInnerUnordered(col, curIdx) = tval(curIdx);
529  tval(curIdx) = Scalar(0.);
530  }
531  }
532 
533  if (nonzeroCol < diagSize && abs(beta) >= pivotThreshold) {
534  m_R.insertBackByOuterInner(col, nonzeroCol) = beta;
535  // The householder coefficient
536  m_hcoeffs(nonzeroCol) = tau;
537  // Record the householder reflections
538  for (Index itq = 0; itq < nzcolQ; ++itq) {
539  Index iQ = Qidx(itq);
540  m_Q.insertBackByOuterInnerUnordered(nonzeroCol, iQ) = tval(iQ);
541  tval(iQ) = Scalar(0.);
542  }
543  nonzeroCol++;
544  if (nonzeroCol < diagSize) m_Q.startVec(nonzeroCol);
545  } else {
546  // Zero pivot found: move implicitly this column to the end
547  for (Index j = nonzeroCol; j < n - 1; j++) std::swap(m_pivotperm.indices()(j), m_pivotperm.indices()[j + 1]);
548 
549  // Recompute the column elimination tree
550  internal::coletree(m_pmat, m_etree, m_firstRowElt, m_pivotperm.indices().data());
551  m_isEtreeOk = false;
552  }
553  }
554 
555  m_hcoeffs.tail(diagSize - nonzeroCol).setZero();
556 
557  // Finalize the column pointers of the sparse matrices R and Q
558  m_Q.finalize();
559  m_Q.makeCompressed();
560  m_R.finalize();
561  m_R.makeCompressed();
562  m_isQSorted = false;
563 
564  m_nonzeropivots = nonzeroCol;
565 
566  if (nonzeroCol < n) {
567  // Permute the triangular factor to put the 'dead' columns to the end
568  QRMatrixType tempR(m_R);
569  m_R = tempR * m_pivotperm;
570 
571  // Update the column permutation
572  m_outputPerm_c = m_outputPerm_c * m_pivotperm;
573  }
574 
575  m_isInitialized = true;
576  m_factorizationIsok = true;
577  m_info = Success;
578 }
579 
580 template <typename SparseQRType, typename Derived>
581 struct SparseQR_QProduct : ReturnByValue<SparseQR_QProduct<SparseQRType, Derived> > {
582  typedef typename SparseQRType::QRMatrixType MatrixType;
583  typedef typename SparseQRType::Scalar Scalar;
584  // Get the references
585  SparseQR_QProduct(const SparseQRType& qr, const Derived& other, bool transpose)
586  : m_qr(qr), m_other(other), m_transpose(transpose) {}
587  inline Index rows() const { return m_qr.matrixQ().rows(); }
588  inline Index cols() const { return m_other.cols(); }
589 
590  // Assign to a vector
591  template <typename DesType>
592  void evalTo(DesType& res) const {
593  Index m = m_qr.rows();
594  Index n = m_qr.cols();
595  Index diagSize = (std::min)(m, n);
596  res = m_other;
597  if (m_transpose) {
598  eigen_assert(m_qr.m_Q.rows() == m_other.rows() && "Non conforming object sizes");
599  // Compute res = Q' * other column by column
600  for (Index j = 0; j < res.cols(); j++) {
601  for (Index k = 0; k < diagSize; k++) {
602  Scalar tau = Scalar(0);
603  tau = m_qr.m_Q.col(k).dot(res.col(j));
604  if (tau == Scalar(0)) continue;
605  tau = tau * m_qr.m_hcoeffs(k);
606  res.col(j) -= tau * m_qr.m_Q.col(k);
607  }
608  }
609  } else {
610  eigen_assert(m_qr.matrixQ().cols() == m_other.rows() && "Non conforming object sizes");
611 
612  res.conservativeResize(rows(), cols());
613 
614  // Compute res = Q * other column by column
615  for (Index j = 0; j < res.cols(); j++) {
616  Index start_k = internal::is_identity<Derived>::value ? numext::mini(j, diagSize - 1) : diagSize - 1;
617  for (Index k = start_k; k >= 0; k--) {
618  Scalar tau = Scalar(0);
619  tau = m_qr.m_Q.col(k).dot(res.col(j));
620  if (tau == Scalar(0)) continue;
621  tau = tau * numext::conj(m_qr.m_hcoeffs(k));
622  res.col(j) -= tau * m_qr.m_Q.col(k);
623  }
624  }
625  }
626  }
627 
628  const SparseQRType& m_qr;
629  const Derived& m_other;
630  bool m_transpose; // TODO this actually means adjoint
631 };
632 
633 template <typename SparseQRType>
634 struct SparseQRMatrixQReturnType : public EigenBase<SparseQRMatrixQReturnType<SparseQRType> > {
635  typedef typename SparseQRType::Scalar Scalar;
636  typedef Matrix<Scalar, Dynamic, Dynamic> DenseMatrix;
637  enum { RowsAtCompileTime = Dynamic, ColsAtCompileTime = Dynamic };
638  explicit SparseQRMatrixQReturnType(const SparseQRType& qr) : m_qr(qr) {}
639  template <typename Derived>
640  SparseQR_QProduct<SparseQRType, Derived> operator*(const MatrixBase<Derived>& other) {
641  return SparseQR_QProduct<SparseQRType, Derived>(m_qr, other.derived(), false);
642  }
643  // To use for operations with the adjoint of Q
644  SparseQRMatrixQTransposeReturnType<SparseQRType> adjoint() const {
645  return SparseQRMatrixQTransposeReturnType<SparseQRType>(m_qr);
646  }
647  inline Index rows() const { return m_qr.rows(); }
648  inline Index cols() const { return m_qr.rows(); }
649  // To use for operations with the transpose of Q FIXME this is the same as adjoint at the moment
650  SparseQRMatrixQTransposeReturnType<SparseQRType> transpose() const {
651  return SparseQRMatrixQTransposeReturnType<SparseQRType>(m_qr);
652  }
653  const SparseQRType& m_qr;
654 };
655 
656 // TODO this actually represents the adjoint of Q
657 template <typename SparseQRType>
658 struct SparseQRMatrixQTransposeReturnType {
659  explicit SparseQRMatrixQTransposeReturnType(const SparseQRType& qr) : m_qr(qr) {}
660  template <typename Derived>
661  SparseQR_QProduct<SparseQRType, Derived> operator*(const MatrixBase<Derived>& other) {
662  return SparseQR_QProduct<SparseQRType, Derived>(m_qr, other.derived(), true);
663  }
664  const SparseQRType& m_qr;
665 };
666 
667 namespace internal {
668 
669 template <typename SparseQRType>
670 struct evaluator_traits<SparseQRMatrixQReturnType<SparseQRType> > {
671  typedef typename SparseQRType::MatrixType MatrixType;
672  typedef typename storage_kind_to_evaluator_kind<typename MatrixType::StorageKind>::Kind Kind;
673  typedef SparseShape Shape;
674 };
675 
676 template <typename DstXprType, typename SparseQRType>
677 struct Assignment<DstXprType, SparseQRMatrixQReturnType<SparseQRType>,
678  internal::assign_op<typename DstXprType::Scalar, typename DstXprType::Scalar>, Sparse2Sparse> {
679  typedef SparseQRMatrixQReturnType<SparseQRType> SrcXprType;
680  typedef typename DstXprType::Scalar Scalar;
681  typedef typename DstXprType::StorageIndex StorageIndex;
682  static void run(DstXprType& dst, const SrcXprType& src, const internal::assign_op<Scalar, Scalar>& /*func*/) {
683  typename DstXprType::PlainObject idMat(src.rows(), src.cols());
684  idMat.setIdentity();
685  // Sort the sparse householder reflectors if needed
686  const_cast<SparseQRType*>(&src.m_qr)->_sort_matrix_Q();
687  dst = SparseQR_QProduct<SparseQRType, DstXprType>(src.m_qr, idMat, false);
688  }
689 };
690 
691 template <typename DstXprType, typename SparseQRType>
692 struct Assignment<DstXprType, SparseQRMatrixQReturnType<SparseQRType>,
693  internal::assign_op<typename DstXprType::Scalar, typename DstXprType::Scalar>, Sparse2Dense> {
694  typedef SparseQRMatrixQReturnType<SparseQRType> SrcXprType;
695  typedef typename DstXprType::Scalar Scalar;
696  typedef typename DstXprType::StorageIndex StorageIndex;
697  static void run(DstXprType& dst, const SrcXprType& src, const internal::assign_op<Scalar, Scalar>& /*func*/) {
698  dst = src.m_qr.matrixQ() * DstXprType::Identity(src.m_qr.rows(), src.m_qr.rows());
699  }
700 };
701 
702 } // end namespace internal
703 
704 } // end namespace Eigen
705 
706 #endif
const Product< MatrixDerived, PermutationDerived, DefaultProduct > operator*(const MatrixBase< MatrixDerived > &matrix, const PermutationBase< PermutationDerived > &permutation)
Definition: PermutationMatrix.h:474
constexpr Derived & derived()
Definition: EigenBase.h:49
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_sqrt_op< typename Derived::Scalar >, const Derived > sqrt(const Eigen::ArrayBase< Derived > &x)
Index rows() const
Definition: SparseMatrix.h:159
const PermutationType & colsPermutation() const
Definition: SparseQR.h:190
Namespace containing all symbols from the Eigen library.
Definition: B01_Experimental.dox:1
void compute(const MatrixType &mat)
Definition: SparseQR.h:128
Index cols() const
Definition: SparseQR.h:141
FixedBlockXpr<...,... >::Type topLeftCorner(NRowsType cRows, NColsType cCols)
Definition: SparseMatrixBase.h:285
Holds information about the various numeric (i.e. scalar) types allowed by Eigen. ...
Definition: NumTraits.h:232
constexpr const Scalar * data() const
Definition: PlainObjectBase.h:247
const Solve< SparseQR, Rhs > solve(const MatrixBase< Rhs > &B) const
Definition: SparseQR.h:244
void factorize(const MatrixType &mat)
Performs the numerical QR factorization of the input matrix.
Definition: SparseQR.h:356
Index size() const
Definition: PermutationMatrix.h:96
Derived & setZero(Index size)
Definition: CwiseNullaryOp.h:567
Sparse left-looking QR factorization with numerical column pivoting.
Definition: SparseQR.h:20
Index rows() const
Definition: SparseMatrixBase.h:182
Base class of any sparse matrices or sparse expressions.
Definition: ForwardDeclarations.h:481
ComputationInfo info() const
Reports whether previous computation was successful.
Definition: SparseQR.h:266
Index rank() const
Definition: SparseQR.h:162
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:82
const QRMatrixType & matrixR() const
Definition: SparseQR.h:156
Derived & setConstant(Index size, const Scalar &val)
Definition: CwiseNullaryOp.h:363
Definition: Constants.h:447
SparseQR(const MatrixType &mat)
Definition: SparseQR.h:117
Definition: Constants.h:440
constexpr Derived & derived()
Definition: EigenBase.h:49
void analyzePattern(const MatrixType &mat)
Preprocessing step of a QR factorization.
Definition: SparseQR.h:314
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_abs_op< typename Derived::Scalar >, const Derived > abs(const Eigen::ArrayBase< Derived > &x)
Index cols() const
Definition: SparseMatrix.h:161
void setPivotThreshold(const RealScalar &threshold)
Definition: SparseQR.h:234
SparseQRMatrixQReturnType< SparseQR > matrixQ() const
Definition: SparseQR.h:185
const int Dynamic
Definition: Constants.h:25
Pseudo expression representing a solving operation.
Definition: Solve.h:62
Index rows() const
Definition: SparseQR.h:137
constexpr Index rows() const noexcept
Definition: EigenBase.h:59
ComputationInfo
Definition: Constants.h:438
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:52
std::string lastErrorMessage() const
Definition: SparseQR.h:198