PartialPivLU.h
Go to the documentation of this file.
1 // This file is part of Eigen, a lightweight C++ template library
2 // for linear algebra.
3 //
4 // Copyright (C) 2006-2009 Benoit Jacob <jacob.benoit.1@gmail.com>
5 // Copyright (C) 2009 Gael Guennebaud <gael.guennebaud@inria.fr>
6 //
7 // This Source Code Form is subject to the terms of the Mozilla
8 // Public License v. 2.0. If a copy of the MPL was not distributed
9 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
10 
11 #ifndef EIGEN_PARTIALLU_H
12 #define EIGEN_PARTIALLU_H
13 
14 #include "./InternalHeaderCheck.h"
15 
16 namespace Eigen {
17 
18 namespace internal {
19 template<typename MatrixType_, typename PermutationIndex_> struct traits<PartialPivLU<MatrixType_, PermutationIndex_> >
20  : traits<MatrixType_>
21 {
22  typedef MatrixXpr XprKind;
23  typedef SolverStorage StorageKind;
24  typedef PermutationIndex_ StorageIndex;
25  typedef traits<MatrixType_> BaseTraits;
26  enum {
27  Flags = BaseTraits::Flags & RowMajorBit,
28  CoeffReadCost = Dynamic
29  };
30 };
31 
32 template<typename T,typename Derived>
33 struct enable_if_ref;
34 // {
35 // typedef Derived type;
36 // };
37 
38 template<typename T,typename Derived>
39 struct enable_if_ref<Ref<T>,Derived> {
40  typedef Derived type;
41 };
42 
43 } // end namespace internal
44 
78 template<typename MatrixType_, typename PermutationIndex_> class PartialPivLU
79  : public SolverBase<PartialPivLU<MatrixType_, PermutationIndex_> >
80 {
81  public:
82 
83  typedef MatrixType_ MatrixType;
85  friend class SolverBase<PartialPivLU>;
86 
88  enum {
89  MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
90  MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
91  };
92  using PermutationIndex = PermutationIndex_;
96 
103  PartialPivLU();
104 
111  explicit PartialPivLU(Index size);
112 
120  template<typename InputType>
121  explicit PartialPivLU(const EigenBase<InputType>& matrix);
122 
130  template<typename InputType>
131  explicit PartialPivLU(EigenBase<InputType>& matrix);
132 
133  template<typename InputType>
135  m_lu = matrix.derived();
136  compute();
137  return *this;
138  }
139 
146  inline const MatrixType& matrixLU() const
147  {
148  eigen_assert(m_isInitialized && "PartialPivLU is not initialized.");
149  return m_lu;
150  }
151 
154  inline const PermutationType& permutationP() const
155  {
156  eigen_assert(m_isInitialized && "PartialPivLU is not initialized.");
157  return m_p;
158  }
159 
160  #ifdef EIGEN_PARSED_BY_DOXYGEN
178  template<typename Rhs>
179  inline const Solve<PartialPivLU, Rhs>
180  solve(const MatrixBase<Rhs>& b) const;
181  #endif
182 
186  inline RealScalar rcond() const
187  {
188  eigen_assert(m_isInitialized && "PartialPivLU is not initialized.");
190  }
191 
199  inline const Inverse<PartialPivLU> inverse() const
200  {
201  eigen_assert(m_isInitialized && "PartialPivLU is not initialized.");
202  return Inverse<PartialPivLU>(*this);
203  }
204 
218  Scalar determinant() const;
219 
221 
222  EIGEN_CONSTEXPR inline Index rows() const EIGEN_NOEXCEPT { return m_lu.rows(); }
223  EIGEN_CONSTEXPR inline Index cols() const EIGEN_NOEXCEPT { return m_lu.cols(); }
224 
225  #ifndef EIGEN_PARSED_BY_DOXYGEN
226  template<typename RhsType, typename DstType>
228  void _solve_impl(const RhsType &rhs, DstType &dst) const {
229  /* The decomposition PA = LU can be rewritten as A = P^{-1} L U.
230  * So we proceed as follows:
231  * Step 1: compute c = Pb.
232  * Step 2: replace c by the solution x to Lx = c.
233  * Step 3: replace c by the solution x to Ux = c.
234  */
235 
236  // Step 1
237  dst = permutationP() * rhs;
238 
239  // Step 2
240  m_lu.template triangularView<UnitLower>().solveInPlace(dst);
241 
242  // Step 3
243  m_lu.template triangularView<Upper>().solveInPlace(dst);
244  }
245 
246  template<bool Conjugate, typename RhsType, typename DstType>
248  void _solve_impl_transposed(const RhsType &rhs, DstType &dst) const {
249  /* The decomposition PA = LU can be rewritten as A^T = U^T L^T P.
250  * So we proceed as follows:
251  * Step 1: compute c as the solution to L^T c = b
252  * Step 2: replace c by the solution x to U^T x = c.
253  * Step 3: update c = P^-1 c.
254  */
255 
256  eigen_assert(rhs.rows() == m_lu.cols());
257 
258  // Step 1
259  dst = m_lu.template triangularView<Upper>().transpose()
260  .template conjugateIf<Conjugate>().solve(rhs);
261  // Step 2
262  m_lu.template triangularView<UnitLower>().transpose()
263  .template conjugateIf<Conjugate>().solveInPlace(dst);
264  // Step 3
265  dst = permutationP().transpose() * dst;
266  }
267  #endif
268 
269  protected:
270 
272 
273  void compute();
274 
278  RealScalar m_l1_norm;
279  signed char m_det_p;
281 };
282 
283 template<typename MatrixType, typename PermutationIndex>
285  : m_lu(),
286  m_p(),
287  m_rowsTranspositions(),
288  m_l1_norm(0),
289  m_det_p(0),
290  m_isInitialized(false)
291 {
292 }
293 
294 template<typename MatrixType, typename PermutationIndex>
296  : m_lu(size, size),
297  m_p(size),
298  m_rowsTranspositions(size),
299  m_l1_norm(0),
300  m_det_p(0),
301  m_isInitialized(false)
302 {
303 }
304 
305 template<typename MatrixType, typename PermutationIndex>
306 template<typename InputType>
308  : m_lu(matrix.rows(),matrix.cols()),
309  m_p(matrix.rows()),
310  m_rowsTranspositions(matrix.rows()),
311  m_l1_norm(0),
312  m_det_p(0),
313  m_isInitialized(false)
314 {
315  compute(matrix.derived());
316 }
317 
318 template<typename MatrixType, typename PermutationIndex>
319 template<typename InputType>
321  : m_lu(matrix.derived()),
322  m_p(matrix.rows()),
323  m_rowsTranspositions(matrix.rows()),
324  m_l1_norm(0),
325  m_det_p(0),
326  m_isInitialized(false)
327 {
328  compute();
329 }
330 
331 namespace internal {
332 
334 template<typename Scalar, int StorageOrder, typename PivIndex, int SizeAtCompileTime=Dynamic>
335 struct partial_lu_impl
336 {
337  static constexpr int UnBlockedBound = 16;
338  static constexpr bool UnBlockedAtCompileTime = SizeAtCompileTime!=Dynamic && SizeAtCompileTime<=UnBlockedBound;
339  static constexpr int ActualSizeAtCompileTime = UnBlockedAtCompileTime ? SizeAtCompileTime : Dynamic;
340  // Remaining rows and columns at compile-time:
341  static constexpr int RRows = SizeAtCompileTime==2 ? 1 : Dynamic;
342  static constexpr int RCols = SizeAtCompileTime==2 ? 1 : Dynamic;
344  typedef Ref<MatrixType> MatrixTypeRef;
346  typedef typename MatrixType::RealScalar RealScalar;
347 
358  static Index unblocked_lu(MatrixTypeRef& lu, PivIndex* row_transpositions, PivIndex& nb_transpositions)
359  {
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);
365  // For small compile-time matrices it is worth processing the last row separately:
366  // speedup: +100% for 2x2, +10% for others.
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)
371  {
372  int rrows = internal::convert_index<int>(rows-k-1);
373  int rcols = internal::convert_index<int>(cols-k-1);
374 
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;
379 
380  row_transpositions[k] = PivIndex(row_of_biggest_in_col);
381 
382  if(!numext::is_exactly_zero(biggest_in_corner))
383  {
384  if(k != row_of_biggest_in_col)
385  {
386  lu.row(k).swap(lu.row(row_of_biggest_in_col));
387  ++nb_transpositions;
388  }
389 
390  lu.col(k).tail(fix<RRows>(rrows)) /= lu.coeff(k,k);
391  }
392  else if(first_zero_pivot==-1)
393  {
394  // the pivot is exactly zero, we record the index of the first pivot which is exactly 0,
395  // and continue the factorization such we still have A = PLU
396  first_zero_pivot = k;
397  }
398 
399  if(k<rows-1)
400  lu.bottomRightCorner(fix<RRows>(rrows),fix<RCols>(rcols)).noalias() -= lu.col(k).tail(fix<RRows>(rrows)) * lu.row(k).tail(fix<RCols>(rcols));
401  }
402 
403  // special handling of the last entry
404  if(UnBlockedAtCompileTime)
405  {
406  Index k = endk;
407  row_transpositions[k] = PivIndex(k);
408  if (numext::is_exactly_zero(Scoring()(lu(k, k))) && first_zero_pivot == -1)
409  first_zero_pivot = k;
410  }
411 
412  return first_zero_pivot;
413  }
414 
430  static Index blocked_lu(Index rows, Index cols, Scalar* lu_data, Index luStride, PivIndex* row_transpositions, PivIndex& nb_transpositions, Index maxBlockSize=256)
431  {
432  MatrixTypeRef lu = MatrixType::Map(lu_data,rows, cols, OuterStride<>(luStride));
433 
434  const Index size = (std::min)(rows,cols);
435 
436  // if the matrix is too small, no blocking:
437  if(UnBlockedAtCompileTime || size<=UnBlockedBound)
438  {
439  return unblocked_lu(lu, row_transpositions, nb_transpositions);
440  }
441 
442  // automatically adjust the number of subdivisions to the size
443  // of the matrix so that there is enough sub blocks:
444  Index blockSize;
445  {
446  blockSize = size/8;
447  blockSize = (blockSize/16)*16;
448  blockSize = (std::min)((std::max)(blockSize,Index(8)), maxBlockSize);
449  }
450 
451  nb_transpositions = 0;
452  Index first_zero_pivot = -1;
453  for(Index k = 0; k < size; k+=blockSize)
454  {
455  Index bs = (std::min)(size-k,blockSize); // actual size of the block
456  Index trows = rows - k - bs; // trailing rows
457  Index tsize = size - k - bs; // trailing size
458 
459  // partition the matrix:
460  // A00 | A01 | A02
461  // lu = A_0 | A_1 | A_2 = A10 | A11 | A12
462  // A20 | A21 | A22
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);
469 
470  PivIndex nb_transpositions_in_panel;
471  // recursively call the blocked LU algorithm on [A11^T A21^T]^T
472  // with a very small blocking size:
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;
477 
478  nb_transpositions += nb_transpositions_in_panel;
479  // update permutations and apply them to A_0
480  for(Index i=k; i<k+bs; ++i)
481  {
482  Index piv = (row_transpositions[i] += internal::convert_index<PivIndex>(k));
483  A_0.row(i).swap(A_0.row(piv));
484  }
485 
486  if(trows)
487  {
488  // apply permutations to A_2
489  for(Index i=k;i<k+bs; ++i)
490  A_2.row(i).swap(A_2.row(row_transpositions[i]));
491 
492  // A12 = A11^-1 A12
493  A11.template triangularView<UnitLower>().solveInPlace(A12);
494 
495  A22.noalias() -= A21 * A12;
496  }
497  }
498  return first_zero_pivot;
499  }
500 };
501 
504 template<typename MatrixType, typename TranspositionType>
505 void partial_lu_inplace(MatrixType& lu, TranspositionType& row_transpositions, typename TranspositionType::StorageIndex& nb_transpositions)
506 {
507  // Special-case of zero matrix.
508  if (lu.rows() == 0 || lu.cols() == 0) {
509  nb_transpositions = 0;
510  return;
511  }
512  eigen_assert(lu.cols() == row_transpositions.size());
513  eigen_assert(row_transpositions.size() < 2 || (&row_transpositions.coeffRef(1)-&row_transpositions.coeffRef(0)) == 1);
514 
515  partial_lu_impl
516  < typename MatrixType::Scalar, MatrixType::Flags&RowMajorBit?RowMajor:ColMajor,
517  typename TranspositionType::StorageIndex,
518  internal::min_size_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);
520 }
521 
522 } // end namespace internal
523 
524 template<typename MatrixType, typename PermutationIndex>
526 {
528 
529  if(m_lu.cols()>0)
530  m_l1_norm = m_lu.cwiseAbs().colwise().sum().maxCoeff();
531  else
532  m_l1_norm = RealScalar(0);
533 
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();
536 
537  m_rowsTranspositions.resize(size);
538 
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;
542 
543  m_p = m_rowsTranspositions;
544 
545  m_isInitialized = true;
546 }
547 
548 template<typename MatrixType, typename PermutationIndex>
550 {
551  eigen_assert(m_isInitialized && "PartialPivLU is not initialized.");
552  return Scalar(m_det_p) * m_lu.diagonal().prod();
553 }
554 
558 template<typename MatrixType, typename PermutationIndex>
560 {
561  eigen_assert(m_isInitialized && "LU is not initialized.");
562  // LU
563  MatrixType res = m_lu.template triangularView<UnitLower>().toDenseMatrix()
564  * m_lu.template triangularView<Upper>();
565 
566  // P^{-1}(LU)
567  res = m_p.inverse() * res;
568 
569  return res;
570 }
571 
572 
574 namespace internal {
575 
576 
577 template<typename DstXprType, typename MatrixType, typename PermutationIndex>
578 struct Assignment<DstXprType, Inverse<PartialPivLU<MatrixType, PermutationIndex> >, internal::assign_op<typename DstXprType::Scalar,typename PartialPivLU<MatrixType, PermutationIndex>::Scalar>, Dense2Dense>
579 {
581  typedef Inverse<LuType> SrcXprType;
582  static void run(DstXprType &dst, const SrcXprType &src, const internal::assign_op<typename DstXprType::Scalar,typename LuType::Scalar> &)
583  {
584  dst = src.nestedExpression().solve(MatrixType::Identity(src.rows(), src.cols()));
585  }
586 };
587 } // end namespace internal
588 
589 
597 template<typename Derived>
598 template<typename PermutationIndex>
599 inline const PartialPivLU<typename MatrixBase<Derived>::PlainObject, PermutationIndex>
601 {
603 }
604 
613 template<typename Derived>
614 template<typename PermutationIndex>
615 inline const PartialPivLU<typename MatrixBase<Derived>::PlainObject, PermutationIndex>
617 {
619 }
620 
621 } // end namespace Eigen
622 
623 #endif // EIGEN_PARTIALLU_H
Array< int, 3, 1 > b
#define EIGEN_GENERIC_PUBLIC_INTERFACE(Derived)
Definition: Macros.h:1149
#define EIGEN_NOEXCEPT
Definition: Macros.h:1260
#define EIGEN_CONSTEXPR
Definition: Macros.h:747
#define EIGEN_DEVICE_FUNC
Definition: Macros.h:883
#define eigen_assert(x)
Definition: Macros.h:902
cout<< "Here is the matrix m:"<< endl<< m<< endl;Matrix< ptrdiff_t, 3, 1 > res
#define EIGEN_STATIC_ASSERT_NON_INTEGER(TYPE)
Definition: StaticAssert.h:81
Matrix< float, 1, Dynamic > MatrixType
NumTraits< Scalar >::Real RealScalar
Definition: DenseBase.h:68
Expression of the inverse of another expression.
Definition: Inverse.h:46
const XprTypeNestedCleaned & nestedExpression() const
Definition: Inverse.h:62
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:52
const PartialPivLU< PlainObject, PermutationIndex > lu() const
const PartialPivLU< PlainObject, PermutationIndex > partialPivLu() const
The matrix class, also used for vectors and row-vectors.
Definition: Matrix.h:182
Base::PlainObject PlainObject
Definition: Matrix.h:194
LU decomposition of a matrix with partial pivoting, and related features.
Definition: PartialPivLU.h:80
SolverBase< PartialPivLU > Base
Definition: PartialPivLU.h:84
MatrixType_ MatrixType
Definition: PartialPivLU.h:83
RealScalar m_l1_norm
Definition: PartialPivLU.h:278
TranspositionType m_rowsTranspositions
Definition: PartialPivLU.h:277
PermutationType m_p
Definition: PartialPivLU.h:276
PartialPivLU & compute(const EigenBase< InputType > &matrix)
Definition: PartialPivLU.h:134
MatrixType::PlainObject PlainObject
Definition: PartialPivLU.h:95
const Solve< PartialPivLU, Rhs > solve(const MatrixBase< Rhs > &b) const
MatrixType reconstructedMatrix() const
Definition: PartialPivLU.h:559
const PermutationType & permutationP() const
Definition: PartialPivLU.h:154
signed char m_det_p
Definition: PartialPivLU.h:279
PermutationIndex_ PermutationIndex
Definition: PartialPivLU.h:92
RealScalar rcond() const
Definition: PartialPivLU.h:186
EIGEN_CONSTEXPR Index rows() const EIGEN_NOEXCEPT
Definition: PartialPivLU.h:222
EIGEN_CONSTEXPR Index cols() const EIGEN_NOEXCEPT
Definition: PartialPivLU.h:223
PartialPivLU()
Default Constructor.
Definition: PartialPivLU.h:284
Transpositions< RowsAtCompileTime, MaxRowsAtCompileTime, PermutationIndex > TranspositionType
Definition: PartialPivLU.h:94
const Inverse< PartialPivLU > inverse() const
Definition: PartialPivLU.h:199
const MatrixType & matrixLU() const
Definition: PartialPivLU.h:146
Scalar determinant() const
Definition: PartialPivLU.h:549
PermutationMatrix< RowsAtCompileTime, MaxRowsAtCompileTime, PermutationIndex > PermutationType
Definition: PartialPivLU.h:93
InverseReturnType transpose() const
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
internal::traits< PartialPivLU< MatrixType_, PermutationIndex_ > >::Scalar Scalar
Definition: SolverBase.h:75
cout<< "Here is the matrix m:"<< endl<< m<< endl;Eigen::FullPivLU< Matrix5x3 > lu(m)
@ ColMajor
Definition: Constants.h:321
@ RowMajor
Definition: Constants.h:323
const unsigned int RowMajorBit
Definition: Constants.h:68
bfloat16() max(const bfloat16 &a, const bfloat16 &b)
Definition: BFloat16.h:690
bfloat16() min(const bfloat16 &a, const bfloat16 &b)
Definition: BFloat16.h:684
void partial_lu_inplace(MatrixType &lu, TranspositionType &row_transpositions, typename TranspositionType::StorageIndex &nb_transpositions)
Definition: PartialPivLU.h:505
constexpr int min_size_prefer_fixed(A a, B b)
Definition: Meta.h:553
Decomposition::RealScalar rcond_estimate_helper(typename Decomposition::RealScalar matrix_norm, const Decomposition &dec)
Reciprocal condition number estimator.
bool is_exactly_zero(const X &x)
Definition: Meta.h:475
: InteropHeaders
Definition: Core:139
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:82
const int Dynamic
Definition: Constants.h:24
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
Holds information about the various numeric (i.e. scalar) types allowed by Eigen.
Definition: NumTraits.h:231