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 // Eigen is free software; you can redistribute it and/or
8 // modify it under the terms of the GNU Lesser General Public
9 // License as published by the Free Software Foundation; either
10 // version 3 of the License, or (at your option) any later version.
11 //
12 // Alternatively, you can redistribute it and/or
13 // modify it under the terms of the GNU General Public License as
14 // published by the Free Software Foundation; either version 2 of
15 // the License, or (at your option) any later version.
16 //
17 // Eigen is distributed in the hope that it will be useful, but WITHOUT ANY
18 // WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
19 // FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License or the
20 // GNU General Public License for more details.
21 //
22 // You should have received a copy of the GNU Lesser General Public
23 // License and a copy of the GNU General Public License along with
24 // Eigen. If not, see <http://www.gnu.org/licenses/>.
25 
26 #ifndef EIGEN_PARTIALLU_H
27 #define EIGEN_PARTIALLU_H
28 
29 namespace Eigen {
30 
62 template<typename _MatrixType> class PartialPivLU
63 {
64  public:
65 
66  typedef _MatrixType MatrixType;
67  enum {
72  MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
73  };
74  typedef typename MatrixType::Scalar Scalar;
76  typedef typename internal::traits<MatrixType>::StorageKind StorageKind;
77  typedef typename MatrixType::Index Index;
80 
81 
88  PartialPivLU();
89 
96  PartialPivLU(Index size);
97 
105  PartialPivLU(const MatrixType& matrix);
106 
107  PartialPivLU& compute(const MatrixType& matrix);
108 
115  inline const MatrixType& matrixLU() const
116  {
117  eigen_assert(m_isInitialized && "PartialPivLU is not initialized.");
118  return m_lu;
119  }
120 
123  inline const PermutationType& permutationP() const
124  {
125  eigen_assert(m_isInitialized && "PartialPivLU is not initialized.");
126  return m_p;
127  }
128 
146  template<typename Rhs>
147  inline const internal::solve_retval<PartialPivLU, Rhs>
148  solve(const MatrixBase<Rhs>& b) const
149  {
150  eigen_assert(m_isInitialized && "PartialPivLU is not initialized.");
151  return internal::solve_retval<PartialPivLU, Rhs>(*this, b.derived());
152  }
153 
161  inline const internal::solve_retval<PartialPivLU,typename MatrixType::IdentityReturnType> inverse() const
162  {
163  eigen_assert(m_isInitialized && "PartialPivLU is not initialized.");
164  return internal::solve_retval<PartialPivLU,typename MatrixType::IdentityReturnType>
165  (*this, MatrixType::Identity(m_lu.rows(), m_lu.cols()));
166  }
167 
181  typename internal::traits<MatrixType>::Scalar determinant() const;
182 
184 
185  inline Index rows() const { return m_lu.rows(); }
186  inline Index cols() const { return m_lu.cols(); }
187 
188  protected:
194 };
195 
196 template<typename MatrixType>
198  : m_lu(),
199  m_p(),
200  m_rowsTranspositions(),
201  m_det_p(0),
202  m_isInitialized(false)
203 {
204 }
205 
206 template<typename MatrixType>
208  : m_lu(size, size),
209  m_p(size),
210  m_rowsTranspositions(size),
211  m_det_p(0),
212  m_isInitialized(false)
213 {
214 }
215 
216 template<typename MatrixType>
218  : m_lu(matrix.rows(), matrix.rows()),
219  m_p(matrix.rows()),
220  m_rowsTranspositions(matrix.rows()),
221  m_det_p(0),
222  m_isInitialized(false)
223 {
224  compute(matrix);
225 }
226 
227 namespace internal {
228 
230 template<typename Scalar, int StorageOrder, typename PivIndex>
231 struct partial_lu_impl
232 {
233  // FIXME add a stride to Map, so that the following mapping becomes easier,
234  // another option would be to create an expression being able to automatically
235  // warp any Map, Matrix, and Block expressions as a unique type, but since that's exactly
236  // a Map + stride, why not adding a stride to Map, and convenient ctors from a Matrix,
237  // and Block.
239  typedef Block<MapLU, Dynamic, Dynamic> MatrixType;
240  typedef Block<MatrixType,Dynamic,Dynamic> BlockType;
241  typedef typename MatrixType::RealScalar RealScalar;
242  typedef typename MatrixType::Index Index;
243 
254  static Index unblocked_lu(MatrixType& lu, PivIndex* row_transpositions, PivIndex& nb_transpositions)
255  {
256  const Index rows = lu.rows();
257  const Index cols = lu.cols();
258  const Index size = (std::min)(rows,cols);
259  nb_transpositions = 0;
260  int first_zero_pivot = -1;
261  for(Index k = 0; k < size; ++k)
262  {
263  Index rrows = rows-k-1;
264  Index rcols = cols-k-1;
265 
266  Index row_of_biggest_in_col;
267  RealScalar biggest_in_corner
268  = lu.col(k).tail(rows-k).cwiseAbs().maxCoeff(&row_of_biggest_in_col);
269  row_of_biggest_in_col += k;
270 
271  row_transpositions[k] = row_of_biggest_in_col;
272 
273  if(biggest_in_corner != RealScalar(0))
274  {
275  if(k != row_of_biggest_in_col)
276  {
277  lu.row(k).swap(lu.row(row_of_biggest_in_col));
278  ++nb_transpositions;
279  }
280 
281  // FIXME shall we introduce a safe quotient expression in cas 1/lu.coeff(k,k)
282  // overflow but not the actual quotient?
283  lu.col(k).tail(rrows) /= lu.coeff(k,k);
284  }
285  else if(first_zero_pivot==-1)
286  {
287  // the pivot is exactly zero, we record the index of the first pivot which is exactly 0,
288  // and continue the factorization such we still have A = PLU
289  first_zero_pivot = k;
290  }
291 
292  if(k<rows-1)
293  lu.bottomRightCorner(rrows,rcols).noalias() -= lu.col(k).tail(rrows) * lu.row(k).tail(rcols);
294  }
295  return first_zero_pivot;
296  }
297 
313  static Index blocked_lu(Index rows, Index cols, Scalar* lu_data, Index luStride, PivIndex* row_transpositions, PivIndex& nb_transpositions, Index maxBlockSize=256)
314  {
315  MapLU lu1(lu_data,StorageOrder==RowMajor?rows:luStride,StorageOrder==RowMajor?luStride:cols);
316  MatrixType lu(lu1,0,0,rows,cols);
317 
318  const Index size = (std::min)(rows,cols);
319 
320  // if the matrix is too small, no blocking:
321  if(size<=16)
322  {
323  return unblocked_lu(lu, row_transpositions, nb_transpositions);
324  }
325 
326  // automatically adjust the number of subdivisions to the size
327  // of the matrix so that there is enough sub blocks:
328  Index blockSize;
329  {
330  blockSize = size/8;
331  blockSize = (blockSize/16)*16;
332  blockSize = (std::min)((std::max)(blockSize,Index(8)), maxBlockSize);
333  }
334 
335  nb_transpositions = 0;
336  int first_zero_pivot = -1;
337  for(Index k = 0; k < size; k+=blockSize)
338  {
339  Index bs = (std::min)(size-k,blockSize); // actual size of the block
340  Index trows = rows - k - bs; // trailing rows
341  Index tsize = size - k - bs; // trailing size
342 
343  // partition the matrix:
344  // A00 | A01 | A02
345  // lu = A_0 | A_1 | A_2 = A10 | A11 | A12
346  // A20 | A21 | A22
347  BlockType A_0(lu,0,0,rows,k);
348  BlockType A_2(lu,0,k+bs,rows,tsize);
349  BlockType A11(lu,k,k,bs,bs);
350  BlockType A12(lu,k,k+bs,bs,tsize);
351  BlockType A21(lu,k+bs,k,trows,bs);
352  BlockType A22(lu,k+bs,k+bs,trows,tsize);
353 
354  PivIndex nb_transpositions_in_panel;
355  // recursively call the blocked LU algorithm on [A11^T A21^T]^T
356  // with a very small blocking size:
357  Index ret = blocked_lu(trows+bs, bs, &lu.coeffRef(k,k), luStride,
358  row_transpositions+k, nb_transpositions_in_panel, 16);
359  if(ret>=0 && first_zero_pivot==-1)
360  first_zero_pivot = k+ret;
361 
362  nb_transpositions += nb_transpositions_in_panel;
363  // update permutations and apply them to A_0
364  for(Index i=k; i<k+bs; ++i)
365  {
366  Index piv = (row_transpositions[i] += k);
367  A_0.row(i).swap(A_0.row(piv));
368  }
369 
370  if(trows)
371  {
372  // apply permutations to A_2
373  for(Index i=k;i<k+bs; ++i)
374  A_2.row(i).swap(A_2.row(row_transpositions[i]));
375 
376  // A12 = A11^-1 A12
377  A11.template triangularView<UnitLower>().solveInPlace(A12);
378 
379  A22.noalias() -= A21 * A12;
380  }
381  }
382  return first_zero_pivot;
383  }
384 };
385 
388 template<typename MatrixType, typename TranspositionType>
389 void partial_lu_inplace(MatrixType& lu, TranspositionType& row_transpositions, typename TranspositionType::Index& nb_transpositions)
390 {
391  eigen_assert(lu.cols() == row_transpositions.size());
392  eigen_assert((&row_transpositions.coeffRef(1)-&row_transpositions.coeffRef(0)) == 1);
393 
394  partial_lu_impl
395  <typename MatrixType::Scalar, MatrixType::Flags&RowMajorBit?RowMajor:ColMajor, typename TranspositionType::Index>
396  ::blocked_lu(lu.rows(), lu.cols(), &lu.coeffRef(0,0), lu.outerStride(), &row_transpositions.coeffRef(0), nb_transpositions);
397 }
398 
399 } // end namespace internal
400 
401 template<typename MatrixType>
403 {
404  m_lu = matrix;
405 
406  eigen_assert(matrix.rows() == matrix.cols() && "PartialPivLU is only for square (and moreover invertible) matrices");
407  const Index size = matrix.rows();
408 
409  m_rowsTranspositions.resize(size);
410 
411  typename TranspositionType::Index nb_transpositions;
412  internal::partial_lu_inplace(m_lu, m_rowsTranspositions, nb_transpositions);
413  m_det_p = (nb_transpositions%2) ? -1 : 1;
414 
415  m_p = m_rowsTranspositions;
416 
417  m_isInitialized = true;
418  return *this;
419 }
420 
421 template<typename MatrixType>
422 typename internal::traits<MatrixType>::Scalar PartialPivLU<MatrixType>::determinant() const
423 {
424  eigen_assert(m_isInitialized && "PartialPivLU is not initialized.");
425  return Scalar(m_det_p) * m_lu.diagonal().prod();
426 }
427 
431 template<typename MatrixType>
433 {
434  eigen_assert(m_isInitialized && "LU is not initialized.");
435  // LU
436  MatrixType res = m_lu.template triangularView<UnitLower>().toDenseMatrix()
437  * m_lu.template triangularView<Upper>();
438 
439  // P^{-1}(LU)
440  res = m_p.inverse() * res;
441 
442  return res;
443 }
444 
445 /***** Implementation of solve() *****************************************************/
446 
447 namespace internal {
448 
449 template<typename _MatrixType, typename Rhs>
450 struct solve_retval<PartialPivLU<_MatrixType>, Rhs>
451  : solve_retval_base<PartialPivLU<_MatrixType>, Rhs>
452 {
454 
455  template<typename Dest> void evalTo(Dest& dst) const
456  {
457  /* The decomposition PA = LU can be rewritten as A = P^{-1} L U.
458  * So we proceed as follows:
459  * Step 1: compute c = Pb.
460  * Step 2: replace c by the solution x to Lx = c.
461  * Step 3: replace c by the solution x to Ux = c.
462  */
463 
464  eigen_assert(rhs().rows() == dec().matrixLU().rows());
465 
466  // Step 1
467  dst = dec().permutationP() * rhs();
468 
469  // Step 2
470  dec().matrixLU().template triangularView<UnitLower>().solveInPlace(dst);
471 
472  // Step 3
473  dec().matrixLU().template triangularView<Upper>().solveInPlace(dst);
474  }
475 };
476 
477 } // end namespace internal
478 
479 /******** MatrixBase methods *******/
480 
487 template<typename Derived>
488 inline const PartialPivLU<typename MatrixBase<Derived>::PlainObject>
490 {
491  return PartialPivLU<PlainObject>(eval());
492 }
493 
494 #if EIGEN2_SUPPORT_STAGE > STAGE20_RESOLVE_API_CONFLICTS
495 
503 template<typename Derived>
506 {
507  return PartialPivLU<PlainObject>(eval());
508 }
509 #endif
510 
511 } // end namespace Eigen
512 
513 #endif // EIGEN_PARTIALLU_H