$darkmode
Eigen  5.0.1-dev
PartialPivLU.h
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 // IWYU pragma: private
15 #include "./InternalHeaderCheck.h"
16 
17 namespace Eigen {
18 
19 namespace internal {
20 template <typename MatrixType_, typename PermutationIndex_>
21 struct traits<PartialPivLU<MatrixType_, PermutationIndex_> > : traits<MatrixType_> {
22  typedef MatrixXpr XprKind;
23  typedef SolverStorage StorageKind;
24  typedef PermutationIndex_ StorageIndex;
25  typedef traits<MatrixType_> BaseTraits;
26  enum { Flags = BaseTraits::Flags & RowMajorBit, CoeffReadCost = Dynamic };
27 };
28 
29 template <typename T, typename Derived>
30 struct enable_if_ref;
31 // {
32 // typedef Derived type;
33 // };
34 
35 template <typename T, typename Derived>
36 struct enable_if_ref<Ref<T>, Derived> {
37  typedef Derived type;
38 };
39 
40 } // end namespace internal
41 
76 template <typename MatrixType_, typename PermutationIndex_>
77 class PartialPivLU : public SolverBase<PartialPivLU<MatrixType_, PermutationIndex_> > {
78  public:
79  typedef MatrixType_ MatrixType;
80  typedef SolverBase<PartialPivLU> Base;
81  friend class SolverBase<PartialPivLU>;
82 
83  EIGEN_GENERIC_PUBLIC_INTERFACE(PartialPivLU)
84  enum {
85  MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
86  MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
87  };
88  using PermutationIndex = PermutationIndex_;
89  typedef PermutationMatrix<RowsAtCompileTime, MaxRowsAtCompileTime, PermutationIndex> PermutationType;
90  typedef Transpositions<RowsAtCompileTime, MaxRowsAtCompileTime, PermutationIndex> TranspositionType;
91  typedef typename MatrixType::PlainObject PlainObject;
92 
100  eigen_assert(m_isInitialized && "PartialPivLU is not initialized.");
101  return Success;
102  }
103 
110  PartialPivLU();
111 
118  explicit PartialPivLU(Index size);
119 
127  template <typename InputType>
128  explicit PartialPivLU(const EigenBase<InputType>& matrix);
129 
137  template <typename InputType>
138  explicit PartialPivLU(EigenBase<InputType>& matrix);
139 
140  template <typename InputType>
141  PartialPivLU& compute(const EigenBase<InputType>& matrix) {
142  m_lu = matrix.derived();
143  compute();
144  return *this;
145  }
146 
153  inline const MatrixType& matrixLU() const {
154  eigen_assert(m_isInitialized && "PartialPivLU is not initialized.");
155  return m_lu;
156  }
157 
160  inline const PermutationType& permutationP() const {
161  eigen_assert(m_isInitialized && "PartialPivLU is not initialized.");
162  return m_p;
163  }
164 
165 #ifdef EIGEN_PARSED_BY_DOXYGEN
166 
183  template <typename Rhs>
184  inline const Solve<PartialPivLU, Rhs> solve(const MatrixBase<Rhs>& b) const;
185 #endif
186 
190  inline RealScalar rcond() const {
191  eigen_assert(m_isInitialized && "PartialPivLU is not initialized.");
192  return internal::rcond_estimate_helper(m_l1_norm, *this);
193  }
194 
202  inline const Inverse<PartialPivLU> inverse() const {
203  eigen_assert(m_isInitialized && "PartialPivLU is not initialized.");
204  return Inverse<PartialPivLU>(*this);
205  }
206 
220  Scalar determinant() const;
221 
222  MatrixType reconstructedMatrix() const;
223 
224  constexpr Index rows() const noexcept { return m_lu.rows(); }
225  constexpr Index cols() const noexcept { return m_lu.cols(); }
226 
227 #ifndef EIGEN_PARSED_BY_DOXYGEN
228  template <typename RhsType, typename DstType>
229  EIGEN_DEVICE_FUNC void _solve_impl(const RhsType& rhs, DstType& dst) const {
230  /* The decomposition PA = LU can be rewritten as A = P^{-1} L U.
231  * So we proceed as follows:
232  * Step 1: compute c = Pb.
233  * Step 2: replace c by the solution x to Lx = c.
234  * Step 3: replace c by the solution x to Ux = c.
235  */
236 
237  // Step 1
238  dst = permutationP() * rhs;
239 
240  // Step 2
241  m_lu.template triangularView<UnitLower>().solveInPlace(dst);
242 
243  // Step 3
244  m_lu.template triangularView<Upper>().solveInPlace(dst);
245  }
246 
247  template <bool Conjugate, typename RhsType, typename DstType>
248  EIGEN_DEVICE_FUNC 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().template conjugateIf<Conjugate>().solve(rhs);
260  // Step 2
261  m_lu.template triangularView<UnitLower>().transpose().template conjugateIf<Conjugate>().solveInPlace(dst);
262  // Step 3
263  dst = permutationP().transpose() * dst;
264  }
265 #endif
266 
267  protected:
268  EIGEN_STATIC_ASSERT_NON_INTEGER(Scalar)
269 
270  void compute();
271 
272  MatrixType m_lu;
273  PermutationType m_p;
274  TranspositionType m_rowsTranspositions;
275  RealScalar m_l1_norm;
276  signed char m_det_p;
277  bool m_isInitialized;
278 };
279 
280 template <typename MatrixType, typename PermutationIndex>
282  : m_lu(), m_p(), m_rowsTranspositions(), m_l1_norm(0), m_det_p(0), m_isInitialized(false) {}
283 
284 template <typename MatrixType, typename PermutationIndex>
286  : m_lu(size, size), m_p(size), m_rowsTranspositions(size), m_l1_norm(0), m_det_p(0), m_isInitialized(false) {}
287 
288 template <typename MatrixType, typename PermutationIndex>
289 template <typename InputType>
291  : m_lu(matrix.rows(), matrix.cols()),
292  m_p(matrix.rows()),
293  m_rowsTranspositions(matrix.rows()),
294  m_l1_norm(0),
295  m_det_p(0),
296  m_isInitialized(false) {
297  compute(matrix.derived());
298 }
299 
300 template <typename MatrixType, typename PermutationIndex>
301 template <typename InputType>
303  : m_lu(matrix.derived()),
304  m_p(matrix.rows()),
305  m_rowsTranspositions(matrix.rows()),
306  m_l1_norm(0),
307  m_det_p(0),
308  m_isInitialized(false) {
309  compute();
310 }
311 
312 namespace internal {
313 
315 template <typename Scalar, int StorageOrder, typename PivIndex, int SizeAtCompileTime = Dynamic>
316 struct partial_lu_impl {
317  static constexpr int UnBlockedBound = 16;
318  static constexpr bool UnBlockedAtCompileTime = SizeAtCompileTime != Dynamic && SizeAtCompileTime <= UnBlockedBound;
319  static constexpr int ActualSizeAtCompileTime = UnBlockedAtCompileTime ? SizeAtCompileTime : Dynamic;
320  // Remaining rows and columns at compile-time:
321  static constexpr int RRows = SizeAtCompileTime == 2 ? 1 : Dynamic;
322  static constexpr int RCols = SizeAtCompileTime == 2 ? 1 : Dynamic;
324  typedef Ref<MatrixType> MatrixTypeRef;
326  typedef typename MatrixType::RealScalar RealScalar;
327 
338  static Index unblocked_lu(MatrixTypeRef& lu, PivIndex* row_transpositions, PivIndex& nb_transpositions) {
339  typedef scalar_score_coeff_op<Scalar> Scoring;
340  typedef typename Scoring::result_type Score;
341  const Index rows = lu.rows();
342  const Index cols = lu.cols();
343  const Index size = (std::min)(rows, cols);
344  // For small compile-time matrices it is worth processing the last row separately:
345  // speedup: +100% for 2x2, +10% for others.
346  const Index endk = UnBlockedAtCompileTime ? size - 1 : size;
347  nb_transpositions = 0;
348  Index first_zero_pivot = -1;
349  for (Index k = 0; k < endk; ++k) {
350  int rrows = internal::convert_index<int>(rows - k - 1);
351  int rcols = internal::convert_index<int>(cols - k - 1);
352 
353  Index row_of_biggest_in_col;
354  Score biggest_in_corner = lu.col(k).tail(rows - k).unaryExpr(Scoring()).maxCoeff(&row_of_biggest_in_col);
355  row_of_biggest_in_col += k;
356 
357  row_transpositions[k] = PivIndex(row_of_biggest_in_col);
358 
359  if (!numext::is_exactly_zero(biggest_in_corner)) {
360  if (k != row_of_biggest_in_col) {
361  lu.row(k).swap(lu.row(row_of_biggest_in_col));
362  ++nb_transpositions;
363  }
364 
365  lu.col(k).tail(fix<RRows>(rrows)) /= lu.coeff(k, k);
366  } else if (first_zero_pivot == -1) {
367  // the pivot is exactly zero, we record the index of the first pivot which is exactly 0,
368  // and continue the factorization such we still have A = PLU
369  first_zero_pivot = k;
370  }
371 
372  if (k < rows - 1)
373  lu.bottomRightCorner(fix<RRows>(rrows), fix<RCols>(rcols)).noalias() -=
374  lu.col(k).tail(fix<RRows>(rrows)) * lu.row(k).tail(fix<RCols>(rcols));
375  }
376 
377  // special handling of the last entry
378  if (UnBlockedAtCompileTime) {
379  Index k = endk;
380  row_transpositions[k] = PivIndex(k);
381  if (numext::is_exactly_zero(Scoring()(lu(k, k))) && first_zero_pivot == -1) first_zero_pivot = k;
382  }
383 
384  return first_zero_pivot;
385  }
386 
402  static Index blocked_lu(Index rows, Index cols, Scalar* lu_data, Index luStride, PivIndex* row_transpositions,
403  PivIndex& nb_transpositions, Index maxBlockSize = 256) {
404  MatrixTypeRef lu = MatrixType::Map(lu_data, rows, cols, OuterStride<>(luStride));
405 
406  const Index size = (std::min)(rows, cols);
407 
408  // if the matrix is too small, no blocking:
409  if (UnBlockedAtCompileTime || size <= UnBlockedBound) {
410  return unblocked_lu(lu, row_transpositions, nb_transpositions);
411  }
412 
413  // automatically adjust the number of subdivisions to the size
414  // of the matrix so that there is enough sub blocks:
415  Index blockSize;
416  {
417  blockSize = size / 8;
418  blockSize = (blockSize / 16) * 16;
419  blockSize = (std::min)((std::max)(blockSize, Index(8)), maxBlockSize);
420  }
421 
422  nb_transpositions = 0;
423  Index first_zero_pivot = -1;
424  for (Index k = 0; k < size; k += blockSize) {
425  Index bs = (std::min)(size - k, blockSize); // actual size of the block
426  Index trows = rows - k - bs; // trailing rows
427  Index tsize = size - k - bs; // trailing size
428 
429  // partition the matrix:
430  // A00 | A01 | A02
431  // lu = A_0 | A_1 | A_2 = A10 | A11 | A12
432  // A20 | A21 | A22
433  BlockType A_0 = lu.block(0, 0, rows, k);
434  BlockType A_2 = lu.block(0, k + bs, rows, tsize);
435  BlockType A11 = lu.block(k, k, bs, bs);
436  BlockType A12 = lu.block(k, k + bs, bs, tsize);
437  BlockType A21 = lu.block(k + bs, k, trows, bs);
438  BlockType A22 = lu.block(k + bs, k + bs, trows, tsize);
439 
440  PivIndex nb_transpositions_in_panel;
441  // recursively call the blocked LU algorithm on [A11^T A21^T]^T
442  // with a very small blocking size:
443  Index ret = blocked_lu(trows + bs, bs, &lu.coeffRef(k, k), luStride, row_transpositions + k,
444  nb_transpositions_in_panel, 16);
445  if (ret >= 0 && first_zero_pivot == -1) first_zero_pivot = k + ret;
446 
447  nb_transpositions += nb_transpositions_in_panel;
448  // update permutations and apply them to A_0
449  for (Index i = k; i < k + bs; ++i) {
450  Index piv = (row_transpositions[i] += internal::convert_index<PivIndex>(k));
451  A_0.row(i).swap(A_0.row(piv));
452  }
453 
454  if (trows) {
455  // apply permutations to A_2
456  for (Index i = k; i < k + bs; ++i) A_2.row(i).swap(A_2.row(row_transpositions[i]));
457 
458  // A12 = A11^-1 A12
459  A11.template triangularView<UnitLower>().solveInPlace(A12);
460 
461  A22.noalias() -= A21 * A12;
462  }
463  }
464  return first_zero_pivot;
465  }
466 };
467 
470 template <typename MatrixType, typename TranspositionType>
471 void partial_lu_inplace(MatrixType& lu, TranspositionType& row_transpositions,
472  typename TranspositionType::StorageIndex& nb_transpositions) {
473  // Special-case of zero matrix.
474  if (lu.rows() == 0 || lu.cols() == 0) {
475  nb_transpositions = 0;
476  return;
477  }
478  eigen_assert(lu.cols() == row_transpositions.size());
479  eigen_assert(row_transpositions.size() < 2 ||
480  (&row_transpositions.coeffRef(1) - &row_transpositions.coeffRef(0)) == 1);
481 
482  partial_lu_impl<typename MatrixType::Scalar, MatrixType::Flags & RowMajorBit ? RowMajor : ColMajor,
483  typename TranspositionType::StorageIndex,
484  internal::min_size_prefer_fixed(MatrixType::RowsAtCompileTime, MatrixType::ColsAtCompileTime)>::
485  blocked_lu(lu.rows(), lu.cols(), &lu.coeffRef(0, 0), lu.outerStride(), &row_transpositions.coeffRef(0),
486  nb_transpositions);
487 }
488 
489 } // end namespace internal
490 
491 template <typename MatrixType, typename PermutationIndex>
492 void PartialPivLU<MatrixType, PermutationIndex>::compute() {
493  eigen_assert(m_lu.rows() < NumTraits<PermutationIndex>::highest());
494 
495  if (m_lu.cols() > 0)
496  m_l1_norm = m_lu.cwiseAbs().colwise().sum().maxCoeff();
497  else
498  m_l1_norm = RealScalar(0);
499 
500  eigen_assert(m_lu.rows() == m_lu.cols() && "PartialPivLU is only for square (and moreover invertible) matrices");
501  const Index size = m_lu.rows();
502 
503  m_rowsTranspositions.resize(size);
504 
505  typename TranspositionType::StorageIndex nb_transpositions;
506  internal::partial_lu_inplace(m_lu, m_rowsTranspositions, nb_transpositions);
507  m_det_p = (nb_transpositions % 2) ? -1 : 1;
508 
509  m_p = m_rowsTranspositions;
510 
511  m_isInitialized = true;
512 }
513 
514 template <typename MatrixType, typename PermutationIndex>
515 typename PartialPivLU<MatrixType, PermutationIndex>::Scalar PartialPivLU<MatrixType, PermutationIndex>::determinant()
516  const {
517  eigen_assert(m_isInitialized && "PartialPivLU is not initialized.");
518  return Scalar(m_det_p) * m_lu.diagonal().prod();
519 }
520 
524 template <typename MatrixType, typename PermutationIndex>
526  eigen_assert(m_isInitialized && "LU is not initialized.");
527  // LU
528  MatrixType res = m_lu.template triangularView<UnitLower>().toDenseMatrix() * m_lu.template triangularView<Upper>();
529 
530  // P^{-1}(LU)
531  res = m_p.inverse() * res;
532 
533  return res;
534 }
535 
536 /***** Implementation details *****************************************************/
537 
538 namespace internal {
539 
540 /***** Implementation of inverse() *****************************************************/
541 template <typename DstXprType, typename MatrixType, typename PermutationIndex>
542 struct Assignment<
543  DstXprType, Inverse<PartialPivLU<MatrixType, PermutationIndex> >,
544  internal::assign_op<typename DstXprType::Scalar, typename PartialPivLU<MatrixType, PermutationIndex>::Scalar>,
545  Dense2Dense> {
547  typedef Inverse<LuType> SrcXprType;
548  static void run(DstXprType& dst, const SrcXprType& src,
549  const internal::assign_op<typename DstXprType::Scalar, typename LuType::Scalar>&) {
550  dst = src.nestedExpression().solve(MatrixType::Identity(src.rows(), src.cols()));
551  }
552 };
553 } // end namespace internal
554 
555 /******** MatrixBase methods *******/
556 
563 template <typename Derived>
564 template <typename PermutationIndex>
565 inline const PartialPivLU<typename MatrixBase<Derived>::PlainObject, PermutationIndex>
568 }
569 
578 template <typename Derived>
579 template <typename PermutationIndex>
582 }
583 
584 } // end namespace Eigen
585 
586 #endif // EIGEN_PARTIALLU_H
const Solve< PartialPivLU, Rhs > solve(const MatrixBase< Rhs > &b) const
PartialPivLU()
Default Constructor.
Definition: PartialPivLU.h:281
constexpr Index size() const noexcept
Definition: EigenBase.h:64
Definition: Constants.h:318
const MatrixType & matrixLU() const
Definition: PartialPivLU.h:153
const PermutationType & permutationP() const
Definition: PartialPivLU.h:160
ComputationInfo info() const
Reports whether the LU factorization was successful.
Definition: PartialPivLU.h:99
RealScalar rcond() const
Definition: PartialPivLU.h:190
Scalar determinant() const
Definition: PartialPivLU.h:515
LU decomposition of a matrix with partial pivoting, and related features.
Definition: ForwardDeclarations.h:421
Namespace containing all symbols from the Eigen library.
Definition: B01_Experimental.dox:1
Eigen::Index Index
The interface type of indices.
Definition: EigenBase.h:43
const unsigned int RowMajorBit
Definition: Constants.h:70
Definition: EigenBase.h:33
Expression of the inverse of another expression.
Definition: Inverse.h:43
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:82
MatrixType reconstructedMatrix() const
Definition: PartialPivLU.h:525
Definition: Constants.h:440
constexpr Derived & derived()
Definition: EigenBase.h:49
const ConstTransposeReturnType transpose() const
Definition: SolverBase.h:128
A matrix or vector expression mapping an existing expression.
Definition: Ref.h:264
const Inverse< PartialPivLU > inverse() const
Definition: PartialPivLU.h:202
InverseReturnType transpose() const
Definition: PermutationMatrix.h:180
Definition: Constants.h:320
const int Dynamic
Definition: Constants.h:25
Pseudo expression representing a solving operation.
Definition: Solve.h:62
The matrix class, also used for vectors and row-vectors.
Definition: Matrix.h:186
ComputationInfo
Definition: Constants.h:438
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:52