11 #ifndef EIGEN_FULLPIVOTINGHOUSEHOLDERQR_H 12 #define EIGEN_FULLPIVOTINGHOUSEHOLDERQR_H 26 template<
typename MatrixType>
63 RowsAtCompileTime = MatrixType::RowsAtCompileTime,
64 ColsAtCompileTime = MatrixType::ColsAtCompileTime,
65 MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
66 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
74 typedef Matrix<StorageIndex, 1,
90 m_rows_transpositions(),
91 m_cols_transpositions(),
94 m_isInitialized(false),
95 m_usePrescribedThreshold(false) {}
105 m_hCoeffs((
std::
min)(rows,cols)),
106 m_rows_transpositions((
std::
min)(rows,cols)),
107 m_cols_transpositions((
std::
min)(rows,cols)),
108 m_cols_permutation(cols),
110 m_isInitialized(false),
111 m_usePrescribedThreshold(false) {}
125 template<
typename InputType>
127 : m_qr(matrix.
rows(), matrix.
cols()),
131 m_cols_permutation(matrix.
cols()),
132 m_temp(matrix.
cols()),
133 m_isInitialized(false),
134 m_usePrescribedThreshold(false)
145 template<
typename InputType>
147 : m_qr(matrix.derived()),
151 m_cols_permutation(matrix.
cols()),
152 m_temp(matrix.
cols()),
153 m_isInitialized(false),
154 m_usePrescribedThreshold(false)
174 template<
typename Rhs>
178 eigen_assert(m_isInitialized &&
"FullPivHouseholderQR is not initialized.");
184 MatrixQReturnType matrixQ(
void)
const;
190 eigen_assert(m_isInitialized &&
"FullPivHouseholderQR is not initialized.");
194 template<
typename InputType>
200 eigen_assert(m_isInitialized &&
"FullPivHouseholderQR is not initialized.");
201 return m_cols_permutation;
207 eigen_assert(m_isInitialized &&
"FullPivHouseholderQR is not initialized.");
208 return m_rows_transpositions;
249 eigen_assert(m_isInitialized &&
"FullPivHouseholderQR is not initialized.");
250 RealScalar premultiplied_threshold =
abs(m_maxpivot) * threshold();
252 for(
Index i = 0;
i < m_nonzero_pivots; ++
i)
253 result += (
abs(m_qr.coeff(
i,
i)) > premultiplied_threshold);
265 eigen_assert(m_isInitialized &&
"FullPivHouseholderQR is not initialized.");
266 return cols() - rank();
278 eigen_assert(m_isInitialized &&
"FullPivHouseholderQR is not initialized.");
279 return rank() ==
cols();
291 eigen_assert(m_isInitialized &&
"FullPivHouseholderQR is not initialized.");
292 return rank() ==
rows();
303 eigen_assert(m_isInitialized &&
"FullPivHouseholderQR is not initialized.");
304 return isInjective() && isSurjective();
314 eigen_assert(m_isInitialized &&
"FullPivHouseholderQR is not initialized.");
325 const HCoeffsType&
hCoeffs()
const {
return m_hCoeffs; }
346 m_usePrescribedThreshold =
true;
347 m_prescribedThreshold = threshold;
361 m_usePrescribedThreshold =
false;
371 eigen_assert(m_isInitialized || m_usePrescribedThreshold);
372 return m_usePrescribedThreshold ? m_prescribedThreshold
387 eigen_assert(m_isInitialized &&
"LU is not initialized.");
388 return m_nonzero_pivots;
396 #ifndef EIGEN_PARSED_BY_DOXYGEN 397 template<
typename RhsType,
typename DstType>
399 void _solve_impl(
const RhsType &rhs, DstType &dst)
const;
409 void computeInPlace();
424 template<
typename MatrixType>
428 eigen_assert(m_isInitialized &&
"FullPivHouseholderQR is not initialized.");
429 eigen_assert(m_qr.rows() == m_qr.cols() &&
"You can't take the determinant of a non-square matrix!");
430 return abs(m_qr.diagonal().prod());
433 template<
typename MatrixType>
436 eigen_assert(m_isInitialized &&
"FullPivHouseholderQR is not initialized.");
437 eigen_assert(m_qr.rows() == m_qr.cols() &&
"You can't take the determinant of a non-square matrix!");
438 return m_qr.diagonal().cwiseAbs().array().log().sum();
447 template<
typename MatrixType>
448 template<
typename InputType>
456 template<
typename MatrixType>
459 check_template_parameters();
467 m_hCoeffs.resize(size);
473 m_rows_transpositions.resize(size);
474 m_cols_transpositions.resize(size);
475 Index number_of_transpositions = 0;
479 m_nonzero_pivots =
size;
484 Index row_of_biggest_in_corner, col_of_biggest_in_corner;
486 typedef typename Scoring::result_type Score;
488 Score score = m_qr.bottomRightCorner(rows-k, cols-k)
489 .unaryExpr(Scoring())
490 .maxCoeff(&row_of_biggest_in_corner, &col_of_biggest_in_corner);
491 row_of_biggest_in_corner += k;
492 col_of_biggest_in_corner += k;
494 if(k==0) biggest = biggest_in_corner;
499 m_nonzero_pivots = k;
502 m_rows_transpositions.coeffRef(
i) =
i;
503 m_cols_transpositions.coeffRef(
i) =
i;
504 m_hCoeffs.coeffRef(
i) =
Scalar(0);
509 m_rows_transpositions.coeffRef(k) = row_of_biggest_in_corner;
510 m_cols_transpositions.coeffRef(k) = col_of_biggest_in_corner;
511 if(k != row_of_biggest_in_corner) {
512 m_qr.row(k).tail(cols-k).swap(m_qr.row(row_of_biggest_in_corner).tail(cols-k));
513 ++number_of_transpositions;
515 if(k != col_of_biggest_in_corner) {
516 m_qr.col(k).swap(m_qr.col(col_of_biggest_in_corner));
517 ++number_of_transpositions;
521 m_qr.col(k).tail(rows-k).makeHouseholderInPlace(m_hCoeffs.coeffRef(k), beta);
522 m_qr.coeffRef(k,k) = beta;
525 if(
abs(beta) > m_maxpivot) m_maxpivot =
abs(beta);
527 m_qr.bottomRightCorner(rows-k, cols-k-1)
528 .applyHouseholderOnTheLeft(m_qr.col(k).tail(rows-k-1), m_hCoeffs.coeffRef(k), &m_temp.coeffRef(k+1));
531 m_cols_permutation.setIdentity(cols);
533 m_cols_permutation.applyTranspositionOnTheRight(k, m_cols_transpositions.coeff(k));
535 m_det_pq = (number_of_transpositions%2) ? -1 : 1;
536 m_isInitialized =
true;
539 #ifndef EIGEN_PARSED_BY_DOXYGEN 540 template<
typename _MatrixType>
541 template<
typename RhsType,
typename DstType>
545 const Index l_rank = rank();
555 typename RhsType::PlainObject
c(rhs);
558 for (
Index k = 0; k < l_rank; ++k)
561 c.row(k).swap(c.row(m_rows_transpositions.coeff(k)));
562 c.bottomRightCorner(remainingSize, rhs.cols())
563 .applyHouseholderOnTheLeft(m_qr.col(k).tail(remainingSize-1),
564 m_hCoeffs.coeff(k), &temp.coeffRef(0));
567 m_qr.topLeftCorner(l_rank, l_rank)
568 .template triangularView<Upper>()
569 .solveInPlace(c.topRows(l_rank));
571 for(
Index i = 0;
i < l_rank; ++
i) dst.row(m_cols_permutation.indices().coeff(
i)) = c.row(
i);
572 for(
Index i = l_rank;
i <
cols(); ++
i) dst.row(m_cols_permutation.indices().coeff(
i)).
setZero();
578 template<
typename DstXprType,
typename MatrixType>
596 :
public ReturnByValue<FullPivHouseholderQRMatrixQReturnType<MatrixType> >
605 const HCoeffsType& hCoeffs,
606 const IntDiagSizeVectorType& rowsTranspositions)
609 m_rowsTranspositions(rowsTranspositions)
612 template <
typename ResultType>
616 WorkVectorType workspace(rows);
617 evalTo(result, workspace);
620 template <
typename ResultType>
631 result.setIdentity(rows, rows);
632 for (
Index k = size-1; k >= 0; k--)
634 result.block(k, k, rows-k, rows-k)
635 .applyHouseholderOnTheLeft(m_qr.col(k).tail(rows-k-1),
conj(m_hCoeffs.coeff(k)), &workspace.
coeffRef(k));
636 result.row(k).swap(result.row(m_rowsTranspositions.coeff(k)));
644 typename MatrixType::Nested
m_qr;
656 template<
typename MatrixType>
659 eigen_assert(m_isInitialized &&
"FullPivHouseholderQR is not initialized.");
667 template<
typename Derived>
676 #endif // EIGEN_FULLPIVOTINGHOUSEHOLDERQR_H Index dimensionOfKernel() const
bool isInvertible() const
FullPivHouseholderQR & setThreshold(Default_t)
EIGEN_DEVICE_FUNC bool isMuchSmallerThan(const Scalar &x, const OtherScalar &y, const typename NumTraits< Scalar >::Real &precision=NumTraits< Scalar >::dummy_precision())
const AutoDiffScalar< DerType > & conj(const AutoDiffScalar< DerType > &x)
EIGEN_DEVICE_FUNC void _solve_impl(const RhsType &rhs, DstType &dst) const
EIGEN_DEVICE_FUNC Index rows() const
Householder rank-revealing QR decomposition of a matrix with full pivoting.
MatrixType::StorageIndex StorageIndex
void evalTo(ResultType &result) const
static void check_template_parameters()
MatrixQReturnType matrixQ(void) const
bool m_usePrescribedThreshold
RealScalar maxPivot() const
Namespace containing all symbols from the Eigen library.
const IntDiagSizeVectorType & rowsTranspositions() const
bool isSurjective() const
MatrixType::RealScalar RealScalar
FullPivHouseholderQRMatrixQReturnType(const MatrixType &qr, const HCoeffsType &hCoeffs, const IntDiagSizeVectorType &rowsTranspositions)
#define EIGEN_SIZE_MIN_PREFER_FIXED(a, b)
void evalTo(ResultType &result, WorkVectorType &workspace) const
HouseholderQR< MatrixXf > qr(A)
FullPivHouseholderQR & compute(const EigenBase< InputType > &matrix)
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Scalar & coeffRef(Index rowId, Index colId)
FullPivHouseholderQR(Index rows, Index cols)
Default Constructor with memory preallocation.
internal::plain_diag_type< MatrixType >::type HCoeffsType
Expression of the inverse of another expression.
const Solve< FullPivHouseholderQR, Rhs > solve(const MatrixBase< Rhs > &b) const
Matrix< typename MatrixType::Scalar, 1, MatrixType::RowsAtCompileTime, RowMajor, 1, MatrixType::MaxRowsAtCompileTime > WorkVectorType
const HCoeffsType & hCoeffs() const
HCoeffsType::Nested m_hCoeffs
FullPivHouseholderQR()
Default Constructor.
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE void resize(Index rows, Index cols)
internal::plain_row_type< MatrixType >::type RowVectorType
const PermutationType & colsPermutation() const
PermutationType m_cols_permutation
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
IntDiagSizeVectorType m_rows_transpositions
EIGEN_DEVICE_FUNC Index cols() const
#define EIGEN_STATIC_ASSERT_NON_INTEGER(TYPE)
internal::FullPivHouseholderQRMatrixQReturnType< MatrixType > MatrixQReturnType
Expression type for return value of FullPivHouseholderQR::matrixQ()
FullPivHouseholderQR(EigenBase< InputType > &matrix)
Constructs a QR factorization from a given matrix.
internal::plain_col_type< MatrixType >::type ColVectorType
FullPivHouseholderQR & setThreshold(const RealScalar &threshold)
NumTraits< Scalar >::Real RealScalar
MatrixType::Scalar Scalar
Inverse< QrType > SrcXprType
IntDiagSizeVectorType::Nested m_rowsTranspositions
PermutationMatrix< ColsAtCompileTime, MaxColsAtCompileTime > PermutationType
FullPivHouseholderQR(const EigenBase< InputType > &matrix)
Constructs a QR factorization from a given matrix.
MatrixType::RealScalar logAbsDeterminant() const
RealScalar threshold() const
FullPivHouseholderQR< MatrixType >::IntDiagSizeVectorType IntDiagSizeVectorType
FullPivHouseholderQR< MatrixType > QrType
Matrix< StorageIndex, 1, EIGEN_SIZE_MIN_PREFER_DYNAMIC(ColsAtCompileTime, RowsAtCompileTime), RowMajor, 1, EIGEN_SIZE_MIN_PREFER_FIXED(MaxColsAtCompileTime, MaxRowsAtCompileTime)> IntDiagSizeVectorType
MatrixType::RealScalar absDeterminant() const
const MatrixType & matrixQR() const
MatrixType::PlainObject PlainObject
IntDiagSizeVectorType m_cols_transpositions
RealScalar m_prescribedThreshold
internal::nested_eval< T, 1 >::type eval(const T &xpr)
Pseudo expression representing a solving operation.
MatrixType::PlainObject ReturnType
EIGEN_DONT_INLINE void compute(Solver &solver, const MatrixType &A)
Map< Matrix< T, Dynamic, Dynamic, ColMajor >, 0, OuterStride<> > matrix(T *data, int rows, int cols, int stride)
The matrix class, also used for vectors and row-vectors.
static void run(DstXprType &dst, const SrcXprType &src, const internal::assign_op< typename DstXprType::Scalar, typename QrType::Scalar > &)
EIGEN_DEVICE_FUNC const XprTypeNestedCleaned & nestedExpression() const
#define EIGEN_SIZE_MIN_PREFER_DYNAMIC(a, b)
EIGEN_DEVICE_FUNC Derived & derived()
Base class for all dense matrices, vectors, and expressions.
internal::plain_diag_type< MatrixType >::type HCoeffsType
ScalarWithExceptions conj(const ScalarWithExceptions &x)
const FullPivHouseholderQR< PlainObject > fullPivHouseholderQr() const
const Inverse< FullPivHouseholderQR > inverse() const
Index nonzeroPivots() const