13 #include "./InternalHeaderCheck.h"
18 template<
typename MatrixType_>
struct traits<FullPivLU<MatrixType_> >
21 typedef MatrixXpr XprKind;
22 typedef SolverStorage StorageKind;
23 typedef int StorageIndex;
66 typedef MatrixType_ MatrixType;
72 MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
73 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
75 typedef typename internal::plain_row_type<MatrixType, StorageIndex>::type IntRowVectorType;
76 typedef typename internal::plain_col_type<MatrixType, StorageIndex>::type IntColVectorType;
79 typedef typename MatrixType::PlainObject PlainObject;
102 template<
typename InputType>
111 template<
typename InputType>
121 template<
typename InputType>
136 eigen_assert(m_isInitialized &&
"LU is not initialized.");
149 eigen_assert(m_isInitialized &&
"LU is not initialized.");
150 return m_nonzero_pivots;
164 eigen_assert(m_isInitialized &&
"LU is not initialized.");
174 eigen_assert(m_isInitialized &&
"LU is not initialized.");
192 inline const internal::kernel_retval<FullPivLU>
kernel()
const
194 eigen_assert(m_isInitialized &&
"LU is not initialized.");
195 return internal::kernel_retval<FullPivLU>(*
this);
217 inline const internal::image_retval<FullPivLU>
218 image(
const MatrixType& originalMatrix)
const
220 eigen_assert(m_isInitialized &&
"LU is not initialized.");
221 return internal::image_retval<FullPivLU>(*
this, originalMatrix);
224 #ifdef EIGEN_PARSED_BY_DOXYGEN
244 template<
typename Rhs>
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 EIGEN_CONSTEXPR
410 inline Index rows() const EIGEN_NOEXCEPT {
return m_lu.rows(); }
411 EIGEN_DEVICE_FUNC EIGEN_CONSTEXPR
412 inline Index cols() const EIGEN_NOEXCEPT {
return m_lu.cols(); }
414 #ifndef EIGEN_PARSED_BY_DOXYGEN
415 template<
typename RhsType,
typename DstType>
416 void _solve_impl(
const RhsType &rhs, DstType &dst)
const;
418 template<
bool Conjugate,
typename RhsType,
typename DstType>
419 void _solve_impl_transposed(
const RhsType &rhs, DstType &dst)
const;
424 EIGEN_STATIC_ASSERT_NON_INTEGER(Scalar)
426 void computeInPlace();
429 PermutationPType m_p;
430 PermutationQType m_q;
431 IntColVectorType m_rowsTranspositions;
432 IntRowVectorType m_colsTranspositions;
433 Index m_nonzero_pivots;
434 RealScalar m_l1_norm;
435 RealScalar m_maxpivot, m_prescribedThreshold;
436 signed char m_det_pq;
437 bool m_isInitialized, m_usePrescribedThreshold;
440 template<
typename MatrixType>
442 : m_isInitialized(false), m_usePrescribedThreshold(false)
446 template<
typename MatrixType>
451 m_rowsTranspositions(rows),
452 m_colsTranspositions(cols),
453 m_isInitialized(false),
454 m_usePrescribedThreshold(false)
458 template<
typename MatrixType>
459 template<
typename InputType>
461 : m_lu(matrix.rows(), matrix.cols()),
464 m_rowsTranspositions(matrix.rows()),
465 m_colsTranspositions(matrix.cols()),
466 m_isInitialized(false),
467 m_usePrescribedThreshold(false)
472 template<
typename MatrixType>
473 template<
typename InputType>
475 : m_lu(matrix.derived()),
478 m_rowsTranspositions(matrix.rows()),
479 m_colsTranspositions(matrix.cols()),
480 m_isInitialized(false),
481 m_usePrescribedThreshold(false)
486 template<
typename MatrixType>
492 m_l1_norm = m_lu.cwiseAbs().colwise().sum().maxCoeff();
494 const Index size = m_lu.diagonalSize();
495 const Index rows = m_lu.rows();
496 const Index cols = m_lu.cols();
500 m_rowsTranspositions.resize(m_lu.rows());
501 m_colsTranspositions.resize(m_lu.cols());
502 Index number_of_transpositions = 0;
504 m_nonzero_pivots = size;
505 m_maxpivot = RealScalar(0);
507 for(
Index k = 0; k < size; ++k)
512 Index row_of_biggest_in_corner, col_of_biggest_in_corner;
513 typedef internal::scalar_score_coeff_op<Scalar> Scoring;
514 typedef typename Scoring::result_type Score;
515 Score biggest_in_corner;
516 biggest_in_corner = m_lu.bottomRightCorner(rows-k, cols-k)
517 .unaryExpr(Scoring())
518 .maxCoeff(&row_of_biggest_in_corner, &col_of_biggest_in_corner);
519 row_of_biggest_in_corner += k;
520 col_of_biggest_in_corner += k;
522 if(numext::is_exactly_zero(biggest_in_corner))
526 m_nonzero_pivots = k;
527 for(
Index i = k; i < size; ++i)
529 m_rowsTranspositions.coeffRef(i) = internal::convert_index<StorageIndex>(i);
530 m_colsTranspositions.coeffRef(i) = internal::convert_index<StorageIndex>(i);
535 RealScalar abs_pivot = internal::abs_knowing_score<Scalar>()(m_lu(row_of_biggest_in_corner, col_of_biggest_in_corner), biggest_in_corner);
536 if(abs_pivot > m_maxpivot) m_maxpivot = abs_pivot;
541 m_rowsTranspositions.coeffRef(k) = internal::convert_index<StorageIndex>(row_of_biggest_in_corner);
542 m_colsTranspositions.coeffRef(k) = internal::convert_index<StorageIndex>(col_of_biggest_in_corner);
543 if(k != row_of_biggest_in_corner) {
544 m_lu.row(k).swap(m_lu.row(row_of_biggest_in_corner));
545 ++number_of_transpositions;
547 if(k != col_of_biggest_in_corner) {
548 m_lu.col(k).swap(m_lu.col(col_of_biggest_in_corner));
549 ++number_of_transpositions;
556 m_lu.col(k).tail(rows-k-1) /= m_lu.coeff(k,k);
558 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);
564 m_p.setIdentity(rows);
565 for(
Index k = size-1; k >= 0; --k)
566 m_p.applyTranspositionOnTheRight(k, m_rowsTranspositions.coeff(k));
568 m_q.setIdentity(cols);
569 for(
Index k = 0; k < size; ++k)
570 m_q.applyTranspositionOnTheRight(k, m_colsTranspositions.coeff(k));
572 m_det_pq = (number_of_transpositions%2) ? -1 : 1;
574 m_isInitialized =
true;
577 template<
typename MatrixType>
580 eigen_assert(m_isInitialized &&
"LU is not initialized.");
581 eigen_assert(m_lu.rows() == m_lu.cols() &&
"You can't take the determinant of a non-square matrix!");
582 return Scalar(m_det_pq) * Scalar(m_lu.diagonal().prod());
588 template<
typename MatrixType>
591 eigen_assert(m_isInitialized &&
"LU is not initialized.");
592 const Index smalldim = (std::min)(m_lu.rows(), m_lu.cols());
594 MatrixType res(m_lu.rows(),m_lu.cols());
596 res = m_lu.leftCols(smalldim)
597 .template triangularView<UnitLower>().toDenseMatrix()
598 * m_lu.topRows(smalldim)
599 .template triangularView<Upper>().toDenseMatrix();
602 res = m_p.inverse() * res;
605 res = res * m_q.inverse();
613 template<
typename MatrixType_>
614 struct kernel_retval<
FullPivLU<MatrixType_> >
615 : kernel_retval_base<FullPivLU<MatrixType_> >
619 enum { MaxSmallDimAtCompileTime = min_size_prefer_fixed(
620 MatrixType::MaxColsAtCompileTime,
621 MatrixType::MaxRowsAtCompileTime)
624 template<
typename Dest>
void evalTo(Dest& dst)
const
627 const Index cols = dec().matrixLU().cols(), dimker = cols - rank();
653 Matrix<Index, Dynamic, 1, 0, MaxSmallDimAtCompileTime, 1> pivots(rank());
654 RealScalar premultiplied_threshold = dec().maxPivot() * dec().threshold();
656 for(
Index i = 0; i < dec().nonzeroPivots(); ++i)
657 if(
abs(dec().matrixLU().coeff(i,i)) > premultiplied_threshold)
658 pivots.coeffRef(p++) = i;
659 eigen_internal_assert(p == rank());
665 Matrix<
typename MatrixType::Scalar,
Dynamic,
Dynamic, MatrixType::Options,
666 MaxSmallDimAtCompileTime, MatrixType::MaxColsAtCompileTime>
667 m(dec().matrixLU().block(0, 0, rank(), cols));
668 for(
Index i = 0; i < rank(); ++i)
670 if(i) m.row(i).head(i).setZero();
671 m.row(i).tail(cols-i) = dec().matrixLU().row(pivots.coeff(i)).tail(cols-i);
673 m.block(0, 0, rank(), rank());
674 m.block(0, 0, rank(), rank()).template triangularView<StrictlyLower>().setZero();
675 for(
Index i = 0; i < rank(); ++i)
676 m.col(i).swap(m.col(pivots.coeff(i)));
681 m.topLeftCorner(rank(), rank())
682 .
template triangularView<Upper>().solveInPlace(
683 m.topRightCorner(rank(), dimker)
687 for(
Index i = rank()-1; i >= 0; --i)
688 m.col(i).swap(m.col(pivots.coeff(i)));
691 for(
Index i = 0; i < rank(); ++i) dst.row(dec().permutationQ().indices().coeff(i)) = -m.row(i).tail(dimker);
692 for(
Index i = rank(); i < cols; ++i) dst.row(dec().permutationQ().indices().coeff(i)).setZero();
693 for(
Index k = 0; k < dimker; ++k) dst.coeffRef(dec().permutationQ().indices().coeff(rank()+k), k) = Scalar(1);
699 template<
typename MatrixType_>
700 struct image_retval<FullPivLU<MatrixType_> >
701 : image_retval_base<FullPivLU<MatrixType_> >
703 EIGEN_MAKE_IMAGE_HELPERS(FullPivLU<MatrixType_>)
705 enum { MaxSmallDimAtCompileTime = min_size_prefer_fixed(
706 MatrixType::MaxColsAtCompileTime,
707 MatrixType::MaxRowsAtCompileTime)
710 template<
typename Dest>
void evalTo(Dest& dst)
const
722 Matrix<Index, Dynamic, 1, 0, MaxSmallDimAtCompileTime, 1> pivots(rank());
723 RealScalar premultiplied_threshold = dec().maxPivot() * dec().threshold();
725 for(
Index i = 0; i < dec().nonzeroPivots(); ++i)
726 if(
abs(dec().matrixLU().coeff(i,i)) > premultiplied_threshold)
727 pivots.coeffRef(p++) = i;
728 eigen_internal_assert(p == rank());
730 for(
Index i = 0; i < rank(); ++i)
731 dst.col(i) = originalMatrix().col(dec().permutationQ().indices().coeff(pivots.coeff(i)));
739 #ifndef EIGEN_PARSED_BY_DOXYGEN
740 template<
typename MatrixType_>
741 template<
typename RhsType,
typename DstType>
742 void FullPivLU<MatrixType_>::_solve_impl(
const RhsType &rhs, DstType &dst)
const
752 const Index rows = this->rows(),
754 nonzero_pivots = this->rank();
755 const Index smalldim = (std::min)(rows, cols);
757 if(nonzero_pivots == 0)
763 typename RhsType::PlainObject c(rhs.rows(), rhs.cols());
766 c = permutationP() * rhs;
769 m_lu.topLeftCorner(smalldim,smalldim)
770 .template triangularView<UnitLower>()
771 .solveInPlace(c.topRows(smalldim));
773 c.bottomRows(rows-cols) -= m_lu.bottomRows(rows-cols) * c.topRows(cols);
776 m_lu.topLeftCorner(nonzero_pivots, nonzero_pivots)
777 .template triangularView<Upper>()
778 .solveInPlace(c.topRows(nonzero_pivots));
781 for(
Index i = 0; i < nonzero_pivots; ++i)
782 dst.row(permutationQ().indices().coeff(i)) = c.row(i);
783 for(
Index i = nonzero_pivots; i < m_lu.cols(); ++i)
784 dst.row(permutationQ().indices().coeff(i)).setZero();
787 template<
typename MatrixType_>
788 template<
bool Conjugate,
typename RhsType,
typename DstType>
789 void FullPivLU<MatrixType_>::_solve_impl_transposed(
const RhsType &rhs, DstType &dst)
const
802 const Index rows = this->rows(), cols = this->cols(),
803 nonzero_pivots = this->rank();
804 const Index smalldim = (std::min)(rows, cols);
806 if(nonzero_pivots == 0)
812 typename RhsType::PlainObject c(rhs.rows(), rhs.cols());
815 c = permutationQ().inverse() * rhs;
818 m_lu.topLeftCorner(nonzero_pivots, nonzero_pivots)
819 .template triangularView<Upper>()
821 .template conjugateIf<Conjugate>()
822 .solveInPlace(c.topRows(nonzero_pivots));
825 m_lu.topLeftCorner(smalldim, smalldim)
826 .template triangularView<UnitLower>()
828 .template conjugateIf<Conjugate>()
829 .solveInPlace(c.topRows(smalldim));
832 PermutationPType invp = permutationP().inverse().eval();
833 for(
Index i = 0; i < smalldim; ++i)
834 dst.row(invp.indices().coeff(i)) = c.row(i);
835 for(
Index i = smalldim; i < rows; ++i)
836 dst.row(invp.indices().coeff(i)).setZero();
845 template<
typename DstXprType,
typename MatrixType>
846 struct Assignment<DstXprType, Inverse<FullPivLU<MatrixType> >, internal::assign_op<typename DstXprType::Scalar,typename FullPivLU<MatrixType>::Scalar>, Dense2Dense>
848 typedef FullPivLU<MatrixType> LuType;
849 typedef Inverse<LuType> SrcXprType;
850 static void run(DstXprType &dst,
const SrcXprType &src,
const internal::assign_op<typename DstXprType::Scalar,typename MatrixType::Scalar> &)
852 dst = src.nestedExpression().
solve(MatrixType::Identity(src.rows(), src.cols()));
865 template<
typename Derived>
866 inline const FullPivLU<typename MatrixBase<Derived>::PlainObject>
LU decomposition of a matrix with complete pivoting, and related features.
Definition: FullPivLU.h:64
RealScalar maxPivot() const
Definition: FullPivLU.h:156
bool isInvertible() const
Definition: FullPivLU.h:387
MatrixType reconstructedMatrix() const
Definition: FullPivLU.h:589
const PermutationQType & permutationQ() const
Definition: FullPivLU.h:172
const Inverse< FullPivLU > inverse() const
Definition: FullPivLU.h:400
Index dimensionOfKernel() const
Definition: FullPivLU.h:349
FullPivLU & compute(const EigenBase< InputType > &matrix)
Definition: FullPivLU.h:122
const MatrixType & matrixLU() const
Definition: FullPivLU.h:134
const PermutationPType & permutationP() const
Definition: FullPivLU.h:162
const Solve< FullPivLU, Rhs > solve(const MatrixBase< Rhs > &b) const
internal::traits< MatrixType >::Scalar determinant() const
Definition: FullPivLU.h:578
Index rank() const
Definition: FullPivLU.h:332
bool isInjective() const
Definition: FullPivLU.h:362
const internal::image_retval< FullPivLU > image(const MatrixType &originalMatrix) const
Definition: FullPivLU.h:218
Index nonzeroPivots() const
Definition: FullPivLU.h:147
FullPivLU & setThreshold(const RealScalar &threshold)
Definition: FullPivLU.h:292
const internal::kernel_retval< FullPivLU > kernel() const
Definition: FullPivLU.h:192
RealScalar threshold() const
Definition: FullPivLU.h:317
RealScalar rcond() const
Definition: FullPivLU.h:252
FullPivLU & setThreshold(Default_t)
Definition: FullPivLU.h:307
FullPivLU()
Default Constructor.
Definition: FullPivLU.h:441
bool isSurjective() const
Definition: FullPivLU.h:375
Expression of the inverse of another expression.
Definition: Inverse.h:46
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:52
const FullPivLU< PlainObject > fullPivLu() const
Definition: FullPivLU.h:867
Pseudo expression representing a solving operation.
Definition: Solve.h:65
A base class for matrix decomposition and solvers.
Definition: SolverBase.h:71
const Solve< Derived, Rhs > solve(const MatrixBase< Rhs > &b) const
Definition: SolverBase.h:107
Namespace containing all symbols from the Eigen library.
Definition: Core:139
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:59
const int Dynamic
Definition: Constants.h:24
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_abs_op< typename Derived::Scalar >, const Derived > abs(const Eigen::ArrayBase< Derived > &x)
Definition: EigenBase.h:32
Derived & derived()
Definition: EigenBase.h:48
Eigen::Index Index
The interface type of indices.
Definition: EigenBase.h:41
Holds information about the various numeric (i.e. scalar) types allowed by Eigen.
Definition: NumTraits.h:231