FullPivLU.h
1 // This file is part of Eigen, a lightweight C++ template library
2 // for linear algebra.
3 //
4 // Copyright (C) 2006-2009 Benoit Jacob <jacob.benoit.1@gmail.com>
5 //
6 // This Source Code Form is subject to the terms of the Mozilla
7 // Public License v. 2.0. If a copy of the MPL was not distributed
8 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
9 
10 #ifndef EIGEN_LU_H
11 #define EIGEN_LU_H
12 
13 namespace Eigen {
14 
45 template<typename _MatrixType> class FullPivLU
46 {
47  public:
48  typedef _MatrixType MatrixType;
49  enum {
50  RowsAtCompileTime = MatrixType::RowsAtCompileTime,
51  ColsAtCompileTime = MatrixType::ColsAtCompileTime,
52  Options = MatrixType::Options,
53  MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
54  MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
55  };
56  typedef typename MatrixType::Scalar Scalar;
57  typedef typename NumTraits<typename MatrixType::Scalar>::Real RealScalar;
58  typedef typename internal::traits<MatrixType>::StorageKind StorageKind;
59  typedef typename MatrixType::Index Index;
60  typedef typename internal::plain_row_type<MatrixType, Index>::type IntRowVectorType;
61  typedef typename internal::plain_col_type<MatrixType, Index>::type IntColVectorType;
64 
71  FullPivLU();
72 
79  FullPivLU(Index rows, Index cols);
80 
86  FullPivLU(const MatrixType& matrix);
87 
95  FullPivLU& compute(const MatrixType& matrix);
96 
103  inline const MatrixType& matrixLU() const
104  {
105  eigen_assert(m_isInitialized && "LU is not initialized.");
106  return m_lu;
107  }
108 
116  inline Index nonzeroPivots() const
117  {
118  eigen_assert(m_isInitialized && "LU is not initialized.");
119  return m_nonzero_pivots;
120  }
121 
125  RealScalar maxPivot() const { return m_maxpivot; }
126 
131  inline const PermutationPType& permutationP() const
132  {
133  eigen_assert(m_isInitialized && "LU is not initialized.");
134  return m_p;
135  }
136 
141  inline const PermutationQType& permutationQ() const
142  {
143  eigen_assert(m_isInitialized && "LU is not initialized.");
144  return m_q;
145  }
146 
161  inline const internal::kernel_retval<FullPivLU> kernel() const
162  {
163  eigen_assert(m_isInitialized && "LU is not initialized.");
164  return internal::kernel_retval<FullPivLU>(*this);
165  }
166 
186  inline const internal::image_retval<FullPivLU>
187  image(const MatrixType& originalMatrix) const
188  {
189  eigen_assert(m_isInitialized && "LU is not initialized.");
190  return internal::image_retval<FullPivLU>(*this, originalMatrix);
191  }
192 
212  template<typename Rhs>
213  inline const internal::solve_retval<FullPivLU, Rhs>
214  solve(const MatrixBase<Rhs>& b) const
215  {
216  eigen_assert(m_isInitialized && "LU is not initialized.");
217  return internal::solve_retval<FullPivLU, Rhs>(*this, b.derived());
218  }
219 
235  typename internal::traits<MatrixType>::Scalar determinant() const;
236 
254  FullPivLU& setThreshold(const RealScalar& threshold)
255  {
256  m_usePrescribedThreshold = true;
257  m_prescribedThreshold = threshold;
258  return *this;
259  }
260 
270  {
271  m_usePrescribedThreshold = false;
272  return *this;
273  }
274 
279  RealScalar threshold() const
280  {
281  eigen_assert(m_isInitialized || m_usePrescribedThreshold);
282  return m_usePrescribedThreshold ? m_prescribedThreshold
283  // this formula comes from experimenting (see "LU precision tuning" thread on the list)
284  // and turns out to be identical to Higham's formula used already in LDLt.
285  : NumTraits<Scalar>::epsilon() * m_lu.diagonalSize();
286  }
287 
294  inline Index rank() const
295  {
296  eigen_assert(m_isInitialized && "LU is not initialized.");
297  RealScalar premultiplied_threshold = internal::abs(m_maxpivot) * threshold();
298  Index result = 0;
299  for(Index i = 0; i < m_nonzero_pivots; ++i)
300  result += (internal::abs(m_lu.coeff(i,i)) > premultiplied_threshold);
301  return result;
302  }
303 
310  inline Index dimensionOfKernel() const
311  {
312  eigen_assert(m_isInitialized && "LU is not initialized.");
313  return cols() - rank();
314  }
315 
323  inline bool isInjective() const
324  {
325  eigen_assert(m_isInitialized && "LU is not initialized.");
326  return rank() == cols();
327  }
328 
336  inline bool isSurjective() const
337  {
338  eigen_assert(m_isInitialized && "LU is not initialized.");
339  return rank() == rows();
340  }
341 
348  inline bool isInvertible() const
349  {
350  eigen_assert(m_isInitialized && "LU is not initialized.");
351  return isInjective() && (m_lu.rows() == m_lu.cols());
352  }
353 
361  inline const internal::solve_retval<FullPivLU,typename MatrixType::IdentityReturnType> inverse() const
362  {
363  eigen_assert(m_isInitialized && "LU is not initialized.");
364  eigen_assert(m_lu.rows() == m_lu.cols() && "You can't take the inverse of a non-square matrix!");
365  return internal::solve_retval<FullPivLU,typename MatrixType::IdentityReturnType>
366  (*this, MatrixType::Identity(m_lu.rows(), m_lu.cols()));
367  }
368 
369  MatrixType reconstructedMatrix() const;
370 
371  inline Index rows() const { return m_lu.rows(); }
372  inline Index cols() const { return m_lu.cols(); }
373 
374  protected:
375  MatrixType m_lu;
376  PermutationPType m_p;
377  PermutationQType m_q;
378  IntColVectorType m_rowsTranspositions;
379  IntRowVectorType m_colsTranspositions;
380  Index m_det_pq, m_nonzero_pivots;
381  RealScalar m_maxpivot, m_prescribedThreshold;
382  bool m_isInitialized, m_usePrescribedThreshold;
383 };
384 
385 template<typename MatrixType>
387  : m_isInitialized(false), m_usePrescribedThreshold(false)
388 {
389 }
390 
391 template<typename MatrixType>
392 FullPivLU<MatrixType>::FullPivLU(Index rows, Index cols)
393  : m_lu(rows, cols),
394  m_p(rows),
395  m_q(cols),
396  m_rowsTranspositions(rows),
397  m_colsTranspositions(cols),
398  m_isInitialized(false),
399  m_usePrescribedThreshold(false)
400 {
401 }
402 
403 template<typename MatrixType>
404 FullPivLU<MatrixType>::FullPivLU(const MatrixType& matrix)
405  : m_lu(matrix.rows(), matrix.cols()),
406  m_p(matrix.rows()),
407  m_q(matrix.cols()),
408  m_rowsTranspositions(matrix.rows()),
409  m_colsTranspositions(matrix.cols()),
410  m_isInitialized(false),
411  m_usePrescribedThreshold(false)
412 {
413  compute(matrix);
414 }
415 
416 template<typename MatrixType>
418 {
419  m_isInitialized = true;
420  m_lu = matrix;
421 
422  const Index size = matrix.diagonalSize();
423  const Index rows = matrix.rows();
424  const Index cols = matrix.cols();
425 
426  // will store the transpositions, before we accumulate them at the end.
427  // can't accumulate on-the-fly because that will be done in reverse order for the rows.
428  m_rowsTranspositions.resize(matrix.rows());
429  m_colsTranspositions.resize(matrix.cols());
430  Index number_of_transpositions = 0; // number of NONTRIVIAL transpositions, i.e. m_rowsTranspositions[i]!=i
431 
432  m_nonzero_pivots = size; // the generic case is that in which all pivots are nonzero (invertible case)
433  m_maxpivot = RealScalar(0);
434 
435  for(Index k = 0; k < size; ++k)
436  {
437  // First, we need to find the pivot.
438 
439  // biggest coefficient in the remaining bottom-right corner (starting at row k, col k)
440  Index row_of_biggest_in_corner, col_of_biggest_in_corner;
441  RealScalar biggest_in_corner;
442  biggest_in_corner = m_lu.bottomRightCorner(rows-k, cols-k)
443  .cwiseAbs()
444  .maxCoeff(&row_of_biggest_in_corner, &col_of_biggest_in_corner);
445  row_of_biggest_in_corner += k; // correct the values! since they were computed in the corner,
446  col_of_biggest_in_corner += k; // need to add k to them.
447 
448  if(biggest_in_corner==RealScalar(0))
449  {
450  // before exiting, make sure to initialize the still uninitialized transpositions
451  // in a sane state without destroying what we already have.
452  m_nonzero_pivots = k;
453  for(Index i = k; i < size; ++i)
454  {
455  m_rowsTranspositions.coeffRef(i) = i;
456  m_colsTranspositions.coeffRef(i) = i;
457  }
458  break;
459  }
460 
461  if(biggest_in_corner > m_maxpivot) m_maxpivot = biggest_in_corner;
462 
463  // Now that we've found the pivot, we need to apply the row/col swaps to
464  // bring it to the location (k,k).
465 
466  m_rowsTranspositions.coeffRef(k) = row_of_biggest_in_corner;
467  m_colsTranspositions.coeffRef(k) = col_of_biggest_in_corner;
468  if(k != row_of_biggest_in_corner) {
469  m_lu.row(k).swap(m_lu.row(row_of_biggest_in_corner));
470  ++number_of_transpositions;
471  }
472  if(k != col_of_biggest_in_corner) {
473  m_lu.col(k).swap(m_lu.col(col_of_biggest_in_corner));
474  ++number_of_transpositions;
475  }
476 
477  // Now that the pivot is at the right location, we update the remaining
478  // bottom-right corner by Gaussian elimination.
479 
480  if(k<rows-1)
481  m_lu.col(k).tail(rows-k-1) /= m_lu.coeff(k,k);
482  if(k<size-1)
483  m_lu.block(k+1,k+1,rows-k-1,cols-k-1).noalias() -= m_lu.col(k).tail(rows-k-1) * m_lu.row(k).tail(cols-k-1);
484  }
485 
486  // the main loop is over, we still have to accumulate the transpositions to find the
487  // permutations P and Q
488 
489  m_p.setIdentity(rows);
490  for(Index k = size-1; k >= 0; --k)
491  m_p.applyTranspositionOnTheRight(k, m_rowsTranspositions.coeff(k));
492 
493  m_q.setIdentity(cols);
494  for(Index k = 0; k < size; ++k)
495  m_q.applyTranspositionOnTheRight(k, m_colsTranspositions.coeff(k));
496 
497  m_det_pq = (number_of_transpositions%2) ? -1 : 1;
498  return *this;
499 }
500 
501 template<typename MatrixType>
502 typename internal::traits<MatrixType>::Scalar FullPivLU<MatrixType>::determinant() const
503 {
504  eigen_assert(m_isInitialized && "LU is not initialized.");
505  eigen_assert(m_lu.rows() == m_lu.cols() && "You can't take the determinant of a non-square matrix!");
506  return Scalar(m_det_pq) * Scalar(m_lu.diagonal().prod());
507 }
508 
512 template<typename MatrixType>
514 {
515  eigen_assert(m_isInitialized && "LU is not initialized.");
516  const Index smalldim = (std::min)(m_lu.rows(), m_lu.cols());
517  // LU
518  MatrixType res(m_lu.rows(),m_lu.cols());
519  // FIXME the .toDenseMatrix() should not be needed...
520  res = m_lu.leftCols(smalldim)
521  .template triangularView<UnitLower>().toDenseMatrix()
522  * m_lu.topRows(smalldim)
523  .template triangularView<Upper>().toDenseMatrix();
524 
525  // P^{-1}(LU)
526  res = m_p.inverse() * res;
527 
528  // (P^{-1}LU)Q^{-1}
529  res = res * m_q.inverse();
530 
531  return res;
532 }
533 
534 /********* Implementation of kernel() **************************************************/
535 
536 namespace internal {
537 template<typename _MatrixType>
538 struct kernel_retval<FullPivLU<_MatrixType> >
539  : kernel_retval_base<FullPivLU<_MatrixType> >
540 {
541  EIGEN_MAKE_KERNEL_HELPERS(FullPivLU<_MatrixType>)
542 
543  enum { MaxSmallDimAtCompileTime = EIGEN_SIZE_MIN_PREFER_FIXED(
544  MatrixType::MaxColsAtCompileTime,
545  MatrixType::MaxRowsAtCompileTime)
546  };
547 
548  template<typename Dest> void evalTo(Dest& dst) const
549  {
550  const Index cols = dec().matrixLU().cols(), dimker = cols - rank();
551  if(dimker == 0)
552  {
553  // The Kernel is just {0}, so it doesn't have a basis properly speaking, but let's
554  // avoid crashing/asserting as that depends on floating point calculations. Let's
555  // just return a single column vector filled with zeros.
556  dst.setZero();
557  return;
558  }
559 
560  /* Let us use the following lemma:
561  *
562  * Lemma: If the matrix A has the LU decomposition PAQ = LU,
563  * then Ker A = Q(Ker U).
564  *
565  * Proof: trivial: just keep in mind that P, Q, L are invertible.
566  */
567 
568  /* Thus, all we need to do is to compute Ker U, and then apply Q.
569  *
570  * U is upper triangular, with eigenvalues sorted so that any zeros appear at the end.
571  * Thus, the diagonal of U ends with exactly
572  * dimKer zero's. Let us use that to construct dimKer linearly
573  * independent vectors in Ker U.
574  */
575 
576  Matrix<Index, Dynamic, 1, 0, MaxSmallDimAtCompileTime, 1> pivots(rank());
577  RealScalar premultiplied_threshold = dec().maxPivot() * dec().threshold();
578  Index p = 0;
579  for(Index i = 0; i < dec().nonzeroPivots(); ++i)
580  if(abs(dec().matrixLU().coeff(i,i)) > premultiplied_threshold)
581  pivots.coeffRef(p++) = i;
582  eigen_internal_assert(p == rank());
583 
584  // we construct a temporaty trapezoid matrix m, by taking the U matrix and
585  // permuting the rows and cols to bring the nonnegligible pivots to the top of
586  // the main diagonal. We need that to be able to apply our triangular solvers.
587  // FIXME when we get triangularView-for-rectangular-matrices, this can be simplified
588  Matrix<typename MatrixType::Scalar, Dynamic, Dynamic, MatrixType::Options,
589  MaxSmallDimAtCompileTime, MatrixType::MaxColsAtCompileTime>
590  m(dec().matrixLU().block(0, 0, rank(), cols));
591  for(Index i = 0; i < rank(); ++i)
592  {
593  if(i) m.row(i).head(i).setZero();
594  m.row(i).tail(cols-i) = dec().matrixLU().row(pivots.coeff(i)).tail(cols-i);
595  }
596  m.block(0, 0, rank(), rank());
597  m.block(0, 0, rank(), rank()).template triangularView<StrictlyLower>().setZero();
598  for(Index i = 0; i < rank(); ++i)
599  m.col(i).swap(m.col(pivots.coeff(i)));
600 
601  // ok, we have our trapezoid matrix, we can apply the triangular solver.
602  // notice that the math behind this suggests that we should apply this to the
603  // negative of the RHS, but for performance we just put the negative sign elsewhere, see below.
604  m.topLeftCorner(rank(), rank())
605  .template triangularView<Upper>().solveInPlace(
606  m.topRightCorner(rank(), dimker)
607  );
608 
609  // now we must undo the column permutation that we had applied!
610  for(Index i = rank()-1; i >= 0; --i)
611  m.col(i).swap(m.col(pivots.coeff(i)));
612 
613  // see the negative sign in the next line, that's what we were talking about above.
614  for(Index i = 0; i < rank(); ++i) dst.row(dec().permutationQ().indices().coeff(i)) = -m.row(i).tail(dimker);
615  for(Index i = rank(); i < cols; ++i) dst.row(dec().permutationQ().indices().coeff(i)).setZero();
616  for(Index k = 0; k < dimker; ++k) dst.coeffRef(dec().permutationQ().indices().coeff(rank()+k), k) = Scalar(1);
617  }
618 };
619 
620 /***** Implementation of image() *****************************************************/
621 
622 template<typename _MatrixType>
623 struct image_retval<FullPivLU<_MatrixType> >
624  : image_retval_base<FullPivLU<_MatrixType> >
625 {
626  EIGEN_MAKE_IMAGE_HELPERS(FullPivLU<_MatrixType>)
627 
628  enum { MaxSmallDimAtCompileTime = EIGEN_SIZE_MIN_PREFER_FIXED(
629  MatrixType::MaxColsAtCompileTime,
630  MatrixType::MaxRowsAtCompileTime)
631  };
632 
633  template<typename Dest> void evalTo(Dest& dst) const
634  {
635  if(rank() == 0)
636  {
637  // The Image is just {0}, so it doesn't have a basis properly speaking, but let's
638  // avoid crashing/asserting as that depends on floating point calculations. Let's
639  // just return a single column vector filled with zeros.
640  dst.setZero();
641  return;
642  }
643 
644  Matrix<Index, Dynamic, 1, 0, MaxSmallDimAtCompileTime, 1> pivots(rank());
645  RealScalar premultiplied_threshold = dec().maxPivot() * dec().threshold();
646  Index p = 0;
647  for(Index i = 0; i < dec().nonzeroPivots(); ++i)
648  if(abs(dec().matrixLU().coeff(i,i)) > premultiplied_threshold)
649  pivots.coeffRef(p++) = i;
650  eigen_internal_assert(p == rank());
651 
652  for(Index i = 0; i < rank(); ++i)
653  dst.col(i) = originalMatrix().col(dec().permutationQ().indices().coeff(pivots.coeff(i)));
654  }
655 };
656 
657 /***** Implementation of solve() *****************************************************/
658 
659 template<typename _MatrixType, typename Rhs>
660 struct solve_retval<FullPivLU<_MatrixType>, Rhs>
661  : solve_retval_base<FullPivLU<_MatrixType>, Rhs>
662 {
663  EIGEN_MAKE_SOLVE_HELPERS(FullPivLU<_MatrixType>,Rhs)
664 
665  template<typename Dest> void evalTo(Dest& dst) const
666  {
667  /* The decomposition PAQ = LU can be rewritten as A = P^{-1} L U Q^{-1}.
668  * So we proceed as follows:
669  * Step 1: compute c = P * rhs.
670  * Step 2: replace c by the solution x to Lx = c. Exists because L is invertible.
671  * Step 3: replace c by the solution x to Ux = c. May or may not exist.
672  * Step 4: result = Q * c;
673  */
674 
675  const Index rows = dec().rows(), cols = dec().cols(),
676  nonzero_pivots = dec().nonzeroPivots();
677  eigen_assert(rhs().rows() == rows);
678  const Index smalldim = (std::min)(rows, cols);
679 
680  if(nonzero_pivots == 0)
681  {
682  dst.setZero();
683  return;
684  }
685 
686  typename Rhs::PlainObject c(rhs().rows(), rhs().cols());
687 
688  // Step 1
689  c = dec().permutationP() * rhs();
690 
691  // Step 2
692  dec().matrixLU()
693  .topLeftCorner(smalldim,smalldim)
694  .template triangularView<UnitLower>()
695  .solveInPlace(c.topRows(smalldim));
696  if(rows>cols)
697  {
698  c.bottomRows(rows-cols)
699  -= dec().matrixLU().bottomRows(rows-cols)
700  * c.topRows(cols);
701  }
702 
703  // Step 3
704  dec().matrixLU()
705  .topLeftCorner(nonzero_pivots, nonzero_pivots)
706  .template triangularView<Upper>()
707  .solveInPlace(c.topRows(nonzero_pivots));
708 
709  // Step 4
710  for(Index i = 0; i < nonzero_pivots; ++i)
711  dst.row(dec().permutationQ().indices().coeff(i)) = c.row(i);
712  for(Index i = nonzero_pivots; i < dec().matrixLU().cols(); ++i)
713  dst.row(dec().permutationQ().indices().coeff(i)).setZero();
714  }
715 };
716 
717 } // end namespace internal
718 
719 /******* MatrixBase methods *****************************************************************/
720 
727 template<typename Derived>
728 inline const FullPivLU<typename MatrixBase<Derived>::PlainObject>
730 {
731  return FullPivLU<PlainObject>(eval());
732 }
733 
734 } // end namespace Eigen
735 
736 #endif // EIGEN_LU_H