Merge of itasc branch. Project files, scons and cmake should be working. Makefile...
[blender.git] / extern / Eigen2 / Eigen / src / QR / HessenbergDecomposition.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_HESSENBERGDECOMPOSITION_H
26 #define EIGEN_HESSENBERGDECOMPOSITION_H
27
28 /** \ingroup QR_Module
29   * \nonstableyet
30   *
31   * \class HessenbergDecomposition
32   *
33   * \brief Reduces a squared matrix to an Hessemberg form
34   *
35   * \param MatrixType the type of the matrix of which we are computing the Hessenberg decomposition
36   *
37   * This class performs an Hessenberg decomposition of a matrix \f$ A \f$ such that:
38   * \f$ A = Q H Q^* \f$ where \f$ Q \f$ is unitary and \f$ H \f$ a Hessenberg matrix.
39   *
40   * \sa class Tridiagonalization, class Qr
41   */
42 template<typename _MatrixType> class HessenbergDecomposition
43 {
44   public:
45
46     typedef _MatrixType MatrixType;
47     typedef typename MatrixType::Scalar Scalar;
48     typedef typename NumTraits<Scalar>::Real RealScalar;
49
50     enum {
51       Size = MatrixType::RowsAtCompileTime,
52       SizeMinusOne = MatrixType::RowsAtCompileTime==Dynamic
53                    ? Dynamic
54                    : MatrixType::RowsAtCompileTime-1
55     };
56
57     typedef Matrix<Scalar, SizeMinusOne, 1> CoeffVectorType;
58     typedef Matrix<RealScalar, Size, 1> DiagonalType;
59     typedef Matrix<RealScalar, SizeMinusOne, 1> SubDiagonalType;
60
61     typedef typename NestByValue<DiagonalCoeffs<MatrixType> >::RealReturnType DiagonalReturnType;
62
63     typedef typename NestByValue<DiagonalCoeffs<
64         NestByValue<Block<MatrixType,SizeMinusOne,SizeMinusOne> > > >::RealReturnType SubDiagonalReturnType;
65
66     /** This constructor initializes a HessenbergDecomposition object for
67       * further use with HessenbergDecomposition::compute()
68       */
69     HessenbergDecomposition(int size = Size==Dynamic ? 2 : Size)
70       : m_matrix(size,size), m_hCoeffs(size-1)
71     {}
72
73     HessenbergDecomposition(const MatrixType& matrix)
74       : m_matrix(matrix),
75         m_hCoeffs(matrix.cols()-1)
76     {
77       _compute(m_matrix, m_hCoeffs);
78     }
79
80     /** Computes or re-compute the Hessenberg decomposition for the matrix \a matrix.
81       *
82       * This method allows to re-use the allocated data.
83       */
84     void compute(const MatrixType& matrix)
85     {
86       m_matrix = matrix;
87       m_hCoeffs.resize(matrix.rows()-1,1);
88       _compute(m_matrix, m_hCoeffs);
89     }
90
91     /** \returns the householder coefficients allowing to
92       * reconstruct the matrix Q from the packed data.
93       *
94       * \sa packedMatrix()
95       */
96     CoeffVectorType householderCoefficients(void) const { return m_hCoeffs; }
97
98     /** \returns the internal result of the decomposition.
99       *
100       * The returned matrix contains the following information:
101       *  - the upper part and lower sub-diagonal represent the Hessenberg matrix H
102       *  - the rest of the lower part contains the Householder vectors that, combined with
103       *    Householder coefficients returned by householderCoefficients(),
104       *    allows to reconstruct the matrix Q as follow:
105       *       Q = H_{N-1} ... H_1 H_0
106       *    where the matrices H are the Householder transformation:
107       *       H_i = (I - h_i * v_i * v_i')
108       *    where h_i == householderCoefficients()[i] and v_i is a Householder vector:
109       *       v_i = [ 0, ..., 0, 1, M(i+2,i), ..., M(N-1,i) ]
110       *
111       * See LAPACK for further details on this packed storage.
112       */
113     const MatrixType& packedMatrix(void) const { return m_matrix; }
114
115     MatrixType matrixQ(void) const;
116     MatrixType matrixH(void) const;
117
118   private:
119
120     static void _compute(MatrixType& matA, CoeffVectorType& hCoeffs);
121
122   protected:
123     MatrixType m_matrix;
124     CoeffVectorType m_hCoeffs;
125 };
126
127 #ifndef EIGEN_HIDE_HEAVY_CODE
128
129 /** \internal
130   * Performs a tridiagonal decomposition of \a matA in place.
131   *
132   * \param matA the input selfadjoint matrix
133   * \param hCoeffs returned Householder coefficients
134   *
135   * The result is written in the lower triangular part of \a matA.
136   *
137   * Implemented from Golub's "Matrix Computations", algorithm 8.3.1.
138   *
139   * \sa packedMatrix()
140   */
141 template<typename MatrixType>
142 void HessenbergDecomposition<MatrixType>::_compute(MatrixType& matA, CoeffVectorType& hCoeffs)
143 {
144   assert(matA.rows()==matA.cols());
145   int n = matA.rows();
146   for (int i = 0; i<n-2; ++i)
147   {
148     // let's consider the vector v = i-th column starting at position i+1
149
150     // start of the householder transformation
151     // squared norm of the vector v skipping the first element
152     RealScalar v1norm2 = matA.col(i).end(n-(i+2)).squaredNorm();
153
154     if (ei_isMuchSmallerThan(v1norm2,static_cast<Scalar>(1)))
155     {
156       hCoeffs.coeffRef(i) = 0.;
157     }
158     else
159     {
160       Scalar v0 = matA.col(i).coeff(i+1);
161       RealScalar beta = ei_sqrt(ei_abs2(v0)+v1norm2);
162       if (ei_real(v0)>=0.)
163         beta = -beta;
164       matA.col(i).end(n-(i+2)) *= (Scalar(1)/(v0-beta));
165       matA.col(i).coeffRef(i+1) = beta;
166       Scalar h = (beta - v0) / beta;
167       // end of the householder transformation
168
169       // Apply similarity transformation to remaining columns,
170       // i.e., A = H' A H where H = I - h v v' and v = matA.col(i).end(n-i-1)
171       matA.col(i).coeffRef(i+1) = 1;
172
173       // first let's do A = H A
174       matA.corner(BottomRight,n-i-1,n-i-1) -= ((ei_conj(h) * matA.col(i).end(n-i-1)) *
175         (matA.col(i).end(n-i-1).adjoint() * matA.corner(BottomRight,n-i-1,n-i-1))).lazy();
176
177       // now let's do A = A H
178       matA.corner(BottomRight,n,n-i-1) -= ((matA.corner(BottomRight,n,n-i-1) * matA.col(i).end(n-i-1))
179                                         * (h * matA.col(i).end(n-i-1).adjoint())).lazy();
180
181       matA.col(i).coeffRef(i+1) = beta;
182       hCoeffs.coeffRef(i) = h;
183     }
184   }
185   if (NumTraits<Scalar>::IsComplex)
186   {
187     // Householder transformation on the remaining single scalar
188     int i = n-2;
189     Scalar v0 = matA.coeff(i+1,i);
190
191     RealScalar beta = ei_sqrt(ei_abs2(v0));
192     if (ei_real(v0)>=0.)
193       beta = -beta;
194     Scalar h = (beta - v0) / beta;
195     hCoeffs.coeffRef(i) = h;
196
197     // A = H* A
198     matA.corner(BottomRight,n-i-1,n-i) -= ei_conj(h) * matA.corner(BottomRight,n-i-1,n-i);
199
200     // A = A H
201     matA.col(n-1) -= h * matA.col(n-1);
202   }
203   else
204   {
205     hCoeffs.coeffRef(n-2) = 0;
206   }
207 }
208
209 /** reconstructs and returns the matrix Q */
210 template<typename MatrixType>
211 typename HessenbergDecomposition<MatrixType>::MatrixType
212 HessenbergDecomposition<MatrixType>::matrixQ(void) const
213 {
214   int n = m_matrix.rows();
215   MatrixType matQ = MatrixType::Identity(n,n);
216   for (int i = n-2; i>=0; i--)
217   {
218     Scalar tmp = m_matrix.coeff(i+1,i);
219     m_matrix.const_cast_derived().coeffRef(i+1,i) = 1;
220
221     matQ.corner(BottomRight,n-i-1,n-i-1) -=
222       ((m_hCoeffs.coeff(i) * m_matrix.col(i).end(n-i-1)) *
223       (m_matrix.col(i).end(n-i-1).adjoint() * matQ.corner(BottomRight,n-i-1,n-i-1)).lazy()).lazy();
224
225     m_matrix.const_cast_derived().coeffRef(i+1,i) = tmp;
226   }
227   return matQ;
228 }
229
230 #endif // EIGEN_HIDE_HEAVY_CODE
231
232 /** constructs and returns the matrix H.
233   * Note that the matrix H is equivalent to the upper part of the packed matrix
234   * (including the lower sub-diagonal). Therefore, it might be often sufficient
235   * to directly use the packed matrix instead of creating a new one.
236   */
237 template<typename MatrixType>
238 typename HessenbergDecomposition<MatrixType>::MatrixType
239 HessenbergDecomposition<MatrixType>::matrixH(void) const
240 {
241   // FIXME should this function (and other similar) rather take a matrix as argument
242   // and fill it (to avoid temporaries)
243   int n = m_matrix.rows();
244   MatrixType matH = m_matrix;
245   if (n>2)
246     matH.corner(BottomLeft,n-2, n-2).template part<LowerTriangular>().setZero();
247   return matH;
248 }
249
250 #endif // EIGEN_HESSENBERGDECOMPOSITION_H