11 #ifndef EIGEN_FULLPIVOTINGHOUSEHOLDERQR_H
12 #define EIGEN_FULLPIVOTINGHOUSEHOLDERQR_H
18 template<
typename MatrixType>
struct FullPivHouseholderQRMatrixQReturnType;
20 template<
typename MatrixType>
21 struct traits<FullPivHouseholderQRMatrixQReturnType<MatrixType> >
23 typedef typename MatrixType::PlainObject ReturnType;
49 template<
typename _MatrixType>
class FullPivHouseholderQR
53 typedef _MatrixType MatrixType;
55 RowsAtCompileTime = MatrixType::RowsAtCompileTime,
56 ColsAtCompileTime = MatrixType::ColsAtCompileTime,
57 Options = MatrixType::Options,
58 MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
59 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
61 typedef typename MatrixType::Scalar Scalar;
62 typedef typename MatrixType::RealScalar RealScalar;
63 typedef typename MatrixType::Index Index;
64 typedef internal::FullPivHouseholderQRMatrixQReturnType<MatrixType> MatrixQReturnType;
65 typedef typename internal::plain_diag_type<MatrixType>::type HCoeffsType;
66 typedef Matrix<Index, 1,
67 EIGEN_SIZE_MIN_PREFER_DYNAMIC(ColsAtCompileTime,RowsAtCompileTime),
RowMajor, 1,
68 EIGEN_SIZE_MIN_PREFER_FIXED(MaxColsAtCompileTime,MaxRowsAtCompileTime)> IntDiagSizeVectorType;
69 typedef PermutationMatrix<ColsAtCompileTime, MaxColsAtCompileTime> PermutationType;
70 typedef typename internal::plain_row_type<MatrixType>::type RowVectorType;
71 typedef typename internal::plain_col_type<MatrixType>::type ColVectorType;
81 m_rows_transpositions(),
82 m_cols_transpositions(),
85 m_isInitialized(false),
86 m_usePrescribedThreshold(false) {}
96 m_hCoeffs((std::min)(rows,cols)),
97 m_rows_transpositions((std::min)(rows,cols)),
98 m_cols_transpositions((std::min)(rows,cols)),
99 m_cols_permutation(cols),
101 m_isInitialized(false),
102 m_usePrescribedThreshold(false) {}
117 : m_qr(matrix.rows(), matrix.cols()),
118 m_hCoeffs((std::min)(matrix.rows(), matrix.cols())),
119 m_rows_transpositions((std::min)(matrix.rows(), matrix.cols())),
120 m_cols_transpositions((std::min)(matrix.rows(), matrix.cols())),
121 m_cols_permutation(matrix.cols()),
122 m_temp(matrix.cols()),
123 m_isInitialized(false),
124 m_usePrescribedThreshold(false)
147 template<
typename Rhs>
148 inline const internal::solve_retval<FullPivHouseholderQR, Rhs>
151 eigen_assert(m_isInitialized &&
"FullPivHouseholderQR is not initialized.");
152 return internal::solve_retval<FullPivHouseholderQR, Rhs>(*
this, b.derived());
157 MatrixQReturnType
matrixQ(
void)
const;
163 eigen_assert(m_isInitialized &&
"FullPivHouseholderQR is not initialized.");
172 eigen_assert(m_isInitialized &&
"FullPivHouseholderQR is not initialized.");
173 return m_cols_permutation;
179 eigen_assert(m_isInitialized &&
"FullPivHouseholderQR is not initialized.");
180 return m_rows_transpositions;
221 eigen_assert(m_isInitialized &&
"FullPivHouseholderQR is not initialized.");
222 RealScalar premultiplied_threshold = abs(m_maxpivot) *
threshold();
224 for(Index i = 0; i < m_nonzero_pivots; ++i)
225 result += (abs(m_qr.coeff(i,i)) > premultiplied_threshold);
237 eigen_assert(m_isInitialized &&
"FullPivHouseholderQR is not initialized.");
238 return cols() -
rank();
250 eigen_assert(m_isInitialized &&
"FullPivHouseholderQR is not initialized.");
251 return rank() == cols();
263 eigen_assert(m_isInitialized &&
"FullPivHouseholderQR is not initialized.");
264 return rank() == rows();
275 eigen_assert(m_isInitialized &&
"FullPivHouseholderQR is not initialized.");
284 internal::solve_retval<FullPivHouseholderQR, typename MatrixType::IdentityReturnType>
287 eigen_assert(m_isInitialized &&
"FullPivHouseholderQR is not initialized.");
288 return internal::solve_retval<FullPivHouseholderQR,typename MatrixType::IdentityReturnType>
289 (*
this, MatrixType::Identity(m_qr.rows(), m_qr.cols()));
292 inline Index rows()
const {
return m_qr.rows(); }
293 inline Index cols()
const {
return m_qr.cols(); }
299 const HCoeffsType&
hCoeffs()
const {
return m_hCoeffs; }
320 m_usePrescribedThreshold =
true;
335 m_usePrescribedThreshold =
false;
345 eigen_assert(m_isInitialized || m_usePrescribedThreshold);
346 return m_usePrescribedThreshold ? m_prescribedThreshold
361 eigen_assert(m_isInitialized &&
"LU is not initialized.");
362 return m_nonzero_pivots;
372 static void check_template_parameters()
374 EIGEN_STATIC_ASSERT_NON_INTEGER(Scalar);
378 HCoeffsType m_hCoeffs;
379 IntDiagSizeVectorType m_rows_transpositions;
380 IntDiagSizeVectorType m_cols_transpositions;
381 PermutationType m_cols_permutation;
382 RowVectorType m_temp;
383 bool m_isInitialized, m_usePrescribedThreshold;
384 RealScalar m_prescribedThreshold, m_maxpivot;
385 Index m_nonzero_pivots;
386 RealScalar m_precision;
390 template<
typename MatrixType>
394 eigen_assert(m_isInitialized &&
"FullPivHouseholderQR is not initialized.");
395 eigen_assert(m_qr.rows() == m_qr.cols() &&
"You can't take the determinant of a non-square matrix!");
396 return abs(m_qr.diagonal().prod());
399 template<
typename MatrixType>
402 eigen_assert(m_isInitialized &&
"FullPivHouseholderQR is not initialized.");
403 eigen_assert(m_qr.rows() == m_qr.cols() &&
"You can't take the determinant of a non-square matrix!");
404 return m_qr.diagonal().cwiseAbs().array().log().sum();
413 template<
typename MatrixType>
416 check_template_parameters();
419 Index rows = matrix.rows();
420 Index cols = matrix.cols();
421 Index size = (std::min)(rows,cols);
424 m_hCoeffs.resize(size);
430 m_rows_transpositions.resize(size);
431 m_cols_transpositions.resize(size);
432 Index number_of_transpositions = 0;
434 RealScalar biggest(0);
436 m_nonzero_pivots = size;
437 m_maxpivot = RealScalar(0);
439 for (Index k = 0; k < size; ++k)
441 Index row_of_biggest_in_corner, col_of_biggest_in_corner;
442 RealScalar biggest_in_corner;
444 biggest_in_corner = m_qr.bottomRightCorner(rows-k, cols-k)
446 .maxCoeff(&row_of_biggest_in_corner, &col_of_biggest_in_corner);
447 row_of_biggest_in_corner += k;
448 col_of_biggest_in_corner += k;
449 if(k==0) biggest = biggest_in_corner;
452 if(internal::isMuchSmallerThan(biggest_in_corner, biggest, m_precision))
454 m_nonzero_pivots = k;
455 for(Index i = k; i < size; i++)
457 m_rows_transpositions.coeffRef(i) = i;
458 m_cols_transpositions.coeffRef(i) = i;
459 m_hCoeffs.coeffRef(i) = Scalar(0);
464 m_rows_transpositions.coeffRef(k) = row_of_biggest_in_corner;
465 m_cols_transpositions.coeffRef(k) = col_of_biggest_in_corner;
466 if(k != row_of_biggest_in_corner) {
467 m_qr.row(k).tail(cols-k).swap(m_qr.row(row_of_biggest_in_corner).tail(cols-k));
468 ++number_of_transpositions;
470 if(k != col_of_biggest_in_corner) {
471 m_qr.col(k).swap(m_qr.col(col_of_biggest_in_corner));
472 ++number_of_transpositions;
476 m_qr.col(k).tail(rows-k).makeHouseholderInPlace(m_hCoeffs.coeffRef(k), beta);
477 m_qr.coeffRef(k,k) = beta;
480 if(abs(beta) > m_maxpivot) m_maxpivot = abs(beta);
482 m_qr.bottomRightCorner(rows-k, cols-k-1)
483 .applyHouseholderOnTheLeft(m_qr.col(k).tail(rows-k-1), m_hCoeffs.coeffRef(k), &m_temp.coeffRef(k+1));
486 m_cols_permutation.setIdentity(cols);
487 for(Index k = 0; k < size; ++k)
488 m_cols_permutation.applyTranspositionOnTheRight(k, m_cols_transpositions.coeff(k));
490 m_det_pq = (number_of_transpositions%2) ? -1 : 1;
491 m_isInitialized =
true;
498 template<
typename _MatrixType,
typename Rhs>
500 : solve_retval_base<FullPivHouseholderQR<_MatrixType>, Rhs>
504 template<typename Dest>
void evalTo(Dest& dst)
const
506 const Index rows = dec().rows(), cols = dec().cols();
507 eigen_assert(rhs().rows() == rows);
517 typename Rhs::PlainObject c(rhs());
519 Matrix<Scalar,1,Rhs::ColsAtCompileTime> temp(rhs().cols());
520 for (Index k = 0; k < dec().rank(); ++k)
522 Index remainingSize = rows-k;
523 c.row(k).swap(c.row(dec().rowsTranspositions().coeff(k)));
524 c.bottomRightCorner(remainingSize, rhs().cols())
525 .applyHouseholderOnTheLeft(dec().matrixQR().col(k).tail(remainingSize-1),
526 dec().hCoeffs().coeff(k), &temp.coeffRef(0));
530 .topLeftCorner(dec().rank(), dec().rank())
531 .template triangularView<Upper>()
532 .solveInPlace(c.topRows(dec().rank()));
534 for(Index i = 0; i < dec().rank(); ++i) dst.row(dec().colsPermutation().indices().coeff(i)) = c.row(i);
535 for(Index i = dec().rank(); i < cols; ++i) dst.row(dec().colsPermutation().indices().coeff(i)).setZero();
545 template<
typename MatrixType>
struct FullPivHouseholderQRMatrixQReturnType
546 :
public ReturnByValue<FullPivHouseholderQRMatrixQReturnType<MatrixType> >
549 typedef typename MatrixType::Index Index;
550 typedef typename FullPivHouseholderQR<MatrixType>::IntDiagSizeVectorType IntDiagSizeVectorType;
551 typedef typename internal::plain_diag_type<MatrixType>::type HCoeffsType;
552 typedef Matrix<
typename MatrixType::Scalar, 1, MatrixType::RowsAtCompileTime,
RowMajor, 1,
553 MatrixType::MaxRowsAtCompileTime> WorkVectorType;
555 FullPivHouseholderQRMatrixQReturnType(
const MatrixType& qr,
556 const HCoeffsType& hCoeffs,
557 const IntDiagSizeVectorType& rowsTranspositions)
560 m_rowsTranspositions(rowsTranspositions)
563 template <
typename ResultType>
564 void evalTo(ResultType& result)
const
566 const Index rows = m_qr.rows();
567 WorkVectorType workspace(rows);
568 evalTo(result, workspace);
571 template <
typename ResultType>
572 void evalTo(ResultType& result, WorkVectorType& workspace)
const
578 const Index rows = m_qr.rows();
579 const Index cols = m_qr.cols();
580 const Index size = (std::min)(rows, cols);
581 workspace.resize(rows);
582 result.setIdentity(rows, rows);
583 for (Index k = size-1; k >= 0; k--)
585 result.block(k, k, rows-k, rows-k)
586 .applyHouseholderOnTheLeft(m_qr.col(k).tail(rows-k-1), conj(m_hCoeffs.coeff(k)), &workspace.coeffRef(k));
587 result.row(k).swap(result.row(m_rowsTranspositions.coeff(k)));
591 Index rows()
const {
return m_qr.rows(); }
592 Index cols()
const {
return m_qr.rows(); }
595 typename MatrixType::Nested m_qr;
596 typename HCoeffsType::Nested m_hCoeffs;
597 typename IntDiagSizeVectorType::Nested m_rowsTranspositions;
602 template<
typename MatrixType>
605 eigen_assert(m_isInitialized &&
"FullPivHouseholderQR is not initialized.");
606 return MatrixQReturnType(m_qr, m_hCoeffs, m_rows_transpositions);
613 template<
typename Derived>
622 #endif // EIGEN_FULLPIVOTINGHOUSEHOLDERQR_H
MatrixType::RealScalar absDeterminant() const
Definition: FullPivHouseholderQR.h:391
const internal::solve_retval< FullPivHouseholderQR, Rhs > solve(const MatrixBase< Rhs > &b) const
Definition: FullPivHouseholderQR.h:149
bool isInjective() const
Definition: FullPivHouseholderQR.h:248
Householder rank-revealing QR decomposition of a matrix with full pivoting.
Definition: ForwardDeclarations.h:223
const PermutationType & colsPermutation() const
Definition: FullPivHouseholderQR.h:170
Index dimensionOfKernel() const
Definition: FullPivHouseholderQR.h:235
Holds information about the various numeric (i.e. scalar) types allowed by Eigen. ...
Definition: NumTraits.h:88
const HCoeffsType & hCoeffs() const
Definition: FullPivHouseholderQR.h:299
FullPivHouseholderQR(Index rows, Index cols)
Default Constructor with memory preallocation.
Definition: FullPivHouseholderQR.h:94
bool isInvertible() const
Definition: FullPivHouseholderQR.h:273
bool isSurjective() const
Definition: FullPivHouseholderQR.h:261
FullPivHouseholderQR & compute(const MatrixType &matrix)
Definition: FullPivHouseholderQR.h:414
RealScalar maxPivot() const
Definition: FullPivHouseholderQR.h:368
const FullPivHouseholderQR< PlainObject > fullPivHouseholderQr() const
Definition: FullPivHouseholderQR.h:615
FullPivHouseholderQR()
Default Constructor.
Definition: FullPivHouseholderQR.h:78
Index rank() const
Definition: FullPivHouseholderQR.h:218
RealScalar threshold() const
Definition: FullPivHouseholderQR.h:343
const MatrixType & matrixQR() const
Definition: FullPivHouseholderQR.h:161
FullPivHouseholderQR & setThreshold(Default_t)
Definition: FullPivHouseholderQR.h:333
FullPivHouseholderQR(const MatrixType &matrix)
Constructs a QR factorization from a given matrix.
Definition: FullPivHouseholderQR.h:116
Index nonzeroPivots() const
Definition: FullPivHouseholderQR.h:359
const internal::solve_retval< FullPivHouseholderQR, typename MatrixType::IdentityReturnType > inverse() const
Definition: FullPivHouseholderQR.h:285
MatrixType::RealScalar logAbsDeterminant() const
Definition: FullPivHouseholderQR.h:400
Definition: Constants.h:266
const IntDiagSizeVectorType & rowsTranspositions() const
Definition: FullPivHouseholderQR.h:177
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:48
MatrixQReturnType matrixQ(void) const
Definition: FullPivHouseholderQR.h:603
FullPivHouseholderQR & setThreshold(const RealScalar &threshold)
Definition: FullPivHouseholderQR.h:318