16 template<
typename _MatrixType>
struct traits<FullPivLU<_MatrixType> >
19 typedef MatrixXpr XprKind;
20 typedef SolverStorage StorageKind;
59 template<
typename _MatrixType>
class FullPivLU
60 :
public SolverBase<FullPivLU<_MatrixType> >
63 typedef _MatrixType MatrixType;
64 typedef SolverBase<FullPivLU> Base;
69 MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
70 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
72 typedef typename internal::plain_row_type<MatrixType, StorageIndex>::type IntRowVectorType;
73 typedef typename internal::plain_col_type<MatrixType, StorageIndex>::type IntColVectorType;
74 typedef PermutationMatrix<ColsAtCompileTime, MaxColsAtCompileTime> PermutationQType;
75 typedef PermutationMatrix<RowsAtCompileTime, MaxRowsAtCompileTime> PermutationPType;
76 typedef typename MatrixType::PlainObject PlainObject;
99 template<
typename InputType>
100 explicit FullPivLU(
const EigenBase<InputType>& matrix);
108 template<
typename InputType>
109 explicit FullPivLU(EigenBase<InputType>& matrix);
118 template<
typename InputType>
133 eigen_assert(m_isInitialized &&
"LU is not initialized.");
146 eigen_assert(m_isInitialized &&
"LU is not initialized.");
147 return m_nonzero_pivots;
161 eigen_assert(m_isInitialized &&
"LU is not initialized.");
171 eigen_assert(m_isInitialized &&
"LU is not initialized.");
189 inline const internal::kernel_retval<FullPivLU>
kernel()
const 191 eigen_assert(m_isInitialized &&
"LU is not initialized.");
192 return internal::kernel_retval<FullPivLU>(*
this);
214 inline const internal::image_retval<FullPivLU>
215 image(
const MatrixType& originalMatrix)
const 217 eigen_assert(m_isInitialized &&
"LU is not initialized.");
218 return internal::image_retval<FullPivLU>(*
this, originalMatrix);
241 template<
typename Rhs>
245 eigen_assert(m_isInitialized &&
"LU is not initialized.");
254 eigen_assert(m_isInitialized &&
"PartialPivLU is not initialized.");
255 return internal::rcond_estimate_helper(m_l1_norm, *
this);
273 typename internal::traits<MatrixType>::Scalar
determinant()
const;
294 m_usePrescribedThreshold =
true;
309 m_usePrescribedThreshold =
false;
319 eigen_assert(m_isInitialized || m_usePrescribedThreshold);
320 return m_usePrescribedThreshold ? m_prescribedThreshold
335 eigen_assert(m_isInitialized &&
"LU is not initialized.");
336 RealScalar premultiplied_threshold =
abs(m_maxpivot) *
threshold();
338 for(
Index i = 0; i < m_nonzero_pivots; ++i)
339 result += (
abs(m_lu.coeff(i,i)) > premultiplied_threshold);
351 eigen_assert(m_isInitialized &&
"LU is not initialized.");
352 return cols() -
rank();
364 eigen_assert(m_isInitialized &&
"LU is not initialized.");
365 return rank() == cols();
377 eigen_assert(m_isInitialized &&
"LU is not initialized.");
378 return rank() == rows();
389 eigen_assert(m_isInitialized &&
"LU is not initialized.");
390 return isInjective() && (m_lu.rows() == m_lu.cols());
402 eigen_assert(m_isInitialized &&
"LU is not initialized.");
403 eigen_assert(m_lu.rows() == m_lu.cols() &&
"You can't take the inverse of a non-square matrix!");
409 EIGEN_DEVICE_FUNC
inline Index rows()
const {
return m_lu.rows(); }
410 EIGEN_DEVICE_FUNC
inline Index cols()
const {
return m_lu.cols(); }
412 #ifndef EIGEN_PARSED_BY_DOXYGEN 413 template<
typename RhsType,
typename DstType>
415 void _solve_impl(
const RhsType &rhs, DstType &dst)
const;
417 template<
bool Conjugate,
typename RhsType,
typename DstType>
419 void _solve_impl_transposed(
const RhsType &rhs, DstType &dst)
const;
424 static void check_template_parameters()
426 EIGEN_STATIC_ASSERT_NON_INTEGER(Scalar);
429 void computeInPlace();
432 PermutationPType m_p;
433 PermutationQType m_q;
434 IntColVectorType m_rowsTranspositions;
435 IntRowVectorType m_colsTranspositions;
436 Index m_nonzero_pivots;
437 RealScalar m_l1_norm;
438 RealScalar m_maxpivot, m_prescribedThreshold;
439 signed char m_det_pq;
440 bool m_isInitialized, m_usePrescribedThreshold;
443 template<
typename MatrixType>
445 : m_isInitialized(false), m_usePrescribedThreshold(false)
449 template<
typename MatrixType>
454 m_rowsTranspositions(rows),
455 m_colsTranspositions(cols),
456 m_isInitialized(false),
457 m_usePrescribedThreshold(false)
461 template<
typename MatrixType>
462 template<
typename InputType>
464 : m_lu(matrix.rows(), matrix.cols()),
467 m_rowsTranspositions(matrix.rows()),
468 m_colsTranspositions(matrix.cols()),
469 m_isInitialized(false),
470 m_usePrescribedThreshold(false)
475 template<
typename MatrixType>
476 template<
typename InputType>
478 : m_lu(matrix.derived()),
481 m_rowsTranspositions(matrix.rows()),
482 m_colsTranspositions(matrix.cols()),
483 m_isInitialized(false),
484 m_usePrescribedThreshold(false)
489 template<
typename MatrixType>
492 check_template_parameters();
497 m_l1_norm = m_lu.cwiseAbs().colwise().sum().maxCoeff();
499 const Index size = m_lu.diagonalSize();
500 const Index rows = m_lu.rows();
501 const Index cols = m_lu.cols();
505 m_rowsTranspositions.resize(m_lu.rows());
506 m_colsTranspositions.resize(m_lu.cols());
507 Index number_of_transpositions = 0;
509 m_nonzero_pivots = size;
510 m_maxpivot = RealScalar(0);
512 for(
Index k = 0; k < size; ++k)
517 Index row_of_biggest_in_corner, col_of_biggest_in_corner;
518 typedef internal::scalar_score_coeff_op<Scalar> Scoring;
519 typedef typename Scoring::result_type Score;
520 Score biggest_in_corner;
521 biggest_in_corner = m_lu.bottomRightCorner(rows-k, cols-k)
522 .unaryExpr(Scoring())
523 .maxCoeff(&row_of_biggest_in_corner, &col_of_biggest_in_corner);
524 row_of_biggest_in_corner += k;
525 col_of_biggest_in_corner += k;
527 if(biggest_in_corner==Score(0))
531 m_nonzero_pivots = k;
532 for(
Index i = k; i < size; ++i)
534 m_rowsTranspositions.coeffRef(i) = i;
535 m_colsTranspositions.coeffRef(i) = i;
540 RealScalar abs_pivot = internal::abs_knowing_score<Scalar>()(m_lu(row_of_biggest_in_corner, col_of_biggest_in_corner), biggest_in_corner);
541 if(abs_pivot > m_maxpivot) m_maxpivot = abs_pivot;
546 m_rowsTranspositions.coeffRef(k) = row_of_biggest_in_corner;
547 m_colsTranspositions.coeffRef(k) = col_of_biggest_in_corner;
548 if(k != row_of_biggest_in_corner) {
549 m_lu.row(k).swap(m_lu.row(row_of_biggest_in_corner));
550 ++number_of_transpositions;
552 if(k != col_of_biggest_in_corner) {
553 m_lu.col(k).swap(m_lu.col(col_of_biggest_in_corner));
554 ++number_of_transpositions;
561 m_lu.col(k).tail(rows-k-1) /= m_lu.coeff(k,k);
563 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);
569 m_p.setIdentity(rows);
570 for(
Index k = size-1; k >= 0; --k)
571 m_p.applyTranspositionOnTheRight(k, m_rowsTranspositions.coeff(k));
573 m_q.setIdentity(cols);
574 for(
Index k = 0; k < size; ++k)
575 m_q.applyTranspositionOnTheRight(k, m_colsTranspositions.coeff(k));
577 m_det_pq = (number_of_transpositions%2) ? -1 : 1;
579 m_isInitialized =
true;
582 template<
typename MatrixType>
585 eigen_assert(m_isInitialized &&
"LU is not initialized.");
586 eigen_assert(m_lu.rows() == m_lu.cols() &&
"You can't take the determinant of a non-square matrix!");
587 return Scalar(m_det_pq) * Scalar(m_lu.diagonal().prod());
593 template<
typename MatrixType>
596 eigen_assert(m_isInitialized &&
"LU is not initialized.");
597 const Index smalldim = (std::min)(m_lu.rows(), m_lu.cols());
599 MatrixType res(m_lu.rows(),m_lu.cols());
601 res = m_lu.leftCols(smalldim)
602 .template triangularView<UnitLower>().toDenseMatrix()
603 * m_lu.topRows(smalldim)
604 .template triangularView<Upper>().toDenseMatrix();
607 res = m_p.inverse() * res;
610 res = res * m_q.inverse();
618 template<
typename _MatrixType>
619 struct kernel_retval<
FullPivLU<_MatrixType> >
620 : kernel_retval_base<FullPivLU<_MatrixType> >
624 enum { MaxSmallDimAtCompileTime = EIGEN_SIZE_MIN_PREFER_FIXED(
625 MatrixType::MaxColsAtCompileTime,
626 MatrixType::MaxRowsAtCompileTime)
629 template<
typename Dest>
void evalTo(Dest& dst)
const 632 const Index cols = dec().matrixLU().cols(), dimker = cols - rank();
658 Matrix<Index, Dynamic, 1, 0, MaxSmallDimAtCompileTime, 1> pivots(rank());
659 RealScalar premultiplied_threshold = dec().maxPivot() * dec().threshold();
661 for(
Index i = 0; i < dec().nonzeroPivots(); ++i)
662 if(
abs(dec().matrixLU().coeff(i,i)) > premultiplied_threshold)
663 pivots.coeffRef(p++) = i;
664 eigen_internal_assert(p == rank());
670 Matrix<
typename MatrixType::Scalar,
Dynamic,
Dynamic, MatrixType::Options,
671 MaxSmallDimAtCompileTime, MatrixType::MaxColsAtCompileTime>
672 m(dec().matrixLU().block(0, 0, rank(), cols));
673 for(
Index i = 0; i < rank(); ++i)
675 if(i) m.row(i).head(i).setZero();
676 m.row(i).tail(cols-i) = dec().matrixLU().row(pivots.coeff(i)).tail(cols-i);
678 m.block(0, 0, rank(), rank());
679 m.block(0, 0, rank(), rank()).template triangularView<StrictlyLower>().setZero();
680 for(
Index i = 0; i < rank(); ++i)
681 m.col(i).swap(m.col(pivots.coeff(i)));
686 m.topLeftCorner(rank(), rank())
687 .template triangularView<Upper>().solveInPlace(
688 m.topRightCorner(rank(), dimker)
692 for(
Index i = rank()-1; i >= 0; --i)
693 m.col(i).swap(m.col(pivots.coeff(i)));
696 for(
Index i = 0; i < rank(); ++i) dst.row(dec().permutationQ().indices().coeff(i)) = -m.row(i).tail(dimker);
697 for(
Index i = rank(); i < cols; ++i) dst.row(dec().permutationQ().indices().coeff(i)).setZero();
698 for(
Index k = 0; k < dimker; ++k) dst.coeffRef(dec().permutationQ().indices().coeff(rank()+k), k) = Scalar(1);
704 template<
typename _MatrixType>
705 struct image_retval<FullPivLU<_MatrixType> >
706 : image_retval_base<FullPivLU<_MatrixType> >
708 EIGEN_MAKE_IMAGE_HELPERS(FullPivLU<_MatrixType>)
710 enum { MaxSmallDimAtCompileTime = EIGEN_SIZE_MIN_PREFER_FIXED(
711 MatrixType::MaxColsAtCompileTime,
712 MatrixType::MaxRowsAtCompileTime)
715 template<
typename Dest>
void evalTo(Dest& dst)
const 727 Matrix<Index, Dynamic, 1, 0, MaxSmallDimAtCompileTime, 1> pivots(rank());
728 RealScalar premultiplied_threshold = dec().maxPivot() * dec().threshold();
730 for(
Index i = 0; i < dec().nonzeroPivots(); ++i)
731 if(
abs(dec().matrixLU().coeff(i,i)) > premultiplied_threshold)
732 pivots.coeffRef(p++) = i;
733 eigen_internal_assert(p == rank());
735 for(
Index i = 0; i < rank(); ++i)
736 dst.col(i) = originalMatrix().col(dec().permutationQ().indices().coeff(pivots.coeff(i)));
744 #ifndef EIGEN_PARSED_BY_DOXYGEN 745 template<
typename _MatrixType>
746 template<
typename RhsType,
typename DstType>
747 void FullPivLU<_MatrixType>::_solve_impl(
const RhsType &rhs, DstType &dst)
const 757 const Index rows = this->rows(),
759 nonzero_pivots = this->rank();
760 eigen_assert(rhs.rows() == rows);
761 const Index smalldim = (std::min)(rows, cols);
763 if(nonzero_pivots == 0)
769 typename RhsType::PlainObject c(rhs.rows(), rhs.cols());
772 c = permutationP() * rhs;
775 m_lu.topLeftCorner(smalldim,smalldim)
776 .template triangularView<UnitLower>()
777 .solveInPlace(c.topRows(smalldim));
779 c.bottomRows(rows-cols) -= m_lu.bottomRows(rows-cols) * c.topRows(cols);
782 m_lu.topLeftCorner(nonzero_pivots, nonzero_pivots)
783 .template triangularView<Upper>()
784 .solveInPlace(c.topRows(nonzero_pivots));
787 for(
Index i = 0; i < nonzero_pivots; ++i)
788 dst.row(permutationQ().indices().coeff(i)) = c.row(i);
789 for(
Index i = nonzero_pivots; i < m_lu.cols(); ++i)
790 dst.row(permutationQ().indices().coeff(i)).setZero();
793 template<
typename _MatrixType>
794 template<
bool Conjugate,
typename RhsType,
typename DstType>
795 void FullPivLU<_MatrixType>::_solve_impl_transposed(
const RhsType &rhs, DstType &dst)
const 808 const Index rows = this->rows(), cols = this->cols(),
809 nonzero_pivots = this->rank();
810 eigen_assert(rhs.rows() == cols);
811 const Index smalldim = (std::min)(rows, cols);
813 if(nonzero_pivots == 0)
819 typename RhsType::PlainObject c(rhs.rows(), rhs.cols());
822 c = permutationQ().inverse() * rhs;
826 m_lu.topLeftCorner(nonzero_pivots, nonzero_pivots)
827 .template triangularView<Upper>()
829 .solveInPlace(c.topRows(nonzero_pivots));
831 m_lu.topLeftCorner(smalldim, smalldim)
832 .template triangularView<UnitLower>()
834 .solveInPlace(c.topRows(smalldim));
837 m_lu.topLeftCorner(nonzero_pivots, nonzero_pivots)
838 .template triangularView<Upper>()
840 .solveInPlace(c.topRows(nonzero_pivots));
842 m_lu.topLeftCorner(smalldim, smalldim)
843 .template triangularView<UnitLower>()
845 .solveInPlace(c.topRows(smalldim));
849 PermutationPType invp = permutationP().inverse().eval();
850 for(
Index i = 0; i < smalldim; ++i)
851 dst.row(invp.indices().coeff(i)) = c.row(i);
852 for(
Index i = smalldim; i < rows; ++i)
853 dst.row(invp.indices().coeff(i)).setZero();
862 template<
typename DstXprType,
typename MatrixType>
863 struct Assignment<DstXprType, Inverse<FullPivLU<MatrixType> >, internal::assign_op<typename DstXprType::Scalar,typename FullPivLU<MatrixType>::Scalar>, Dense2Dense>
865 typedef FullPivLU<MatrixType> LuType;
866 typedef Inverse<LuType> SrcXprType;
867 static void run(DstXprType &dst,
const SrcXprType &src,
const internal::assign_op<typename DstXprType::Scalar,typename MatrixType::Scalar> &)
869 dst = src.nestedExpression().solve(MatrixType::Identity(src.rows(), src.cols()));
882 template<
typename Derived>
883 inline const FullPivLU<typename MatrixBase<Derived>::PlainObject>
const internal::image_retval< FullPivLU > image(const MatrixType &originalMatrix) const
Definition: FullPivLU.h:215
FullPivLU & setThreshold(const RealScalar &threshold)
Definition: FullPivLU.h:292
FullPivLU & compute(const EigenBase< InputType > &matrix)
Definition: FullPivLU.h:119
const PermutationQType & permutationQ() const
Definition: FullPivLU.h:169
RealScalar rcond() const
Definition: FullPivLU.h:252
Namespace containing all symbols from the Eigen library.
Definition: Core:306
Holds information about the various numeric (i.e. scalar) types allowed by Eigen.
Definition: NumTraits.h:150
Derived & derived()
Definition: EigenBase.h:45
Eigen::Index Index
The interface type of indices.
Definition: EigenBase.h:38
const FullPivLU< PlainObject > fullPivLu() const
Definition: FullPivLU.h:884
Index nonzeroPivots() const
Definition: FullPivLU.h:144
Definition: EigenBase.h:29
Index rank() const
Definition: FullPivLU.h:332
Expression of the inverse of another expression.
Definition: Inverse.h:43
RealScalar threshold() const
Definition: FullPivLU.h:317
internal::traits< MatrixType >::Scalar determinant() const
Definition: FullPivLU.h:583
const MatrixType & matrixLU() const
Definition: FullPivLU.h:131
Derived & derived()
Definition: EigenBase.h:45
FullPivLU & setThreshold(Default_t)
Definition: FullPivLU.h:307
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:33
const Inverse< FullPivLU > inverse() const
Definition: FullPivLU.h:400
bool isInjective() const
Definition: FullPivLU.h:362
bool isSurjective() const
Definition: FullPivLU.h:375
MatrixType reconstructedMatrix() const
Definition: FullPivLU.h:594
bool isInvertible() const
Definition: FullPivLU.h:387
const PermutationPType & permutationP() const
Definition: FullPivLU.h:159
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_abs_op< typename Derived::Scalar >, const Derived > abs(const Eigen::ArrayBase< Derived > &x)
const internal::kernel_retval< FullPivLU > kernel() const
Definition: FullPivLU.h:189
const Solve< FullPivLU, Rhs > solve(const MatrixBase< Rhs > &b) const
Definition: FullPivLU.h:243
LU decomposition of a matrix with complete pivoting, and related features.
Definition: ForwardDeclarations.h:249
FullPivLU()
Default Constructor.
Definition: FullPivLU.h:444
const int Dynamic
Definition: Constants.h:21
Pseudo expression representing a solving operation.
Definition: Solve.h:62
Index dimensionOfKernel() const
Definition: FullPivLU.h:349
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:48
RealScalar maxPivot() const
Definition: FullPivLU.h:153