Merge of itasc branch. Project files, scons and cmake should be working. Makefile...
[blender.git] / extern / Eigen2 / Eigen / src / QR / QR.h
1 // This file is part of Eigen, a lightweight C++ template library
2 // for linear algebra.
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_QR_H
26 #define EIGEN_QR_H
27
28 /** \ingroup QR_Module
29   * \nonstableyet
30   *
31   * \class QR
32   *
33   * \brief QR decomposition of a matrix
34   *
35   * \param MatrixType the type of the matrix of which we are computing the QR decomposition
36   *
37   * This class performs a QR decomposition using Householder transformations. The result is
38   * stored in a compact way compatible with LAPACK.
39   *
40   * \sa MatrixBase::qr()
41   */
42 template<typename MatrixType> class QR
43 {
44   public:
45
46     typedef typename MatrixType::Scalar Scalar;
47     typedef typename MatrixType::RealScalar RealScalar;
48     typedef Block<MatrixType, MatrixType::ColsAtCompileTime, MatrixType::ColsAtCompileTime> MatrixRBlockType;
49     typedef Matrix<Scalar, MatrixType::ColsAtCompileTime, MatrixType::ColsAtCompileTime> MatrixTypeR;
50     typedef Matrix<Scalar, MatrixType::ColsAtCompileTime, 1> VectorType;
51
52     /** 
53     * \brief Default Constructor.
54     *
55     * The default constructor is useful in cases in which the user intends to
56     * perform decompositions via QR::compute(const MatrixType&).
57     */
58     QR() : m_qr(), m_hCoeffs(), m_isInitialized(false) {}
59
60     QR(const MatrixType& matrix)
61       : m_qr(matrix.rows(), matrix.cols()),
62         m_hCoeffs(matrix.cols()),
63         m_isInitialized(false)
64     {
65       compute(matrix);
66     }
67     
68     /** \deprecated use isInjective()
69       * \returns whether or not the matrix is of full rank
70       *
71       * \note Since the rank is computed only once, i.e. the first time it is needed, this
72       *       method almost does not perform any further computation.
73       */
74     EIGEN_DEPRECATED bool isFullRank() const 
75     { 
76       ei_assert(m_isInitialized && "QR is not initialized.");
77       return rank() == m_qr.cols(); 
78     }
79     
80     /** \returns the rank of the matrix of which *this is the QR decomposition.
81       *
82       * \note Since the rank is computed only once, i.e. the first time it is needed, this
83       *       method almost does not perform any further computation.
84       */
85     int rank() const;
86     
87     /** \returns the dimension of the kernel of the matrix of which *this is the QR decomposition.
88       *
89       * \note Since the rank is computed only once, i.e. the first time it is needed, this
90       *       method almost does not perform any further computation.
91       */
92     inline int dimensionOfKernel() const
93     {
94       ei_assert(m_isInitialized && "QR is not initialized.");
95       return m_qr.cols() - rank();
96     }
97     
98     /** \returns true if the matrix of which *this is the QR decomposition represents an injective
99       *          linear map, i.e. has trivial kernel; false otherwise.
100       *
101       * \note Since the rank is computed only once, i.e. the first time it is needed, this
102       *       method almost does not perform any further computation.
103       */
104     inline bool isInjective() const
105     {
106       ei_assert(m_isInitialized && "QR is not initialized.");
107       return rank() == m_qr.cols();
108     }
109     
110     /** \returns true if the matrix of which *this is the QR decomposition represents a surjective
111       *          linear map; false otherwise.
112       *
113       * \note Since the rank is computed only once, i.e. the first time it is needed, this
114       *       method almost does not perform any further computation.
115       */
116     inline bool isSurjective() const
117     {
118       ei_assert(m_isInitialized && "QR is not initialized.");
119       return rank() == m_qr.rows();
120     }
121
122     /** \returns true if the matrix of which *this is the QR decomposition is invertible.
123       *
124       * \note Since the rank is computed only once, i.e. the first time it is needed, this
125       *       method almost does not perform any further computation.
126       */
127     inline bool isInvertible() const
128     {
129       ei_assert(m_isInitialized && "QR is not initialized.");
130       return isInjective() && isSurjective();
131     }
132     
133     /** \returns a read-only expression of the matrix R of the actual the QR decomposition */
134     const Part<NestByValue<MatrixRBlockType>, UpperTriangular>
135     matrixR(void) const
136     {
137       ei_assert(m_isInitialized && "QR is not initialized.");
138       int cols = m_qr.cols();
139       return MatrixRBlockType(m_qr, 0, 0, cols, cols).nestByValue().template part<UpperTriangular>();
140     }
141
142     /** This method finds a solution x to the equation Ax=b, where A is the matrix of which
143       * *this is the QR decomposition, if any exists.
144       *
145       * \param b the right-hand-side of the equation to solve.
146       *
147       * \param result a pointer to the vector/matrix in which to store the solution, if any exists.
148       *          Resized if necessary, so that result->rows()==A.cols() and result->cols()==b.cols().
149       *          If no solution exists, *result is left with undefined coefficients.
150       *
151       * \returns true if any solution exists, false if no solution exists.
152       *
153       * \note If there exist more than one solution, this method will arbitrarily choose one.
154       *       If you need a complete analysis of the space of solutions, take the one solution obtained
155       *       by this method and add to it elements of the kernel, as determined by kernel().
156       *
157       * \note The case where b is a matrix is not yet implemented. Also, this
158       *       code is space inefficient.
159       *
160       * Example: \include QR_solve.cpp
161       * Output: \verbinclude QR_solve.out
162       *
163       * \sa MatrixBase::solveTriangular(), kernel(), computeKernel(), inverse(), computeInverse()
164       */
165     template<typename OtherDerived, typename ResultType>
166     bool solve(const MatrixBase<OtherDerived>& b, ResultType *result) const;
167
168     MatrixType matrixQ(void) const;
169
170     void compute(const MatrixType& matrix);
171
172   protected:
173     MatrixType m_qr;
174     VectorType m_hCoeffs;
175     mutable int m_rank;
176     mutable bool m_rankIsUptodate;
177     bool m_isInitialized;
178 };
179
180 /** \returns the rank of the matrix of which *this is the QR decomposition. */
181 template<typename MatrixType>
182 int QR<MatrixType>::rank() const
183 {
184   ei_assert(m_isInitialized && "QR is not initialized.");
185   if (!m_rankIsUptodate)
186   {
187     RealScalar maxCoeff = m_qr.diagonal().cwise().abs().maxCoeff();
188     int n = m_qr.cols();
189     m_rank = 0;
190     while(m_rank<n && !ei_isMuchSmallerThan(m_qr.diagonal().coeff(m_rank), maxCoeff))
191       ++m_rank;
192     m_rankIsUptodate = true;
193   }
194   return m_rank;
195 }
196
197 #ifndef EIGEN_HIDE_HEAVY_CODE
198
199 template<typename MatrixType>
200 void QR<MatrixType>::compute(const MatrixType& matrix)
201
202   m_rankIsUptodate = false;
203   m_qr = matrix;
204   m_hCoeffs.resize(matrix.cols());
205
206   int rows = matrix.rows();
207   int cols = matrix.cols();
208   RealScalar eps2 = precision<RealScalar>()*precision<RealScalar>();
209
210   for (int k = 0; k < cols; ++k)
211   {
212     int remainingSize = rows-k;
213
214     RealScalar beta;
215     Scalar v0 = m_qr.col(k).coeff(k);
216
217     if (remainingSize==1)
218     {
219       if (NumTraits<Scalar>::IsComplex)
220       {
221         // Householder transformation on the remaining single scalar
222         beta = ei_abs(v0);
223         if (ei_real(v0)>0)
224           beta = -beta;
225         m_qr.coeffRef(k,k) = beta;
226         m_hCoeffs.coeffRef(k) = (beta - v0) / beta;
227       }
228       else
229       {
230         m_hCoeffs.coeffRef(k) = 0;
231       }
232     }
233     else if ((beta=m_qr.col(k).end(remainingSize-1).squaredNorm())>eps2)
234     // FIXME what about ei_imag(v0) ??
235     {
236       // form k-th Householder vector
237       beta = ei_sqrt(ei_abs2(v0)+beta);
238       if (ei_real(v0)>=0.)
239         beta = -beta;
240       m_qr.col(k).end(remainingSize-1) /= v0-beta;
241       m_qr.coeffRef(k,k) = beta;
242       Scalar h = m_hCoeffs.coeffRef(k) = (beta - v0) / beta;
243
244       // apply the Householder transformation (I - h v v') to remaining columns, i.e.,
245       // R <- (I - h v v') * R   where v = [1,m_qr(k+1,k), m_qr(k+2,k), ...]
246       int remainingCols = cols - k -1;
247       if (remainingCols>0)
248       {
249         m_qr.coeffRef(k,k) = Scalar(1);
250         m_qr.corner(BottomRight, remainingSize, remainingCols) -= ei_conj(h) * m_qr.col(k).end(remainingSize)
251             * (m_qr.col(k).end(remainingSize).adjoint() * m_qr.corner(BottomRight, remainingSize, remainingCols));
252         m_qr.coeffRef(k,k) = beta;
253       }
254     }
255     else
256     {
257       m_hCoeffs.coeffRef(k) = 0;
258     }
259   }
260   m_isInitialized = true;
261 }
262
263 template<typename MatrixType>
264 template<typename OtherDerived, typename ResultType>
265 bool QR<MatrixType>::solve(
266   const MatrixBase<OtherDerived>& b,
267   ResultType *result
268 ) const
269 {
270   ei_assert(m_isInitialized && "QR is not initialized.");
271   const int rows = m_qr.rows();
272   ei_assert(b.rows() == rows);
273   result->resize(rows, b.cols());
274   
275   // TODO(keir): There is almost certainly a faster way to multiply by
276   // Q^T without explicitly forming matrixQ(). Investigate.
277   *result = matrixQ().transpose()*b;
278   
279   if(!isSurjective())
280   {
281     // is result is in the image of R ?
282     RealScalar biggest_in_res = result->corner(TopLeft, m_rank, result->cols()).cwise().abs().maxCoeff();
283     for(int col = 0; col < result->cols(); ++col)
284       for(int row = m_rank; row < result->rows(); ++row)
285         if(!ei_isMuchSmallerThan(result->coeff(row,col), biggest_in_res))
286           return false;
287   }
288   m_qr.corner(TopLeft, m_rank, m_rank)
289       .template marked<UpperTriangular>()
290       .solveTriangularInPlace(result->corner(TopLeft, m_rank, result->cols()));
291
292   return true;
293 }
294
295 /** \returns the matrix Q */
296 template<typename MatrixType>
297 MatrixType QR<MatrixType>::matrixQ() const
298 {
299   ei_assert(m_isInitialized && "QR is not initialized.");
300   // compute the product Q_0 Q_1 ... Q_n-1,
301   // where Q_k is the k-th Householder transformation I - h_k v_k v_k'
302   // and v_k is the k-th Householder vector [1,m_qr(k+1,k), m_qr(k+2,k), ...]
303   int rows = m_qr.rows();
304   int cols = m_qr.cols();
305   MatrixType res = MatrixType::Identity(rows, cols);
306   for (int k = cols-1; k >= 0; k--)
307   {
308     // to make easier the computation of the transformation, let's temporarily
309     // overwrite m_qr(k,k) such that the end of m_qr.col(k) is exactly our Householder vector.
310     Scalar beta = m_qr.coeff(k,k);
311     m_qr.const_cast_derived().coeffRef(k,k) = 1;
312     int endLength = rows-k;
313     res.corner(BottomRight,endLength, cols-k) -= ((m_hCoeffs.coeff(k) * m_qr.col(k).end(endLength))
314       * (m_qr.col(k).end(endLength).adjoint() * res.corner(BottomRight,endLength, cols-k)).lazy()).lazy();
315     m_qr.const_cast_derived().coeffRef(k,k) = beta;
316   }
317   return res;
318 }
319
320 #endif // EIGEN_HIDE_HEAVY_CODE
321
322 /** \return the QR decomposition of \c *this.
323   *
324   * \sa class QR
325   */
326 template<typename Derived>
327 const QR<typename MatrixBase<Derived>::PlainMatrixType>
328 MatrixBase<Derived>::qr() const
329 {
330   return QR<PlainMatrixType>(eval());
331 }
332
333
334 #endif // EIGEN_QR_H