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
|
// This file is part of Eigen, a lightweight C++ template library
// for linear algebra.
//
// Copyright (C) 2008-2009 Gael Guennebaud <g.gael@free.fr>
// Copyright (C) 2009 Benoit Jacob <jacob.benoit.1@gmail.com>
//
// 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_QR_H
#define EIGEN_QR_H
/** \ingroup QR_Module
* \nonstableyet
*
* \class HouseholderQR
*
* \brief Householder QR decomposition of a matrix
*
* \param MatrixType the type of the matrix of which we are computing the QR decomposition
*
* This class performs a QR decomposition using Householder transformations. The result is
* stored in a compact way compatible with LAPACK.
*
* Note that no pivoting is performed. This is \b not a rank-revealing decomposition.
* If you want that feature, use FullPivotingHouseholderQR or ColPivotingHouseholderQR instead.
*
* This Householder QR decomposition is faster, but less numerically stable and less feature-full than
* FullPivotingHouseholderQR or ColPivotingHouseholderQR.
*
* \sa MatrixBase::householderQr()
*/
template<typename MatrixType> class HouseholderQR
{
public:
enum {
RowsAtCompileTime = MatrixType::RowsAtCompileTime,
ColsAtCompileTime = MatrixType::ColsAtCompileTime,
Options = MatrixType::Options,
DiagSizeAtCompileTime = EIGEN_ENUM_MIN(ColsAtCompileTime,RowsAtCompileTime)
};
typedef typename MatrixType::Scalar Scalar;
typedef typename MatrixType::RealScalar RealScalar;
typedef Matrix<Scalar, RowsAtCompileTime, RowsAtCompileTime, AutoAlign | (ei_traits<MatrixType>::Flags&RowMajorBit ? RowMajor : ColMajor)> MatrixQType;
typedef Matrix<Scalar, DiagSizeAtCompileTime, 1> HCoeffsType;
typedef Matrix<Scalar, 1, ColsAtCompileTime> RowVectorType;
typedef typename HouseholderSequence<MatrixType,HCoeffsType>::ConjugateReturnType HouseholderSequenceType;
/**
* \brief Default Constructor.
*
* The default constructor is useful in cases in which the user intends to
* perform decompositions via HouseholderQR::compute(const MatrixType&).
*/
HouseholderQR() : m_qr(), m_hCoeffs(), m_isInitialized(false) {}
HouseholderQR(const MatrixType& matrix)
: m_qr(matrix.rows(), matrix.cols()),
m_hCoeffs(std::min(matrix.rows(),matrix.cols())),
m_isInitialized(false)
{
compute(matrix);
}
/** This method finds a solution x to the equation Ax=b, where A is the matrix of which
* *this is the QR decomposition, if any exists.
*
* \param b the right-hand-side of the equation to solve.
*
* \param result a pointer to the vector/matrix in which to store the solution, if any exists.
* Resized if necessary, so that result->rows()==A.cols() and result->cols()==b.cols().
* If no solution exists, *result is left with undefined coefficients.
*
* \note The case where b is a matrix is not yet implemented. Also, this
* code is space inefficient.
*
* Example: \include HouseholderQR_solve.cpp
* Output: \verbinclude HouseholderQR_solve.out
*/
template<typename OtherDerived, typename ResultType>
void solve(const MatrixBase<OtherDerived>& b, ResultType *result) const;
MatrixQType matrixQ() const;
HouseholderSequenceType matrixQAsHouseholderSequence() const
{
return HouseholderSequenceType(m_qr, m_hCoeffs.conjugate());
}
/** \returns a reference to the matrix where the Householder QR decomposition is stored
* in a LAPACK-compatible way.
*/
const MatrixType& matrixQR() const
{
ei_assert(m_isInitialized && "HouseholderQR is not initialized.");
return m_qr;
}
HouseholderQR& compute(const MatrixType& matrix);
/** \returns the absolute value of the determinant of the matrix of which
* *this is the QR decomposition. It has only linear complexity
* (that is, O(n) where n is the dimension of the square matrix)
* as the QR decomposition has already been computed.
*
* \note This is only for square matrices.
*
* \warning a determinant can be very big or small, so for matrices
* of large enough dimension, there is a risk of overflow/underflow.
* One way to work around that is to use logAbsDeterminant() instead.
*
* \sa logAbsDeterminant(), MatrixBase::determinant()
*/
typename MatrixType::RealScalar absDeterminant() const;
/** \returns the natural log of the absolute value of the determinant of the matrix of which
* *this is the QR decomposition. It has only linear complexity
* (that is, O(n) where n is the dimension of the square matrix)
* as the QR decomposition has already been computed.
*
* \note This is only for square matrices.
*
* \note This method is useful to work around the risk of overflow/underflow that's inherent
* to determinant computation.
*
* \sa absDeterminant(), MatrixBase::determinant()
*/
typename MatrixType::RealScalar logAbsDeterminant() const;
protected:
MatrixType m_qr;
HCoeffsType m_hCoeffs;
bool m_isInitialized;
};
#ifndef EIGEN_HIDE_HEAVY_CODE
template<typename MatrixType>
typename MatrixType::RealScalar HouseholderQR<MatrixType>::absDeterminant() const
{
ei_assert(m_isInitialized && "HouseholderQR is not initialized.");
ei_assert(m_qr.rows() == m_qr.cols() && "You can't take the determinant of a non-square matrix!");
return ei_abs(m_qr.diagonal().prod());
}
template<typename MatrixType>
typename MatrixType::RealScalar HouseholderQR<MatrixType>::logAbsDeterminant() const
{
ei_assert(m_isInitialized && "HouseholderQR is not initialized.");
ei_assert(m_qr.rows() == m_qr.cols() && "You can't take the determinant of a non-square matrix!");
return m_qr.diagonal().cwise().abs().cwise().log().sum();
}
template<typename MatrixType>
HouseholderQR<MatrixType>& HouseholderQR<MatrixType>::compute(const MatrixType& matrix)
{
int rows = matrix.rows();
int cols = matrix.cols();
int size = std::min(rows,cols);
m_qr = matrix;
m_hCoeffs.resize(size);
RowVectorType temp(cols);
for (int k = 0; k < size; ++k)
{
int remainingRows = rows - k;
int remainingCols = cols - k - 1;
RealScalar beta;
m_qr.col(k).end(remainingRows).makeHouseholderInPlace(&m_hCoeffs.coeffRef(k), &beta);
m_qr.coeffRef(k,k) = beta;
// apply H to remaining part of m_qr from the left
m_qr.corner(BottomRight, remainingRows, remainingCols)
.applyHouseholderOnTheLeft(m_qr.col(k).end(remainingRows-1), m_hCoeffs.coeffRef(k), &temp.coeffRef(k+1));
}
m_isInitialized = true;
return *this;
}
template<typename MatrixType>
template<typename OtherDerived, typename ResultType>
void HouseholderQR<MatrixType>::solve(
const MatrixBase<OtherDerived>& b,
ResultType *result
) const
{
ei_assert(m_isInitialized && "HouseholderQR is not initialized.");
result->resize(m_qr.cols(), b.cols());
const int rows = m_qr.rows();
const int rank = std::min(m_qr.rows(), m_qr.cols());
ei_assert(b.rows() == rows);
typename OtherDerived::PlainMatrixType c(b);
// Note that the matrix Q = H_0^* H_1^*... so its inverse is Q^* = (H_0 H_1 ...)^T
c.applyOnTheLeft(makeHouseholderSequence(m_qr.corner(TopLeft,rows,rank), m_hCoeffs.start(rank)).transpose());
m_qr.corner(TopLeft, rank, rank)
.template triangularView<UpperTriangular>()
.solveInPlace(c.corner(TopLeft, rank, c.cols()));
result->corner(TopLeft, rank, c.cols()) = c.corner(TopLeft,rank, c.cols());
result->corner(BottomLeft, result->rows()-rank, c.cols()).setZero();
}
/** \returns the matrix Q */
template<typename MatrixType>
typename HouseholderQR<MatrixType>::MatrixQType HouseholderQR<MatrixType>::matrixQ() const
{
ei_assert(m_isInitialized && "HouseholderQR is not initialized.");
return matrixQAsHouseholderSequence();
}
#endif // EIGEN_HIDE_HEAVY_CODE
/** \return the Householder QR decomposition of \c *this.
*
* \sa class HouseholderQR
*/
template<typename Derived>
const HouseholderQR<typename MatrixBase<Derived>::PlainMatrixType>
MatrixBase<Derived>::householderQr() const
{
return HouseholderQR<PlainMatrixType>(eval());
}
#endif // EIGEN_QR_H
|