45 template<
typename _MatrixType>
class FullPivLU
56 typedef typename MatrixType::Scalar
Scalar;
59 typedef typename MatrixType::Index
Index;
187 image(
const MatrixType& originalMatrix)
const 212 template<
typename Rhs>
301 result += (
abs(
m_lu.coeff(i,i)) > premultiplied_threshold);
365 eigen_assert(
m_lu.rows() ==
m_lu.cols() &&
"You can't take the inverse of a non-square matrix!");
367 (*
this, MatrixType::Identity(
m_lu.rows(),
m_lu.cols()));
386 template<
typename MatrixType>
392 template<
typename MatrixType>
404 template<
typename MatrixType>
417 template<
typename MatrixType>
426 const Index size = matrix.diagonalSize();
434 Index number_of_transpositions = 0;
439 for(
Index k = 0; k < size; ++k)
444 Index row_of_biggest_in_corner, col_of_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;
457 for(
Index i = k; i < size; ++i)
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);
494 for(
Index k = size-1; k >= 0; --k)
498 for(
Index k = 0; k < size; ++k)
501 m_det_pq = (number_of_transpositions%2) ? -1 : 1;
505 template<
typename MatrixType>
509 eigen_assert(
m_lu.rows() ==
m_lu.cols() &&
"You can't take the determinant of a non-square matrix!");
516 template<
typename MatrixType>
524 res =
m_lu.leftCols(smalldim)
525 .template triangularView<UnitLower>().toDenseMatrix()
526 *
m_lu.topRows(smalldim)
527 .template triangularView<Upper>().toDenseMatrix();
541 template<
typename _MatrixType>
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();
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)
594 MaxSmallDimAtCompileTime, MatrixType::MaxColsAtCompileTime>
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);
602 m.block(0, 0,
rank(),
rank()).template triangularView<StrictlyLower>().setZero();
604 m.col(i).swap(m.col(pivots.
coeff(i)));
610 .template triangularView<Upper>().solveInPlace(
611 m.topRightCorner(rank(), dimker)
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);
627 template<
typename _MatrixType>
634 MatrixType::MaxColsAtCompileTime,
635 MatrixType::MaxRowsAtCompileTime)
638 template<
typename Dest>
void evalTo(Dest& dst)
const 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)
659 dst.col(i) = originalMatrix().col(dec().
permutationQ().indices().coeff(pivots.
coeff(i)));
665 template<
typename _MatrixType,
typename Rhs>
671 template<typename Dest>
void evalTo(Dest& dst)
const 682 nonzero_pivots = dec().nonzeroPivots();
684 const Index smalldim = (std::min)(rows,
cols);
686 if(nonzero_pivots == 0)
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>
PermutationMatrix< ColsAtCompileTime, MaxColsAtCompileTime > PermutationQType
bool m_usePrescribedThreshold
internal::traits< MatrixType >::StorageKind StorageKind
const internal::image_retval< FullPivLU > image(const MatrixType &originalMatrix) const
FullPivLU & setThreshold(const RealScalar &threshold)
RealScalar m_prescribedThreshold
internal::traits< MatrixType >::Scalar determinant() const
Derived & applyTranspositionOnTheRight(Index i, Index j)
iterative scaling algorithm to equilibrate rows and column norms in matrices
Block< Derived > block(Index startRow, Index startCol, Index blockRows, Index blockCols)
Holds information about the various numeric (i.e. scalar) types allowed by Eigen. ...
#define eigen_internal_assert(x)
#define EIGEN_SIZE_MIN_PREFER_FIXED(a, b)
RealScalar threshold() const
void evalTo(Dest &dst) const
Index nonzeroPivots() const
EIGEN_STRONG_INLINE const CwiseUnaryOp< internal::scalar_abs_op< Scalar >, const Derived > abs() const
bool isSurjective() const
MatrixType::Scalar Scalar
EIGEN_STRONG_INLINE const Scalar & coeff(Index rowId, Index colId) const
IntRowVectorType m_colsTranspositions
FullPivLU & setThreshold(Default_t)
const PermutationQType & permutationQ() const
EIGEN_STRONG_INLINE Scalar & coeffRef(Index rowId, Index colId)
SegmentReturnType tail(Index vecSize)
const internal::solve_retval< FullPivLU, typename MatrixType::IdentityReturnType > inverse() const
Provides a generic way to set and pass user-specified options.
NumTraits< typename MatrixType::Scalar >::Real RealScalar
FullPivLU & compute(const MatrixType &matrix)
internal::plain_col_type< MatrixType, Index >::type IntColVectorType
PermutationMatrix< RowsAtCompileTime, MaxRowsAtCompileTime > PermutationPType
MatrixType reconstructedMatrix() const
Transpose< PermutationBase > inverse() const
void rhs(const real_t *x, real_t *f)
const PermutationPType & permutationP() const
#define EIGEN_MAKE_KERNEL_HELPERS(DecompositionType)
LU decomposition of a matrix with complete pivoting, and related features.
bool isInvertible() const
const FullPivLU< PlainObject > fullPivLu() const
FullPivLU()
Default Constructor.
IntColVectorType m_rowsTranspositions
#define EIGEN_MAKE_IMAGE_HELPERS(DecompositionType)
internal::plain_row_type< MatrixType, Index >::type IntRowVectorType
#define EIGEN_MAKE_SOLVE_HELPERS(DecompositionType, Rhs)
const MatrixType & matrixLU() const
void evalTo(Dest &dst) const
The matrix class, also used for vectors and row-vectors.
Index dimensionOfKernel() const
RealScalar maxPivot() const
Base class for all dense matrices, vectors, and expressions.
const internal::solve_retval< FullPivLU, Rhs > solve(const MatrixBase< Rhs > &b) const
const internal::kernel_retval< FullPivLU > kernel() const