soc-2008-mxcurioni: merged changes to revision 23516
[blender.git] / extern / Eigen2 / Eigen / src / Cholesky / LDLT.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) 2008 Gael Guennebaud <g.gael@free.fr>
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_LDLT_H
26 #define EIGEN_LDLT_H
27
28 /** \ingroup cholesky_Module
29   *
30   * \class LDLT
31   *
32   * \brief Robust Cholesky decomposition of a matrix and associated features
33   *
34   * \param MatrixType the type of the matrix of which we are computing the LDL^T Cholesky decomposition
35   *
36   * This class performs a Cholesky decomposition without square root of a symmetric, positive definite
37   * matrix A such that A = L D L^* = U^* D U, where L is lower triangular with a unit diagonal
38   * and D is a diagonal matrix.
39   *
40   * Compared to a standard Cholesky decomposition, avoiding the square roots allows for faster and more
41   * stable computation.
42   *
43   * Note that during the decomposition, only the upper triangular part of A is considered. Therefore,
44   * the strict lower part does not have to store correct values.
45   *
46   * \sa MatrixBase::ldlt(), class LLT
47   */
48 template<typename MatrixType> class LDLT
49 {
50   public:
51
52     typedef typename MatrixType::Scalar Scalar;
53     typedef typename NumTraits<typename MatrixType::Scalar>::Real RealScalar;
54     typedef Matrix<Scalar, MatrixType::ColsAtCompileTime, 1> VectorType;
55
56     LDLT(const MatrixType& matrix)
57       : m_matrix(matrix.rows(), matrix.cols())
58     {
59       compute(matrix);
60     }
61
62     /** \returns the lower triangular matrix L */
63     inline Part<MatrixType, UnitLowerTriangular> matrixL(void) const { return m_matrix; }
64
65     /** \returns the coefficients of the diagonal matrix D */
66     inline DiagonalCoeffs<MatrixType> vectorD(void) const { return m_matrix.diagonal(); }
67
68     /** \returns true if the matrix is positive definite */
69     inline bool isPositiveDefinite(void) const { return m_isPositiveDefinite; }
70
71     template<typename RhsDerived, typename ResultType>
72     bool solve(const MatrixBase<RhsDerived> &b, ResultType *result) const;
73
74     template<typename Derived>
75     bool solveInPlace(MatrixBase<Derived> &bAndX) const;
76
77     void compute(const MatrixType& matrix);
78
79   protected:
80     /** \internal
81       * Used to compute and store the cholesky decomposition A = L D L^* = U^* D U.
82       * The strict upper part is used during the decomposition, the strict lower
83       * part correspond to the coefficients of L (its diagonal is equal to 1 and
84       * is not stored), and the diagonal entries correspond to D.
85       */
86     MatrixType m_matrix;
87
88     bool m_isPositiveDefinite;
89 };
90
91 /** Compute / recompute the LLT decomposition A = L D L^* = U^* D U of \a matrix
92   */
93 template<typename MatrixType>
94 void LDLT<MatrixType>::compute(const MatrixType& a)
95 {
96   assert(a.rows()==a.cols());
97   const int size = a.rows();
98   m_matrix.resize(size, size);
99   m_isPositiveDefinite = true;
100   const RealScalar eps = ei_sqrt(precision<Scalar>());
101
102   if (size<=1)
103   {
104     m_matrix = a;
105     return;
106   }
107
108   // Let's preallocate a temporay vector to evaluate the matrix-vector product into it.
109   // Unlike the standard LLT decomposition, here we cannot evaluate it to the destination
110   // matrix because it a sub-row which is not compatible suitable for efficient packet evaluation.
111   // (at least if we assume the matrix is col-major)
112   Matrix<Scalar,MatrixType::RowsAtCompileTime,1> _temporary(size);
113
114   // Note that, in this algorithm the rows of the strict upper part of m_matrix is used to store
115   // column vector, thus the strange .conjugate() and .transpose()...
116
117   m_matrix.row(0) = a.row(0).conjugate();
118   m_matrix.col(0).end(size-1) = m_matrix.row(0).end(size-1) / m_matrix.coeff(0,0);
119   for (int j = 1; j < size; ++j)
120   {
121     RealScalar tmp = ei_real(a.coeff(j,j) - (m_matrix.row(j).start(j) * m_matrix.col(j).start(j).conjugate()).coeff(0,0));
122     m_matrix.coeffRef(j,j) = tmp;
123
124     if (tmp < eps)
125     {
126       m_isPositiveDefinite = false;
127       return;
128     }
129
130     int endSize = size-j-1;
131     if (endSize>0)
132     {
133       _temporary.end(endSize) = ( m_matrix.block(j+1,0, endSize, j)
134                                   * m_matrix.col(j).start(j).conjugate() ).lazy();
135
136       m_matrix.row(j).end(endSize) = a.row(j).end(endSize).conjugate()
137                                    - _temporary.end(endSize).transpose();
138
139       m_matrix.col(j).end(endSize) = m_matrix.row(j).end(endSize) / tmp;
140     }
141   }
142 }
143
144 /** Computes the solution x of \f$ A x = b \f$ using the current decomposition of A.
145   * The result is stored in \a result
146   *
147   * \returns true in case of success, false otherwise.
148   *
149   * In other words, it computes \f$ b = A^{-1} b \f$ with
150   * \f$ {L^{*}}^{-1} D^{-1} L^{-1} b \f$ from right to left.
151   *
152   * \sa LDLT::solveInPlace(), MatrixBase::ldlt()
153   */
154 template<typename MatrixType>
155 template<typename RhsDerived, typename ResultType>
156 bool LDLT<MatrixType>
157 ::solve(const MatrixBase<RhsDerived> &b, ResultType *result) const
158 {
159   const int size = m_matrix.rows();
160   ei_assert(size==b.rows() && "LLT::solve(): invalid number of rows of the right hand side matrix b");
161   *result = b;
162   return solveInPlace(*result);
163 }
164
165 /** This is the \em in-place version of solve().
166   *
167   * \param bAndX represents both the right-hand side matrix b and result x.
168   *
169   * This version avoids a copy when the right hand side matrix b is not
170   * needed anymore.
171   *
172   * \sa LDLT::solve(), MatrixBase::ldlt()
173   */
174 template<typename MatrixType>
175 template<typename Derived>
176 bool LDLT<MatrixType>::solveInPlace(MatrixBase<Derived> &bAndX) const
177 {
178   const int size = m_matrix.rows();
179   ei_assert(size==bAndX.rows());
180   if (!m_isPositiveDefinite)
181     return false;
182   matrixL().solveTriangularInPlace(bAndX);
183   bAndX = (m_matrix.cwise().inverse().template part<Diagonal>() * bAndX).lazy();
184   m_matrix.adjoint().template part<UnitUpperTriangular>().solveTriangularInPlace(bAndX);
185   return true;
186 }
187
188 /** \cholesky_module
189   * \returns the Cholesky decomposition without square root of \c *this
190   */
191 template<typename Derived>
192 inline const LDLT<typename MatrixBase<Derived>::PlainMatrixType>
193 MatrixBase<Derived>::ldlt() const
194 {
195   return derived();
196 }
197
198 #endif // EIGEN_LDLT_H