$darkmode
Eigen  5.0.1-dev
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 // IWYU pragma: private
14 #include "./InternalHeaderCheck.h"
15 
16 namespace Eigen {
17 
18 namespace internal {
19 template <typename MatrixType_, typename PermutationIndex_>
20 struct traits<FullPivLU<MatrixType_, PermutationIndex_> > : traits<MatrixType_> {
21  typedef MatrixXpr XprKind;
22  typedef SolverStorage StorageKind;
23  typedef PermutationIndex_ StorageIndex;
24  enum { Flags = 0 };
25 };
26 
27 } // end namespace internal
28 
62 template <typename MatrixType_, typename PermutationIndex_>
63 class FullPivLU : public SolverBase<FullPivLU<MatrixType_, PermutationIndex_> > {
64  public:
65  typedef MatrixType_ MatrixType;
66  typedef SolverBase<FullPivLU> Base;
67  friend class SolverBase<FullPivLU>;
68 
69  EIGEN_GENERIC_PUBLIC_INTERFACE(FullPivLU)
70  enum {
71  MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
72  MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
73  };
74  using PermutationIndex = PermutationIndex_;
75  typedef typename internal::plain_row_type<MatrixType, PermutationIndex>::type IntRowVectorType;
76  typedef typename internal::plain_col_type<MatrixType, PermutationIndex>::type IntColVectorType;
77  typedef PermutationMatrix<ColsAtCompileTime, MaxColsAtCompileTime, PermutationIndex> PermutationQType;
78  typedef PermutationMatrix<RowsAtCompileTime, MaxRowsAtCompileTime, PermutationIndex> PermutationPType;
79  typedef typename MatrixType::PlainObject PlainObject;
80 
88  eigen_assert(m_isInitialized && "FullPivLU is not initialized.");
89  return Success;
90  }
91 
98  FullPivLU();
99 
106  FullPivLU(Index rows, Index cols);
107 
113  template <typename InputType>
114  explicit FullPivLU(const EigenBase<InputType>& matrix);
115 
123  template <typename InputType>
124  explicit FullPivLU(EigenBase<InputType>& matrix);
125 
133  template <typename InputType>
135  m_lu = matrix.derived();
136  computeInPlace();
137  return *this;
138  }
139 
146  inline const MatrixType& matrixLU() const {
147  eigen_assert(m_isInitialized && "LU is not initialized.");
148  return m_lu;
149  }
150 
158  inline Index nonzeroPivots() const {
159  eigen_assert(m_isInitialized && "LU is not initialized.");
160  return m_nonzero_pivots;
161  }
162 
166  RealScalar maxPivot() const { return m_maxpivot; }
167 
172  EIGEN_DEVICE_FUNC inline const PermutationPType& permutationP() const {
173  eigen_assert(m_isInitialized && "LU is not initialized.");
174  return m_p;
175  }
176 
181  inline const PermutationQType& permutationQ() const {
182  eigen_assert(m_isInitialized && "LU is not initialized.");
183  return m_q;
184  }
185 
200  inline const internal::kernel_retval<FullPivLU> kernel() const {
201  eigen_assert(m_isInitialized && "LU is not initialized.");
202  return internal::kernel_retval<FullPivLU>(*this);
203  }
204 
224  inline const internal::image_retval<FullPivLU> image(const MatrixType& originalMatrix) const {
225  eigen_assert(m_isInitialized && "LU is not initialized.");
226  return internal::image_retval<FullPivLU>(*this, originalMatrix);
227  }
228 
229 #ifdef EIGEN_PARSED_BY_DOXYGEN
230 
249  template <typename Rhs>
250  inline const Solve<FullPivLU, Rhs> solve(const MatrixBase<Rhs>& b) const;
251 #endif
252 
256  inline RealScalar rcond() const {
257  eigen_assert(m_isInitialized && "FullPivLU is not initialized.");
258  if (!isInvertible()) {
259  return RealScalar(0);
260  }
261  return internal::rcond_estimate_helper(m_l1_norm, *this);
262  }
263 
279  typename internal::traits<MatrixType>::Scalar determinant() const;
280 
298  FullPivLU& setThreshold(const RealScalar& threshold) {
299  m_usePrescribedThreshold = true;
300  m_prescribedThreshold = threshold;
301  return *this;
302  }
303 
312  FullPivLU& setThreshold(Default_t) {
313  m_usePrescribedThreshold = false;
314  return *this;
315  }
316 
321  RealScalar threshold() const {
322  eigen_assert(m_isInitialized || m_usePrescribedThreshold);
323  return m_usePrescribedThreshold ? m_prescribedThreshold
324  // this formula comes from experimenting (see "LU precision tuning" thread on the
325  // list) and turns out to be identical to Higham's formula used already in LDLt.
326  : NumTraits<Scalar>::epsilon() * RealScalar(m_lu.diagonalSize());
327  }
328 
335  inline Index rank() const {
336  using std::abs;
337  eigen_assert(m_isInitialized && "LU is not initialized.");
338  RealScalar premultiplied_threshold = abs(m_maxpivot) * threshold();
339  Index result = 0;
340  for (Index i = 0; i < m_nonzero_pivots; ++i) result += (abs(m_lu.coeff(i, i)) > premultiplied_threshold);
341  return result;
342  }
343 
350  inline Index dimensionOfKernel() const {
351  eigen_assert(m_isInitialized && "LU is not initialized.");
352  return cols() - rank();
353  }
354 
362  inline bool isInjective() const {
363  eigen_assert(m_isInitialized && "LU is not initialized.");
364  return rank() == cols();
365  }
366 
374  inline bool isSurjective() const {
375  eigen_assert(m_isInitialized && "LU is not initialized.");
376  return rank() == rows();
377  }
378 
385  inline bool isInvertible() const {
386  eigen_assert(m_isInitialized && "LU is not initialized.");
387  return isInjective() && (m_lu.rows() == m_lu.cols());
388  }
389 
397  inline const Inverse<FullPivLU> inverse() const {
398  eigen_assert(m_isInitialized && "LU is not initialized.");
399  eigen_assert(m_lu.rows() == m_lu.cols() && "You can't take the inverse of a non-square matrix!");
400  return Inverse<FullPivLU>(*this);
401  }
402 
403  MatrixType reconstructedMatrix() const;
404 
405  EIGEN_DEVICE_FUNC constexpr Index rows() const noexcept { return m_lu.rows(); }
406  EIGEN_DEVICE_FUNC constexpr Index cols() const noexcept { return m_lu.cols(); }
407 
408 #ifndef EIGEN_PARSED_BY_DOXYGEN
409  template <typename RhsType, typename DstType>
410  void _solve_impl(const RhsType& rhs, DstType& dst) const;
411 
412  template <bool Conjugate, typename RhsType, typename DstType>
413  void _solve_impl_transposed(const RhsType& rhs, DstType& dst) const;
414 #endif
415 
416  protected:
417  EIGEN_STATIC_ASSERT_NON_INTEGER(Scalar)
418 
419  void computeInPlace();
420 
421  MatrixType m_lu;
422  PermutationPType m_p;
423  PermutationQType m_q;
424  IntColVectorType m_rowsTranspositions;
425  IntRowVectorType m_colsTranspositions;
426  Index m_nonzero_pivots;
427  RealScalar m_l1_norm;
428  RealScalar m_maxpivot, m_prescribedThreshold;
429  signed char m_det_pq;
430  bool m_isInitialized, m_usePrescribedThreshold;
431 };
432 
433 template <typename MatrixType, typename PermutationIndex>
434 FullPivLU<MatrixType, PermutationIndex>::FullPivLU() : m_isInitialized(false), m_usePrescribedThreshold(false) {}
435 
436 template <typename MatrixType, typename PermutationIndex>
438  : m_lu(rows, cols),
439  m_p(rows),
440  m_q(cols),
441  m_rowsTranspositions(rows),
442  m_colsTranspositions(cols),
443  m_isInitialized(false),
444  m_usePrescribedThreshold(false) {}
445 
446 template <typename MatrixType, typename PermutationIndex>
447 template <typename InputType>
449  : m_lu(matrix.rows(), matrix.cols()),
450  m_p(matrix.rows()),
451  m_q(matrix.cols()),
452  m_rowsTranspositions(matrix.rows()),
453  m_colsTranspositions(matrix.cols()),
454  m_isInitialized(false),
455  m_usePrescribedThreshold(false) {
456  compute(matrix.derived());
457 }
458 
459 template <typename MatrixType, typename PermutationIndex>
460 template <typename InputType>
462  : m_lu(matrix.derived()),
463  m_p(matrix.rows()),
464  m_q(matrix.cols()),
465  m_rowsTranspositions(matrix.rows()),
466  m_colsTranspositions(matrix.cols()),
467  m_isInitialized(false),
468  m_usePrescribedThreshold(false) {
469  computeInPlace();
470 }
471 
472 template <typename MatrixType, typename PermutationIndex>
474  eigen_assert(m_lu.rows() <= NumTraits<PermutationIndex>::highest() &&
475  m_lu.cols() <= NumTraits<PermutationIndex>::highest());
476 
477  m_l1_norm = m_lu.cwiseAbs().colwise().sum().maxCoeff();
478 
479  const Index size = m_lu.diagonalSize();
480  const Index rows = m_lu.rows();
481  const Index cols = m_lu.cols();
482 
483  // will store the transpositions, before we accumulate them at the end.
484  // can't accumulate on-the-fly because that will be done in reverse order for the rows.
485  m_rowsTranspositions.resize(m_lu.rows());
486  m_colsTranspositions.resize(m_lu.cols());
487  Index number_of_transpositions = 0; // number of NONTRIVIAL transpositions, i.e. m_rowsTranspositions[i]!=i
488 
489  m_nonzero_pivots = size; // the generic case is that in which all pivots are nonzero (invertible case)
490  m_maxpivot = RealScalar(0);
491 
492  for (Index k = 0; k < size; ++k) {
493  // First, we need to find the pivot.
494 
495  // biggest coefficient in the remaining bottom-right corner (starting at row k, col k)
496  Index row_of_biggest_in_corner, col_of_biggest_in_corner;
497  typedef internal::scalar_score_coeff_op<Scalar> Scoring;
498  typedef typename Scoring::result_type Score;
499  Score biggest_in_corner;
500  biggest_in_corner = m_lu.bottomRightCorner(rows - k, cols - k)
501  .unaryExpr(Scoring())
502  .maxCoeff(&row_of_biggest_in_corner, &col_of_biggest_in_corner);
503  row_of_biggest_in_corner += k; // correct the values! since they were computed in the corner,
504  col_of_biggest_in_corner += k; // need to add k to them.
505 
506  if (numext::is_exactly_zero(biggest_in_corner)) {
507  // before exiting, make sure to initialize the still uninitialized transpositions
508  // in a sane state without destroying what we already have.
509  m_nonzero_pivots = k;
510  for (Index i = k; i < size; ++i) {
511  m_rowsTranspositions.coeffRef(i) = internal::convert_index<StorageIndex>(i);
512  m_colsTranspositions.coeffRef(i) = internal::convert_index<StorageIndex>(i);
513  }
514  break;
515  }
516 
517  RealScalar abs_pivot = internal::abs_knowing_score<Scalar>()(
518  m_lu(row_of_biggest_in_corner, col_of_biggest_in_corner), biggest_in_corner);
519  if (abs_pivot > m_maxpivot) m_maxpivot = abs_pivot;
520 
521  // Now that we've found the pivot, we need to apply the row/col swaps to
522  // bring it to the location (k,k).
523 
524  m_rowsTranspositions.coeffRef(k) = internal::convert_index<StorageIndex>(row_of_biggest_in_corner);
525  m_colsTranspositions.coeffRef(k) = internal::convert_index<StorageIndex>(col_of_biggest_in_corner);
526  if (k != row_of_biggest_in_corner) {
527  m_lu.row(k).swap(m_lu.row(row_of_biggest_in_corner));
528  ++number_of_transpositions;
529  }
530  if (k != col_of_biggest_in_corner) {
531  m_lu.col(k).swap(m_lu.col(col_of_biggest_in_corner));
532  ++number_of_transpositions;
533  }
534 
535  // Now that the pivot is at the right location, we update the remaining
536  // bottom-right corner by Gaussian elimination.
537 
538  if (k < rows - 1) m_lu.col(k).tail(rows - k - 1) /= m_lu.coeff(k, k);
539  if (k < size - 1)
540  m_lu.block(k + 1, k + 1, rows - k - 1, cols - k - 1).noalias() -=
541  m_lu.col(k).tail(rows - k - 1) * m_lu.row(k).tail(cols - k - 1);
542  }
543 
544  // the main loop is over, we still have to accumulate the transpositions to find the
545  // permutations P and Q
546 
547  m_p.setIdentity(rows);
548  for (Index k = size - 1; k >= 0; --k) m_p.applyTranspositionOnTheRight(k, m_rowsTranspositions.coeff(k));
549 
550  m_q.setIdentity(cols);
551  for (Index k = 0; k < size; ++k) m_q.applyTranspositionOnTheRight(k, m_colsTranspositions.coeff(k));
552 
553  m_det_pq = (number_of_transpositions % 2) ? -1 : 1;
554 
555  m_isInitialized = true;
556 }
557 
558 template <typename MatrixType, typename PermutationIndex>
559 typename internal::traits<MatrixType>::Scalar FullPivLU<MatrixType, PermutationIndex>::determinant() const {
560  eigen_assert(m_isInitialized && "LU is not initialized.");
561  eigen_assert(m_lu.rows() == m_lu.cols() && "You can't take the determinant of a non-square matrix!");
562  return Scalar(m_det_pq) * Scalar(m_lu.diagonal().prod());
563 }
564 
568 template <typename MatrixType, typename PermutationIndex>
570  eigen_assert(m_isInitialized && "LU is not initialized.");
571  const Index smalldim = (std::min)(m_lu.rows(), m_lu.cols());
572  // LU
573  MatrixType res(m_lu.rows(), m_lu.cols());
574  // FIXME the .toDenseMatrix() should not be needed...
575  res = m_lu.leftCols(smalldim).template triangularView<UnitLower>().toDenseMatrix() *
576  m_lu.topRows(smalldim).template triangularView<Upper>().toDenseMatrix();
577 
578  // P^{-1}(LU)
579  res = m_p.inverse() * res;
580 
581  // (P^{-1}LU)Q^{-1}
582  res = res * m_q.inverse();
583 
584  return res;
585 }
586 
587 /********* Implementation of kernel() **************************************************/
588 
589 namespace internal {
590 template <typename MatrixType_, typename PermutationIndex_>
591 struct kernel_retval<FullPivLU<MatrixType_, PermutationIndex_> >
592  : kernel_retval_base<FullPivLU<MatrixType_, PermutationIndex_> > {
593  using DecompositionType = FullPivLU<MatrixType_, PermutationIndex_>;
594  EIGEN_MAKE_KERNEL_HELPERS(DecompositionType)
595 
596  enum {
597  MaxSmallDimAtCompileTime = min_size_prefer_fixed(MatrixType::MaxColsAtCompileTime, MatrixType::MaxRowsAtCompileTime)
598  };
599 
600  template <typename Dest>
601  void evalTo(Dest& dst) const {
602  using std::abs;
603  const Index cols = dec().matrixLU().cols(), dimker = cols - rank();
604  if (dimker == 0) {
605  // The Kernel is just {0}, so it doesn't have a basis properly speaking, but let's
606  // avoid crashing/asserting as that depends on floating point calculations. Let's
607  // just return a single column vector filled with zeros.
608  dst.setZero();
609  return;
610  }
611 
612  /* Let us use the following lemma:
613  *
614  * Lemma: If the matrix A has the LU decomposition PAQ = LU,
615  * then Ker A = Q(Ker U).
616  *
617  * Proof: trivial: just keep in mind that P, Q, L are invertible.
618  */
619 
620  /* Thus, all we need to do is to compute Ker U, and then apply Q.
621  *
622  * U is upper triangular, with eigenvalues sorted so that any zeros appear at the end.
623  * Thus, the diagonal of U ends with exactly
624  * dimKer zero's. Let us use that to construct dimKer linearly
625  * independent vectors in Ker U.
626  */
627 
628  Matrix<Index, Dynamic, 1, 0, MaxSmallDimAtCompileTime, 1> pivots(rank());
629  RealScalar premultiplied_threshold = dec().maxPivot() * dec().threshold();
630  Index p = 0;
631  for (Index i = 0; i < dec().nonzeroPivots(); ++i)
632  if (abs(dec().matrixLU().coeff(i, i)) > premultiplied_threshold) pivots.coeffRef(p++) = i;
633  eigen_internal_assert(p == rank());
634 
635  // we construct a temporaty trapezoid matrix m, by taking the U matrix and
636  // permuting the rows and cols to bring the nonnegligible pivots to the top of
637  // the main diagonal. We need that to be able to apply our triangular solvers.
638  // FIXME when we get triangularView-for-rectangular-matrices, this can be simplified
639  Matrix<typename MatrixType::Scalar, Dynamic, Dynamic, traits<MatrixType>::Options, MaxSmallDimAtCompileTime,
640  MatrixType::MaxColsAtCompileTime>
641  m(dec().matrixLU().block(0, 0, rank(), cols));
642  for (Index i = 0; i < rank(); ++i) {
643  if (i) m.row(i).head(i).setZero();
644  m.row(i).tail(cols - i) = dec().matrixLU().row(pivots.coeff(i)).tail(cols - i);
645  }
646  m.block(0, 0, rank(), rank());
647  m.block(0, 0, rank(), rank()).template triangularView<StrictlyLower>().setZero();
648  for (Index i = 0; i < rank(); ++i) m.col(i).swap(m.col(pivots.coeff(i)));
649 
650  // ok, we have our trapezoid matrix, we can apply the triangular solver.
651  // notice that the math behind this suggests that we should apply this to the
652  // negative of the RHS, but for performance we just put the negative sign elsewhere, see below.
653  m.topLeftCorner(rank(), rank()).template triangularView<Upper>().solveInPlace(m.topRightCorner(rank(), dimker));
654 
655  // now we must undo the column permutation that we had applied!
656  for (Index i = rank() - 1; i >= 0; --i) m.col(i).swap(m.col(pivots.coeff(i)));
657 
658  // see the negative sign in the next line, that's what we were talking about above.
659  for (Index i = 0; i < rank(); ++i) dst.row(dec().permutationQ().indices().coeff(i)) = -m.row(i).tail(dimker);
660  for (Index i = rank(); i < cols; ++i) dst.row(dec().permutationQ().indices().coeff(i)).setZero();
661  for (Index k = 0; k < dimker; ++k) dst.coeffRef(dec().permutationQ().indices().coeff(rank() + k), k) = Scalar(1);
662  }
663 };
664 
665 /***** Implementation of image() *****************************************************/
666 
667 template <typename MatrixType_, typename PermutationIndex_>
668 struct image_retval<FullPivLU<MatrixType_, PermutationIndex_> >
669  : image_retval_base<FullPivLU<MatrixType_, PermutationIndex_> > {
670  using DecompositionType = FullPivLU<MatrixType_, PermutationIndex_>;
671  EIGEN_MAKE_IMAGE_HELPERS(DecompositionType)
672 
673  enum {
674  MaxSmallDimAtCompileTime = min_size_prefer_fixed(MatrixType::MaxColsAtCompileTime, MatrixType::MaxRowsAtCompileTime)
675  };
676 
677  template <typename Dest>
678  void evalTo(Dest& dst) const {
679  using std::abs;
680  if (rank() == 0) {
681  // The Image is just {0}, so it doesn't have a basis properly speaking, but let's
682  // avoid crashing/asserting as that depends on floating point calculations. Let's
683  // just return a single column vector filled with zeros.
684  dst.setZero();
685  return;
686  }
687 
688  Matrix<Index, Dynamic, 1, 0, MaxSmallDimAtCompileTime, 1> pivots(rank());
689  RealScalar premultiplied_threshold = dec().maxPivot() * dec().threshold();
690  Index p = 0;
691  for (Index i = 0; i < dec().nonzeroPivots(); ++i)
692  if (abs(dec().matrixLU().coeff(i, i)) > premultiplied_threshold) pivots.coeffRef(p++) = i;
693  eigen_internal_assert(p == rank());
694 
695  for (Index i = 0; i < rank(); ++i)
696  dst.col(i) = originalMatrix().col(dec().permutationQ().indices().coeff(pivots.coeff(i)));
697  }
698 };
699 
700 /***** Implementation of solve() *****************************************************/
701 
702 } // end namespace internal
703 
704 #ifndef EIGEN_PARSED_BY_DOXYGEN
705 template <typename MatrixType_, typename PermutationIndex_>
706 template <typename RhsType, typename DstType>
707 void FullPivLU<MatrixType_, PermutationIndex_>::_solve_impl(const RhsType& rhs, DstType& dst) const {
708  /* The decomposition PAQ = LU can be rewritten as A = P^{-1} L U Q^{-1}.
709  * So we proceed as follows:
710  * Step 1: compute c = P * rhs.
711  * Step 2: replace c by the solution x to Lx = c. Exists because L is invertible.
712  * Step 3: replace c by the solution x to Ux = c. May or may not exist.
713  * Step 4: result = Q * c;
714  */
715 
716  const Index rows = this->rows(), cols = this->cols(), nonzero_pivots = this->rank();
717  const Index smalldim = (std::min)(rows, cols);
718 
719  if (nonzero_pivots == 0) {
720  dst.setZero();
721  return;
722  }
723 
724  typename RhsType::PlainObject c(rhs.rows(), rhs.cols());
725 
726  // Step 1
727  c = permutationP() * rhs;
728 
729  // Step 2
730  m_lu.topLeftCorner(smalldim, smalldim).template triangularView<UnitLower>().solveInPlace(c.topRows(smalldim));
731  if (rows > cols) c.bottomRows(rows - cols).noalias() -= m_lu.bottomRows(rows - cols) * c.topRows(cols);
732 
733  // Step 3
734  m_lu.topLeftCorner(nonzero_pivots, nonzero_pivots)
735  .template triangularView<Upper>()
736  .solveInPlace(c.topRows(nonzero_pivots));
737 
738  // Step 4
739  for (Index i = 0; i < nonzero_pivots; ++i) dst.row(permutationQ().indices().coeff(i)) = c.row(i);
740  for (Index i = nonzero_pivots; i < m_lu.cols(); ++i) dst.row(permutationQ().indices().coeff(i)).setZero();
741 }
742 
743 template <typename MatrixType_, typename PermutationIndex_>
744 template <bool Conjugate, typename RhsType, typename DstType>
745 void FullPivLU<MatrixType_, PermutationIndex_>::_solve_impl_transposed(const RhsType& rhs, DstType& dst) const {
746  /* The decomposition PAQ = LU can be rewritten as A = P^{-1} L U Q^{-1},
747  * and since permutations are real and unitary, we can write this
748  * as A^T = Q U^T L^T P,
749  * So we proceed as follows:
750  * Step 1: compute c = Q^T rhs.
751  * Step 2: replace c by the solution x to U^T x = c. May or may not exist.
752  * Step 3: replace c by the solution x to L^T x = c.
753  * Step 4: result = P^T c.
754  * If Conjugate is true, replace "^T" by "^*" above.
755  */
756 
757  const Index rows = this->rows(), cols = this->cols(), nonzero_pivots = this->rank();
758  const Index smalldim = (std::min)(rows, cols);
759 
760  if (nonzero_pivots == 0) {
761  dst.setZero();
762  return;
763  }
764 
765  typename RhsType::PlainObject c(rhs.rows(), rhs.cols());
766 
767  // Step 1
768  c = permutationQ().inverse() * rhs;
769 
770  // Step 2
771  m_lu.topLeftCorner(nonzero_pivots, nonzero_pivots)
772  .template triangularView<Upper>()
773  .transpose()
774  .template conjugateIf<Conjugate>()
775  .solveInPlace(c.topRows(nonzero_pivots));
776 
777  // Step 3
778  m_lu.topLeftCorner(smalldim, smalldim)
779  .template triangularView<UnitLower>()
780  .transpose()
781  .template conjugateIf<Conjugate>()
782  .solveInPlace(c.topRows(smalldim));
783 
784  // Step 4
785  PermutationPType invp = permutationP().inverse().eval();
786  for (Index i = 0; i < smalldim; ++i) dst.row(invp.indices().coeff(i)) = c.row(i);
787  for (Index i = smalldim; i < rows; ++i) dst.row(invp.indices().coeff(i)).setZero();
788 }
789 
790 #endif
791 
792 namespace internal {
793 
794 /***** Implementation of inverse() *****************************************************/
795 template <typename DstXprType, typename MatrixType, typename PermutationIndex>
796 struct Assignment<
797  DstXprType, Inverse<FullPivLU<MatrixType, PermutationIndex> >,
798  internal::assign_op<typename DstXprType::Scalar, typename FullPivLU<MatrixType, PermutationIndex>::Scalar>,
799  Dense2Dense> {
800  typedef FullPivLU<MatrixType, PermutationIndex> LuType;
801  typedef Inverse<LuType> SrcXprType;
802  static void run(DstXprType& dst, const SrcXprType& src,
803  const internal::assign_op<typename DstXprType::Scalar, typename MatrixType::Scalar>&) {
804  dst = src.nestedExpression().solve(MatrixType::Identity(src.rows(), src.cols()));
805  }
806 };
807 } // end namespace internal
808 
809 /******* MatrixBase methods *****************************************************************/
810 
817 template <typename Derived>
818 template <typename PermutationIndex>
820  const {
822 }
823 
824 } // end namespace Eigen
825 
826 #endif // EIGEN_LU_H
RealScalar rcond() const
Definition: FullPivLU.h:256
FullPivLU & setThreshold(Default_t)
Definition: FullPivLU.h:312
const Inverse< FullPivLU > inverse() const
Definition: FullPivLU.h:397
const MatrixType & matrixLU() const
Definition: FullPivLU.h:146
FullPivLU & compute(const EigenBase< InputType > &matrix)
Definition: FullPivLU.h:134
MatrixType reconstructedMatrix() const
Definition: FullPivLU.h:569
Namespace containing all symbols from the Eigen library.
Definition: B01_Experimental.dox:1
ComputationInfo info() const
Reports whether the LU factorization was successful.
Definition: FullPivLU.h:87
Holds information about the various numeric (i.e. scalar) types allowed by Eigen. ...
Definition: NumTraits.h:232
const PermutationPType & permutationP() const
Definition: FullPivLU.h:172
Eigen::Index Index
The interface type of indices.
Definition: EigenBase.h:43
Index dimensionOfKernel() const
Definition: FullPivLU.h:350
bool isInjective() const
Definition: FullPivLU.h:362
Definition: EigenBase.h:33
bool isInvertible() const
Definition: FullPivLU.h:385
Expression of the inverse of another expression.
Definition: Inverse.h:43
Index rank() const
Definition: FullPivLU.h:335
const internal::image_retval< FullPivLU > image(const MatrixType &originalMatrix) const
Definition: FullPivLU.h:224
RealScalar threshold() const
Definition: FullPivLU.h:321
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:82
const PermutationQType & permutationQ() const
Definition: FullPivLU.h:181
FullPivLU & setThreshold(const RealScalar &threshold)
Definition: FullPivLU.h:298
Definition: Constants.h:440
constexpr Derived & derived()
Definition: EigenBase.h:49
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_abs_op< typename Derived::Scalar >, const Derived > abs(const Eigen::ArrayBase< Derived > &x)
FullPivLU()
Default Constructor.
Definition: FullPivLU.h:434
Index nonzeroPivots() const
Definition: FullPivLU.h:158
LU decomposition of a matrix with complete pivoting, and related features.
Definition: ForwardDeclarations.h:419
RealScalar maxPivot() const
Definition: FullPivLU.h:166
const Solve< FullPivLU, Rhs > solve(const MatrixBase< Rhs > &b) const
internal::traits< MatrixType >::Scalar determinant() const
Definition: FullPivLU.h:559
Pseudo expression representing a solving operation.
Definition: Solve.h:62
ComputationInfo
Definition: Constants.h:438
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:52
bool isSurjective() const
Definition: FullPivLU.h:374
const internal::kernel_retval< FullPivLU > kernel() const
Definition: FullPivLU.h:200