SparseQR.h
Go to the documentation of this file.
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 namespace Eigen {
15 
16 template<typename MatrixType, typename OrderingType> class SparseQR;
17 template<typename SparseQRType> struct SparseQRMatrixQReturnType;
18 template<typename SparseQRType> struct SparseQRMatrixQTransposeReturnType;
19 template<typename SparseQRType, typename Derived> struct SparseQR_QProduct;
20 namespace internal {
21  template <typename SparseQRType> struct traits<SparseQRMatrixQReturnType<SparseQRType> >
22  {
24  typedef typename ReturnType::StorageIndex StorageIndex;
25  typedef typename ReturnType::StorageKind StorageKind;
26  enum {
27  RowsAtCompileTime = Dynamic,
28  ColsAtCompileTime = Dynamic
29  };
30  };
31  template <typename SparseQRType> struct traits<SparseQRMatrixQTransposeReturnType<SparseQRType> >
32  {
34  };
35  template <typename SparseQRType, typename Derived> struct traits<SparseQR_QProduct<SparseQRType, Derived> >
36  {
37  typedef typename Derived::PlainObject ReturnType;
38  };
39 } // End namespace internal
40 
71 template<typename _MatrixType, typename _OrderingType>
72 class SparseQR : public SparseSolverBase<SparseQR<_MatrixType,_OrderingType> >
73 {
74  protected:
76  using Base::m_isInitialized;
77  public:
78  using Base::_solve_impl;
79  typedef _MatrixType MatrixType;
80  typedef _OrderingType OrderingType;
81  typedef typename MatrixType::Scalar Scalar;
83  typedef typename MatrixType::StorageIndex StorageIndex;
88 
89  enum {
90  ColsAtCompileTime = MatrixType::ColsAtCompileTime,
91  MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
92  };
93 
94  public:
95  SparseQR () : m_analysisIsok(false), m_lastError(""), m_useDefaultThreshold(true),m_isQSorted(false),m_isEtreeOk(false)
96  { }
97 
104  explicit SparseQR(const MatrixType& mat) : m_analysisIsok(false), m_lastError(""), m_useDefaultThreshold(true),m_isQSorted(false),m_isEtreeOk(false)
105  {
106  compute(mat);
107  }
108 
115  void compute(const MatrixType& mat)
116  {
117  analyzePattern(mat);
118  factorize(mat);
119  }
120  void analyzePattern(const MatrixType& mat);
121  void factorize(const MatrixType& mat);
122 
125  inline Index rows() const { return m_pmat.rows(); }
126 
129  inline Index cols() const { return m_pmat.cols();}
130 
144  const QRMatrixType& matrixR() const { return m_R; }
145 
150  Index rank() const
151  {
152  eigen_assert(m_isInitialized && "The factorization should be called first, use compute()");
153  return m_nonzeropivots;
154  }
155 
175  { return SparseQRMatrixQReturnType<SparseQR>(*this); }
176 
180  const PermutationType& colsPermutation() const
181  {
182  eigen_assert(m_isInitialized && "Decomposition is not initialized.");
183  return m_outputPerm_c;
184  }
185 
189  std::string lastErrorMessage() const { return m_lastError; }
190 
192  template<typename Rhs, typename Dest>
193  bool _solve_impl(const MatrixBase<Rhs> &B, MatrixBase<Dest> &dest) const
194  {
195  eigen_assert(m_isInitialized && "The factorization should be called first, use compute()");
196  eigen_assert(this->rows() == B.rows() && "SparseQR::solve() : invalid number of rows in the right hand side matrix");
197 
198  Index rank = this->rank();
199 
200  // Compute Q^* * b;
201  typename Dest::PlainObject y, b;
202  y = this->matrixQ().adjoint() * B;
203  b = y;
204 
205  // Solve with the triangular matrix R
206  y.resize((std::max<Index>)(cols(),y.rows()),y.cols());
207  y.topRows(rank) = this->matrixR().topLeftCorner(rank, rank).template triangularView<Upper>().solve(b.topRows(rank));
208  y.bottomRows(y.rows()-rank).setZero();
209 
210  // Apply the column permutation
211  if (m_perm_c.size()) dest = colsPermutation() * y.topRows(cols());
212  else dest = y.topRows(cols());
213 
214  m_info = Success;
215  return true;
216  }
217 
223  void setPivotThreshold(const RealScalar& threshold)
224  {
225  m_useDefaultThreshold = false;
226  m_threshold = threshold;
227  }
228 
233  template<typename Rhs>
234  inline const Solve<SparseQR, Rhs> solve(const MatrixBase<Rhs>& B) const
235  {
236  eigen_assert(m_isInitialized && "The factorization should be called first, use compute()");
237  eigen_assert(this->rows() == B.rows() && "SparseQR::solve() : invalid number of rows in the right hand side matrix");
238  return Solve<SparseQR, Rhs>(*this, B.derived());
239  }
240  template<typename Rhs>
242  {
243  eigen_assert(m_isInitialized && "The factorization should be called first, use compute()");
244  eigen_assert(this->rows() == B.rows() && "SparseQR::solve() : invalid number of rows in the right hand side matrix");
245  return Solve<SparseQR, Rhs>(*this, B.derived());
246  }
247 
257  {
258  eigen_assert(m_isInitialized && "Decomposition is not initialized.");
259  return m_info;
260  }
261 
262 
264  inline void _sort_matrix_Q()
265  {
266  if(this->m_isQSorted) return;
267  // The matrix Q is sorted during the transposition
269  this->m_Q = mQrm;
270  this->m_isQSorted = true;
271  }
272 
273 
274  protected:
278  std::string m_lastError;
279  QRMatrixType m_pmat; // Temporary matrix
280  QRMatrixType m_R; // The triangular factor matrix
281  QRMatrixType m_Q; // The orthogonal reflectors
282  ScalarVector m_hcoeffs; // The Householder coefficients
283  PermutationType m_perm_c; // Fill-reducing Column permutation
284  PermutationType m_pivotperm; // The permutation for rank revealing
285  PermutationType m_outputPerm_c; // The final column permutation
286  RealScalar m_threshold; // Threshold to determine null Householder reflections
287  bool m_useDefaultThreshold; // Use default threshold
288  Index m_nonzeropivots; // Number of non zero pivots found
289  IndexVector m_etree; // Column elimination tree
290  IndexVector m_firstRowElt; // First element in each row
291  bool m_isQSorted; // whether Q is sorted or not
292  bool m_isEtreeOk; // whether the elimination tree match the initial input matrix
293 
294  template <typename, typename > friend struct SparseQR_QProduct;
295 
296 };
297 
307 template <typename MatrixType, typename OrderingType>
309 {
310  eigen_assert(mat.isCompressed() && "SparseQR requires a sparse matrix in compressed mode. Call .makeCompressed() before passing it to SparseQR");
311  // Copy to a column major matrix if the input is rowmajor
313  // Compute the column fill reducing ordering
314  OrderingType ord;
315  ord(matCpy, m_perm_c);
316  Index n = mat.cols();
317  Index m = mat.rows();
318  Index diagSize = (std::min)(m,n);
319 
320  if (!m_perm_c.size())
321  {
322  m_perm_c.resize(n);
323  m_perm_c.indices().setLinSpaced(n, 0,StorageIndex(n-1));
324  }
325 
326  // Compute the column elimination tree of the permuted matrix
327  m_outputPerm_c = m_perm_c.inverse();
328  internal::coletree(matCpy, m_etree, m_firstRowElt, m_outputPerm_c.indices().data());
329  m_isEtreeOk = true;
330 
331  m_R.resize(m, n);
332  m_Q.resize(m, diagSize);
333 
334  // Allocate space for nonzero elements : rough estimation
335  m_R.reserve(2*mat.nonZeros()); //FIXME Get a more accurate estimation through symbolic factorization with the etree
336  m_Q.reserve(2*mat.nonZeros());
337  m_hcoeffs.resize(diagSize);
338  m_analysisIsok = true;
339 }
340 
348 template <typename MatrixType, typename OrderingType>
350 {
351  using std::abs;
352 
353  eigen_assert(m_analysisIsok && "analyzePattern() should be called before this step");
354  StorageIndex m = StorageIndex(mat.rows());
355  StorageIndex n = StorageIndex(mat.cols());
356  StorageIndex diagSize = (std::min)(m,n);
357  IndexVector mark((std::max)(m,n)); mark.setConstant(-1); // Record the visited nodes
358  IndexVector Ridx(n), Qidx(m); // Store temporarily the row indexes for the current column of R and Q
359  Index nzcolR, nzcolQ; // Number of nonzero for the current column of R and Q
360  ScalarVector tval(m); // The dense vector used to compute the current column
361  RealScalar pivotThreshold = m_threshold;
362 
363  m_R.setZero();
364  m_Q.setZero();
365  m_pmat = mat;
366  if(!m_isEtreeOk)
367  {
368  m_outputPerm_c = m_perm_c.inverse();
369  internal::coletree(m_pmat, m_etree, m_firstRowElt, m_outputPerm_c.indices().data());
370  m_isEtreeOk = true;
371  }
372 
373  m_pmat.uncompress(); // To have the innerNonZeroPtr allocated
374 
375  // Apply the fill-in reducing permutation lazily:
376  {
377  // If the input is row major, copy the original column indices,
378  // otherwise directly use the input matrix
379  //
380  IndexVector originalOuterIndicesCpy;
381  const StorageIndex *originalOuterIndices = mat.outerIndexPtr();
382  if(MatrixType::IsRowMajor)
383  {
384  originalOuterIndicesCpy = IndexVector::Map(m_pmat.outerIndexPtr(),n+1);
385  originalOuterIndices = originalOuterIndicesCpy.data();
386  }
387 
388  for (int i = 0; i < n; i++)
389  {
390  Index p = m_perm_c.size() ? m_perm_c.indices()(i) : i;
391  m_pmat.outerIndexPtr()[p] = originalOuterIndices[i];
392  m_pmat.innerNonZeroPtr()[p] = originalOuterIndices[i+1] - originalOuterIndices[i];
393  }
394  }
395 
396  /* Compute the default threshold as in MatLab, see:
397  * Tim Davis, "Algorithm 915, SuiteSparseQR: Multifrontal Multithreaded Rank-Revealing
398  * Sparse QR Factorization, ACM Trans. on Math. Soft. 38(1), 2011, Page 8:3
399  */
400  if(m_useDefaultThreshold)
401  {
402  RealScalar max2Norm = 0.0;
403  for (int j = 0; j < n; j++) max2Norm = numext::maxi(max2Norm, m_pmat.col(j).norm());
404  if(max2Norm==RealScalar(0))
405  max2Norm = RealScalar(1);
406  pivotThreshold = 20 * (m + n) * max2Norm * NumTraits<RealScalar>::epsilon();
407  }
408 
409  // Initialize the numerical permutation
410  m_pivotperm.setIdentity(n);
411 
412  StorageIndex nonzeroCol = 0; // Record the number of valid pivots
413  m_Q.startVec(0);
414 
415  // Left looking rank-revealing QR factorization: compute a column of R and Q at a time
416  for (StorageIndex col = 0; col < n; ++col)
417  {
418  mark.setConstant(-1);
419  m_R.startVec(col);
420  mark(nonzeroCol) = col;
421  Qidx(0) = nonzeroCol;
422  nzcolR = 0; nzcolQ = 1;
423  bool found_diag = nonzeroCol>=m;
424  tval.setZero();
425 
426  // Symbolic factorization: find the nonzero locations of the column k of the factors R and Q, i.e.,
427  // all the nodes (with indexes lower than rank) reachable through the column elimination tree (etree) rooted at node k.
428  // Note: if the diagonal entry does not exist, then its contribution must be explicitly added,
429  // thus the trick with found_diag that permits to do one more iteration on the diagonal element if this one has not been found.
430  for (typename QRMatrixType::InnerIterator itp(m_pmat, col); itp || !found_diag; ++itp)
431  {
432  StorageIndex curIdx = nonzeroCol;
433  if(itp) curIdx = StorageIndex(itp.row());
434  if(curIdx == nonzeroCol) found_diag = true;
435 
436  // Get the nonzeros indexes of the current column of R
437  StorageIndex st = m_firstRowElt(curIdx); // The traversal of the etree starts here
438  if (st < 0 )
439  {
440  m_lastError = "Empty row found during numerical factorization";
441  m_info = InvalidInput;
442  return;
443  }
444 
445  // Traverse the etree
446  Index bi = nzcolR;
447  for (; mark(st) != col; st = m_etree(st))
448  {
449  Ridx(nzcolR) = st; // Add this row to the list,
450  mark(st) = col; // and mark this row as visited
451  nzcolR++;
452  }
453 
454  // Reverse the list to get the topological ordering
455  Index nt = nzcolR-bi;
456  for(Index i = 0; i < nt/2; i++) std::swap(Ridx(bi+i), Ridx(nzcolR-i-1));
457 
458  // Copy the current (curIdx,pcol) value of the input matrix
459  if(itp) tval(curIdx) = itp.value();
460  else tval(curIdx) = Scalar(0);
461 
462  // Compute the pattern of Q(:,k)
463  if(curIdx > nonzeroCol && mark(curIdx) != col )
464  {
465  Qidx(nzcolQ) = curIdx; // Add this row to the pattern of Q,
466  mark(curIdx) = col; // and mark it as visited
467  nzcolQ++;
468  }
469  }
470 
471  // Browse all the indexes of R(:,col) in reverse order
472  for (Index i = nzcolR-1; i >= 0; i--)
473  {
474  Index curIdx = Ridx(i);
475 
476  // Apply the curIdx-th householder vector to the current column (temporarily stored into tval)
477  Scalar tdot(0);
478 
479  // First compute q' * tval
480  tdot = m_Q.col(curIdx).dot(tval);
481 
482  tdot *= m_hcoeffs(curIdx);
483 
484  // Then update tval = tval - q * tau
485  // FIXME: tval -= tdot * m_Q.col(curIdx) should amount to the same (need to check/add support for efficient "dense ?= sparse")
486  for (typename QRMatrixType::InnerIterator itq(m_Q, curIdx); itq; ++itq)
487  tval(itq.row()) -= itq.value() * tdot;
488 
489  // Detect fill-in for the current column of Q
490  if(m_etree(Ridx(i)) == nonzeroCol)
491  {
492  for (typename QRMatrixType::InnerIterator itq(m_Q, curIdx); itq; ++itq)
493  {
494  StorageIndex iQ = StorageIndex(itq.row());
495  if (mark(iQ) != col)
496  {
497  Qidx(nzcolQ++) = iQ; // Add this row to the pattern of Q,
498  mark(iQ) = col; // and mark it as visited
499  }
500  }
501  }
502  } // End update current column
503 
504  Scalar tau = RealScalar(0);
505  RealScalar beta = 0;
506 
507  if(nonzeroCol < diagSize)
508  {
509  // Compute the Householder reflection that eliminate the current column
510  // FIXME this step should call the Householder module.
511  Scalar c0 = nzcolQ ? tval(Qidx(0)) : Scalar(0);
512 
513  // First, the squared norm of Q((col+1):m, col)
514  RealScalar sqrNorm = 0.;
515  for (Index itq = 1; itq < nzcolQ; ++itq) sqrNorm += numext::abs2(tval(Qidx(itq)));
516  if(sqrNorm == RealScalar(0) && numext::imag(c0) == RealScalar(0))
517  {
518  beta = numext::real(c0);
519  tval(Qidx(0)) = 1;
520  }
521  else
522  {
523  using std::sqrt;
524  beta = sqrt(numext::abs2(c0) + sqrNorm);
525  if(numext::real(c0) >= RealScalar(0))
526  beta = -beta;
527  tval(Qidx(0)) = 1;
528  for (Index itq = 1; itq < nzcolQ; ++itq)
529  tval(Qidx(itq)) /= (c0 - beta);
530  tau = numext::conj((beta-c0) / beta);
531 
532  }
533  }
534 
535  // Insert values in R
536  for (Index i = nzcolR-1; i >= 0; i--)
537  {
538  Index curIdx = Ridx(i);
539  if(curIdx < nonzeroCol)
540  {
541  m_R.insertBackByOuterInnerUnordered(col, curIdx) = tval(curIdx);
542  tval(curIdx) = Scalar(0.);
543  }
544  }
545 
546  if(nonzeroCol < diagSize && abs(beta) >= pivotThreshold)
547  {
548  m_R.insertBackByOuterInner(col, nonzeroCol) = beta;
549  // The householder coefficient
550  m_hcoeffs(nonzeroCol) = tau;
551  // Record the householder reflections
552  for (Index itq = 0; itq < nzcolQ; ++itq)
553  {
554  Index iQ = Qidx(itq);
555  m_Q.insertBackByOuterInnerUnordered(nonzeroCol,iQ) = tval(iQ);
556  tval(iQ) = Scalar(0.);
557  }
558  nonzeroCol++;
559  if(nonzeroCol<diagSize)
560  m_Q.startVec(nonzeroCol);
561  }
562  else
563  {
564  // Zero pivot found: move implicitly this column to the end
565  for (Index j = nonzeroCol; j < n-1; j++)
566  std::swap(m_pivotperm.indices()(j), m_pivotperm.indices()[j+1]);
567 
568  // Recompute the column elimination tree
569  internal::coletree(m_pmat, m_etree, m_firstRowElt, m_pivotperm.indices().data());
570  m_isEtreeOk = false;
571  }
572  }
573 
574  m_hcoeffs.tail(diagSize-nonzeroCol).setZero();
575 
576  // Finalize the column pointers of the sparse matrices R and Q
577  m_Q.finalize();
578  m_Q.makeCompressed();
579  m_R.finalize();
580  m_R.makeCompressed();
581  m_isQSorted = false;
582 
583  m_nonzeropivots = nonzeroCol;
584 
585  if(nonzeroCol<n)
586  {
587  // Permute the triangular factor to put the 'dead' columns to the end
588  QRMatrixType tempR(m_R);
589  m_R = tempR * m_pivotperm;
590 
591  // Update the column permutation
592  m_outputPerm_c = m_outputPerm_c * m_pivotperm;
593  }
594 
595  m_isInitialized = true;
596  m_factorizationIsok = true;
597  m_info = Success;
598 }
599 
600 template <typename SparseQRType, typename Derived>
601 struct SparseQR_QProduct : ReturnByValue<SparseQR_QProduct<SparseQRType, Derived> >
602 {
603  typedef typename SparseQRType::QRMatrixType MatrixType;
604  typedef typename SparseQRType::Scalar Scalar;
605  // Get the references
606  SparseQR_QProduct(const SparseQRType& qr, const Derived& other, bool transpose) :
607  m_qr(qr),m_other(other),m_transpose(transpose) {}
608  inline Index rows() const { return m_qr.matrixQ().rows(); }
609  inline Index cols() const { return m_other.cols(); }
610 
611  // Assign to a vector
612  template<typename DesType>
613  void evalTo(DesType& res) const
614  {
615  Index m = m_qr.rows();
616  Index n = m_qr.cols();
617  Index diagSize = (std::min)(m,n);
618  res = m_other;
619  if (m_transpose)
620  {
621  eigen_assert(m_qr.m_Q.rows() == m_other.rows() && "Non conforming object sizes");
622  //Compute res = Q' * other column by column
623  for(Index j = 0; j < res.cols(); j++){
624  for (Index k = 0; k < diagSize; k++)
625  {
626  Scalar tau = Scalar(0);
627  tau = m_qr.m_Q.col(k).dot(res.col(j));
628  if(tau==Scalar(0)) continue;
629  tau = tau * m_qr.m_hcoeffs(k);
630  res.col(j) -= tau * m_qr.m_Q.col(k);
631  }
632  }
633  }
634  else
635  {
636  eigen_assert(m_qr.matrixQ().cols() == m_other.rows() && "Non conforming object sizes");
637 
638  res.conservativeResize(rows(), cols());
639 
640  // Compute res = Q * other column by column
641  for(Index j = 0; j < res.cols(); j++)
642  {
643  for (Index k = diagSize-1; k >=0; k--)
644  {
645  Scalar tau = Scalar(0);
646  tau = m_qr.m_Q.col(k).dot(res.col(j));
647  if(tau==Scalar(0)) continue;
648  tau = tau * numext::conj(m_qr.m_hcoeffs(k));
649  res.col(j) -= tau * m_qr.m_Q.col(k);
650  }
651  }
652  }
653  }
654 
655  const SparseQRType& m_qr;
656  const Derived& m_other;
657  bool m_transpose; // TODO this actually means adjoint
658 };
659 
660 template<typename SparseQRType>
661 struct SparseQRMatrixQReturnType : public EigenBase<SparseQRMatrixQReturnType<SparseQRType> >
662 {
663  typedef typename SparseQRType::Scalar Scalar;
665  enum {
666  RowsAtCompileTime = Dynamic,
667  ColsAtCompileTime = Dynamic
668  };
669  explicit SparseQRMatrixQReturnType(const SparseQRType& qr) : m_qr(qr) {}
670  template<typename Derived>
672  {
673  return SparseQR_QProduct<SparseQRType,Derived>(m_qr,other.derived(),false);
674  }
675  // To use for operations with the adjoint of Q
677  {
679  }
680  inline Index rows() const { return m_qr.rows(); }
681  inline Index cols() const { return m_qr.rows(); }
682  // To use for operations with the transpose of Q FIXME this is the same as adjoint at the moment
684  {
686  }
687  const SparseQRType& m_qr;
688 };
689 
690 // TODO this actually represents the adjoint of Q
691 template<typename SparseQRType>
693 {
694  explicit SparseQRMatrixQTransposeReturnType(const SparseQRType& qr) : m_qr(qr) {}
695  template<typename Derived>
697  {
698  return SparseQR_QProduct<SparseQRType,Derived>(m_qr,other.derived(), true);
699  }
700  const SparseQRType& m_qr;
701 };
702 
703 namespace internal {
704 
705 template<typename SparseQRType>
707 {
711 };
712 
713 template< typename DstXprType, typename SparseQRType>
714 struct Assignment<DstXprType, SparseQRMatrixQReturnType<SparseQRType>, internal::assign_op<typename DstXprType::Scalar,typename DstXprType::Scalar>, Sparse2Sparse>
715 {
717  typedef typename DstXprType::Scalar Scalar;
718  typedef typename DstXprType::StorageIndex StorageIndex;
719  static void run(DstXprType &dst, const SrcXprType &src, const internal::assign_op<Scalar,Scalar> &/*func*/)
720  {
721  typename DstXprType::PlainObject idMat(src.rows(), src.cols());
722  idMat.setIdentity();
723  // Sort the sparse householder reflectors if needed
724  const_cast<SparseQRType *>(&src.m_qr)->_sort_matrix_Q();
725  dst = SparseQR_QProduct<SparseQRType, DstXprType>(src.m_qr, idMat, false);
726  }
727 };
728 
729 template< typename DstXprType, typename SparseQRType>
730 struct Assignment<DstXprType, SparseQRMatrixQReturnType<SparseQRType>, internal::assign_op<typename DstXprType::Scalar,typename DstXprType::Scalar>, Sparse2Dense>
731 {
733  typedef typename DstXprType::Scalar Scalar;
734  typedef typename DstXprType::StorageIndex StorageIndex;
735  static void run(DstXprType &dst, const SrcXprType &src, const internal::assign_op<Scalar,Scalar> &/*func*/)
736  {
737  dst = src.m_qr.matrixQ() * DstXprType::Identity(src.m_qr.rows(), src.m_qr.rows());
738  }
739 };
740 
741 } // end namespace internal
742 
743 } // end namespace Eigen
744 
745 #endif
Matrix3f m
_OrderingType OrderingType
Definition: SparseQR.h:80
MatrixType::StorageIndex StorageIndex
Definition: SparseQR.h:83
SCALAR Scalar
Definition: bench_gemm.cpp:33
PermutationMatrix< Dynamic, Dynamic, StorageIndex > PermutationType
Definition: SparseQR.h:87
#define max(a, b)
Definition: datatypes.h:20
Matrix< Scalar, Dynamic, 1 > ScalarVector
Definition: SparseQR.h:86
EIGEN_DEVICE_FUNC Derived & setZero(Index size)
QRMatrixType m_Q
Definition: SparseQR.h:281
float real
Definition: datatypes.h:10
_MatrixType MatrixType
Definition: SparseQR.h:79
Scalar * b
Definition: benchVecAdd.cpp:17
ScalarVector m_hcoeffs
Definition: SparseQR.h:282
int coletree(const MatrixType &mat, IndexVector &parent, IndexVector &firstRowElt, typename MatrixType::StorageIndex *perm=0)
#define min(a, b)
Definition: datatypes.h:19
Index cols() const
Definition: SparseQR.h:609
Index rows() const
Definition: SparseQR.h:608
PermutationType m_perm_c
Definition: SparseQR.h:283
SparseQRMatrixQReturnType(const SparseQRType &qr)
Definition: SparseQR.h:669
Matrix< StorageIndex, Dynamic, 1 > IndexVector
Definition: SparseQR.h:85
const SparseQRType & m_qr
Definition: SparseQR.h:655
A base class for sparse solvers.
int n
QRMatrixType m_pmat
Definition: SparseQR.h:279
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE const Scalar * data() const
Index rank() const
Definition: SparseQR.h:150
storage_kind_to_evaluator_kind< typename MatrixType::StorageKind >::Kind Kind
Definition: SparseQR.h:709
SparseQR_QProduct(const SparseQRType &qr, const Derived &other, bool transpose)
Definition: SparseQR.h:606
EIGEN_DEVICE_FUNC const SqrtReturnType sqrt() const
Namespace containing all symbols from the Eigen library.
Definition: jet.h:637
MatrixXf MatrixType
Holds information about the various numeric (i.e. scalar) types allowed by Eigen. ...
Definition: NumTraits.h:150
EIGEN_DEVICE_FUNC EIGEN_ALWAYS_INLINE T maxi(const T &x, const T &y)
Eigen::Index Index
The interface type of indices.
Definition: EigenBase.h:38
SparseQRType::Scalar Scalar
Definition: SparseQR.h:663
EIGEN_DEVICE_FUNC RowsBlockXpr topRows(Index n)
Definition: DenseBase.h:433
ComputationInfo info() const
Reports whether previous computation was successful.
Definition: SparseQR.h:256
const QRMatrixType & matrixR() const
Definition: SparseQR.h:144
SparseQRType::QRMatrixType MatrixType
Definition: SparseQR.h:603
HouseholderQR< MatrixXf > qr(A)
void factorize(const MatrixType &mat)
Performs the numerical QR factorization of the input matrix.
Definition: SparseQR.h:349
const Derived & m_other
Definition: SparseQR.h:656
void _sort_matrix_Q()
Definition: SparseQR.h:264
SparseQR_QProduct< SparseQRType, Derived > operator*(const MatrixBase< Derived > &other)
Definition: SparseQR.h:696
cout<< "Here is the matrix m:"<< endl<< m<< endl;Matrix< ptrdiff_t, 3, 1 > res
Sparse left-looking rank-revealing QR factorization.
Definition: SparseQR.h:16
IndexVector m_firstRowElt
Definition: SparseQR.h:290
Base class of any sparse matrices or sparse expressions.
SparseSolverBase< SparseQR< _MatrixType, _OrderingType > > Base
Definition: SparseQR.h:75
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:33
Matrix< Scalar, Dynamic, Dynamic > DenseMatrix
Definition: SparseQR.h:664
#define eigen_assert(x)
Definition: Macros.h:579
Index cols() const
Definition: SparseQR.h:129
PermutationType m_pivotperm
Definition: SparseQR.h:284
SparseQRMatrixQTransposeReturnType< SparseQRType > adjoint() const
Definition: SparseQR.h:676
ComputationInfo m_info
Definition: SparseQR.h:277
EIGEN_DEVICE_FUNC Derived & setConstant(Index size, const Scalar &val)
QRMatrixType m_R
Definition: SparseQR.h:280
bool m_isEtreeOk
Definition: SparseQR.h:292
const Solve< SparseQR, Rhs > solve(const SparseMatrixBase< Rhs > &B) const
Definition: SparseQR.h:241
NumTraits< Scalar >::Real RealScalar
Definition: bench_gemm.cpp:34
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE const Abs2ReturnType abs2() const
MatrixType::Scalar Scalar
Definition: SparseQR.h:81
void evalTo(DesType &res) const
Definition: SparseQR.h:613
bool _solve_impl(const MatrixBase< Rhs > &B, MatrixBase< Dest > &dest) const
Definition: SparseQR.h:193
bool m_isQSorted
Definition: SparseQR.h:291
void analyzePattern(const MatrixType &mat)
Preprocessing step of a QR factorization.
Definition: SparseQR.h:308
SparseQRType::Scalar Scalar
Definition: SparseQR.h:604
PermutationType m_outputPerm_c
Definition: SparseQR.h:285
Index rows() const
Definition: SparseQR.h:125
const Derived & derived() const
const PermutationType & colsPermutation() const
Definition: SparseQR.h:180
const Solve< SparseQR, Rhs > solve(const MatrixBase< Rhs > &B) const
Definition: SparseQR.h:234
float * p
SparseQRMatrixQReturnType< SparseQR > matrixQ() const
Definition: SparseQR.h:174
const SparseQRType & m_qr
Definition: SparseQR.h:687
m col(1)
EIGEN_DEVICE_FUNC const ImagReturnType imag() const
SparseQRMatrixQTransposeReturnType(const SparseQRType &qr)
Definition: SparseQR.h:694
std::string m_lastError
Definition: SparseQR.h:278
const int Dynamic
Definition: Constants.h:21
Pseudo expression representing a solving operation.
Definition: Solve.h:62
MatrixType::RealScalar RealScalar
Definition: SparseQR.h:82
EIGEN_DONT_INLINE void compute(Solver &solver, const MatrixType &A)
bool m_analysisIsok
Definition: SparseQR.h:275
#define abs(x)
Definition: datatypes.h:17
ComputationInfo
Definition: Constants.h:430
SparseQR(const MatrixType &mat)
Definition: SparseQR.h:104
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:48
IndexVector m_etree
Definition: SparseQR.h:289
void swap(mpfr::mpreal &x, mpfr::mpreal &y)
Definition: mpreal.h:2986
std::ptrdiff_t j
RealScalar m_threshold
Definition: SparseQR.h:286
SparseQRMatrixQTransposeReturnType< SparseQRType > transpose() const
Definition: SparseQR.h:683
Index m_nonzeropivots
Definition: SparseQR.h:288
void setPivotThreshold(const RealScalar &threshold)
Definition: SparseQR.h:223
std::string lastErrorMessage() const
Definition: SparseQR.h:189
bool m_useDefaultThreshold
Definition: SparseQR.h:287
ScalarWithExceptions conj(const ScalarWithExceptions &x)
Definition: exceptions.cpp:74
SparseMatrix< Scalar, ColMajor, StorageIndex > QRMatrixType
Definition: SparseQR.h:84
SparseQR_QProduct< SparseQRType, Derived > operator*(const MatrixBase< Derived > &other)
Definition: SparseQR.h:671
v setZero(3)
bool m_factorizationIsok
Definition: SparseQR.h:276
void compute(const MatrixType &mat)
Definition: SparseQR.h:115


gtsam
Author(s):
autogenerated on Sat May 8 2021 02:44:36