11 #ifndef EIGEN_PARTIALLU_H
12 #define EIGEN_PARTIALLU_H
14 #include "./InternalHeaderCheck.h"
19 template<
typename MatrixType_>
struct traits<PartialPivLU<MatrixType_> >
22 typedef MatrixXpr XprKind;
23 typedef SolverStorage StorageKind;
24 typedef int StorageIndex;
25 typedef traits<MatrixType_> BaseTraits;
32 template<
typename T,
typename Derived>
38 template<
typename T,
typename Derived>
39 struct enable_if_ref<Ref<T>,Derived> {
79 :
public SolverBase<PartialPivLU<MatrixType_> >
83 typedef MatrixType_ MatrixType;
89 MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
90 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
94 typedef typename MatrixType::PlainObject PlainObject;
119 template<
typename InputType>
129 template<
typename InputType>
132 template<
typename InputType>
147 eigen_assert(m_isInitialized &&
"PartialPivLU is not initialized.");
155 eigen_assert(m_isInitialized &&
"PartialPivLU is not initialized.");
159 #ifdef EIGEN_PARSED_BY_DOXYGEN
177 template<
typename Rhs>
187 eigen_assert(m_isInitialized &&
"PartialPivLU is not initialized.");
188 return internal::rcond_estimate_helper(m_l1_norm, *
this);
200 eigen_assert(m_isInitialized &&
"PartialPivLU is not initialized.");
221 EIGEN_CONSTEXPR
inline Index rows() const EIGEN_NOEXCEPT {
return m_lu.rows(); }
222 EIGEN_CONSTEXPR
inline Index cols() const EIGEN_NOEXCEPT {
return m_lu.cols(); }
224 #ifndef EIGEN_PARSED_BY_DOXYGEN
225 template<
typename RhsType,
typename DstType>
227 void _solve_impl(
const RhsType &rhs, DstType &dst)
const {
239 m_lu.template triangularView<UnitLower>().solveInPlace(dst);
242 m_lu.template triangularView<Upper>().solveInPlace(dst);
245 template<
bool Conjugate,
typename RhsType,
typename DstType>
247 void _solve_impl_transposed(
const RhsType &rhs, DstType &dst)
const {
255 eigen_assert(rhs.rows() == m_lu.cols());
258 dst = m_lu.template triangularView<Upper>().transpose()
259 .template conjugateIf<Conjugate>().solve(rhs);
261 m_lu.template triangularView<UnitLower>().transpose()
262 .template conjugateIf<Conjugate>().solveInPlace(dst);
270 EIGEN_STATIC_ASSERT_NON_INTEGER(Scalar)
276 TranspositionType m_rowsTranspositions;
277 RealScalar m_l1_norm;
279 bool m_isInitialized;
282 template<
typename MatrixType>
286 m_rowsTranspositions(),
289 m_isInitialized(false)
293 template<
typename MatrixType>
297 m_rowsTranspositions(size),
300 m_isInitialized(false)
304 template<
typename MatrixType>
305 template<
typename InputType>
307 : m_lu(matrix.rows(),matrix.cols()),
309 m_rowsTranspositions(matrix.rows()),
312 m_isInitialized(false)
317 template<
typename MatrixType>
318 template<
typename InputType>
320 : m_lu(matrix.derived()),
322 m_rowsTranspositions(matrix.rows()),
325 m_isInitialized(false)
333 template<
typename Scalar,
int StorageOrder,
typename PivIndex,
int SizeAtCompileTime=Dynamic>
334 struct partial_lu_impl
336 static constexpr
int UnBlockedBound = 16;
337 static constexpr
bool UnBlockedAtCompileTime = SizeAtCompileTime!=
Dynamic && SizeAtCompileTime<=UnBlockedBound;
338 static constexpr
int ActualSizeAtCompileTime = UnBlockedAtCompileTime ? SizeAtCompileTime :
Dynamic;
340 static constexpr
int RRows = SizeAtCompileTime==2 ? 1 :
Dynamic;
341 static constexpr
int RCols = SizeAtCompileTime==2 ? 1 :
Dynamic;
345 typedef typename MatrixType::RealScalar RealScalar;
357 static Index unblocked_lu(MatrixTypeRef& lu, PivIndex* row_transpositions, PivIndex& nb_transpositions)
359 typedef scalar_score_coeff_op<Scalar> Scoring;
360 typedef typename Scoring::result_type Score;
361 const Index rows = lu.rows();
362 const Index cols = lu.cols();
363 const Index size = (std::min)(rows,cols);
366 const Index endk = UnBlockedAtCompileTime ? size-1 : size;
367 nb_transpositions = 0;
368 Index first_zero_pivot = -1;
369 for(
Index k = 0; k < endk; ++k)
371 int rrows = internal::convert_index<int>(rows-k-1);
372 int rcols = internal::convert_index<int>(cols-k-1);
374 Index row_of_biggest_in_col;
375 Score biggest_in_corner
376 = lu.col(k).tail(rows-k).unaryExpr(Scoring()).maxCoeff(&row_of_biggest_in_col);
377 row_of_biggest_in_col += k;
379 row_transpositions[k] = PivIndex(row_of_biggest_in_col);
381 if(!numext::is_exactly_zero(biggest_in_corner))
383 if(k != row_of_biggest_in_col)
385 lu.row(k).swap(lu.row(row_of_biggest_in_col));
389 lu.col(k).tail(fix<RRows>(rrows)) /= lu.coeff(k,k);
391 else if(first_zero_pivot==-1)
395 first_zero_pivot = k;
399 lu.bottomRightCorner(fix<RRows>(rrows),fix<RCols>(rcols)).noalias() -= lu.col(k).tail(fix<RRows>(rrows)) * lu.row(k).tail(fix<RCols>(rcols));
403 if(UnBlockedAtCompileTime)
406 row_transpositions[k] = PivIndex(k);
407 if (numext::is_exactly_zero(Scoring()(lu(k, k))) && first_zero_pivot == -1)
408 first_zero_pivot = k;
411 return first_zero_pivot;
429 static Index blocked_lu(
Index rows,
Index cols, Scalar* lu_data,
Index luStride, PivIndex* row_transpositions, PivIndex& nb_transpositions,
Index maxBlockSize=256)
431 MatrixTypeRef lu =
MatrixType::Map(lu_data,rows, cols, OuterStride<>(luStride));
433 const Index size = (std::min)(rows,cols);
436 if(UnBlockedAtCompileTime || size<=UnBlockedBound)
438 return unblocked_lu(lu, row_transpositions, nb_transpositions);
446 blockSize = (blockSize/16)*16;
447 blockSize = (std::min)((std::max)(blockSize,
Index(8)), maxBlockSize);
450 nb_transpositions = 0;
451 Index first_zero_pivot = -1;
452 for(
Index k = 0; k < size; k+=blockSize)
454 Index bs = (std::min)(size-k,blockSize);
455 Index trows = rows - k - bs;
456 Index tsize = size - k - bs;
462 BlockType A_0 = lu.block(0,0,rows,k);
463 BlockType A_2 = lu.block(0,k+bs,rows,tsize);
464 BlockType A11 = lu.block(k,k,bs,bs);
465 BlockType A12 = lu.block(k,k+bs,bs,tsize);
466 BlockType A21 = lu.block(k+bs,k,trows,bs);
467 BlockType A22 = lu.block(k+bs,k+bs,trows,tsize);
469 PivIndex nb_transpositions_in_panel;
472 Index ret = blocked_lu(trows+bs, bs, &lu.coeffRef(k,k), luStride,
473 row_transpositions+k, nb_transpositions_in_panel, 16);
474 if(ret>=0 && first_zero_pivot==-1)
475 first_zero_pivot = k+ret;
477 nb_transpositions += nb_transpositions_in_panel;
479 for(
Index i=k; i<k+bs; ++i)
481 Index piv = (row_transpositions[i] += internal::convert_index<PivIndex>(k));
482 A_0.row(i).swap(A_0.row(piv));
488 for(
Index i=k;i<k+bs; ++i)
489 A_2.row(i).swap(A_2.row(row_transpositions[i]));
492 A11.template triangularView<UnitLower>().solveInPlace(A12);
494 A22.noalias() -= A21 * A12;
497 return first_zero_pivot;
503 template<
typename MatrixType,
typename TranspositionType>
504 void partial_lu_inplace(MatrixType& lu, TranspositionType& row_transpositions,
typename TranspositionType::StorageIndex& nb_transpositions)
507 if (lu.rows() == 0 || lu.cols() == 0) {
508 nb_transpositions = 0;
511 eigen_assert(lu.cols() == row_transpositions.size());
512 eigen_assert(row_transpositions.size() < 2 || (&row_transpositions.coeffRef(1)-&row_transpositions.coeffRef(0)) == 1);
516 typename TranspositionType::StorageIndex,
517 internal::min_size_prefer_fixed(MatrixType::RowsAtCompileTime, MatrixType::ColsAtCompileTime)>
518 ::blocked_lu(lu.rows(), lu.cols(), &lu.coeffRef(0,0), lu.outerStride(), &row_transpositions.coeffRef(0), nb_transpositions);
523 template<
typename MatrixType>
524 void PartialPivLU<MatrixType>::compute()
527 eigen_assert(m_lu.rows()<NumTraits<int>::highest());
530 m_l1_norm = m_lu.cwiseAbs().colwise().sum().maxCoeff();
532 m_l1_norm = RealScalar(0);
534 eigen_assert(m_lu.rows() == m_lu.cols() &&
"PartialPivLU is only for square (and moreover invertible) matrices");
535 const Index size = m_lu.rows();
537 m_rowsTranspositions.resize(size);
539 typename TranspositionType::StorageIndex nb_transpositions;
540 internal::partial_lu_inplace(m_lu, m_rowsTranspositions, nb_transpositions);
541 m_det_p = (nb_transpositions%2) ? -1 : 1;
543 m_p = m_rowsTranspositions;
545 m_isInitialized =
true;
548 template<
typename MatrixType>
551 eigen_assert(m_isInitialized &&
"PartialPivLU is not initialized.");
552 return Scalar(m_det_p) * m_lu.diagonal().prod();
558 template<
typename MatrixType>
561 eigen_assert(m_isInitialized &&
"LU is not initialized.");
563 MatrixType res = m_lu.template triangularView<UnitLower>().toDenseMatrix()
564 * m_lu.template triangularView<Upper>();
567 res = m_p.inverse() * res;
577 template<
typename DstXprType,
typename MatrixType>
578 struct Assignment<DstXprType,
Inverse<
PartialPivLU<MatrixType> >, internal::assign_op<typename DstXprType::Scalar,typename PartialPivLU<MatrixType>::Scalar>, Dense2Dense>
582 static void run(DstXprType &dst,
const SrcXprType &src,
const internal::assign_op<typename DstXprType::Scalar,typename LuType::Scalar> &)
584 dst = src.nestedExpression().solve(MatrixType::Identity(src.rows(), src.cols()));
597 template<
typename Derived>
598 inline const PartialPivLU<typename MatrixBase<Derived>::PlainObject>
612 template<
typename Derived>
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 PartialPivLU< PlainObject > partialPivLu() const
Definition: PartialPivLU.h:599
const PartialPivLU< PlainObject > lu() const
Definition: PartialPivLU.h:614
The matrix class, also used for vectors and row-vectors.
Definition: Matrix.h:182
LU decomposition of a matrix with partial pivoting, and related features.
Definition: PartialPivLU.h:80
const Solve< PartialPivLU, Rhs > solve(const MatrixBase< Rhs > &b) const
RealScalar rcond() const
Definition: PartialPivLU.h:185
Scalar determinant() const
Definition: PartialPivLU.h:549
PartialPivLU()
Default Constructor.
Definition: PartialPivLU.h:283
const PermutationType & permutationP() const
Definition: PartialPivLU.h:153
const Inverse< PartialPivLU > inverse() const
Definition: PartialPivLU.h:198
MatrixType reconstructedMatrix() const
Definition: PartialPivLU.h:559
const MatrixType & matrixLU() const
Definition: PartialPivLU.h:145
InverseReturnType transpose() const
Definition: PermutationMatrix.h:193
static ConstMapType Map(const Scalar *data)
Definition: PlainObjectBase.h:642
A matrix or vector expression mapping an existing expression.
Definition: Ref.h:285
Pseudo expression representing a solving operation.
Definition: Solve.h:65
A base class for matrix decomposition and solvers.
Definition: SolverBase.h:71
@ ColMajor
Definition: Constants.h:321
@ RowMajor
Definition: Constants.h:323
const unsigned int RowMajorBit
Definition: Constants.h:68
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
Definition: EigenBase.h:32
Derived & derived()
Definition: EigenBase.h:48
Eigen::Index Index
The interface type of indices.
Definition: EigenBase.h:41
EIGEN_CONSTEXPR Index size() const EIGEN_NOEXCEPT
Definition: EigenBase.h:69