45 template<
typename _MatrixType>
class FullPivLU
48 typedef _MatrixType MatrixType;
50 RowsAtCompileTime = MatrixType::RowsAtCompileTime,
51 ColsAtCompileTime = MatrixType::ColsAtCompileTime,
52 Options = MatrixType::Options,
53 MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
54 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
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;
62 typedef PermutationMatrix<ColsAtCompileTime, MaxColsAtCompileTime> PermutationQType;
63 typedef PermutationMatrix<RowsAtCompileTime, MaxRowsAtCompileTime> PermutationPType;
105 eigen_assert(m_isInitialized &&
"LU is not initialized.");
118 eigen_assert(m_isInitialized &&
"LU is not initialized.");
119 return m_nonzero_pivots;
133 eigen_assert(m_isInitialized &&
"LU is not initialized.");
143 eigen_assert(m_isInitialized &&
"LU is not initialized.");
161 inline const internal::kernel_retval<FullPivLU>
kernel()
const
163 eigen_assert(m_isInitialized &&
"LU is not initialized.");
164 return internal::kernel_retval<FullPivLU>(*this);
186 inline const internal::image_retval<FullPivLU>
187 image(
const MatrixType& originalMatrix)
const
189 eigen_assert(m_isInitialized &&
"LU is not initialized.");
190 return internal::image_retval<FullPivLU>(*
this, originalMatrix);
212 template<
typename Rhs>
213 inline const internal::solve_retval<FullPivLU, Rhs>
216 eigen_assert(m_isInitialized &&
"LU is not initialized.");
217 return internal::solve_retval<FullPivLU, Rhs>(*
this, b.derived());
235 typename internal::traits<MatrixType>::Scalar
determinant()
const;
256 m_usePrescribedThreshold =
true;
271 m_usePrescribedThreshold =
false;
281 eigen_assert(m_isInitialized || m_usePrescribedThreshold);
282 return m_usePrescribedThreshold ? m_prescribedThreshold
297 eigen_assert(m_isInitialized &&
"LU is not initialized.");
298 RealScalar premultiplied_threshold = abs(m_maxpivot) *
threshold();
300 for(Index i = 0; i < m_nonzero_pivots; ++i)
301 result += (abs(m_lu.coeff(i,i)) > premultiplied_threshold);
313 eigen_assert(m_isInitialized &&
"LU is not initialized.");
314 return cols() -
rank();
326 eigen_assert(m_isInitialized &&
"LU is not initialized.");
327 return rank() == cols();
339 eigen_assert(m_isInitialized &&
"LU is not initialized.");
340 return rank() == rows();
351 eigen_assert(m_isInitialized &&
"LU is not initialized.");
352 return isInjective() && (m_lu.rows() == m_lu.cols());
362 inline const internal::solve_retval<FullPivLU,typename MatrixType::IdentityReturnType>
inverse()
const
364 eigen_assert(m_isInitialized &&
"LU is not initialized.");
365 eigen_assert(m_lu.rows() == m_lu.cols() &&
"You can't take the inverse of a non-square matrix!");
366 return internal::solve_retval<FullPivLU,typename MatrixType::IdentityReturnType>
367 (*
this, MatrixType::Identity(m_lu.rows(), m_lu.cols()));
372 inline Index rows()
const {
return m_lu.rows(); }
373 inline Index cols()
const {
return m_lu.cols(); }
377 PermutationPType m_p;
378 PermutationQType m_q;
379 IntColVectorType m_rowsTranspositions;
380 IntRowVectorType m_colsTranspositions;
381 Index m_det_pq, m_nonzero_pivots;
382 RealScalar m_maxpivot, m_prescribedThreshold;
383 bool m_isInitialized, m_usePrescribedThreshold;
386 template<
typename MatrixType>
388 : m_isInitialized(false), m_usePrescribedThreshold(false)
392 template<
typename MatrixType>
397 m_rowsTranspositions(rows),
398 m_colsTranspositions(cols),
399 m_isInitialized(false),
400 m_usePrescribedThreshold(false)
404 template<
typename MatrixType>
406 : m_lu(matrix.rows(), matrix.cols()),
409 m_rowsTranspositions(matrix.rows()),
410 m_colsTranspositions(matrix.cols()),
411 m_isInitialized(false),
412 m_usePrescribedThreshold(false)
417 template<
typename MatrixType>
423 m_isInitialized =
true;
426 const Index size = matrix.diagonalSize();
427 const Index rows = matrix.rows();
428 const Index cols = matrix.cols();
432 m_rowsTranspositions.resize(matrix.rows());
433 m_colsTranspositions.resize(matrix.cols());
434 Index number_of_transpositions = 0;
436 m_nonzero_pivots = size;
437 m_maxpivot = RealScalar(0);
439 for(Index k = 0; k < size; ++k)
444 Index row_of_biggest_in_corner, col_of_biggest_in_corner;
445 RealScalar biggest_in_corner;
446 biggest_in_corner = m_lu.bottomRightCorner(rows-k, cols-k)
448 .maxCoeff(&row_of_biggest_in_corner, &col_of_biggest_in_corner);
449 row_of_biggest_in_corner += k;
450 col_of_biggest_in_corner += k;
452 if(biggest_in_corner==RealScalar(0))
456 m_nonzero_pivots = k;
457 for(Index i = k; i < size; ++i)
459 m_rowsTranspositions.coeffRef(i) = i;
460 m_colsTranspositions.coeffRef(i) = i;
465 if(biggest_in_corner > m_maxpivot) m_maxpivot = biggest_in_corner;
470 m_rowsTranspositions.coeffRef(k) = row_of_biggest_in_corner;
471 m_colsTranspositions.coeffRef(k) = col_of_biggest_in_corner;
472 if(k != row_of_biggest_in_corner) {
473 m_lu.row(k).swap(m_lu.row(row_of_biggest_in_corner));
474 ++number_of_transpositions;
476 if(k != col_of_biggest_in_corner) {
477 m_lu.col(k).swap(m_lu.col(col_of_biggest_in_corner));
478 ++number_of_transpositions;
485 m_lu.col(k).tail(rows-k-1) /= m_lu.coeff(k,k);
487 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);
493 m_p.setIdentity(rows);
494 for(Index k = size-1; k >= 0; --k)
495 m_p.applyTranspositionOnTheRight(k, m_rowsTranspositions.coeff(k));
497 m_q.setIdentity(cols);
498 for(Index k = 0; k < size; ++k)
499 m_q.applyTranspositionOnTheRight(k, m_colsTranspositions.coeff(k));
501 m_det_pq = (number_of_transpositions%2) ? -1 : 1;
505 template<
typename MatrixType>
508 eigen_assert(m_isInitialized &&
"LU is not initialized.");
509 eigen_assert(m_lu.rows() == m_lu.cols() &&
"You can't take the determinant of a non-square matrix!");
510 return Scalar(m_det_pq) * Scalar(m_lu.diagonal().prod());
516 template<
typename MatrixType>
519 eigen_assert(m_isInitialized &&
"LU is not initialized.");
520 const Index smalldim = (std::min)(m_lu.rows(), m_lu.cols());
522 MatrixType res(m_lu.rows(),m_lu.cols());
524 res = m_lu.leftCols(smalldim)
525 .template triangularView<UnitLower>().toDenseMatrix()
526 * m_lu.topRows(smalldim)
527 .template triangularView<Upper>().toDenseMatrix();
530 res = m_p.inverse() * res;
533 res = res * m_q.inverse();
541 template<
typename _MatrixType>
542 struct kernel_retval<
FullPivLU<_MatrixType> >
543 : kernel_retval_base<FullPivLU<_MatrixType> >
547 enum { MaxSmallDimAtCompileTime = EIGEN_SIZE_MIN_PREFER_FIXED(
548 MatrixType::MaxColsAtCompileTime,
549 MatrixType::MaxRowsAtCompileTime)
552 template<
typename Dest>
void evalTo(Dest& dst)
const
555 const Index cols = dec().matrixLU().cols(), dimker = cols - rank();
581 Matrix<Index, Dynamic, 1, 0, MaxSmallDimAtCompileTime, 1> pivots(rank());
582 RealScalar premultiplied_threshold = dec().maxPivot() * dec().threshold();
584 for(Index i = 0; i < dec().nonzeroPivots(); ++i)
585 if(abs(dec().matrixLU().coeff(i,i)) > premultiplied_threshold)
586 pivots.coeffRef(p++) = i;
587 eigen_internal_assert(p == rank());
593 Matrix<
typename MatrixType::Scalar,
Dynamic,
Dynamic, MatrixType::Options,
594 MaxSmallDimAtCompileTime, MatrixType::MaxColsAtCompileTime>
595 m(dec().matrixLU().block(0, 0, rank(), cols));
596 for(Index i = 0; i < rank(); ++i)
598 if(i) m.row(i).head(i).setZero();
599 m.row(i).tail(cols-i) = dec().matrixLU().row(pivots.coeff(i)).tail(cols-i);
601 m.block(0, 0, rank(), rank());
602 m.block(0, 0, rank(), rank()).template triangularView<StrictlyLower>().setZero();
603 for(Index i = 0; i < rank(); ++i)
604 m.col(i).swap(m.col(pivots.coeff(i)));
609 m.topLeftCorner(rank(), rank())
610 .template triangularView<Upper>().solveInPlace(
611 m.topRightCorner(rank(), dimker)
615 for(Index i = rank()-1; i >= 0; --i)
616 m.col(i).swap(m.col(pivots.coeff(i)));
619 for(Index i = 0; i < rank(); ++i) dst.row(dec().permutationQ().indices().coeff(i)) = -m.row(i).tail(dimker);
620 for(Index i = rank(); i < cols; ++i) dst.row(dec().permutationQ().indices().coeff(i)).setZero();
621 for(Index k = 0; k < dimker; ++k) dst.coeffRef(dec().permutationQ().indices().coeff(rank()+k), k) = Scalar(1);
627 template<
typename _MatrixType>
628 struct image_retval<FullPivLU<_MatrixType> >
629 : image_retval_base<FullPivLU<_MatrixType> >
631 EIGEN_MAKE_IMAGE_HELPERS(FullPivLU<_MatrixType>)
633 enum { MaxSmallDimAtCompileTime = EIGEN_SIZE_MIN_PREFER_FIXED(
634 MatrixType::MaxColsAtCompileTime,
635 MatrixType::MaxRowsAtCompileTime)
638 template<
typename Dest>
void evalTo(Dest& dst)
const
650 Matrix<Index, Dynamic, 1, 0, MaxSmallDimAtCompileTime, 1> pivots(rank());
651 RealScalar premultiplied_threshold = dec().maxPivot() * dec().threshold();
653 for(Index i = 0; i < dec().nonzeroPivots(); ++i)
654 if(abs(dec().matrixLU().coeff(i,i)) > premultiplied_threshold)
655 pivots.coeffRef(p++) = i;
656 eigen_internal_assert(p == rank());
658 for(Index i = 0; i < rank(); ++i)
659 dst.col(i) = originalMatrix().col(dec().permutationQ().indices().coeff(pivots.coeff(i)));
665 template<
typename _MatrixType,
typename Rhs>
666 struct solve_retval<FullPivLU<_MatrixType>, Rhs>
667 : solve_retval_base<FullPivLU<_MatrixType>, Rhs>
669 EIGEN_MAKE_SOLVE_HELPERS(FullPivLU<_MatrixType>,Rhs)
671 template<typename Dest>
void evalTo(Dest& dst)
const
681 const Index rows = dec().rows(), cols = dec().cols(),
682 nonzero_pivots = dec().nonzeroPivots();
683 eigen_assert(rhs().rows() == rows);
684 const Index smalldim = (std::min)(rows, cols);
686 if(nonzero_pivots == 0)
692 typename Rhs::PlainObject c(rhs().rows(), rhs().cols());
695 c = dec().permutationP() * rhs();
699 .topLeftCorner(smalldim,smalldim)
700 .template triangularView<UnitLower>()
701 .solveInPlace(c.topRows(smalldim));
704 c.bottomRows(rows-cols)
705 -= dec().matrixLU().bottomRows(rows-cols)
711 .topLeftCorner(nonzero_pivots, nonzero_pivots)
712 .template triangularView<Upper>()
713 .solveInPlace(c.topRows(nonzero_pivots));
716 for(Index i = 0; i < nonzero_pivots; ++i)
717 dst.row(dec().permutationQ().indices().coeff(i)) = c.row(i);
718 for(Index i = nonzero_pivots; i < dec().matrixLU().cols(); ++i)
719 dst.row(dec().permutationQ().indices().coeff(i)).setZero();
733 template<
typename Derived>
734 inline const FullPivLU<typename MatrixBase<Derived>::PlainObject>
bool isInvertible() const
Definition: FullPivLU.h:349
RealScalar threshold() const
Definition: FullPivLU.h:279
const MatrixType & matrixLU() const
Definition: FullPivLU.h:103
const FullPivLU< PlainObject > fullPivLu() const
Definition: FullPivLU.h:735
internal::traits< MatrixType >::Scalar determinant() const
Definition: FullPivLU.h:506
MatrixType reconstructedMatrix() const
Definition: FullPivLU.h:517
Holds information about the various numeric (i.e. scalar) types allowed by Eigen. ...
Definition: NumTraits.h:88
const int Dynamic
Definition: Constants.h:21
bool isInjective() const
Definition: FullPivLU.h:324
FullPivLU()
Default Constructor.
Definition: FullPivLU.h:387
FullPivLU & compute(const MatrixType &matrix)
Definition: FullPivLU.h:418
const PermutationPType & permutationP() const
Definition: FullPivLU.h:131
Index rank() const
Definition: FullPivLU.h:294
const internal::image_retval< FullPivLU > image(const MatrixType &originalMatrix) const
Definition: FullPivLU.h:187
FullPivLU & setThreshold(Default_t)
Definition: FullPivLU.h:269
Index nonzeroPivots() const
Definition: FullPivLU.h:116
Index dimensionOfKernel() const
Definition: FullPivLU.h:311
LU decomposition of a matrix with complete pivoting, and related features.
Definition: ForwardDeclarations.h:216
FullPivLU & setThreshold(const RealScalar &threshold)
Definition: FullPivLU.h:254
const internal::kernel_retval< FullPivLU > kernel() const
Definition: FullPivLU.h:161
RealScalar maxPivot() const
Definition: FullPivLU.h:125
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:48
const internal::solve_retval< FullPivLU, Rhs > solve(const MatrixBase< Rhs > &b) const
Definition: FullPivLU.h:214
bool isSurjective() const
Definition: FullPivLU.h:337
const internal::solve_retval< FullPivLU, typename MatrixType::IdentityReturnType > inverse() const
Definition: FullPivLU.h:362
const PermutationQType & permutationQ() const
Definition: FullPivLU.h:141