• Home
  • Line#
  • Scopes#
  • Navigate#
  • Raw
  • Download
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   {
23     typedef typename SparseQRType::MatrixType ReturnType;
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   {
33     typedef typename SparseQRType::MatrixType ReturnType;
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 
41 /**
42   * \ingroup SparseQR_Module
43   * \class SparseQR
44   * \brief Sparse left-looking rank-revealing QR factorization
45   *
46   * This class implements a left-looking rank-revealing QR decomposition
47   * of sparse matrices. When a column has a norm less than a given tolerance
48   * it is implicitly permuted to the end. The QR factorization thus obtained is
49   * given by A*P = Q*R where R is upper triangular or trapezoidal.
50   *
51   * P is the column permutation which is the product of the fill-reducing and the
52   * rank-revealing permutations. Use colsPermutation() to get it.
53   *
54   * Q is the orthogonal matrix represented as products of Householder reflectors.
55   * Use matrixQ() to get an expression and matrixQ().adjoint() to get the adjoint.
56   * You can then apply it to a vector.
57   *
58   * R is the sparse triangular or trapezoidal matrix. The later occurs when A is rank-deficient.
59   * matrixR().topLeftCorner(rank(), rank()) always returns a triangular factor of full rank.
60   *
61   * \tparam _MatrixType The type of the sparse matrix A, must be a column-major SparseMatrix<>
62   * \tparam _OrderingType The fill-reducing ordering method. See the \link OrderingMethods_Module
63   *  OrderingMethods \endlink module for the list of built-in and external ordering methods.
64   *
65   * \implsparsesolverconcept
66   *
67   * \warning The input sparse matrix A must be in compressed mode (see SparseMatrix::makeCompressed()).
68   * \warning For complex matrices matrixQ().transpose() will actually return the adjoint matrix.
69   *
70   */
71 template<typename _MatrixType, typename _OrderingType>
72 class SparseQR : public SparseSolverBase<SparseQR<_MatrixType,_OrderingType> >
73 {
74   protected:
75     typedef SparseSolverBase<SparseQR<_MatrixType,_OrderingType> > Base;
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;
82     typedef typename MatrixType::RealScalar RealScalar;
83     typedef typename MatrixType::StorageIndex StorageIndex;
84     typedef SparseMatrix<Scalar,ColMajor,StorageIndex> QRMatrixType;
85     typedef Matrix<StorageIndex, Dynamic, 1> IndexVector;
86     typedef Matrix<Scalar, Dynamic, 1> ScalarVector;
87     typedef PermutationMatrix<Dynamic, Dynamic, StorageIndex> PermutationType;
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 
98     /** Construct a QR factorization of the matrix \a mat.
99       *
100       * \warning The matrix \a mat must be in compressed mode (see SparseMatrix::makeCompressed()).
101       *
102       * \sa compute()
103       */
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 
109     /** Computes the QR factorization of the sparse matrix \a mat.
110       *
111       * \warning The matrix \a mat must be in compressed mode (see SparseMatrix::makeCompressed()).
112       *
113       * \sa analyzePattern(), factorize()
114       */
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 
123     /** \returns the number of rows of the represented matrix.
124       */
125     inline Index rows() const { return m_pmat.rows(); }
126 
127     /** \returns the number of columns of the represented matrix.
128       */
129     inline Index cols() const { return m_pmat.cols();}
130 
131     /** \returns a const reference to the \b sparse upper triangular matrix R of the QR factorization.
132       * \warning The entries of the returned matrix are not sorted. This means that using it in algorithms
133       *          expecting sorted entries will fail. This include random coefficient accesses (SpaseMatrix::coeff()),
134       *          and coefficient-wise operations. Matrix products and triangular solves are fine though.
135       *
136       * To sort the entries, you can assign it to a row-major matrix, and if a column-major matrix
137       * is required, you can copy it again:
138       * \code
139       * SparseMatrix<double>          R  = qr.matrixR();  // column-major, not sorted!
140       * SparseMatrix<double,RowMajor> Rr = qr.matrixR();  // row-major, sorted
141       * SparseMatrix<double>          Rc = Rr;            // column-major, sorted
142       * \endcode
143       */
144     const QRMatrixType& matrixR() const { return m_R; }
145 
146     /** \returns the number of non linearly dependent columns as determined by the pivoting threshold.
147       *
148       * \sa setPivotThreshold()
149       */
150     Index rank() const
151     {
152       eigen_assert(m_isInitialized && "The factorization should be called first, use compute()");
153       return m_nonzeropivots;
154     }
155 
156     /** \returns an expression of the matrix Q as products of sparse Householder reflectors.
157     * The common usage of this function is to apply it to a dense matrix or vector
158     * \code
159     * VectorXd B1, B2;
160     * // Initialize B1
161     * B2 = matrixQ() * B1;
162     * \endcode
163     *
164     * To get a plain SparseMatrix representation of Q:
165     * \code
166     * SparseMatrix<double> Q;
167     * Q = SparseQR<SparseMatrix<double> >(A).matrixQ();
168     * \endcode
169     * Internally, this call simply performs a sparse product between the matrix Q
170     * and a sparse identity matrix. However, due to the fact that the sparse
171     * reflectors are stored unsorted, two transpositions are needed to sort
172     * them before performing the product.
173     */
174     SparseQRMatrixQReturnType<SparseQR> matrixQ() const
175     { return SparseQRMatrixQReturnType<SparseQR>(*this); }
176 
177     /** \returns a const reference to the column permutation P that was applied to A such that A*P = Q*R
178       * It is the combination of the fill-in reducing permutation and numerical column pivoting.
179       */
180     const PermutationType& colsPermutation() const
181     {
182       eigen_assert(m_isInitialized && "Decomposition is not initialized.");
183       return m_outputPerm_c;
184     }
185 
186     /** \returns A string describing the type of error.
187       * This method is provided to ease debugging, not to handle errors.
188       */
189     std::string lastErrorMessage() const { return m_lastError; }
190 
191     /** \internal */
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 
218     /** Sets the threshold that is used to determine linearly dependent columns during the factorization.
219       *
220       * In practice, if during the factorization the norm of the column that has to be eliminated is below
221       * this threshold, then the entire column is treated as zero, and it is moved at the end.
222       */
223     void setPivotThreshold(const RealScalar& threshold)
224     {
225       m_useDefaultThreshold = false;
226       m_threshold = threshold;
227     }
228 
229     /** \returns the solution X of \f$ A X = B \f$ using the current decomposition of A.
230       *
231       * \sa compute()
232       */
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>
241     inline const Solve<SparseQR, Rhs> solve(const SparseMatrixBase<Rhs>& B) const
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 
248     /** \brief Reports whether previous computation was successful.
249       *
250       * \returns \c Success if computation was successful,
251       *          \c NumericalIssue if the QR factorization reports a numerical problem
252       *          \c InvalidInput if the input matrix is invalid
253       *
254       * \sa iparm()
255       */
256     ComputationInfo info() const
257     {
258       eigen_assert(m_isInitialized && "Decomposition is not initialized.");
259       return m_info;
260     }
261 
262 
263     /** \internal */
264     inline void _sort_matrix_Q()
265     {
266       if(this->m_isQSorted) return;
267       // The matrix Q is sorted during the transposition
268       SparseMatrix<Scalar, RowMajor, Index> mQrm(this->m_Q);
269       this->m_Q = mQrm;
270       this->m_isQSorted = true;
271     }
272 
273 
274   protected:
275     bool m_analysisIsok;
276     bool m_factorizationIsok;
277     mutable ComputationInfo m_info;
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 
298 /** \brief Preprocessing step of a QR factorization
299   *
300   * \warning The matrix \a mat must be in compressed mode (see SparseMatrix::makeCompressed()).
301   *
302   * In this step, the fill-reducing permutation is computed and applied to the columns of A
303   * and the column elimination tree is computed as well. Only the sparsity pattern of \a mat is exploited.
304   *
305   * \note In this step it is assumed that there is no empty row in the matrix \a mat.
306   */
307 template <typename MatrixType, typename OrderingType>
308 void SparseQR<MatrixType,OrderingType>::analyzePattern(const MatrixType& mat)
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
312   typename internal::conditional<MatrixType::IsRowMajor,QRMatrixType,const MatrixType&>::type matCpy(mat);
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 
341 /** \brief Performs the numerical QR factorization of the input matrix
342   *
343   * The function SparseQR::analyzePattern(const MatrixType&) must have been called beforehand with
344   * a matrix having the same sparsity pattern than \a mat.
345   *
346   * \param mat The sparse column-major matrix
347   */
348 template <typename MatrixType, typename OrderingType>
349 void SparseQR<MatrixType,OrderingType>::factorize(const MatrixType& mat)
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;
664   typedef Matrix<Scalar,Dynamic,Dynamic> DenseMatrix;
665   enum {
666     RowsAtCompileTime = Dynamic,
667     ColsAtCompileTime = Dynamic
668   };
669   explicit SparseQRMatrixQReturnType(const SparseQRType& qr) : m_qr(qr) {}
670   template<typename Derived>
671   SparseQR_QProduct<SparseQRType, Derived> operator*(const MatrixBase<Derived>& other)
672   {
673     return SparseQR_QProduct<SparseQRType,Derived>(m_qr,other.derived(),false);
674   }
675   // To use for operations with the adjoint of Q
676   SparseQRMatrixQTransposeReturnType<SparseQRType> adjoint() const
677   {
678     return SparseQRMatrixQTransposeReturnType<SparseQRType>(m_qr);
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
683   SparseQRMatrixQTransposeReturnType<SparseQRType> transpose() const
684   {
685     return SparseQRMatrixQTransposeReturnType<SparseQRType>(m_qr);
686   }
687   const SparseQRType& m_qr;
688 };
689 
690 // TODO this actually represents the adjoint of Q
691 template<typename SparseQRType>
692 struct SparseQRMatrixQTransposeReturnType
693 {
694   explicit SparseQRMatrixQTransposeReturnType(const SparseQRType& qr) : m_qr(qr) {}
695   template<typename Derived>
696   SparseQR_QProduct<SparseQRType,Derived> operator*(const MatrixBase<Derived>& other)
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>
706 struct evaluator_traits<SparseQRMatrixQReturnType<SparseQRType> >
707 {
708   typedef typename SparseQRType::MatrixType MatrixType;
709   typedef typename storage_kind_to_evaluator_kind<typename MatrixType::StorageKind>::Kind Kind;
710   typedef SparseShape Shape;
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 {
716   typedef SparseQRMatrixQReturnType<SparseQRType> SrcXprType;
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 {
732   typedef SparseQRMatrixQReturnType<SparseQRType> SrcXprType;
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
746