Eigen  3.2.10
 
Loading...
Searching...
No Matches
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
14namespace Eigen {
15
47template<typename _MatrixType> class PartialPivLU
48{
49 public:
50
51 typedef _MatrixType MatrixType;
52 enum {
53 RowsAtCompileTime = MatrixType::RowsAtCompileTime,
54 ColsAtCompileTime = MatrixType::ColsAtCompileTime,
55 Options = MatrixType::Options,
56 MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
57 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
58 };
59 typedef typename MatrixType::Scalar Scalar;
60 typedef typename NumTraits<typename MatrixType::Scalar>::Real RealScalar;
61 typedef typename internal::traits<MatrixType>::StorageKind StorageKind;
62 typedef typename MatrixType::Index Index;
65
66
74
81 PartialPivLU(Index size);
82
90 PartialPivLU(const MatrixType& matrix);
91
92 PartialPivLU& compute(const MatrixType& matrix);
93
100 inline const MatrixType& matrixLU() const
101 {
102 eigen_assert(m_isInitialized && "PartialPivLU is not initialized.");
103 return m_lu;
104 }
105
108 inline const PermutationType& permutationP() const
109 {
110 eigen_assert(m_isInitialized && "PartialPivLU is not initialized.");
111 return m_p;
112 }
113
131 template<typename Rhs>
132 inline const internal::solve_retval<PartialPivLU, Rhs>
133 solve(const MatrixBase<Rhs>& b) const
134 {
135 eigen_assert(m_isInitialized && "PartialPivLU is not initialized.");
136 return internal::solve_retval<PartialPivLU, Rhs>(*this, b.derived());
137 }
138
146 inline const internal::solve_retval<PartialPivLU,typename MatrixType::IdentityReturnType> inverse() const
147 {
148 eigen_assert(m_isInitialized && "PartialPivLU is not initialized.");
149 return internal::solve_retval<PartialPivLU,typename MatrixType::IdentityReturnType>
150 (*this, MatrixType::Identity(m_lu.rows(), m_lu.cols()));
151 }
152
166 typename internal::traits<MatrixType>::Scalar determinant() const;
167
169
170 inline Index rows() const { return m_lu.rows(); }
171 inline Index cols() const { return m_lu.cols(); }
172
173 protected:
174
175 static void check_template_parameters()
176 {
177 EIGEN_STATIC_ASSERT_NON_INTEGER(Scalar);
178 }
179
180 MatrixType m_lu;
181 PermutationType m_p;
182 TranspositionType m_rowsTranspositions;
183 Index m_det_p;
184 bool m_isInitialized;
185};
186
187template<typename MatrixType>
189 : m_lu(),
190 m_p(),
191 m_rowsTranspositions(),
192 m_det_p(0),
193 m_isInitialized(false)
194{
195}
196
197template<typename MatrixType>
199 : m_lu(size, size),
200 m_p(size),
201 m_rowsTranspositions(size),
202 m_det_p(0),
203 m_isInitialized(false)
204{
205}
206
207template<typename MatrixType>
209 : m_lu(matrix.rows(), matrix.rows()),
210 m_p(matrix.rows()),
211 m_rowsTranspositions(matrix.rows()),
212 m_det_p(0),
213 m_isInitialized(false)
214{
215 compute(matrix);
216}
217
218namespace internal {
219
221template<typename Scalar, int StorageOrder, typename PivIndex>
222struct partial_lu_impl
223{
224 // FIXME add a stride to Map, so that the following mapping becomes easier,
225 // another option would be to create an expression being able to automatically
226 // warp any Map, Matrix, and Block expressions as a unique type, but since that's exactly
227 // a Map + stride, why not adding a stride to Map, and convenient ctors from a Matrix,
228 // and Block.
231 typedef Block<MatrixType,Dynamic,Dynamic> BlockType;
232 typedef typename MatrixType::RealScalar RealScalar;
233 typedef typename MatrixType::Index Index;
234
245 static Index unblocked_lu(MatrixType& lu, PivIndex* row_transpositions, PivIndex& nb_transpositions)
246 {
247 const Index rows = lu.rows();
248 const Index cols = lu.cols();
249 const Index size = (std::min)(rows,cols);
250 nb_transpositions = 0;
251 Index first_zero_pivot = -1;
252 for(Index k = 0; k < size; ++k)
253 {
254 Index rrows = rows-k-1;
255 Index rcols = cols-k-1;
256
257 Index row_of_biggest_in_col;
258 RealScalar biggest_in_corner
259 = lu.col(k).tail(rows-k).cwiseAbs().maxCoeff(&row_of_biggest_in_col);
260 row_of_biggest_in_col += k;
261
262 row_transpositions[k] = PivIndex(row_of_biggest_in_col);
263
264 if(biggest_in_corner != RealScalar(0))
265 {
266 if(k != row_of_biggest_in_col)
267 {
268 lu.row(k).swap(lu.row(row_of_biggest_in_col));
269 ++nb_transpositions;
270 }
271
272 // FIXME shall we introduce a safe quotient expression in cas 1/lu.coeff(k,k)
273 // overflow but not the actual quotient?
274 lu.col(k).tail(rrows) /= lu.coeff(k,k);
275 }
276 else if(first_zero_pivot==-1)
277 {
278 // the pivot is exactly zero, we record the index of the first pivot which is exactly 0,
279 // and continue the factorization such we still have A = PLU
280 first_zero_pivot = k;
281 }
282
283 if(k<rows-1)
284 lu.bottomRightCorner(rrows,rcols).noalias() -= lu.col(k).tail(rrows) * lu.row(k).tail(rcols);
285 }
286 return first_zero_pivot;
287 }
288
304 static Index blocked_lu(Index rows, Index cols, Scalar* lu_data, Index luStride, PivIndex* row_transpositions, PivIndex& nb_transpositions, Index maxBlockSize=256)
305 {
306 MapLU lu1(lu_data,StorageOrder==RowMajor?rows:luStride,StorageOrder==RowMajor?luStride:cols);
307 MatrixType lu(lu1,0,0,rows,cols);
308
309 const Index size = (std::min)(rows,cols);
310
311 // if the matrix is too small, no blocking:
312 if(size<=16)
313 {
314 return unblocked_lu(lu, row_transpositions, nb_transpositions);
315 }
316
317 // automatically adjust the number of subdivisions to the size
318 // of the matrix so that there is enough sub blocks:
319 Index blockSize;
320 {
321 blockSize = size/8;
322 blockSize = (blockSize/16)*16;
323 blockSize = (std::min)((std::max)(blockSize,Index(8)), maxBlockSize);
324 }
325
326 nb_transpositions = 0;
327 Index first_zero_pivot = -1;
328 for(Index k = 0; k < size; k+=blockSize)
329 {
330 Index bs = (std::min)(size-k,blockSize); // actual size of the block
331 Index trows = rows - k - bs; // trailing rows
332 Index tsize = size - k - bs; // trailing size
334 // partition the matrix:
335 // A00 | A01 | A02
336 // lu = A_0 | A_1 | A_2 = A10 | A11 | A12
337 // A20 | A21 | A22
338 BlockType A_0(lu,0,0,rows,k);
339 BlockType A_2(lu,0,k+bs,rows,tsize);
340 BlockType A11(lu,k,k,bs,bs);
341 BlockType A12(lu,k,k+bs,bs,tsize);
342 BlockType A21(lu,k+bs,k,trows,bs);
343 BlockType A22(lu,k+bs,k+bs,trows,tsize);
344
345 PivIndex nb_transpositions_in_panel;
346 // recursively call the blocked LU algorithm on [A11^T A21^T]^T
347 // with a very small blocking size:
348 Index ret = blocked_lu(trows+bs, bs, &lu.coeffRef(k,k), luStride,
349 row_transpositions+k, nb_transpositions_in_panel, 16);
350 if(ret>=0 && first_zero_pivot==-1)
351 first_zero_pivot = k+ret;
352
353 nb_transpositions += nb_transpositions_in_panel;
354 // update permutations and apply them to A_0
355 for(Index i=k; i<k+bs; ++i)
356 {
357 Index piv = (row_transpositions[i] += k);
358 A_0.row(i).swap(A_0.row(piv));
359 }
360
361 if(trows)
362 {
363 // apply permutations to A_2
364 for(Index i=k;i<k+bs; ++i)
365 A_2.row(i).swap(A_2.row(row_transpositions[i]));
366
367 // A12 = A11^-1 A12
368 A11.template triangularView<UnitLower>().solveInPlace(A12);
369
370 A22.noalias() -= A21 * A12;
371 }
372 }
373 return first_zero_pivot;
374 }
375};
376
379template<typename MatrixType, typename TranspositionType>
380void partial_lu_inplace(MatrixType& lu, TranspositionType& row_transpositions, typename TranspositionType::Index& nb_transpositions)
381{
382 eigen_assert(lu.cols() == row_transpositions.size());
383 eigen_assert((&row_transpositions.coeffRef(1)-&row_transpositions.coeffRef(0)) == 1);
384
385 partial_lu_impl
386 <typename MatrixType::Scalar, MatrixType::Flags&RowMajorBit?RowMajor:ColMajor, typename TranspositionType::Index>
387 ::blocked_lu(lu.rows(), lu.cols(), &lu.coeffRef(0,0), lu.outerStride(), &row_transpositions.coeffRef(0), nb_transpositions);
388}
389
390} // end namespace internal
391
392template<typename MatrixType>
393PartialPivLU<MatrixType>& PartialPivLU<MatrixType>::compute(const MatrixType& matrix)
394{
395 check_template_parameters();
396
397 // the row permutation is stored as int indices, so just to be sure:
398 eigen_assert(matrix.rows()<NumTraits<int>::highest());
399
400 m_lu = matrix;
401
402 eigen_assert(matrix.rows() == matrix.cols() && "PartialPivLU is only for square (and moreover invertible) matrices");
403 const Index size = matrix.rows();
404
405 m_rowsTranspositions.resize(size);
406
407 typename TranspositionType::Index nb_transpositions;
408 internal::partial_lu_inplace(m_lu, m_rowsTranspositions, nb_transpositions);
409 m_det_p = (nb_transpositions%2) ? -1 : 1;
410
411 m_p = m_rowsTranspositions;
412
413 m_isInitialized = true;
414 return *this;
415}
416
417template<typename MatrixType>
418typename internal::traits<MatrixType>::Scalar PartialPivLU<MatrixType>::determinant() const
419{
420 eigen_assert(m_isInitialized && "PartialPivLU is not initialized.");
421 return Scalar(m_det_p) * m_lu.diagonal().prod();
422}
423
427template<typename MatrixType>
429{
430 eigen_assert(m_isInitialized && "LU is not initialized.");
431 // LU
432 MatrixType res = m_lu.template triangularView<UnitLower>().toDenseMatrix()
433 * m_lu.template triangularView<Upper>();
434
435 // P^{-1}(LU)
436 res = m_p.inverse() * res;
437
438 return res;
439}
440
441/***** Implementation of solve() *****************************************************/
442
443namespace internal {
444
445template<typename _MatrixType, typename Rhs>
446struct solve_retval<PartialPivLU<_MatrixType>, Rhs>
447 : solve_retval_base<PartialPivLU<_MatrixType>, Rhs>
448{
449 EIGEN_MAKE_SOLVE_HELPERS(PartialPivLU<_MatrixType>,Rhs)
450
451 template<typename Dest> void evalTo(Dest& dst) const
452 {
453 /* The decomposition PA = LU can be rewritten as A = P^{-1} L U.
454 * So we proceed as follows:
455 * Step 1: compute c = Pb.
456 * Step 2: replace c by the solution x to Lx = c.
457 * Step 3: replace c by the solution x to Ux = c.
458 */
459
460 eigen_assert(rhs().rows() == dec().matrixLU().rows());
461
462 // Step 1
463 dst = dec().permutationP() * rhs();
464
465 // Step 2
466 dec().matrixLU().template triangularView<UnitLower>().solveInPlace(dst);
467
468 // Step 3
469 dec().matrixLU().template triangularView<Upper>().solveInPlace(dst);
470 }
471};
472
473} // end namespace internal
474
475/******** MatrixBase methods *******/
476
483template<typename Derived>
489
490#if EIGEN2_SUPPORT_STAGE > STAGE20_RESOLVE_API_CONFLICTS
499template<typename Derived>
505#endif
506
507} // end namespace Eigen
508
509#endif // EIGEN_PARTIALLU_H
Expression of a fixed-size or dynamic-size block.
Definition Block.h:104
EvalReturnType eval() const
Definition DenseBase.h:360
internal::traits< MatrixWrapper< ExpressionType > >::Index Index
Definition DenseBase.h:60
A matrix or vector expression mapping an existing array of data.
Definition Map.h:106
Base class for all dense matrices, vectors, and expressions.
Definition MatrixBase.h:50
const PartialPivLU< PlainObject > partialPivLu() const
Definition PartialPivLU.h:485
const PartialPivLU< PlainObject > lu() const
LU decomposition of a matrix with partial pivoting, and related features.
Definition PartialPivLU.h:48
internal::traits< MatrixType >::Scalar determinant() const
Definition PartialPivLU.h:418
const PermutationType & permutationP() const
Definition PartialPivLU.h:108
const internal::solve_retval< PartialPivLU, typename MatrixType::IdentityReturnType > inverse() const
Definition PartialPivLU.h:146
MatrixType reconstructedMatrix() const
Definition PartialPivLU.h:428
PartialPivLU()
Default Constructor.
Definition PartialPivLU.h:188
const MatrixType & matrixLU() const
Definition PartialPivLU.h:100
const internal::solve_retval< PartialPivLU, Rhs > solve(const MatrixBase< Rhs > &b) const
Definition PartialPivLU.h:133
Permutation matrix.
Definition PermutationMatrix.h:313
Represents a sequence of transpositions (row/column interchange)
Definition Transpositions.h:157
@ RowMajor
Definition Constants.h:266
@ ColMajor
Definition Constants.h:264
const unsigned int RowMajorBit
Definition Constants.h:53
Definition LDLT.h:18