11 #ifndef EIGEN_PARTIALLU_H
12 #define EIGEN_PARTIALLU_H
17 template<
typename _MatrixType>
struct traits<PartialPivLU<_MatrixType> >
20 typedef MatrixXpr XprKind;
21 typedef SolverStorage StorageKind;
22 typedef traits<_MatrixType> BaseTraits;
25 CoeffReadCost = Dynamic
62 template<
typename _MatrixType>
class PartialPivLU
63 :
public SolverBase<PartialPivLU<_MatrixType> >
67 typedef _MatrixType MatrixType;
68 typedef SolverBase<PartialPivLU> Base;
72 MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
73 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
75 typedef PermutationMatrix<RowsAtCompileTime, MaxRowsAtCompileTime> PermutationType;
76 typedef Transpositions<RowsAtCompileTime, MaxRowsAtCompileTime> TranspositionType;
77 typedef typename MatrixType::PlainObject PlainObject;
103 template<
typename InputType>
104 explicit PartialPivLU(
const EigenBase<InputType>& matrix);
106 template<
typename InputType>
107 PartialPivLU& compute(
const EigenBase<InputType>& matrix);
117 eigen_assert(m_isInitialized &&
"PartialPivLU is not initialized.");
125 eigen_assert(m_isInitialized &&
"PartialPivLU is not initialized.");
147 template<
typename Rhs>
151 eigen_assert(m_isInitialized &&
"PartialPivLU is not initialized.");
164 eigen_assert(m_isInitialized &&
"PartialPivLU is not initialized.");
181 typename internal::traits<MatrixType>::Scalar
determinant()
const;
185 inline Index rows()
const {
return m_lu.rows(); }
186 inline Index cols()
const {
return m_lu.cols(); }
188 #ifndef EIGEN_PARSED_BY_DOXYGEN
189 template<
typename RhsType,
typename DstType>
191 void _solve_impl(
const RhsType &rhs, DstType &dst)
const {
199 eigen_assert(rhs.rows() == m_lu.rows());
205 m_lu.template triangularView<UnitLower>().solveInPlace(dst);
208 m_lu.template triangularView<Upper>().solveInPlace(dst);
211 template<
bool Conjugate,
typename RhsType,
typename DstType>
213 void _solve_impl_transposed(
const RhsType &rhs, DstType &dst)
const {
221 eigen_assert(rhs.rows() == m_lu.cols());
225 dst = m_lu.template triangularView<Upper>().
adjoint().solve(rhs);
227 m_lu.template triangularView<UnitLower>().
adjoint().solveInPlace(dst);
230 dst = m_lu.template triangularView<Upper>().
transpose().solve(rhs);
232 m_lu.template triangularView<UnitLower>().
transpose().solveInPlace(dst);
241 static void check_template_parameters()
243 EIGEN_STATIC_ASSERT_NON_INTEGER(Scalar);
248 TranspositionType m_rowsTranspositions;
250 bool m_isInitialized;
253 template<
typename MatrixType>
257 m_rowsTranspositions(),
259 m_isInitialized(false)
263 template<
typename MatrixType>
267 m_rowsTranspositions(size),
269 m_isInitialized(false)
273 template<
typename MatrixType>
274 template<
typename InputType>
276 : m_lu(matrix.rows(), matrix.rows()),
278 m_rowsTranspositions(matrix.rows()),
280 m_isInitialized(false)
288 template<
typename Scalar,
int StorageOrder,
typename PivIndex>
289 struct partial_lu_impl
299 typedef typename MatrixType::RealScalar RealScalar;
311 static Index unblocked_lu(MatrixType& lu, PivIndex* row_transpositions, PivIndex& nb_transpositions)
313 typedef scalar_score_coeff_op<Scalar> Scoring;
314 typedef typename Scoring::result_type Score;
315 const Index rows = lu.rows();
316 const Index cols = lu.cols();
317 const Index size = (std::min)(rows,cols);
318 nb_transpositions = 0;
319 Index first_zero_pivot = -1;
320 for(Index k = 0; k < size; ++k)
322 Index rrows = rows-k-1;
323 Index rcols = cols-k-1;
325 Index row_of_biggest_in_col;
326 Score biggest_in_corner
327 = lu.col(k).tail(rows-k).unaryExpr(Scoring()).maxCoeff(&row_of_biggest_in_col);
328 row_of_biggest_in_col += k;
330 row_transpositions[k] = PivIndex(row_of_biggest_in_col);
332 if(biggest_in_corner != Score(0))
334 if(k != row_of_biggest_in_col)
336 lu.row(k).swap(lu.row(row_of_biggest_in_col));
342 lu.col(k).tail(rrows) /= lu.coeff(k,k);
344 else if(first_zero_pivot==-1)
348 first_zero_pivot = k;
352 lu.bottomRightCorner(rrows,rcols).noalias() -= lu.col(k).tail(rrows) * lu.row(k).tail(rcols);
354 return first_zero_pivot;
372 static Index blocked_lu(Index rows, Index cols, Scalar* lu_data, Index luStride, PivIndex* row_transpositions, PivIndex& nb_transpositions, Index maxBlockSize=256)
374 MapLU lu1(lu_data,StorageOrder==
RowMajor?rows:luStride,StorageOrder==
RowMajor?luStride:cols);
375 MatrixType lu(lu1,0,0,rows,cols);
377 const Index size = (std::min)(rows,cols);
382 return unblocked_lu(lu, row_transpositions, nb_transpositions);
390 blockSize = (blockSize/16)*16;
391 blockSize = (std::min)((std::max)(blockSize,Index(8)), maxBlockSize);
394 nb_transpositions = 0;
395 Index first_zero_pivot = -1;
396 for(Index k = 0; k < size; k+=blockSize)
398 Index bs = (std::min)(size-k,blockSize);
399 Index trows = rows - k - bs;
400 Index tsize = size - k - bs;
406 BlockType A_0(lu,0,0,rows,k);
407 BlockType A_2(lu,0,k+bs,rows,tsize);
408 BlockType A11(lu,k,k,bs,bs);
409 BlockType A12(lu,k,k+bs,bs,tsize);
410 BlockType A21(lu,k+bs,k,trows,bs);
411 BlockType A22(lu,k+bs,k+bs,trows,tsize);
413 PivIndex nb_transpositions_in_panel;
416 Index ret = blocked_lu(trows+bs, bs, &lu.coeffRef(k,k), luStride,
417 row_transpositions+k, nb_transpositions_in_panel, 16);
418 if(ret>=0 && first_zero_pivot==-1)
419 first_zero_pivot = k+ret;
421 nb_transpositions += nb_transpositions_in_panel;
423 for(Index i=k; i<k+bs; ++i)
425 Index piv = (row_transpositions[i] += k);
426 A_0.row(i).swap(A_0.row(piv));
432 for(Index i=k;i<k+bs; ++i)
433 A_2.row(i).swap(A_2.row(row_transpositions[i]));
436 A11.template triangularView<UnitLower>().solveInPlace(A12);
438 A22.noalias() -= A21 * A12;
441 return first_zero_pivot;
447 template<
typename MatrixType,
typename TranspositionType>
448 void partial_lu_inplace(MatrixType& lu, TranspositionType& row_transpositions,
typename TranspositionType::StorageIndex& nb_transpositions)
450 eigen_assert(lu.cols() == row_transpositions.size());
451 eigen_assert((&row_transpositions.coeffRef(1)-&row_transpositions.coeffRef(0)) == 1);
455 ::blocked_lu(lu.rows(), lu.cols(), &lu.coeffRef(0,0), lu.outerStride(), &row_transpositions.coeffRef(0), nb_transpositions);
460 template<
typename MatrixType>
461 template<
typename InputType>
462 PartialPivLU<MatrixType>& PartialPivLU<MatrixType>::compute(
const EigenBase<InputType>& matrix)
464 check_template_parameters();
467 eigen_assert(matrix.rows()<NumTraits<int>::highest());
469 m_lu = matrix.derived();
471 eigen_assert(matrix.rows() == matrix.cols() &&
"PartialPivLU is only for square (and moreover invertible) matrices");
472 const Index size = matrix.rows();
474 m_rowsTranspositions.resize(size);
476 typename TranspositionType::StorageIndex nb_transpositions;
477 internal::partial_lu_inplace(m_lu, m_rowsTranspositions, nb_transpositions);
478 m_det_p = (nb_transpositions%2) ? -1 : 1;
480 m_p = m_rowsTranspositions;
482 m_isInitialized =
true;
486 template<
typename MatrixType>
489 eigen_assert(m_isInitialized &&
"PartialPivLU is not initialized.");
490 return Scalar(m_det_p) * m_lu.diagonal().prod();
496 template<
typename MatrixType>
499 eigen_assert(m_isInitialized &&
"LU is not initialized.");
501 MatrixType res = m_lu.template triangularView<UnitLower>().toDenseMatrix()
502 * m_lu.template triangularView<Upper>();
505 res = m_p.inverse() * res;
515 template<
typename DstXprType,
typename MatrixType,
typename Scalar>
520 static void run(DstXprType &dst,
const SrcXprType &src,
const internal::assign_op<Scalar> &)
522 dst = src.nestedExpression().solve(MatrixType::Identity(src.rows(), src.cols()));
536 template<
typename Derived>
537 inline const PartialPivLU<typename MatrixBase<Derived>::PlainObject>
553 template<
typename Derived>
563 #endif // EIGEN_PARTIALLU_H
const Inverse< PartialPivLU > inverse() const
Definition: PartialPivLU.h:162
PartialPivLU()
Default Constructor.
Definition: PartialPivLU.h:254
A matrix or vector expression mapping an existing array of data.
Definition: Map.h:89
AdjointReturnType adjoint() const
Definition: SolverBase.h:109
LU decomposition of a matrix with partial pivoting, and related features.
Definition: ForwardDeclarations.h:248
Derived & derived()
Definition: EigenBase.h:44
Eigen::Index Index
The interface type of indices.
Definition: EigenBase.h:37
const unsigned int RowMajorBit
Definition: Constants.h:61
const PartialPivLU< PlainObject > partialPivLu() const
Definition: PartialPivLU.h:538
Definition: EigenBase.h:28
Expression of the inverse of another expression.
Definition: Inverse.h:43
Definition: Constants.h:320
Definition: Constants.h:322
const PermutationType & permutationP() const
Definition: PartialPivLU.h:123
InverseReturnType transpose() const
Definition: PermutationMatrix.h:203
ConstTransposeReturnType transpose() const
Definition: SolverBase.h:90
Definition: Eigen_Colamd.h:54
Expression of a fixed-size or dynamic-size block.
Definition: Block.h:104
internal::traits< MatrixType >::Scalar determinant() const
Definition: PartialPivLU.h:487
MatrixType reconstructedMatrix() const
Definition: PartialPivLU.h:497
Pseudo expression representing a solving operation.
Definition: Solve.h:62
const MatrixType & matrixLU() const
Definition: PartialPivLU.h:115
const Solve< PartialPivLU, Rhs > solve(const MatrixBase< Rhs > &b) const
Definition: PartialPivLU.h:149
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:48
const PartialPivLU< PlainObject > lu() const
Definition: PartialPivLU.h:555
Index size() const
Definition: EigenBase.h:65