16 template<
typename _MatrixType>
struct traits<FullPivLU<_MatrixType> >
55 template<
typename _MatrixType>
class FullPivLU
58 typedef _MatrixType MatrixType;
60 RowsAtCompileTime = MatrixType::RowsAtCompileTime,
61 ColsAtCompileTime = MatrixType::ColsAtCompileTime,
62 Options = MatrixType::Options,
63 MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
64 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
66 typedef typename MatrixType::Scalar Scalar;
67 typedef typename NumTraits<typename MatrixType::Scalar>::Real RealScalar;
68 typedef typename internal::traits<MatrixType>::StorageKind StorageKind;
70 typedef typename MatrixType::StorageIndex StorageIndex;
71 typedef typename internal::plain_row_type<MatrixType, StorageIndex>::type IntRowVectorType;
72 typedef typename internal::plain_col_type<MatrixType, StorageIndex>::type IntColVectorType;
73 typedef PermutationMatrix<ColsAtCompileTime, MaxColsAtCompileTime> PermutationQType;
74 typedef PermutationMatrix<RowsAtCompileTime, MaxRowsAtCompileTime> PermutationPType;
75 typedef typename MatrixType::PlainObject PlainObject;
98 explicit FullPivLU(
const MatrixType& matrix);
117 eigen_assert(m_isInitialized &&
"LU is not initialized.");
130 eigen_assert(m_isInitialized &&
"LU is not initialized.");
131 return m_nonzero_pivots;
145 eigen_assert(m_isInitialized &&
"LU is not initialized.");
155 eigen_assert(m_isInitialized &&
"LU is not initialized.");
173 inline const internal::kernel_retval<FullPivLU>
kernel()
const
175 eigen_assert(m_isInitialized &&
"LU is not initialized.");
176 return internal::kernel_retval<FullPivLU>(*this);
198 inline const internal::image_retval<FullPivLU>
199 image(
const MatrixType& originalMatrix)
const
201 eigen_assert(m_isInitialized &&
"LU is not initialized.");
202 return internal::image_retval<FullPivLU>(*
this, originalMatrix);
224 template<
typename Rhs>
228 eigen_assert(m_isInitialized &&
"LU is not initialized.");
247 typename internal::traits<MatrixType>::Scalar
determinant()
const;
268 m_usePrescribedThreshold =
true;
283 m_usePrescribedThreshold =
false;
293 eigen_assert(m_isInitialized || m_usePrescribedThreshold);
294 return m_usePrescribedThreshold ? m_prescribedThreshold
309 eigen_assert(m_isInitialized &&
"LU is not initialized.");
310 RealScalar premultiplied_threshold = abs(m_maxpivot) *
threshold();
312 for(Index i = 0; i < m_nonzero_pivots; ++i)
313 result += (abs(m_lu.coeff(i,i)) > premultiplied_threshold);
325 eigen_assert(m_isInitialized &&
"LU is not initialized.");
326 return cols() -
rank();
338 eigen_assert(m_isInitialized &&
"LU is not initialized.");
339 return rank() == cols();
351 eigen_assert(m_isInitialized &&
"LU is not initialized.");
352 return rank() == rows();
363 eigen_assert(m_isInitialized &&
"LU is not initialized.");
364 return isInjective() && (m_lu.rows() == m_lu.cols());
376 eigen_assert(m_isInitialized &&
"LU is not initialized.");
377 eigen_assert(m_lu.rows() == m_lu.cols() &&
"You can't take the inverse of a non-square matrix!");
383 inline Index rows()
const {
return m_lu.rows(); }
384 inline Index cols()
const {
return m_lu.cols(); }
386 #ifndef EIGEN_PARSED_BY_DOXYGEN
387 template<
typename RhsType,
typename DstType>
389 void _solve_impl(
const RhsType &rhs, DstType &dst)
const;
394 static void check_template_parameters()
396 EIGEN_STATIC_ASSERT_NON_INTEGER(Scalar);
400 PermutationPType m_p;
401 PermutationQType m_q;
402 IntColVectorType m_rowsTranspositions;
403 IntRowVectorType m_colsTranspositions;
404 Index m_det_pq, m_nonzero_pivots;
405 RealScalar m_maxpivot, m_prescribedThreshold;
406 bool m_isInitialized, m_usePrescribedThreshold;
409 template<
typename MatrixType>
411 : m_isInitialized(false), m_usePrescribedThreshold(false)
415 template<
typename MatrixType>
420 m_rowsTranspositions(rows),
421 m_colsTranspositions(cols),
422 m_isInitialized(false),
423 m_usePrescribedThreshold(false)
427 template<
typename MatrixType>
429 : m_lu(matrix.rows(), matrix.cols()),
432 m_rowsTranspositions(matrix.rows()),
433 m_colsTranspositions(matrix.cols()),
434 m_isInitialized(false),
435 m_usePrescribedThreshold(false)
440 template<
typename MatrixType>
443 check_template_parameters();
448 m_isInitialized =
true;
451 const Index size = matrix.diagonalSize();
452 const Index rows = matrix.rows();
453 const Index cols = matrix.cols();
457 m_rowsTranspositions.resize(matrix.rows());
458 m_colsTranspositions.resize(matrix.cols());
459 Index number_of_transpositions = 0;
461 m_nonzero_pivots = size;
462 m_maxpivot = RealScalar(0);
464 for(Index k = 0; k < size; ++k)
469 Index row_of_biggest_in_corner, col_of_biggest_in_corner;
470 typedef internal::scalar_score_coeff_op<Scalar> Scoring;
471 typedef typename Scoring::result_type Score;
472 Score biggest_in_corner;
473 biggest_in_corner = m_lu.bottomRightCorner(rows-k, cols-k)
474 .unaryExpr(Scoring())
475 .maxCoeff(&row_of_biggest_in_corner, &col_of_biggest_in_corner);
476 row_of_biggest_in_corner += k;
477 col_of_biggest_in_corner += k;
479 if(biggest_in_corner==Score(0))
483 m_nonzero_pivots = k;
484 for(Index i = k; i < size; ++i)
486 m_rowsTranspositions.coeffRef(i) = i;
487 m_colsTranspositions.coeffRef(i) = i;
492 RealScalar abs_pivot = internal::abs_knowing_score<Scalar>()(m_lu(row_of_biggest_in_corner, col_of_biggest_in_corner), biggest_in_corner);
493 if(abs_pivot > m_maxpivot) m_maxpivot = abs_pivot;
498 m_rowsTranspositions.coeffRef(k) = row_of_biggest_in_corner;
499 m_colsTranspositions.coeffRef(k) = col_of_biggest_in_corner;
500 if(k != row_of_biggest_in_corner) {
501 m_lu.row(k).swap(m_lu.row(row_of_biggest_in_corner));
502 ++number_of_transpositions;
504 if(k != col_of_biggest_in_corner) {
505 m_lu.col(k).swap(m_lu.col(col_of_biggest_in_corner));
506 ++number_of_transpositions;
513 m_lu.col(k).tail(rows-k-1) /= m_lu.coeff(k,k);
515 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);
521 m_p.setIdentity(rows);
522 for(Index k = size-1; k >= 0; --k)
523 m_p.applyTranspositionOnTheRight(k, m_rowsTranspositions.coeff(k));
525 m_q.setIdentity(cols);
526 for(Index k = 0; k < size; ++k)
527 m_q.applyTranspositionOnTheRight(k, m_colsTranspositions.coeff(k));
529 m_det_pq = (number_of_transpositions%2) ? -1 : 1;
533 template<
typename MatrixType>
536 eigen_assert(m_isInitialized &&
"LU is not initialized.");
537 eigen_assert(m_lu.rows() == m_lu.cols() &&
"You can't take the determinant of a non-square matrix!");
538 return Scalar(m_det_pq) * Scalar(m_lu.diagonal().prod());
544 template<
typename MatrixType>
547 eigen_assert(m_isInitialized &&
"LU is not initialized.");
548 const Index smalldim = (std::min)(m_lu.rows(), m_lu.cols());
550 MatrixType res(m_lu.rows(),m_lu.cols());
552 res = m_lu.leftCols(smalldim)
553 .template triangularView<UnitLower>().toDenseMatrix()
554 * m_lu.topRows(smalldim)
555 .template triangularView<Upper>().toDenseMatrix();
558 res = m_p.inverse() * res;
561 res = res * m_q.inverse();
569 template<
typename _MatrixType>
570 struct kernel_retval<
FullPivLU<_MatrixType> >
571 : kernel_retval_base<FullPivLU<_MatrixType> >
575 enum { MaxSmallDimAtCompileTime = EIGEN_SIZE_MIN_PREFER_FIXED(
576 MatrixType::MaxColsAtCompileTime,
577 MatrixType::MaxRowsAtCompileTime)
580 template<
typename Dest>
void evalTo(Dest& dst)
const
583 const Index cols = dec().matrixLU().cols(), dimker = cols - rank();
609 Matrix<Index, Dynamic, 1, 0, MaxSmallDimAtCompileTime, 1> pivots(rank());
610 RealScalar premultiplied_threshold = dec().maxPivot() * dec().threshold();
612 for(Index i = 0; i < dec().nonzeroPivots(); ++i)
613 if(abs(dec().matrixLU().coeff(i,i)) > premultiplied_threshold)
614 pivots.coeffRef(p++) = i;
615 eigen_internal_assert(p == rank());
621 Matrix<
typename MatrixType::Scalar, Dynamic, Dynamic, MatrixType::Options,
622 MaxSmallDimAtCompileTime, MatrixType::MaxColsAtCompileTime>
623 m(dec().matrixLU().block(0, 0, rank(), cols));
624 for(Index i = 0; i < rank(); ++i)
626 if(i) m.row(i).head(i).setZero();
627 m.row(i).tail(cols-i) = dec().matrixLU().row(pivots.coeff(i)).tail(cols-i);
629 m.block(0, 0, rank(), rank());
630 m.block(0, 0, rank(), rank()).template triangularView<StrictlyLower>().setZero();
631 for(Index i = 0; i < rank(); ++i)
632 m.col(i).swap(m.col(pivots.coeff(i)));
637 m.topLeftCorner(rank(), rank())
638 .template triangularView<Upper>().solveInPlace(
639 m.topRightCorner(rank(), dimker)
643 for(Index i = rank()-1; i >= 0; --i)
644 m.col(i).swap(m.col(pivots.coeff(i)));
647 for(Index i = 0; i < rank(); ++i) dst.row(dec().permutationQ().indices().coeff(i)) = -m.row(i).tail(dimker);
648 for(Index i = rank(); i < cols; ++i) dst.row(dec().permutationQ().indices().coeff(i)).setZero();
649 for(Index k = 0; k < dimker; ++k) dst.coeffRef(dec().permutationQ().indices().coeff(rank()+k), k) = Scalar(1);
655 template<
typename _MatrixType>
656 struct image_retval<FullPivLU<_MatrixType> >
657 : image_retval_base<FullPivLU<_MatrixType> >
659 EIGEN_MAKE_IMAGE_HELPERS(FullPivLU<_MatrixType>)
661 enum { MaxSmallDimAtCompileTime = EIGEN_SIZE_MIN_PREFER_FIXED(
662 MatrixType::MaxColsAtCompileTime,
663 MatrixType::MaxRowsAtCompileTime)
666 template<
typename Dest>
void evalTo(Dest& dst)
const
678 Matrix<Index, Dynamic, 1, 0, MaxSmallDimAtCompileTime, 1> pivots(rank());
679 RealScalar premultiplied_threshold = dec().maxPivot() * dec().threshold();
681 for(Index i = 0; i < dec().nonzeroPivots(); ++i)
682 if(abs(dec().matrixLU().coeff(i,i)) > premultiplied_threshold)
683 pivots.coeffRef(p++) = i;
684 eigen_internal_assert(p == rank());
686 for(Index i = 0; i < rank(); ++i)
687 dst.col(i) = originalMatrix().col(dec().permutationQ().indices().coeff(pivots.coeff(i)));
695 #ifndef EIGEN_PARSED_BY_DOXYGEN
696 template<
typename _MatrixType>
697 template<
typename RhsType,
typename DstType>
698 void FullPivLU<_MatrixType>::_solve_impl(
const RhsType &rhs, DstType &dst)
const
708 const Index rows = this->rows(),
710 nonzero_pivots = this->nonzeroPivots();
711 eigen_assert(rhs.rows() == rows);
712 const Index smalldim = (std::min)(rows, cols);
714 if(nonzero_pivots == 0)
720 typename RhsType::PlainObject c(rhs.rows(), rhs.cols());
723 c = permutationP() * rhs;
726 m_lu.topLeftCorner(smalldim,smalldim)
727 .template triangularView<UnitLower>()
728 .solveInPlace(c.topRows(smalldim));
730 c.bottomRows(rows-cols) -= m_lu.bottomRows(rows-cols) * c.topRows(cols);
733 m_lu.topLeftCorner(nonzero_pivots, nonzero_pivots)
734 .template triangularView<Upper>()
735 .solveInPlace(c.topRows(nonzero_pivots));
738 for(Index i = 0; i < nonzero_pivots; ++i)
739 dst.row(permutationQ().indices().coeff(i)) = c.row(i);
740 for(Index i = nonzero_pivots; i < m_lu.cols(); ++i)
741 dst.row(permutationQ().indices().coeff(i)).setZero();
749 template<
typename DstXprType,
typename MatrixType,
typename Scalar>
750 struct Assignment<DstXprType, Inverse<FullPivLU<MatrixType> >,
internal::assign_op<Scalar>, Dense2Dense, Scalar>
752 typedef FullPivLU<MatrixType> LuType;
753 typedef Inverse<LuType> SrcXprType;
754 static void run(DstXprType &dst,
const SrcXprType &src,
const internal::assign_op<Scalar> &)
756 dst = src.nestedExpression().solve(MatrixType::Identity(src.rows(), src.cols()));
770 template<
typename Derived>
771 inline const FullPivLU<typename MatrixBase<Derived>::PlainObject>
bool isInjective() const
Definition: FullPivLU.h:336
const internal::image_retval< FullPivLU > image(const MatrixType &originalMatrix) const
Definition: FullPivLU.h:199
FullPivLU & setThreshold(const RealScalar &threshold)
Definition: FullPivLU.h:266
internal::traits< MatrixType >::Scalar determinant() const
Definition: FullPivLU.h:534
Holds information about the various numeric (i.e. scalar) types allowed by Eigen. ...
Definition: NumTraits.h:107
RealScalar threshold() const
Definition: FullPivLU.h:291
Index nonzeroPivots() const
Definition: FullPivLU.h:128
Expression of the inverse of another expression.
Definition: Inverse.h:45
bool isSurjective() const
Definition: FullPivLU.h:349
FullPivLU & setThreshold(Default_t)
Definition: FullPivLU.h:281
const PermutationQType & permutationQ() const
Definition: FullPivLU.h:153
FullPivLU & compute(const MatrixType &matrix)
Definition: FullPivLU.h:441
MatrixType reconstructedMatrix() const
Definition: FullPivLU.h:545
const PermutationPType & permutationP() const
Definition: FullPivLU.h:143
Definition: Eigen_Colamd.h:54
LU decomposition of a matrix with complete pivoting, and related features.
Definition: ForwardDeclarations.h:246
bool isInvertible() const
Definition: FullPivLU.h:361
const FullPivLU< PlainObject > fullPivLu() const
Definition: FullPivLU.h:772
FullPivLU()
Default Constructor.
Definition: FullPivLU.h:410
Index rank() const
Definition: FullPivLU.h:306
const MatrixType & matrixLU() const
Definition: FullPivLU.h:115
Pseudo expression representing a solving operation.
Definition: Solve.h:63
const Inverse< FullPivLU > inverse() const
Definition: FullPivLU.h:374
Index dimensionOfKernel() const
Definition: FullPivLU.h:323
RealScalar maxPivot() const
Definition: FullPivLU.h:137
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:48
const Solve< FullPivLU, Rhs > solve(const MatrixBase< Rhs > &b) const
Definition: FullPivLU.h:226
const internal::kernel_retval< FullPivLU > kernel() const
Definition: FullPivLU.h:173