Merge of itasc branch. Project files, scons and cmake should be working. Makefile...
[blender.git] / extern / Eigen2 / Eigen / src / LU / LU.h
1 // This file is part of Eigen, a lightweight C++ template library
2 // for linear algebra. Eigen itself is part of the KDE project.
3 //
4 // Copyright (C) 2006-2008 Benoit Jacob <jacob.benoit.1@gmail.com>
5 //
6 // Eigen is free software; you can redistribute it and/or
7 // modify it under the terms of the GNU Lesser General Public
8 // License as published by the Free Software Foundation; either
9 // version 3 of the License, or (at your option) any later version.
10 //
11 // Alternatively, you can redistribute it and/or
12 // modify it under the terms of the GNU General Public License as
13 // published by the Free Software Foundation; either version 2 of
14 // the License, or (at your option) any later version.
15 //
16 // Eigen is distributed in the hope that it will be useful, but WITHOUT ANY
17 // WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
18 // FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License or the
19 // GNU General Public License for more details.
20 //
21 // You should have received a copy of the GNU Lesser General Public
22 // License and a copy of the GNU General Public License along with
23 // Eigen. If not, see <http://www.gnu.org/licenses/>.
24
25 #ifndef EIGEN_LU_H
26 #define EIGEN_LU_H
27
28 /** \ingroup LU_Module
29   *
30   * \class LU
31   *
32   * \brief LU decomposition of a matrix with complete pivoting, and related features
33   *
34   * \param MatrixType the type of the matrix of which we are computing the LU decomposition
35   *
36   * This class represents a LU decomposition of any matrix, with complete pivoting: the matrix A
37   * is decomposed as A = PLUQ where L is unit-lower-triangular, U is upper-triangular, and P and Q
38   * are permutation matrices. This is a rank-revealing LU decomposition. The eigenvalues (diagonal
39   * coefficients) of U are sorted in such a way that any zeros are at the end, so that the rank
40   * of A is the index of the first zero on the diagonal of U (with indices starting at 0) if any.
41   *
42   * This decomposition provides the generic approach to solving systems of linear equations, computing
43   * the rank, invertibility, inverse, kernel, and determinant.
44   *
45   * This LU decomposition is very stable and well tested with large matrices. Even exact rank computation
46   * works at sizes larger than 1000x1000. However there are use cases where the SVD decomposition is inherently
47   * more stable when dealing with numerically damaged input. For example, computing the kernel is more stable with
48   * SVD because the SVD can determine which singular values are negligible while LU has to work at the level of matrix
49   * coefficients that are less meaningful in this respect.
50   *
51   * The data of the LU decomposition can be directly accessed through the methods matrixLU(),
52   * permutationP(), permutationQ().
53   *
54   * As an exemple, here is how the original matrix can be retrieved:
55   * \include class_LU.cpp
56   * Output: \verbinclude class_LU.out
57   *
58   * \sa MatrixBase::lu(), MatrixBase::determinant(), MatrixBase::inverse(), MatrixBase::computeInverse()
59   */
60 template<typename MatrixType> class LU
61 {
62   public:
63
64     typedef typename MatrixType::Scalar Scalar;
65     typedef typename NumTraits<typename MatrixType::Scalar>::Real RealScalar;
66     typedef Matrix<int, 1, MatrixType::ColsAtCompileTime> IntRowVectorType;
67     typedef Matrix<int, MatrixType::RowsAtCompileTime, 1> IntColVectorType;
68     typedef Matrix<Scalar, 1, MatrixType::ColsAtCompileTime> RowVectorType;
69     typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, 1> ColVectorType;
70
71     enum { MaxSmallDimAtCompileTime = EIGEN_ENUM_MIN(
72              MatrixType::MaxColsAtCompileTime,
73              MatrixType::MaxRowsAtCompileTime)
74     };
75
76     typedef Matrix<typename MatrixType::Scalar,
77                   MatrixType::ColsAtCompileTime, // the number of rows in the "kernel matrix" is the number of cols of the original matrix
78                                                  // so that the product "matrix * kernel = zero" makes sense
79                   Dynamic,                       // we don't know at compile-time the dimension of the kernel
80                   MatrixType::Options,
81                   MatrixType::MaxColsAtCompileTime, // see explanation for 2nd template parameter
82                   MatrixType::MaxColsAtCompileTime // the kernel is a subspace of the domain space, whose dimension is the number
83                                                    // of columns of the original matrix
84     > KernelResultType;
85
86     typedef Matrix<typename MatrixType::Scalar,
87                    MatrixType::RowsAtCompileTime, // the image is a subspace of the destination space, whose dimension is the number
88                                                   // of rows of the original matrix
89                    Dynamic,                       // we don't know at compile time the dimension of the image (the rank)
90                    MatrixType::Options,
91                    MatrixType::MaxRowsAtCompileTime, // the image matrix will consist of columns from the original matrix,
92                    MatrixType::MaxColsAtCompileTime  // so it has the same number of rows and at most as many columns.
93     > ImageResultType;
94
95     /** Constructor.
96       *
97       * \param matrix the matrix of which to compute the LU decomposition.
98       */
99     LU(const MatrixType& matrix);
100
101     /** \returns the LU decomposition matrix: the upper-triangular part is U, the
102       * unit-lower-triangular part is L (at least for square matrices; in the non-square
103       * case, special care is needed, see the documentation of class LU).
104       *
105       * \sa matrixL(), matrixU()
106       */
107     inline const MatrixType& matrixLU() const
108     {
109       return m_lu;
110     }
111
112     /** \returns a vector of integers, whose size is the number of rows of the matrix being decomposed,
113       * representing the P permutation i.e. the permutation of the rows. For its precise meaning,
114       * see the examples given in the documentation of class LU.
115       *
116       * \sa permutationQ()
117       */
118     inline const IntColVectorType& permutationP() const
119     {
120       return m_p;
121     }
122
123     /** \returns a vector of integers, whose size is the number of columns of the matrix being
124       * decomposed, representing the Q permutation i.e. the permutation of the columns.
125       * For its precise meaning, see the examples given in the documentation of class LU.
126       *
127       * \sa permutationP()
128       */
129     inline const IntRowVectorType& permutationQ() const
130     {
131       return m_q;
132     }
133
134     /** Computes a basis of the kernel of the matrix, also called the null-space of the matrix.
135       *
136       * \note This method is only allowed on non-invertible matrices, as determined by
137       * isInvertible(). Calling it on an invertible matrix will make an assertion fail.
138       *
139       * \param result a pointer to the matrix in which to store the kernel. The columns of this
140       * matrix will be set to form a basis of the kernel (it will be resized
141       * if necessary).
142       *
143       * Example: \include LU_computeKernel.cpp
144       * Output: \verbinclude LU_computeKernel.out
145       *
146       * \sa kernel(), computeImage(), image()
147       */
148     template<typename KernelMatrixType>
149     void computeKernel(KernelMatrixType *result) const;
150
151     /** Computes a basis of the image of the matrix, also called the column-space or range of he matrix.
152       *
153       * \note Calling this method on the zero matrix will make an assertion fail.
154       *
155       * \param result a pointer to the matrix in which to store the image. The columns of this
156       * matrix will be set to form a basis of the image (it will be resized
157       * if necessary).
158       *
159       * Example: \include LU_computeImage.cpp
160       * Output: \verbinclude LU_computeImage.out
161       *
162       * \sa image(), computeKernel(), kernel()
163       */
164     template<typename ImageMatrixType>
165     void computeImage(ImageMatrixType *result) const;
166
167     /** \returns the kernel of the matrix, also called its null-space. The columns of the returned matrix
168       * will form a basis of the kernel.
169       *
170       * \note: this method is only allowed on non-invertible matrices, as determined by
171       * isInvertible(). Calling it on an invertible matrix will make an assertion fail.
172       *
173       * \note: this method returns a matrix by value, which induces some inefficiency.
174       *        If you prefer to avoid this overhead, use computeKernel() instead.
175       *
176       * Example: \include LU_kernel.cpp
177       * Output: \verbinclude LU_kernel.out
178       *
179       * \sa computeKernel(), image()
180       */
181     const KernelResultType kernel() const;
182
183     /** \returns the image of the matrix, also called its column-space. The columns of the returned matrix
184       * will form a basis of the kernel.
185       *
186       * \note: Calling this method on the zero matrix will make an assertion fail.
187       *
188       * \note: this method returns a matrix by value, which induces some inefficiency.
189       *        If you prefer to avoid this overhead, use computeImage() instead.
190       *
191       * Example: \include LU_image.cpp
192       * Output: \verbinclude LU_image.out
193       *
194       * \sa computeImage(), kernel()
195       */
196     const ImageResultType image() const;
197
198     /** This method finds a solution x to the equation Ax=b, where A is the matrix of which
199       * *this is the LU decomposition, if any exists.
200       *
201       * \param b the right-hand-side of the equation to solve. Can be a vector or a matrix,
202       *          the only requirement in order for the equation to make sense is that
203       *          b.rows()==A.rows(), where A is the matrix of which *this is the LU decomposition.
204       * \param result a pointer to the vector or matrix in which to store the solution, if any exists.
205       *          Resized if necessary, so that result->rows()==A.cols() and result->cols()==b.cols().
206       *          If no solution exists, *result is left with undefined coefficients.
207       *
208       * \returns true if any solution exists, false if no solution exists.
209       *
210       * \note If there exist more than one solution, this method will arbitrarily choose one.
211       *       If you need a complete analysis of the space of solutions, take the one solution obtained
212       *       by this method and add to it elements of the kernel, as determined by kernel().
213       *
214       * Example: \include LU_solve.cpp
215       * Output: \verbinclude LU_solve.out
216       *
217       * \sa MatrixBase::solveTriangular(), kernel(), computeKernel(), inverse(), computeInverse()
218       */
219     template<typename OtherDerived, typename ResultType>
220     bool solve(const MatrixBase<OtherDerived>& b, ResultType *result) const;
221
222     /** \returns the determinant of the matrix of which
223       * *this is the LU decomposition. It has only linear complexity
224       * (that is, O(n) where n is the dimension of the square matrix)
225       * as the LU decomposition has already been computed.
226       *
227       * \note This is only for square matrices.
228       *
229       * \note For fixed-size matrices of size up to 4, MatrixBase::determinant() offers
230       *       optimized paths.
231       *
232       * \warning a determinant can be very big or small, so for matrices
233       * of large enough dimension, there is a risk of overflow/underflow.
234       *
235       * \sa MatrixBase::determinant()
236       */
237     typename ei_traits<MatrixType>::Scalar determinant() const;
238
239     /** \returns the rank of the matrix of which *this is the LU decomposition.
240       *
241       * \note This is computed at the time of the construction of the LU decomposition. This
242       *       method does not perform any further computation.
243       */
244     inline int rank() const
245     {
246       return m_rank;
247     }
248
249     /** \returns the dimension of the kernel of the matrix of which *this is the LU decomposition.
250       *
251       * \note Since the rank is computed at the time of the construction of the LU decomposition, this
252       *       method almost does not perform any further computation.
253       */
254     inline int dimensionOfKernel() const
255     {
256       return m_lu.cols() - m_rank;
257     }
258
259     /** \returns true if the matrix of which *this is the LU decomposition represents an injective
260       *          linear map, i.e. has trivial kernel; false otherwise.
261       *
262       * \note Since the rank is computed at the time of the construction of the LU decomposition, this
263       *       method almost does not perform any further computation.
264       */
265     inline bool isInjective() const
266     {
267       return m_rank == m_lu.cols();
268     }
269
270     /** \returns true if the matrix of which *this is the LU decomposition represents a surjective
271       *          linear map; false otherwise.
272       *
273       * \note Since the rank is computed at the time of the construction of the LU decomposition, this
274       *       method almost does not perform any further computation.
275       */
276     inline bool isSurjective() const
277     {
278       return m_rank == m_lu.rows();
279     }
280
281     /** \returns true if the matrix of which *this is the LU decomposition is invertible.
282       *
283       * \note Since the rank is computed at the time of the construction of the LU decomposition, this
284       *       method almost does not perform any further computation.
285       */
286     inline bool isInvertible() const
287     {
288       return isInjective() && isSurjective();
289     }
290
291     /** Computes the inverse of the matrix of which *this is the LU decomposition.
292       *
293       * \param result a pointer to the matrix into which to store the inverse. Resized if needed.
294       *
295       * \note If this matrix is not invertible, *result is left with undefined coefficients.
296       *       Use isInvertible() to first determine whether this matrix is invertible.
297       *
298       * \sa MatrixBase::computeInverse(), inverse()
299       */
300     inline void computeInverse(MatrixType *result) const
301     {
302       solve(MatrixType::Identity(m_lu.rows(), m_lu.cols()), result);
303     }
304
305     /** \returns the inverse of the matrix of which *this is the LU decomposition.
306       *
307       * \note If this matrix is not invertible, the returned matrix has undefined coefficients.
308       *       Use isInvertible() to first determine whether this matrix is invertible.
309       *
310       * \sa computeInverse(), MatrixBase::inverse()
311       */
312     inline MatrixType inverse() const
313     {
314       MatrixType result;
315       computeInverse(&result);
316       return result;
317     }
318
319   protected:
320     const MatrixType& m_originalMatrix;
321     MatrixType m_lu;
322     IntColVectorType m_p;
323     IntRowVectorType m_q;
324     int m_det_pq;
325     int m_rank;
326     RealScalar m_precision;
327 };
328
329 template<typename MatrixType>
330 LU<MatrixType>::LU(const MatrixType& matrix)
331   : m_originalMatrix(matrix),
332     m_lu(matrix),
333     m_p(matrix.rows()),
334     m_q(matrix.cols())
335 {
336   const int size = matrix.diagonal().size();
337   const int rows = matrix.rows();
338   const int cols = matrix.cols();
339   
340   // this formula comes from experimenting (see "LU precision tuning" thread on the list)
341   // and turns out to be identical to Higham's formula used already in LDLt.
342   m_precision = machine_epsilon<Scalar>() * size;
343
344   IntColVectorType rows_transpositions(matrix.rows());
345   IntRowVectorType cols_transpositions(matrix.cols());
346   int number_of_transpositions = 0;
347
348   RealScalar biggest = RealScalar(0);
349   m_rank = size;
350   for(int k = 0; k < size; ++k)
351   {
352     int row_of_biggest_in_corner, col_of_biggest_in_corner;
353     RealScalar biggest_in_corner;
354
355     biggest_in_corner = m_lu.corner(Eigen::BottomRight, rows-k, cols-k)
356                         .cwise().abs()
357                         .maxCoeff(&row_of_biggest_in_corner, &col_of_biggest_in_corner);
358     row_of_biggest_in_corner += k;
359     col_of_biggest_in_corner += k;
360     if(k==0) biggest = biggest_in_corner;
361
362     // if the corner is negligible, then we have less than full rank, and we can finish early
363     if(ei_isMuchSmallerThan(biggest_in_corner, biggest, m_precision))
364     {
365       m_rank = k;
366       for(int i = k; i < size; i++)
367       {
368         rows_transpositions.coeffRef(i) = i;
369         cols_transpositions.coeffRef(i) = i;
370       }
371       break;
372     }
373
374     rows_transpositions.coeffRef(k) = row_of_biggest_in_corner;
375     cols_transpositions.coeffRef(k) = col_of_biggest_in_corner;
376     if(k != row_of_biggest_in_corner) {
377       m_lu.row(k).swap(m_lu.row(row_of_biggest_in_corner));
378       ++number_of_transpositions;
379     }
380     if(k != col_of_biggest_in_corner) {
381       m_lu.col(k).swap(m_lu.col(col_of_biggest_in_corner));
382       ++number_of_transpositions;
383     }
384     if(k<rows-1)
385       m_lu.col(k).end(rows-k-1) /= m_lu.coeff(k,k);
386     if(k<size-1)
387       for(int col = k + 1; col < cols; ++col)
388         m_lu.col(col).end(rows-k-1) -= m_lu.col(k).end(rows-k-1) * m_lu.coeff(k,col);
389   }
390
391   for(int k = 0; k < matrix.rows(); ++k) m_p.coeffRef(k) = k;
392   for(int k = size-1; k >= 0; --k)
393     std::swap(m_p.coeffRef(k), m_p.coeffRef(rows_transpositions.coeff(k)));
394
395   for(int k = 0; k < matrix.cols(); ++k) m_q.coeffRef(k) = k;
396   for(int k = 0; k < size; ++k)
397     std::swap(m_q.coeffRef(k), m_q.coeffRef(cols_transpositions.coeff(k)));
398
399   m_det_pq = (number_of_transpositions%2) ? -1 : 1;
400 }
401
402 template<typename MatrixType>
403 typename ei_traits<MatrixType>::Scalar LU<MatrixType>::determinant() const
404 {
405   return Scalar(m_det_pq) * m_lu.diagonal().redux(ei_scalar_product_op<Scalar>());
406 }
407
408 template<typename MatrixType>
409 template<typename KernelMatrixType>
410 void LU<MatrixType>::computeKernel(KernelMatrixType *result) const
411 {
412   ei_assert(!isInvertible());
413   const int dimker = dimensionOfKernel(), cols = m_lu.cols();
414   result->resize(cols, dimker);
415
416   /* Let us use the following lemma:
417     *
418     * Lemma: If the matrix A has the LU decomposition PAQ = LU,
419     * then Ker A = Q(Ker U).
420     *
421     * Proof: trivial: just keep in mind that P, Q, L are invertible.
422     */
423
424   /* Thus, all we need to do is to compute Ker U, and then apply Q.
425     *
426     * U is upper triangular, with eigenvalues sorted so that any zeros appear at the end.
427     * Thus, the diagonal of U ends with exactly
428     * m_dimKer zero's. Let us use that to construct m_dimKer linearly
429     * independent vectors in Ker U.
430     */
431
432   Matrix<Scalar, Dynamic, Dynamic, MatrixType::Options,
433          MatrixType::MaxColsAtCompileTime, MatrixType::MaxColsAtCompileTime>
434     y(-m_lu.corner(TopRight, m_rank, dimker));
435
436   m_lu.corner(TopLeft, m_rank, m_rank)
437       .template marked<UpperTriangular>()
438       .solveTriangularInPlace(y);
439
440   for(int i = 0; i < m_rank; ++i) result->row(m_q.coeff(i)) = y.row(i);
441   for(int i = m_rank; i < cols; ++i) result->row(m_q.coeff(i)).setZero();
442   for(int k = 0; k < dimker; ++k) result->coeffRef(m_q.coeff(m_rank+k), k) = Scalar(1);
443 }
444
445 template<typename MatrixType>
446 const typename LU<MatrixType>::KernelResultType
447 LU<MatrixType>::kernel() const
448 {
449   KernelResultType result(m_lu.cols(), dimensionOfKernel());
450   computeKernel(&result);
451   return result;
452 }
453
454 template<typename MatrixType>
455 template<typename ImageMatrixType>
456 void LU<MatrixType>::computeImage(ImageMatrixType *result) const
457 {
458   ei_assert(m_rank > 0);
459   result->resize(m_originalMatrix.rows(), m_rank);
460   for(int i = 0; i < m_rank; ++i)
461     result->col(i) = m_originalMatrix.col(m_q.coeff(i));
462 }
463
464 template<typename MatrixType>
465 const typename LU<MatrixType>::ImageResultType
466 LU<MatrixType>::image() const
467 {
468   ImageResultType result(m_originalMatrix.rows(), m_rank);
469   computeImage(&result);
470   return result;
471 }
472
473 template<typename MatrixType>
474 template<typename OtherDerived, typename ResultType>
475 bool LU<MatrixType>::solve(
476   const MatrixBase<OtherDerived>& b,
477   ResultType *result
478 ) const
479 {
480   /* The decomposition PAQ = LU can be rewritten as A = P^{-1} L U Q^{-1}.
481    * So we proceed as follows:
482    * Step 1: compute c = Pb.
483    * Step 2: replace c by the solution x to Lx = c. Exists because L is invertible.
484    * Step 3: replace c by the solution x to Ux = c. Check if a solution really exists.
485    * Step 4: result = Qc;
486    */
487
488   const int rows = m_lu.rows(), cols = m_lu.cols();
489   ei_assert(b.rows() == rows);
490   const int smalldim = std::min(rows, cols);
491
492   typename OtherDerived::PlainMatrixType c(b.rows(), b.cols());
493
494   // Step 1
495   for(int i = 0; i < rows; ++i) c.row(m_p.coeff(i)) = b.row(i);
496
497   // Step 2
498   m_lu.corner(Eigen::TopLeft,smalldim,smalldim).template marked<UnitLowerTriangular>()
499     .solveTriangularInPlace(
500       c.corner(Eigen::TopLeft, smalldim, c.cols()));
501   if(rows>cols)
502   {
503     c.corner(Eigen::BottomLeft, rows-cols, c.cols())
504       -= m_lu.corner(Eigen::BottomLeft, rows-cols, cols) * c.corner(Eigen::TopLeft, cols, c.cols());
505   }
506
507   // Step 3
508   if(!isSurjective())
509   {
510     // is c is in the image of U ?
511     RealScalar biggest_in_c = m_rank>0 ? c.corner(TopLeft, m_rank, c.cols()).cwise().abs().maxCoeff() : 0;
512     for(int col = 0; col < c.cols(); ++col)
513       for(int row = m_rank; row < c.rows(); ++row)
514         if(!ei_isMuchSmallerThan(c.coeff(row,col), biggest_in_c, m_precision))
515           return false;
516   }
517   m_lu.corner(TopLeft, m_rank, m_rank)
518       .template marked<UpperTriangular>()
519       .solveTriangularInPlace(c.corner(TopLeft, m_rank, c.cols()));
520
521   // Step 4
522   result->resize(m_lu.cols(), b.cols());
523   for(int i = 0; i < m_rank; ++i) result->row(m_q.coeff(i)) = c.row(i);
524   for(int i = m_rank; i < m_lu.cols(); ++i) result->row(m_q.coeff(i)).setZero();
525   return true;
526 }
527
528 /** \lu_module
529   *
530   * \return the LU decomposition of \c *this.
531   *
532   * \sa class LU
533   */
534 template<typename Derived>
535 inline const LU<typename MatrixBase<Derived>::PlainMatrixType>
536 MatrixBase<Derived>::lu() const
537 {
538   return LU<PlainMatrixType>(eval());
539 }
540
541 #endif // EIGEN_LU_H