11 #ifndef EIGEN_COLPIVOTINGHOUSEHOLDERQR_H
12 #define EIGEN_COLPIVOTINGHOUSEHOLDERQR_H
37 template<
typename _MatrixType>
class ColPivHouseholderQR
41 typedef _MatrixType MatrixType;
43 RowsAtCompileTime = MatrixType::RowsAtCompileTime,
44 ColsAtCompileTime = MatrixType::ColsAtCompileTime,
45 Options = MatrixType::Options,
46 MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
47 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
49 typedef typename MatrixType::Scalar Scalar;
50 typedef typename MatrixType::RealScalar RealScalar;
51 typedef typename MatrixType::Index Index;
52 typedef Matrix<Scalar, RowsAtCompileTime, RowsAtCompileTime, Options, MaxRowsAtCompileTime, MaxRowsAtCompileTime> MatrixQType;
53 typedef typename internal::plain_diag_type<MatrixType>::type HCoeffsType;
54 typedef PermutationMatrix<ColsAtCompileTime, MaxColsAtCompileTime> PermutationType;
55 typedef typename internal::plain_row_type<MatrixType, Index>::type IntRowVectorType;
56 typedef typename internal::plain_row_type<MatrixType>::type RowVectorType;
57 typedef typename internal::plain_row_type<MatrixType, RealScalar>::type RealRowVectorType;
58 typedef HouseholderSequence<MatrixType,typename internal::remove_all<typename HCoeffsType::ConjugateReturnType>::type> HouseholderSequenceType;
62 typedef typename PermutationType::Index PermIndexType;
76 m_colsTranspositions(),
79 m_isInitialized(false),
80 m_usePrescribedThreshold(false) {}
90 m_hCoeffs((std::min)(rows,cols)),
91 m_colsPermutation(PermIndexType(cols)),
92 m_colsTranspositions(cols),
95 m_isInitialized(false),
96 m_usePrescribedThreshold(false) {}
111 : m_qr(matrix.rows(), matrix.cols()),
112 m_hCoeffs((std::min)(matrix.rows(),matrix.cols())),
113 m_colsPermutation(PermIndexType(matrix.cols())),
114 m_colsTranspositions(matrix.cols()),
115 m_temp(matrix.cols()),
116 m_colSqNorms(matrix.cols()),
117 m_isInitialized(false),
118 m_usePrescribedThreshold(false)
140 template<
typename Rhs>
141 inline const internal::solve_retval<ColPivHouseholderQR, Rhs>
144 eigen_assert(m_isInitialized &&
"ColPivHouseholderQR is not initialized.");
145 return internal::solve_retval<ColPivHouseholderQR, Rhs>(*
this, b.derived());
149 HouseholderSequenceType matrixQ(
void)
const
158 eigen_assert(m_isInitialized &&
"ColPivHouseholderQR is not initialized.");
173 eigen_assert(m_isInitialized &&
"ColPivHouseholderQR is not initialized.");
182 eigen_assert(m_isInitialized &&
"ColPivHouseholderQR is not initialized.");
183 return m_colsPermutation;
224 eigen_assert(m_isInitialized &&
"ColPivHouseholderQR is not initialized.");
225 RealScalar premultiplied_threshold = abs(m_maxpivot) *
threshold();
227 for(Index i = 0; i < m_nonzero_pivots; ++i)
228 result += (abs(m_qr.coeff(i,i)) > premultiplied_threshold);
240 eigen_assert(m_isInitialized &&
"ColPivHouseholderQR is not initialized.");
241 return cols() -
rank();
253 eigen_assert(m_isInitialized &&
"ColPivHouseholderQR is not initialized.");
254 return rank() == cols();
266 eigen_assert(m_isInitialized &&
"ColPivHouseholderQR is not initialized.");
267 return rank() == rows();
278 eigen_assert(m_isInitialized &&
"ColPivHouseholderQR is not initialized.");
288 internal::solve_retval<ColPivHouseholderQR, typename MatrixType::IdentityReturnType>
291 eigen_assert(m_isInitialized &&
"ColPivHouseholderQR is not initialized.");
292 return internal::solve_retval<ColPivHouseholderQR,typename MatrixType::IdentityReturnType>
293 (*
this, MatrixType::Identity(m_qr.rows(), m_qr.cols()));
296 inline Index rows()
const {
return m_qr.rows(); }
297 inline Index cols()
const {
return m_qr.cols(); }
303 const HCoeffsType&
hCoeffs()
const {
return m_hCoeffs; }
324 m_usePrescribedThreshold =
true;
339 m_usePrescribedThreshold =
false;
349 eigen_assert(m_isInitialized || m_usePrescribedThreshold);
350 return m_usePrescribedThreshold ? m_prescribedThreshold
365 eigen_assert(m_isInitialized &&
"ColPivHouseholderQR is not initialized.");
366 return m_nonzero_pivots;
382 eigen_assert(m_isInitialized &&
"Decomposition is not initialized.");
388 static void check_template_parameters()
390 EIGEN_STATIC_ASSERT_NON_INTEGER(Scalar);
394 HCoeffsType m_hCoeffs;
395 PermutationType m_colsPermutation;
396 IntRowVectorType m_colsTranspositions;
397 RowVectorType m_temp;
398 RealRowVectorType m_colSqNorms;
399 bool m_isInitialized, m_usePrescribedThreshold;
400 RealScalar m_prescribedThreshold, m_maxpivot;
401 Index m_nonzero_pivots;
405 template<
typename MatrixType>
409 eigen_assert(m_isInitialized &&
"ColPivHouseholderQR is not initialized.");
410 eigen_assert(m_qr.rows() == m_qr.cols() &&
"You can't take the determinant of a non-square matrix!");
411 return abs(m_qr.diagonal().prod());
414 template<
typename MatrixType>
417 eigen_assert(m_isInitialized &&
"ColPivHouseholderQR is not initialized.");
418 eigen_assert(m_qr.rows() == m_qr.cols() &&
"You can't take the determinant of a non-square matrix!");
419 return m_qr.diagonal().cwiseAbs().array().log().sum();
428 template<
typename MatrixType>
431 check_template_parameters();
434 Index rows = matrix.rows();
435 Index cols = matrix.cols();
436 Index size = matrix.diagonalSize();
442 m_hCoeffs.resize(size);
446 m_colsTranspositions.resize(matrix.cols());
447 Index number_of_transpositions = 0;
449 m_colSqNorms.resize(cols);
450 for(Index k = 0; k < cols; ++k)
451 m_colSqNorms.coeffRef(k) = m_qr.col(k).squaredNorm();
455 m_nonzero_pivots = size;
456 m_maxpivot = RealScalar(0);
458 for(Index k = 0; k < size; ++k)
461 Index biggest_col_index;
462 RealScalar biggest_col_sq_norm = m_colSqNorms.tail(cols-k).maxCoeff(&biggest_col_index);
463 biggest_col_index += k;
469 biggest_col_sq_norm = m_qr.col(biggest_col_index).tail(rows-k).squaredNorm();
472 m_colSqNorms.coeffRef(biggest_col_index) = biggest_col_sq_norm;
476 if(m_nonzero_pivots==size && biggest_col_sq_norm < threshold_helper * RealScalar(rows-k))
477 m_nonzero_pivots = k;
480 m_colsTranspositions.coeffRef(k) = biggest_col_index;
481 if(k != biggest_col_index) {
482 m_qr.col(k).swap(m_qr.col(biggest_col_index));
483 std::swap(m_colSqNorms.coeffRef(k), m_colSqNorms.coeffRef(biggest_col_index));
484 ++number_of_transpositions;
489 m_qr.col(k).tail(rows-k).makeHouseholderInPlace(m_hCoeffs.coeffRef(k), beta);
492 m_qr.coeffRef(k,k) = beta;
495 if(abs(beta) > m_maxpivot) m_maxpivot = abs(beta);
498 m_qr.bottomRightCorner(rows-k, cols-k-1)
499 .applyHouseholderOnTheLeft(m_qr.col(k).tail(rows-k-1), m_hCoeffs.coeffRef(k), &m_temp.coeffRef(k+1));
502 m_colSqNorms.tail(cols-k-1) -= m_qr.row(k).tail(cols-k-1).cwiseAbs2();
505 m_colsPermutation.setIdentity(PermIndexType(cols));
506 for(PermIndexType k = 0; k < size; ++k)
507 m_colsPermutation.applyTranspositionOnTheRight(k, PermIndexType(m_colsTranspositions.coeff(k)));
509 m_det_pq = (number_of_transpositions%2) ? -1 : 1;
510 m_isInitialized =
true;
517 template<
typename _MatrixType,
typename Rhs>
519 : solve_retval_base<ColPivHouseholderQR<_MatrixType>, Rhs>
523 template<typename Dest>
void evalTo(Dest& dst)
const
525 eigen_assert(rhs().rows() == dec().rows());
527 const Index cols = dec().cols(),
528 nonzero_pivots = dec().nonzeroPivots();
530 if(nonzero_pivots == 0)
536 typename Rhs::PlainObject c(rhs());
540 .setLength(dec().nonzeroPivots())
545 .topLeftCorner(nonzero_pivots, nonzero_pivots)
546 .template triangularView<Upper>()
547 .solveInPlace(c.topRows(nonzero_pivots));
549 for(Index i = 0; i < nonzero_pivots; ++i) dst.row(dec().colsPermutation().indices().coeff(i)) = c.row(i);
550 for(Index i = nonzero_pivots; i < cols; ++i) dst.row(dec().colsPermutation().indices().coeff(i)).setZero();
559 template<
typename MatrixType>
563 eigen_assert(m_isInitialized &&
"ColPivHouseholderQR is not initialized.");
571 template<
typename Derived>
580 #endif // EIGEN_COLPIVOTINGHOUSEHOLDERQR_H
MatrixType::RealScalar absDeterminant() const
Definition: ColPivHouseholderQR.h:406
MatrixType::RealScalar logAbsDeterminant() const
Definition: ColPivHouseholderQR.h:415
RealScalar maxPivot() const
Definition: ColPivHouseholderQR.h:372
RealScalar threshold() const
Definition: ColPivHouseholderQR.h:347
Holds information about the various numeric (i.e. scalar) types allowed by Eigen. ...
Definition: NumTraits.h:88
const HCoeffsType & hCoeffs() const
Definition: ColPivHouseholderQR.h:303
const MatrixType & matrixQR() const
Definition: ColPivHouseholderQR.h:156
HouseholderSequence< VectorsType, CoeffsType > householderSequence(const VectorsType &v, const CoeffsType &h)
Convenience function for constructing a Householder sequence.
Definition: HouseholderSequence.h:422
ColPivHouseholderQR()
Default Constructor.
Definition: ColPivHouseholderQR.h:72
ComputationInfo info() const
Reports whether the QR factorization was succesful.
Definition: ColPivHouseholderQR.h:380
HouseholderSequenceType householderQ(void) const
Definition: ColPivHouseholderQR.h:561
ColPivHouseholderQR(const MatrixType &matrix)
Constructs a QR factorization from a given matrix.
Definition: ColPivHouseholderQR.h:110
Sequence of Householder reflections acting on subspaces with decreasing size.
Definition: ForwardDeclarations.h:227
ColPivHouseholderQR & compute(const MatrixType &matrix)
Definition: ColPivHouseholderQR.h:429
const PermutationType & colsPermutation() const
Definition: ColPivHouseholderQR.h:180
bool isSurjective() const
Definition: ColPivHouseholderQR.h:264
Householder rank-revealing QR decomposition of a matrix with column-pivoting.
Definition: ForwardDeclarations.h:222
bool isInjective() const
Definition: ColPivHouseholderQR.h:251
bool isInvertible() const
Definition: ColPivHouseholderQR.h:276
ColPivHouseholderQR & setThreshold(Default_t)
Definition: ColPivHouseholderQR.h:337
ColPivHouseholderQR & setThreshold(const RealScalar &threshold)
Definition: ColPivHouseholderQR.h:322
const internal::solve_retval< ColPivHouseholderQR, Rhs > solve(const MatrixBase< Rhs > &b) const
Definition: ColPivHouseholderQR.h:142
const MatrixType & matrixR() const
Definition: ColPivHouseholderQR.h:171
Index nonzeroPivots() const
Definition: ColPivHouseholderQR.h:363
Definition: Constants.h:376
const internal::solve_retval< ColPivHouseholderQR, typename MatrixType::IdentityReturnType > inverse() const
Definition: ColPivHouseholderQR.h:289
Index dimensionOfKernel() const
Definition: ColPivHouseholderQR.h:238
ColPivHouseholderQR(Index rows, Index cols)
Default Constructor with memory preallocation.
Definition: ColPivHouseholderQR.h:88
ComputationInfo
Definition: Constants.h:374
Index rank() const
Definition: ColPivHouseholderQR.h:221
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:48
const ColPivHouseholderQR< PlainObject > colPivHouseholderQr() const
Definition: ColPivHouseholderQR.h:573