1
// This file is part of Eigen, a lightweight C++ template library
4
// Copyright (C) 2008-2009 Gael Guennebaud <gael.guennebaud@inria.fr>
5
// Copyright (C) 2009 Benoit Jacob <jacob.benoit.1@gmail.com>
7
// Eigen is free software; you can redistribute it and/or
8
// modify it under the terms of the GNU Lesser General Public
9
// License as published by the Free Software Foundation; either
10
// version 3 of the License, or (at your option) any later version.
12
// Alternatively, you can redistribute it and/or
13
// modify it under the terms of the GNU General Public License as
14
// published by the Free Software Foundation; either version 2 of
15
// the License, or (at your option) any later version.
17
// Eigen is distributed in the hope that it will be useful, but WITHOUT ANY
18
// WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
19
// FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License or the
20
// GNU General Public License for more details.
22
// You should have received a copy of the GNU Lesser General Public
23
// License and a copy of the GNU General Public License along with
24
// Eigen. If not, see <http://www.gnu.org/licenses/>.
26
#ifndef EIGEN_FULLPIVOTINGHOUSEHOLDERQR_H
27
#define EIGEN_FULLPIVOTINGHOUSEHOLDERQR_H
29
/** \ingroup QR_Module
31
* \class FullPivHouseholderQR
33
* \brief Householder rank-revealing QR decomposition of a matrix with full pivoting
35
* \param MatrixType the type of the matrix of which we are computing the QR decomposition
37
* This class performs a rank-revealing QR decomposition of a matrix \b A into matrices \b P, \b Q and \b R
40
* \mathbf{A} \, \mathbf{P} = \mathbf{Q} \, \mathbf{R}
42
* by using Householder transformations. Here, \b P is a permutation matrix, \b Q a unitary matrix and \b R an
43
* upper triangular matrix.
45
* This decomposition performs a very prudent full pivoting in order to be rank-revealing and achieve optimal
46
* numerical stability. The trade-off is that it is slower than HouseholderQR and ColPivHouseholderQR.
48
* \sa MatrixBase::fullPivHouseholderQr()
50
template<typename _MatrixType> class FullPivHouseholderQR
54
typedef _MatrixType MatrixType;
56
RowsAtCompileTime = MatrixType::RowsAtCompileTime,
57
ColsAtCompileTime = MatrixType::ColsAtCompileTime,
58
Options = MatrixType::Options,
59
MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
60
MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
62
typedef typename MatrixType::Scalar Scalar;
63
typedef typename MatrixType::RealScalar RealScalar;
64
typedef typename MatrixType::Index Index;
65
typedef Matrix<Scalar, RowsAtCompileTime, RowsAtCompileTime, Options, MaxRowsAtCompileTime, MaxRowsAtCompileTime> MatrixQType;
66
typedef typename internal::plain_diag_type<MatrixType>::type HCoeffsType;
67
typedef Matrix<Index, 1, ColsAtCompileTime, RowMajor, 1, MaxColsAtCompileTime> IntRowVectorType;
68
typedef PermutationMatrix<ColsAtCompileTime, MaxColsAtCompileTime> PermutationType;
69
typedef typename internal::plain_col_type<MatrixType, Index>::type IntColVectorType;
70
typedef typename internal::plain_row_type<MatrixType>::type RowVectorType;
71
typedef typename internal::plain_col_type<MatrixType>::type ColVectorType;
73
/** \brief Default Constructor.
75
* The default constructor is useful in cases in which the user intends to
76
* perform decompositions via FullPivHouseholderQR::compute(const MatrixType&).
78
FullPivHouseholderQR()
81
m_rows_transpositions(),
82
m_cols_transpositions(),
85
m_isInitialized(false),
86
m_usePrescribedThreshold(false) {}
88
/** \brief Default Constructor with memory preallocation
90
* Like the default constructor but with preallocation of the internal data
91
* according to the specified problem \a size.
92
* \sa FullPivHouseholderQR()
94
FullPivHouseholderQR(Index rows, Index cols)
96
m_hCoeffs(std::min(rows,cols)),
97
m_rows_transpositions(rows),
98
m_cols_transpositions(cols),
99
m_cols_permutation(cols),
100
m_temp(std::min(rows,cols)),
101
m_isInitialized(false),
102
m_usePrescribedThreshold(false) {}
104
FullPivHouseholderQR(const MatrixType& matrix)
105
: m_qr(matrix.rows(), matrix.cols()),
106
m_hCoeffs(std::min(matrix.rows(), matrix.cols())),
107
m_rows_transpositions(matrix.rows()),
108
m_cols_transpositions(matrix.cols()),
109
m_cols_permutation(matrix.cols()),
110
m_temp(std::min(matrix.rows(), matrix.cols())),
111
m_isInitialized(false),
112
m_usePrescribedThreshold(false)
117
/** This method finds a solution x to the equation Ax=b, where A is the matrix of which
118
* *this is the QR decomposition, if any exists.
120
* \param b the right-hand-side of the equation to solve.
122
* \returns a solution.
124
* \note The case where b is a matrix is not yet implemented. Also, this
125
* code is space inefficient.
127
* \note_about_checking_solutions
129
* \note_about_arbitrary_choice_of_solution
131
* Example: \include FullPivHouseholderQR_solve.cpp
132
* Output: \verbinclude FullPivHouseholderQR_solve.out
134
template<typename Rhs>
135
inline const internal::solve_retval<FullPivHouseholderQR, Rhs>
136
solve(const MatrixBase<Rhs>& b) const
138
eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized.");
139
return internal::solve_retval<FullPivHouseholderQR, Rhs>(*this, b.derived());
142
MatrixQType matrixQ(void) const;
144
/** \returns a reference to the matrix where the Householder QR decomposition is stored
146
const MatrixType& matrixQR() const
148
eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized.");
152
FullPivHouseholderQR& compute(const MatrixType& matrix);
154
const PermutationType& colsPermutation() const
156
eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized.");
157
return m_cols_permutation;
160
const IntColVectorType& rowsTranspositions() const
162
eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized.");
163
return m_rows_transpositions;
166
/** \returns the absolute value of the determinant of the matrix of which
167
* *this is the QR decomposition. It has only linear complexity
168
* (that is, O(n) where n is the dimension of the square matrix)
169
* as the QR decomposition has already been computed.
171
* \note This is only for square matrices.
173
* \warning a determinant can be very big or small, so for matrices
174
* of large enough dimension, there is a risk of overflow/underflow.
175
* One way to work around that is to use logAbsDeterminant() instead.
177
* \sa logAbsDeterminant(), MatrixBase::determinant()
179
typename MatrixType::RealScalar absDeterminant() const;
181
/** \returns the natural log of the absolute value of the determinant of the matrix of which
182
* *this is the QR decomposition. It has only linear complexity
183
* (that is, O(n) where n is the dimension of the square matrix)
184
* as the QR decomposition has already been computed.
186
* \note This is only for square matrices.
188
* \note This method is useful to work around the risk of overflow/underflow that's inherent
189
* to determinant computation.
191
* \sa absDeterminant(), MatrixBase::determinant()
193
typename MatrixType::RealScalar logAbsDeterminant() const;
195
/** \returns the rank of the matrix of which *this is the QR decomposition.
197
* \note This method has to determine which pivots should be considered nonzero.
198
* For that, it uses the threshold value that you can control by calling
199
* setThreshold(const RealScalar&).
201
inline Index rank() const
203
eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized.");
204
RealScalar premultiplied_threshold = internal::abs(m_maxpivot) * threshold();
206
for(Index i = 0; i < m_nonzero_pivots; ++i)
207
result += (internal::abs(m_qr.coeff(i,i)) > premultiplied_threshold);
211
/** \returns the dimension of the kernel of the matrix of which *this is the QR decomposition.
213
* \note This method has to determine which pivots should be considered nonzero.
214
* For that, it uses the threshold value that you can control by calling
215
* setThreshold(const RealScalar&).
217
inline Index dimensionOfKernel() const
219
eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized.");
220
return cols() - rank();
223
/** \returns true if the matrix of which *this is the QR decomposition represents an injective
224
* linear map, i.e. has trivial kernel; false otherwise.
226
* \note This method has to determine which pivots should be considered nonzero.
227
* For that, it uses the threshold value that you can control by calling
228
* setThreshold(const RealScalar&).
230
inline bool isInjective() const
232
eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized.");
233
return rank() == cols();
236
/** \returns true if the matrix of which *this is the QR decomposition represents a surjective
237
* linear map; false otherwise.
239
* \note This method has to determine which pivots should be considered nonzero.
240
* For that, it uses the threshold value that you can control by calling
241
* setThreshold(const RealScalar&).
243
inline bool isSurjective() const
245
eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized.");
246
return rank() == rows();
249
/** \returns true if the matrix of which *this is the QR decomposition is invertible.
251
* \note This method has to determine which pivots should be considered nonzero.
252
* For that, it uses the threshold value that you can control by calling
253
* setThreshold(const RealScalar&).
255
inline bool isInvertible() const
257
eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized.");
258
return isInjective() && isSurjective();
261
/** \returns the inverse of the matrix of which *this is the QR decomposition.
263
* \note If this matrix is not invertible, the returned matrix has undefined coefficients.
264
* Use isInvertible() to first determine whether this matrix is invertible.
266
internal::solve_retval<FullPivHouseholderQR, typename MatrixType::IdentityReturnType>
269
eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized.");
270
return internal::solve_retval<FullPivHouseholderQR,typename MatrixType::IdentityReturnType>
271
(*this, MatrixType::Identity(m_qr.rows(), m_qr.cols()));
274
inline Index rows() const { return m_qr.rows(); }
275
inline Index cols() const { return m_qr.cols(); }
276
const HCoeffsType& hCoeffs() const { return m_hCoeffs; }
278
/** Allows to prescribe a threshold to be used by certain methods, such as rank(),
279
* who need to determine when pivots are to be considered nonzero. This is not used for the
280
* QR decomposition itself.
282
* When it needs to get the threshold value, Eigen calls threshold(). By default, this
283
* uses a formula to automatically determine a reasonable threshold.
284
* Once you have called the present method setThreshold(const RealScalar&),
285
* your value is used instead.
287
* \param threshold The new value to use as the threshold.
289
* A pivot will be considered nonzero if its absolute value is strictly greater than
290
* \f$ \vert pivot \vert \leqslant threshold \times \vert maxpivot \vert \f$
291
* where maxpivot is the biggest pivot.
293
* If you want to come back to the default behavior, call setThreshold(Default_t)
295
FullPivHouseholderQR& setThreshold(const RealScalar& threshold)
297
m_usePrescribedThreshold = true;
298
m_prescribedThreshold = threshold;
302
/** Allows to come back to the default behavior, letting Eigen use its default formula for
303
* determining the threshold.
305
* You should pass the special object Eigen::Default as parameter here.
306
* \code qr.setThreshold(Eigen::Default); \endcode
308
* See the documentation of setThreshold(const RealScalar&).
310
FullPivHouseholderQR& setThreshold(Default_t)
312
m_usePrescribedThreshold = false;
316
/** Returns the threshold that will be used by certain methods such as rank().
318
* See the documentation of setThreshold(const RealScalar&).
320
RealScalar threshold() const
322
eigen_assert(m_isInitialized || m_usePrescribedThreshold);
323
return m_usePrescribedThreshold ? m_prescribedThreshold
324
// this formula comes from experimenting (see "LU precision tuning" thread on the list)
325
// and turns out to be identical to Higham's formula used already in LDLt.
326
: NumTraits<Scalar>::epsilon() * m_qr.diagonalSize();
329
/** \returns the number of nonzero pivots in the QR decomposition.
330
* Here nonzero is meant in the exact sense, not in a fuzzy sense.
331
* So that notion isn't really intrinsically interesting, but it is
332
* still useful when implementing algorithms.
336
inline Index nonzeroPivots() const
338
eigen_assert(m_isInitialized && "LU is not initialized.");
339
return m_nonzero_pivots;
342
/** \returns the absolute value of the biggest pivot, i.e. the biggest
343
* diagonal coefficient of U.
345
RealScalar maxPivot() const { return m_maxpivot; }
349
HCoeffsType m_hCoeffs;
350
IntColVectorType m_rows_transpositions;
351
IntRowVectorType m_cols_transpositions;
352
PermutationType m_cols_permutation;
353
RowVectorType m_temp;
354
bool m_isInitialized, m_usePrescribedThreshold;
355
RealScalar m_prescribedThreshold, m_maxpivot;
356
Index m_nonzero_pivots;
357
RealScalar m_precision;
361
template<typename MatrixType>
362
typename MatrixType::RealScalar FullPivHouseholderQR<MatrixType>::absDeterminant() const
364
eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized.");
365
eigen_assert(m_qr.rows() == m_qr.cols() && "You can't take the determinant of a non-square matrix!");
366
return internal::abs(m_qr.diagonal().prod());
369
template<typename MatrixType>
370
typename MatrixType::RealScalar FullPivHouseholderQR<MatrixType>::logAbsDeterminant() const
372
eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized.");
373
eigen_assert(m_qr.rows() == m_qr.cols() && "You can't take the determinant of a non-square matrix!");
374
return m_qr.diagonal().cwiseAbs().array().log().sum();
377
template<typename MatrixType>
378
FullPivHouseholderQR<MatrixType>& FullPivHouseholderQR<MatrixType>::compute(const MatrixType& matrix)
380
Index rows = matrix.rows();
381
Index cols = matrix.cols();
382
Index size = std::min(rows,cols);
385
m_hCoeffs.resize(size);
389
m_precision = NumTraits<Scalar>::epsilon() * size;
391
m_rows_transpositions.resize(matrix.rows());
392
m_cols_transpositions.resize(matrix.cols());
393
Index number_of_transpositions = 0;
395
RealScalar biggest(0);
397
m_nonzero_pivots = size; // the generic case is that in which all pivots are nonzero (invertible case)
398
m_maxpivot = RealScalar(0);
400
for (Index k = 0; k < size; ++k)
402
Index row_of_biggest_in_corner, col_of_biggest_in_corner;
403
RealScalar biggest_in_corner;
405
biggest_in_corner = m_qr.bottomRightCorner(rows-k, cols-k)
407
.maxCoeff(&row_of_biggest_in_corner, &col_of_biggest_in_corner);
408
row_of_biggest_in_corner += k;
409
col_of_biggest_in_corner += k;
410
if(k==0) biggest = biggest_in_corner;
412
// if the corner is negligible, then we have less than full rank, and we can finish early
413
if(internal::isMuchSmallerThan(biggest_in_corner, biggest, m_precision))
415
m_nonzero_pivots = k;
416
for(Index i = k; i < size; i++)
418
m_rows_transpositions.coeffRef(i) = i;
419
m_cols_transpositions.coeffRef(i) = i;
420
m_hCoeffs.coeffRef(i) = Scalar(0);
425
m_rows_transpositions.coeffRef(k) = row_of_biggest_in_corner;
426
m_cols_transpositions.coeffRef(k) = col_of_biggest_in_corner;
427
if(k != row_of_biggest_in_corner) {
428
m_qr.row(k).tail(cols-k).swap(m_qr.row(row_of_biggest_in_corner).tail(cols-k));
429
++number_of_transpositions;
431
if(k != col_of_biggest_in_corner) {
432
m_qr.col(k).swap(m_qr.col(col_of_biggest_in_corner));
433
++number_of_transpositions;
437
m_qr.col(k).tail(rows-k).makeHouseholderInPlace(m_hCoeffs.coeffRef(k), beta);
438
m_qr.coeffRef(k,k) = beta;
440
// remember the maximum absolute value of diagonal coefficients
441
if(internal::abs(beta) > m_maxpivot) m_maxpivot = internal::abs(beta);
443
m_qr.bottomRightCorner(rows-k, cols-k-1)
444
.applyHouseholderOnTheLeft(m_qr.col(k).tail(rows-k-1), m_hCoeffs.coeffRef(k), &m_temp.coeffRef(k+1));
447
m_cols_permutation.setIdentity(cols);
448
for(Index k = 0; k < size; ++k)
449
m_cols_permutation.applyTranspositionOnTheRight(k, m_cols_transpositions.coeff(k));
451
m_det_pq = (number_of_transpositions%2) ? -1 : 1;
452
m_isInitialized = true;
459
template<typename _MatrixType, typename Rhs>
460
struct solve_retval<FullPivHouseholderQR<_MatrixType>, Rhs>
461
: solve_retval_base<FullPivHouseholderQR<_MatrixType>, Rhs>
463
EIGEN_MAKE_SOLVE_HELPERS(FullPivHouseholderQR<_MatrixType>,Rhs)
465
template<typename Dest> void evalTo(Dest& dst) const
467
const Index rows = dec().rows(), cols = dec().cols();
468
eigen_assert(rhs().rows() == rows);
470
// FIXME introduce nonzeroPivots() and use it here. and more generally,
471
// make the same improvements in this dec as in FullPivLU.
478
typename Rhs::PlainObject c(rhs());
480
Matrix<Scalar,1,Rhs::ColsAtCompileTime> temp(rhs().cols());
481
for (Index k = 0; k < dec().rank(); ++k)
483
Index remainingSize = rows-k;
484
c.row(k).swap(c.row(dec().rowsTranspositions().coeff(k)));
485
c.bottomRightCorner(remainingSize, rhs().cols())
486
.applyHouseholderOnTheLeft(dec().matrixQR().col(k).tail(remainingSize-1),
487
dec().hCoeffs().coeff(k), &temp.coeffRef(0));
490
if(!dec().isSurjective())
492
// is c is in the image of R ?
493
RealScalar biggest_in_upper_part_of_c = c.topRows( dec().rank() ).cwiseAbs().maxCoeff();
494
RealScalar biggest_in_lower_part_of_c = c.bottomRows(rows-dec().rank()).cwiseAbs().maxCoeff();
496
const RealScalar m_precision = NumTraits<Scalar>::epsilon() * std::min(rows,cols);
497
// this internal:: prefix is needed by at least gcc 3.4 and ICC
498
if(!internal::isMuchSmallerThan(biggest_in_lower_part_of_c, biggest_in_upper_part_of_c, m_precision))
502
.topLeftCorner(dec().rank(), dec().rank())
503
.template triangularView<Upper>()
504
.solveInPlace(c.topRows(dec().rank()));
506
for(Index i = 0; i < dec().rank(); ++i) dst.row(dec().colsPermutation().indices().coeff(i)) = c.row(i);
507
for(Index i = dec().rank(); i < cols; ++i) dst.row(dec().colsPermutation().indices().coeff(i)).setZero();
511
} // end namespace internal
513
/** \returns the matrix Q */
514
template<typename MatrixType>
515
typename FullPivHouseholderQR<MatrixType>::MatrixQType FullPivHouseholderQR<MatrixType>::matrixQ() const
517
eigen_assert(m_isInitialized && "FullPivHouseholderQR is not initialized.");
518
// compute the product H'_0 H'_1 ... H'_n-1,
519
// where H_k is the k-th Householder transformation I - h_k v_k v_k'
520
// and v_k is the k-th Householder vector [1,m_qr(k+1,k), m_qr(k+2,k), ...]
521
Index rows = m_qr.rows();
522
Index cols = m_qr.cols();
523
Index size = std::min(rows,cols);
524
MatrixQType res = MatrixQType::Identity(rows, rows);
525
Matrix<Scalar,1,MatrixType::RowsAtCompileTime> temp(rows);
526
for (Index k = size-1; k >= 0; k--)
528
res.block(k, k, rows-k, rows-k)
529
.applyHouseholderOnTheLeft(m_qr.col(k).tail(rows-k-1), internal::conj(m_hCoeffs.coeff(k)), &temp.coeffRef(k));
530
res.row(k).swap(res.row(m_rows_transpositions.coeff(k)));
535
/** \return the full-pivoting Householder QR decomposition of \c *this.
537
* \sa class FullPivHouseholderQR
539
template<typename Derived>
540
const FullPivHouseholderQR<typename MatrixBase<Derived>::PlainObject>
541
MatrixBase<Derived>::fullPivHouseholderQr() const
543
return FullPivHouseholderQR<PlainObject>(eval());
546
#endif // EIGEN_FULLPIVOTINGHOUSEHOLDERQR_H