1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
|
// This file is part of Eigen, a lightweight C++ template library
// for linear algebra.
//
// Copyright (C) 2008 Gael Guennebaud <g.gael@free.fr>
//
// Eigen is free software; you can redistribute it and/or
// modify it under the terms of the GNU Lesser General Public
// License as published by the Free Software Foundation; either
// version 3 of the License, or (at your option) any later version.
//
// Alternatively, you can redistribute it and/or
// modify it under the terms of the GNU General Public License as
// published by the Free Software Foundation; either version 2 of
// the License, or (at your option) any later version.
//
// Eigen is distributed in the hope that it will be useful, but WITHOUT ANY
// WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
// FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License or the
// GNU General Public License for more details.
//
// You should have received a copy of the GNU Lesser General Public
// License and a copy of the GNU General Public License along with
// Eigen. If not, see <http://www.gnu.org/licenses/>.
#ifndef EIGEN_TRIDIAGONALIZATION_H
#define EIGEN_TRIDIAGONALIZATION_H
/** \eigenvalues_module \ingroup Eigenvalues_Module
* \nonstableyet
*
* \class Tridiagonalization
*
* \brief Trigiagonal decomposition of a selfadjoint matrix
*
* \param MatrixType the type of the matrix of which we are performing the tridiagonalization
*
* This class performs a tridiagonal decomposition of a selfadjoint matrix \f$ A \f$ such that:
* \f$ A = Q T Q^* \f$ where \f$ Q \f$ is unitary and \f$ T \f$ a real symmetric tridiagonal matrix.
*
* \sa MatrixBase::tridiagonalize()
*/
template<typename _MatrixType> class Tridiagonalization
{
public:
typedef _MatrixType MatrixType;
typedef typename MatrixType::Scalar Scalar;
typedef typename NumTraits<Scalar>::Real RealScalar;
typedef typename ei_packet_traits<Scalar>::type Packet;
enum {
Size = MatrixType::RowsAtCompileTime,
SizeMinusOne = Size == Dynamic ? Dynamic : Size - 1,
Options = MatrixType::Options,
MaxSize = MatrixType::MaxRowsAtCompileTime,
MaxSizeMinusOne = MaxSize == Dynamic ? Dynamic : MaxSize - 1,
PacketSize = ei_packet_traits<Scalar>::size
};
typedef Matrix<Scalar, SizeMinusOne, 1, Options & ~RowMajor, MaxSizeMinusOne, 1> CoeffVectorType;
typedef typename ei_plain_col_type<MatrixType, RealScalar>::type DiagonalType;
typedef Matrix<RealScalar, SizeMinusOne, 1, Options & ~RowMajor, MaxSizeMinusOne, 1> SubDiagonalType;
typedef typename ei_plain_row_type<MatrixType>::type RowVectorType;
typedef typename ei_meta_if<NumTraits<Scalar>::IsComplex,
typename Diagonal<MatrixType,0>::RealReturnType,
Diagonal<MatrixType,0>
>::ret DiagonalReturnType;
typedef typename ei_meta_if<NumTraits<Scalar>::IsComplex,
typename Diagonal<
Block<MatrixType,SizeMinusOne,SizeMinusOne>,0 >::RealReturnType,
Diagonal<
Block<MatrixType,SizeMinusOne,SizeMinusOne>,0 >
>::ret SubDiagonalReturnType;
/** This constructor initializes a Tridiagonalization object for
* further use with Tridiagonalization::compute()
*/
Tridiagonalization(int size = Size==Dynamic ? 2 : Size)
: m_matrix(size,size), m_hCoeffs(size-1)
{}
Tridiagonalization(const MatrixType& matrix)
: m_matrix(matrix), m_hCoeffs(matrix.cols()-1)
{
_compute(m_matrix, m_hCoeffs);
}
/** Computes or re-compute the tridiagonalization for the matrix \a matrix.
*
* This method allows to re-use the allocated data.
*/
void compute(const MatrixType& matrix)
{
m_matrix = matrix;
m_hCoeffs.resize(matrix.rows()-1, 1);
_compute(m_matrix, m_hCoeffs);
}
/** \returns the householder coefficients allowing to
* reconstruct the matrix Q from the packed data.
*
* \sa packedMatrix()
*/
inline CoeffVectorType householderCoefficients(void) const { return m_hCoeffs; }
/** \returns the internal result of the decomposition.
*
* The returned matrix contains the following information:
* - the strict upper part is equal to the input matrix A
* - the diagonal and lower sub-diagonal represent the tridiagonal symmetric matrix (real).
* - the rest of the lower part contains the Householder vectors that, combined with
* Householder coefficients returned by householderCoefficients(),
* allows to reconstruct the matrix Q as follow:
* Q = H_{N-1} ... H_1 H_0
* where the matrices H are the Householder transformations:
* H_i = (I - h_i * v_i * v_i')
* where h_i == householderCoefficients()[i] and v_i is a Householder vector:
* v_i = [ 0, ..., 0, 1, M(i+2,i), ..., M(N-1,i) ]
*
* See LAPACK for further details on this packed storage.
*/
inline const MatrixType& packedMatrix(void) const { return m_matrix; }
MatrixType matrixQ() const;
template<typename QDerived> void matrixQInPlace(MatrixBase<QDerived>* q) const;
MatrixType matrixT() const;
const DiagonalReturnType diagonal(void) const;
const SubDiagonalReturnType subDiagonal(void) const;
static void decomposeInPlace(MatrixType& mat, DiagonalType& diag, SubDiagonalType& subdiag, bool extractQ = true);
static void _compute(MatrixType& matA, CoeffVectorType& hCoeffs);
protected:
static void _decomposeInPlace3x3(MatrixType& mat, DiagonalType& diag, SubDiagonalType& subdiag, bool extractQ = true);
MatrixType m_matrix;
CoeffVectorType m_hCoeffs;
};
/** \returns an expression of the diagonal vector */
template<typename MatrixType>
const typename Tridiagonalization<MatrixType>::DiagonalReturnType
Tridiagonalization<MatrixType>::diagonal(void) const
{
return m_matrix.diagonal();
}
/** \returns an expression of the sub-diagonal vector */
template<typename MatrixType>
const typename Tridiagonalization<MatrixType>::SubDiagonalReturnType
Tridiagonalization<MatrixType>::subDiagonal(void) const
{
int n = m_matrix.rows();
return Block<MatrixType,SizeMinusOne,SizeMinusOne>(m_matrix, 1, 0, n-1,n-1).diagonal();
}
/** constructs and returns the tridiagonal matrix T.
* Note that the matrix T is equivalent to the diagonal and sub-diagonal of the packed matrix.
* Therefore, it might be often sufficient to directly use the packed matrix, or the vector
* expressions returned by diagonal() and subDiagonal() instead of creating a new matrix.
*/
template<typename MatrixType>
typename Tridiagonalization<MatrixType>::MatrixType
Tridiagonalization<MatrixType>::matrixT(void) const
{
// FIXME should this function (and other similar ones) rather take a matrix as argument
// and fill it ? (to avoid temporaries)
int n = m_matrix.rows();
MatrixType matT = m_matrix;
matT.corner(TopRight,n-1, n-1).diagonal() = subDiagonal().template cast<Scalar>().conjugate();
if (n>2)
{
matT.corner(TopRight,n-2, n-2).template triangularView<Upper>().setZero();
matT.corner(BottomLeft,n-2, n-2).template triangularView<Lower>().setZero();
}
return matT;
}
#ifndef EIGEN_HIDE_HEAVY_CODE
/** \internal
* Performs a tridiagonal decomposition of \a matA in place.
*
* \param matA the input selfadjoint matrix
* \param hCoeffs returned Householder coefficients
*
* The result is written in the lower triangular part of \a matA.
*
* Implemented from Golub's "Matrix Computations", algorithm 8.3.1.
*
* \sa packedMatrix()
*/
template<typename MatrixType>
void Tridiagonalization<MatrixType>::_compute(MatrixType& matA, CoeffVectorType& hCoeffs)
{
assert(matA.rows()==matA.cols());
int n = matA.rows();
for (int i = 0; i<n-1; ++i)
{
int remainingSize = n-i-1;
RealScalar beta;
Scalar h;
matA.col(i).tail(remainingSize).makeHouseholderInPlace(h, beta);
// Apply similarity transformation to remaining columns,
// i.e., A = H A H' where H = I - h v v' and v = matA.col(i).tail(n-i-1)
matA.col(i).coeffRef(i+1) = 1;
hCoeffs.tail(n-i-1) = (matA.corner(BottomRight,remainingSize,remainingSize).template selfadjointView<Lower>()
* (ei_conj(h) * matA.col(i).tail(remainingSize)));
hCoeffs.tail(n-i-1) += (ei_conj(h)*Scalar(-0.5)*(hCoeffs.tail(remainingSize).dot(matA.col(i).tail(remainingSize)))) * matA.col(i).tail(n-i-1);
matA.corner(BottomRight, remainingSize, remainingSize).template selfadjointView<Lower>()
.rankUpdate(matA.col(i).tail(remainingSize), hCoeffs.tail(remainingSize), -1);
matA.col(i).coeffRef(i+1) = beta;
hCoeffs.coeffRef(i) = h;
}
}
/** reconstructs and returns the matrix Q */
template<typename MatrixType>
typename Tridiagonalization<MatrixType>::MatrixType
Tridiagonalization<MatrixType>::matrixQ(void) const
{
MatrixType matQ;
matrixQInPlace(&matQ);
return matQ;
}
template<typename MatrixType>
template<typename QDerived>
void Tridiagonalization<MatrixType>::matrixQInPlace(MatrixBase<QDerived>* q) const
{
QDerived& matQ = q->derived();
int n = m_matrix.rows();
matQ = MatrixType::Identity(n,n);
RowVectorType aux(n);
for (int i = n-2; i>=0; i--)
{
matQ.corner(BottomRight,n-i-1,n-i-1)
.applyHouseholderOnTheLeft(m_matrix.col(i).tail(n-i-2), ei_conj(m_hCoeffs.coeff(i)), &aux.coeffRef(0,0));
}
}
/** Performs a full decomposition in place */
template<typename MatrixType>
void Tridiagonalization<MatrixType>::decomposeInPlace(MatrixType& mat, DiagonalType& diag, SubDiagonalType& subdiag, bool extractQ)
{
int n = mat.rows();
ei_assert(mat.cols()==n && diag.size()==n && subdiag.size()==n-1);
if (n==3 && (!NumTraits<Scalar>::IsComplex) )
{
_decomposeInPlace3x3(mat, diag, subdiag, extractQ);
}
else
{
Tridiagonalization tridiag(mat);
diag = tridiag.diagonal();
subdiag = tridiag.subDiagonal();
if (extractQ)
tridiag.matrixQInPlace(&mat);
}
}
/** \internal
* Optimized path for 3x3 matrices.
* Especially useful for plane fitting.
*/
template<typename MatrixType>
void Tridiagonalization<MatrixType>::_decomposeInPlace3x3(MatrixType& mat, DiagonalType& diag, SubDiagonalType& subdiag, bool extractQ)
{
diag[0] = ei_real(mat(0,0));
RealScalar v1norm2 = ei_abs2(mat(0,2));
if (ei_isMuchSmallerThan(v1norm2, RealScalar(1)))
{
diag[1] = ei_real(mat(1,1));
diag[2] = ei_real(mat(2,2));
subdiag[0] = ei_real(mat(0,1));
subdiag[1] = ei_real(mat(1,2));
if (extractQ)
mat.setIdentity();
}
else
{
RealScalar beta = ei_sqrt(ei_abs2(mat(0,1))+v1norm2);
RealScalar invBeta = RealScalar(1)/beta;
Scalar m01 = mat(0,1) * invBeta;
Scalar m02 = mat(0,2) * invBeta;
Scalar q = RealScalar(2)*m01*mat(1,2) + m02*(mat(2,2) - mat(1,1));
diag[1] = ei_real(mat(1,1) + m02*q);
diag[2] = ei_real(mat(2,2) - m02*q);
subdiag[0] = beta;
subdiag[1] = ei_real(mat(1,2) - m01 * q);
if (extractQ)
{
mat(0,0) = 1;
mat(0,1) = 0;
mat(0,2) = 0;
mat(1,0) = 0;
mat(1,1) = m01;
mat(1,2) = m02;
mat(2,0) = 0;
mat(2,1) = m02;
mat(2,2) = -m01;
}
}
}
#endif // EIGEN_HIDE_HEAVY_CODE
#endif // EIGEN_TRIDIAGONALIZATION_H
|