11#ifndef EIGEN_PARTIALLU_H
12#define EIGEN_PARTIALLU_H
17template<
typename _MatrixType>
struct traits<PartialPivLU<_MatrixType> >
20 typedef MatrixXpr XprKind;
21 typedef SolverStorage StorageKind;
22 typedef int StorageIndex;
23 typedef traits<_MatrixType> BaseTraits;
30template<
typename T,
typename Derived>
36template<
typename T,
typename Derived>
37struct enable_if_ref<Ref<T>,Derived> {
77 :
public SolverBase<PartialPivLU<_MatrixType> >
81 typedef _MatrixType MatrixType;
87 MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
88 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
92 typedef typename MatrixType::PlainObject PlainObject;
117 template<
typename InputType>
127 template<
typename InputType>
130 template<
typename InputType>
145 eigen_assert(m_isInitialized &&
"PartialPivLU is not initialized.");
153 eigen_assert(m_isInitialized &&
"PartialPivLU is not initialized.");
157 #ifdef EIGEN_PARSED_BY_DOXYGEN
175 template<
typename Rhs>
185 eigen_assert(m_isInitialized &&
"PartialPivLU is not initialized.");
186 return internal::rcond_estimate_helper(m_l1_norm, *
this);
198 eigen_assert(m_isInitialized &&
"PartialPivLU is not initialized.");
219 EIGEN_CONSTEXPR
inline Index rows() const EIGEN_NOEXCEPT {
return m_lu.rows(); }
220 EIGEN_CONSTEXPR
inline Index cols() const EIGEN_NOEXCEPT {
return m_lu.cols(); }
222 #ifndef EIGEN_PARSED_BY_DOXYGEN
223 template<
typename RhsType,
typename DstType>
225 void _solve_impl(
const RhsType &rhs, DstType &dst)
const {
237 m_lu.template triangularView<UnitLower>().solveInPlace(dst);
240 m_lu.template triangularView<Upper>().solveInPlace(dst);
243 template<
bool Conjugate,
typename RhsType,
typename DstType>
245 void _solve_impl_transposed(
const RhsType &rhs, DstType &dst)
const {
253 eigen_assert(rhs.rows() == m_lu.cols());
256 dst = m_lu.template triangularView<Upper>().transpose()
257 .template conjugateIf<Conjugate>().solve(rhs);
259 m_lu.template triangularView<UnitLower>().transpose()
260 .template conjugateIf<Conjugate>().solveInPlace(dst);
268 static void check_template_parameters()
270 EIGEN_STATIC_ASSERT_NON_INTEGER(Scalar);
277 TranspositionType m_rowsTranspositions;
278 RealScalar m_l1_norm;
280 bool m_isInitialized;
283template<
typename MatrixType>
287 m_rowsTranspositions(),
290 m_isInitialized(false)
294template<
typename MatrixType>
298 m_rowsTranspositions(size),
301 m_isInitialized(false)
305template<
typename MatrixType>
306template<
typename InputType>
308 : m_lu(matrix.rows(),matrix.cols()),
310 m_rowsTranspositions(matrix.rows()),
313 m_isInitialized(false)
318template<
typename MatrixType>
319template<
typename InputType>
321 : m_lu(matrix.derived()),
323 m_rowsTranspositions(matrix.rows()),
326 m_isInitialized(false)
334template<
typename Scalar,
int StorageOrder,
typename PivIndex,
int SizeAtCompileTime=Dynamic>
335struct partial_lu_impl
337 static const int UnBlockedBound = 16;
338 static const bool UnBlockedAtCompileTime = SizeAtCompileTime!=
Dynamic && SizeAtCompileTime<=UnBlockedBound;
339 static const int ActualSizeAtCompileTime = UnBlockedAtCompileTime ? SizeAtCompileTime :
Dynamic;
341 static const int RRows = SizeAtCompileTime==2 ? 1 :
Dynamic;
342 static const int RCols = SizeAtCompileTime==2 ? 1 :
Dynamic;
346 typedef typename MatrixType::RealScalar RealScalar;
358 static Index unblocked_lu(MatrixTypeRef& lu, PivIndex* row_transpositions, PivIndex& nb_transpositions)
360 typedef scalar_score_coeff_op<Scalar> Scoring;
361 typedef typename Scoring::result_type Score;
362 const Index rows = lu.rows();
363 const Index cols = lu.cols();
364 const Index size = (std::min)(rows,cols);
367 const Index endk = UnBlockedAtCompileTime ? size-1 : size;
368 nb_transpositions = 0;
369 Index first_zero_pivot = -1;
370 for(
Index k = 0; k < endk; ++k)
372 int rrows = internal::convert_index<int>(rows-k-1);
373 int rcols = internal::convert_index<int>(cols-k-1);
375 Index row_of_biggest_in_col;
376 Score biggest_in_corner
377 = lu.col(k).tail(rows-k).unaryExpr(Scoring()).maxCoeff(&row_of_biggest_in_col);
378 row_of_biggest_in_col += k;
380 row_transpositions[k] = PivIndex(row_of_biggest_in_col);
382 if(biggest_in_corner != Score(0))
384 if(k != row_of_biggest_in_col)
386 lu.row(k).swap(lu.row(row_of_biggest_in_col));
390 lu.col(k).tail(fix<RRows>(rrows)) /= lu.coeff(k,k);
392 else if(first_zero_pivot==-1)
396 first_zero_pivot = k;
400 lu.bottomRightCorner(fix<RRows>(rrows),fix<RCols>(rcols)).noalias() -= lu.col(k).tail(fix<RRows>(rrows)) * lu.row(k).tail(fix<RCols>(rcols));
404 if(UnBlockedAtCompileTime)
407 row_transpositions[k] = PivIndex(k);
408 if (Scoring()(lu(k, k)) == Score(0) && first_zero_pivot == -1)
409 first_zero_pivot = k;
412 return first_zero_pivot;
430 static Index blocked_lu(
Index rows,
Index cols, Scalar* lu_data,
Index luStride, PivIndex* row_transpositions, PivIndex& nb_transpositions,
Index maxBlockSize=256)
434 const Index size = (std::min)(rows,cols);
437 if(UnBlockedAtCompileTime || size<=UnBlockedBound)
439 return unblocked_lu(lu, row_transpositions, nb_transpositions);
447 blockSize = (blockSize/16)*16;
448 blockSize = (std::min)((std::max)(blockSize,
Index(8)), maxBlockSize);
451 nb_transpositions = 0;
452 Index first_zero_pivot = -1;
453 for(
Index k = 0; k < size; k+=blockSize)
455 Index bs = (std::min)(size-k,blockSize);
456 Index trows = rows - k - bs;
457 Index tsize = size - k - bs;
463 BlockType A_0 = lu.block(0,0,rows,k);
464 BlockType A_2 = lu.block(0,k+bs,rows,tsize);
465 BlockType A11 = lu.block(k,k,bs,bs);
466 BlockType A12 = lu.block(k,k+bs,bs,tsize);
467 BlockType A21 = lu.block(k+bs,k,trows,bs);
468 BlockType A22 = lu.block(k+bs,k+bs,trows,tsize);
470 PivIndex nb_transpositions_in_panel;
473 Index ret = blocked_lu(trows+bs, bs, &lu.coeffRef(k,k), luStride,
474 row_transpositions+k, nb_transpositions_in_panel, 16);
475 if(ret>=0 && first_zero_pivot==-1)
476 first_zero_pivot = k+ret;
478 nb_transpositions += nb_transpositions_in_panel;
480 for(
Index i=k; i<k+bs; ++i)
482 Index piv = (row_transpositions[i] += internal::convert_index<PivIndex>(k));
483 A_0.row(i).swap(A_0.row(piv));
489 for(
Index i=k;i<k+bs; ++i)
490 A_2.row(i).swap(A_2.row(row_transpositions[i]));
493 A11.template triangularView<UnitLower>().solveInPlace(A12);
495 A22.noalias() -= A21 * A12;
498 return first_zero_pivot;
504template<
typename MatrixType,
typename TranspositionType>
505void partial_lu_inplace(MatrixType& lu, TranspositionType& row_transpositions,
typename TranspositionType::StorageIndex& nb_transpositions)
508 if (lu.rows() == 0 || lu.cols() == 0) {
509 nb_transpositions = 0;
512 eigen_assert(lu.cols() == row_transpositions.size());
513 eigen_assert(row_transpositions.size() < 2 || (&row_transpositions.coeffRef(1)-&row_transpositions.coeffRef(0)) == 1);
517 typename TranspositionType::StorageIndex,
518 EIGEN_SIZE_MIN_PREFER_FIXED(MatrixType::RowsAtCompileTime,MatrixType::ColsAtCompileTime)>
519 ::blocked_lu(lu.rows(), lu.cols(), &lu.coeffRef(0,0), lu.outerStride(), &row_transpositions.coeffRef(0), nb_transpositions);
524template<
typename MatrixType>
525void PartialPivLU<MatrixType>::compute()
527 check_template_parameters();
530 eigen_assert(m_lu.rows()<NumTraits<int>::highest());
533 m_l1_norm = m_lu.cwiseAbs().colwise().sum().maxCoeff();
535 m_l1_norm = RealScalar(0);
537 eigen_assert(m_lu.rows() == m_lu.cols() &&
"PartialPivLU is only for square (and moreover invertible) matrices");
538 const Index size = m_lu.rows();
540 m_rowsTranspositions.resize(size);
542 typename TranspositionType::StorageIndex nb_transpositions;
543 internal::partial_lu_inplace(m_lu, m_rowsTranspositions, nb_transpositions);
544 m_det_p = (nb_transpositions%2) ? -1 : 1;
546 m_p = m_rowsTranspositions;
548 m_isInitialized =
true;
551template<
typename MatrixType>
554 eigen_assert(m_isInitialized &&
"PartialPivLU is not initialized.");
555 return Scalar(m_det_p) * m_lu.diagonal().prod();
561template<
typename MatrixType>
564 eigen_assert(m_isInitialized &&
"LU is not initialized.");
566 MatrixType res = m_lu.template triangularView<UnitLower>().toDenseMatrix()
567 * m_lu.template triangularView<Upper>();
570 res = m_p.inverse() * res;
580template<
typename DstXprType,
typename MatrixType>
581struct Assignment<DstXprType,
Inverse<
PartialPivLU<MatrixType> >, internal::assign_op<typename DstXprType::Scalar,typename PartialPivLU<MatrixType>::Scalar>, Dense2Dense>
585 static void run(DstXprType &dst,
const SrcXprType &src,
const internal::assign_op<typename DstXprType::Scalar,typename LuType::Scalar> &)
587 dst = src.nestedExpression().solve(MatrixType::Identity(src.rows(), src.cols()));
600template<
typename Derived>
601inline const PartialPivLU<typename MatrixBase<Derived>::PlainObject>
615template<
typename Derived>
Expression of the inverse of another expression.
Definition: Inverse.h:44
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:50
The matrix class, also used for vectors and row-vectors.
Definition: Matrix.h:180
Convenience specialization of Stride to specify only an outer stride See class Map for some examples.
Definition: Stride.h:107
LU decomposition of a matrix with partial pivoting, and related features.
Definition: PartialPivLU.h:78
const Inverse< PartialPivLU > inverse() const
Definition: PartialPivLU.h:196
const Solve< PartialPivLU, Rhs > solve(const MatrixBase< Rhs > &b) const
const PermutationType & permutationP() const
Definition: PartialPivLU.h:151
RealScalar rcond() const
Definition: PartialPivLU.h:183
Scalar determinant() const
Definition: PartialPivLU.h:552
PartialPivLU()
Default Constructor.
Definition: PartialPivLU.h:284
const MatrixType & matrixLU() const
Definition: PartialPivLU.h:143
MatrixType reconstructedMatrix() const
Definition: PartialPivLU.h:562
InverseReturnType transpose() const
Definition: PermutationMatrix.h:191
Permutation matrix.
Definition: PermutationMatrix.h:298
static ConstMapType Map(const Scalar *data)
Definition: PlainObjectBase.h:644
A matrix or vector expression mapping an existing expression.
Definition: Ref.h:283
Pseudo expression representing a solving operation.
Definition: Solve.h:63
A base class for matrix decomposition and solvers.
Definition: SolverBase.h:69
Represents a sequence of transpositions (row/column interchange)
Definition: Transpositions.h:156
@ ColMajor
Definition: Constants.h:319
@ RowMajor
Definition: Constants.h:321
const unsigned int RowMajorBit
Definition: Constants.h:66
Namespace containing all symbols from the Eigen library.
Definition: Core:141
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:74
const int Dynamic
Definition: Constants.h:22
Definition: EigenBase.h:30
Eigen::Index Index
The interface type of indices.
Definition: EigenBase.h:39
Derived & derived()
Definition: EigenBase.h:46
EIGEN_CONSTEXPR Index size() const EIGEN_NOEXCEPT
Definition: EigenBase.h:67