16 #include "./InternalHeaderCheck.h" 21 template <
typename MatrixType_>
22 struct traits<HouseholderQR<MatrixType_>> : traits<MatrixType_> {
23 typedef MatrixXpr XprKind;
24 typedef SolverStorage StorageKind;
25 typedef int StorageIndex;
58 template <
typename MatrixType_>
59 class HouseholderQR :
public SolverBase<HouseholderQR<MatrixType_>> {
61 typedef MatrixType_ MatrixType;
62 typedef SolverBase<HouseholderQR> Base;
67 MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
68 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
71 MaxRowsAtCompileTime, MaxRowsAtCompileTime>
73 typedef typename internal::plain_diag_type<MatrixType>::type HCoeffsType;
74 typedef typename internal::plain_row_type<MatrixType>::type RowVectorType;
75 typedef HouseholderSequence<MatrixType, internal::remove_all_t<typename HCoeffsType::ConjugateReturnType>>
76 HouseholderSequenceType;
85 eigen_assert(m_isInitialized &&
"HouseHolderQR is not initialized.");
95 HouseholderQR() : m_qr(), m_hCoeffs(), m_temp(), m_isInitialized(false) {}
104 : m_qr(rows, cols), m_hCoeffs((
std::min)(rows, cols)), m_temp(cols), m_isInitialized(false) {}
118 template <
typename InputType>
120 : m_qr(matrix.rows(), matrix.cols()),
121 m_hCoeffs((
std::min)(matrix.rows(), matrix.cols())),
122 m_temp(matrix.cols()),
123 m_isInitialized(false) {
134 template <
typename InputType>
137 m_hCoeffs((
std::min)(matrix.rows(), matrix.cols())),
138 m_temp(matrix.cols()),
139 m_isInitialized(false) {
143 #ifdef EIGEN_PARSED_BY_DOXYGEN 158 template <
typename Rhs>
172 eigen_assert(m_isInitialized &&
"HouseholderQR is not initialized.");
180 eigen_assert(m_isInitialized &&
"HouseholderQR is not initialized.");
184 template <
typename InputType>
259 inline Index rows()
const {
return m_qr.rows(); }
260 inline Index cols()
const {
return m_qr.cols(); }
266 const HCoeffsType&
hCoeffs()
const {
return m_hCoeffs; }
268 #ifndef EIGEN_PARSED_BY_DOXYGEN 269 template <
typename RhsType,
typename DstType>
270 void _solve_impl(
const RhsType& rhs, DstType& dst)
const;
272 template <
bool Conjugate,
typename RhsType,
typename DstType>
273 void _solve_impl_transposed(
const RhsType& rhs, DstType& dst)
const;
277 EIGEN_STATIC_ASSERT_NON_INTEGER(Scalar)
282 HCoeffsType m_hCoeffs;
283 RowVectorType m_temp;
284 bool m_isInitialized;
290 template <
typename HCoeffs,
typename Scalar,
bool IsComplex>
291 struct householder_determinant {
292 static void run(
const HCoeffs& hCoeffs, Scalar& out_det) {
294 Index size = hCoeffs.rows();
295 for (
Index i = 0; i < size; i++) {
299 if (hCoeffs(i) != Scalar(0)) out_det *= -numext::conj(hCoeffs(i)) / hCoeffs(i);
305 template <
typename HCoeffs,
typename Scalar>
306 struct householder_determinant<HCoeffs, Scalar, false> {
307 static void run(
const HCoeffs& hCoeffs, Scalar& out_det) {
308 bool negated =
false;
309 Index size = hCoeffs.rows();
310 for (
Index i = 0; i < size; i++) {
312 if (hCoeffs(i) != Scalar(0)) negated ^=
true;
314 out_det = negated ? Scalar(-1) : Scalar(1);
320 template <
typename MatrixType>
322 eigen_assert(m_isInitialized &&
"HouseholderQR is not initialized.");
323 eigen_assert(m_qr.rows() == m_qr.cols() &&
"You can't take the determinant of a non-square matrix!");
325 internal::householder_determinant<HCoeffsType, Scalar, NumTraits<Scalar>::IsComplex>::run(m_hCoeffs, detQ);
326 return m_qr.diagonal().prod() * detQ;
329 template <
typename MatrixType>
332 eigen_assert(m_isInitialized &&
"HouseholderQR is not initialized.");
333 eigen_assert(m_qr.rows() == m_qr.cols() &&
"You can't take the determinant of a non-square matrix!");
334 return abs(m_qr.diagonal().prod());
337 template <
typename MatrixType>
339 eigen_assert(m_isInitialized &&
"HouseholderQR is not initialized.");
340 eigen_assert(m_qr.rows() == m_qr.cols() &&
"You can't take the determinant of a non-square matrix!");
341 return m_qr.diagonal().cwiseAbs().array().log().sum();
344 template <
typename MatrixType>
346 eigen_assert(m_isInitialized &&
"HouseholderQR is not initialized.");
347 eigen_assert(m_qr.rows() == m_qr.cols() &&
"You can't take the determinant of a non-square matrix!");
349 internal::householder_determinant<HCoeffsType, Scalar, NumTraits<Scalar>::IsComplex>::run(m_hCoeffs, detQ);
350 return detQ * m_qr.diagonal().array().sign().prod();
356 template <
typename MatrixQR,
typename HCoeffs>
357 void householder_qr_inplace_unblocked(MatrixQR& mat, HCoeffs& hCoeffs,
typename MatrixQR::Scalar* tempData = 0) {
358 typedef typename MatrixQR::Scalar Scalar;
359 typedef typename MatrixQR::RealScalar RealScalar;
360 Index rows = mat.rows();
361 Index cols = mat.cols();
362 Index size = (std::min)(rows, cols);
364 eigen_assert(hCoeffs.size() == size);
370 tempData = tempVector.data();
373 for (Index k = 0; k < size; ++k) {
374 Index remainingRows = rows - k;
375 Index remainingCols = cols - k - 1;
378 mat.col(k).tail(remainingRows).makeHouseholderInPlace(hCoeffs.coeffRef(k), beta);
379 mat.coeffRef(k, k) = beta;
382 mat.bottomRightCorner(remainingRows, remainingCols)
383 .applyHouseholderOnTheLeft(mat.col(k).tail(remainingRows - 1), hCoeffs.coeffRef(k), tempData + k + 1);
396 template <
typename MatrixQR,
typename HCoeffs,
typename VectorQR>
397 void householder_qr_inplace_update(MatrixQR& mat, HCoeffs& hCoeffs,
const VectorQR& newColumn,
398 typename MatrixQR::Index k,
typename MatrixQR::Scalar* tempData) {
399 typedef typename MatrixQR::Index
Index;
400 typedef typename MatrixQR::RealScalar RealScalar;
401 Index rows = mat.rows();
403 eigen_assert(k < mat.cols());
404 eigen_assert(k < rows);
405 eigen_assert(hCoeffs.size() == mat.cols());
406 eigen_assert(newColumn.size() == rows);
407 eigen_assert(tempData);
410 mat.col(k) = newColumn;
412 for (Index i = 0; i < k; ++i) {
413 Index remainingRows = rows - i;
416 .applyHouseholderOnTheLeft(mat.col(i).tail(remainingRows - 1), hCoeffs.coeffRef(i), tempData + i + 1);
420 mat.col(k).tail(rows - k).makeHouseholderInPlace(hCoeffs.coeffRef(k), beta);
421 mat.coeffRef(k, k) = beta;
425 template <
typename MatrixQR,
typename HCoeffs,
typename MatrixQRScalar =
typename MatrixQR::Scalar,
426 bool InnerStrideIsOne = (MatrixQR::InnerStrideAtCompileTime == 1 && HCoeffs::InnerStrideAtCompileTime == 1)>
427 struct householder_qr_inplace_blocked {
429 static void run(MatrixQR& mat, HCoeffs& hCoeffs, Index maxBlockSize = 32,
typename MatrixQR::Scalar* tempData = 0) {
430 typedef typename MatrixQR::Scalar Scalar;
431 typedef Block<MatrixQR, Dynamic, Dynamic> BlockType;
433 Index rows = mat.rows();
434 Index cols = mat.cols();
435 Index size = (std::min)(rows, cols);
437 typedef Matrix<Scalar, Dynamic, 1, ColMajor, MatrixQR::MaxColsAtCompileTime, 1> TempType;
440 tempVector.resize(cols);
441 tempData = tempVector.data();
444 Index blockSize = (std::min)(maxBlockSize, size);
447 for (k = 0; k < size; k += blockSize) {
448 Index bs = (std::min)(size - k, blockSize);
449 Index tcols = cols - k - bs;
450 Index brows = rows - k;
460 BlockType A11_21 = mat.block(k, k, brows, bs);
461 Block<HCoeffs, Dynamic, 1> hCoeffsSegment = hCoeffs.segment(k, bs);
463 householder_qr_inplace_unblocked(A11_21, hCoeffsSegment, tempData);
466 BlockType A21_22 = mat.block(k, k + bs, brows, tcols);
467 apply_block_householder_on_the_left(A21_22, A11_21, hCoeffsSegment,
false);
475 #ifndef EIGEN_PARSED_BY_DOXYGEN 476 template <
typename MatrixType_>
477 template <
typename RhsType,
typename DstType>
478 void HouseholderQR<MatrixType_>::_solve_impl(
const RhsType& rhs, DstType& dst)
const {
479 const Index rank = (std::min)(rows(), cols());
481 typename RhsType::PlainObject c(rhs);
483 c.applyOnTheLeft(householderQ().setLength(rank).adjoint());
485 m_qr.topLeftCorner(rank, rank).template triangularView<Upper>().solveInPlace(c.topRows(rank));
487 dst.topRows(rank) = c.topRows(rank);
488 dst.bottomRows(cols() - rank).setZero();
491 template <
typename MatrixType_>
492 template <
bool Conjugate,
typename RhsType,
typename DstType>
493 void HouseholderQR<MatrixType_>::_solve_impl_transposed(
const RhsType& rhs, DstType& dst)
const {
494 const Index rank = (std::min)(rows(), cols());
496 typename RhsType::PlainObject c(rhs);
498 m_qr.topLeftCorner(rank, rank)
499 .template triangularView<Upper>()
501 .template conjugateIf<Conjugate>()
502 .solveInPlace(c.topRows(rank));
504 dst.topRows(rank) = c.topRows(rank);
505 dst.bottomRows(rows() - rank).setZero();
507 dst.applyOnTheLeft(householderQ().setLength(rank).
template conjugateIf<!Conjugate>());
517 template <
typename MatrixType>
519 Index rows = m_qr.rows();
520 Index cols = m_qr.cols();
521 Index size = (std::min)(rows, cols);
523 m_hCoeffs.resize(size);
527 internal::householder_qr_inplace_blocked<MatrixType, HCoeffsType>::run(m_qr, m_hCoeffs, 48, m_temp.data());
529 m_isInitialized =
true;
536 template <
typename Derived>
HouseholderSequenceType householderQ() const
Definition: HouseholderQR.h:171
const MatrixType & matrixQR() const
Definition: HouseholderQR.h:179
Definition: Constants.h:318
constexpr void resize(Index rows, Index cols)
Definition: PlainObjectBase.h:268
MatrixType::Scalar signDeterminant() const
Definition: HouseholderQR.h:345
MatrixType::RealScalar absDeterminant() const
Definition: HouseholderQR.h:330
Namespace containing all symbols from the Eigen library.
Definition: B01_Experimental.dox:1
MatrixType::Scalar determinant() const
Definition: HouseholderQR.h:321
Definition: BFloat16.h:231
Eigen::Index Index
The interface type of indices.
Definition: EigenBase.h:43
const unsigned int RowMajorBit
Definition: Constants.h:70
Definition: EigenBase.h:33
Sequence of Householder reflections acting on subspaces with decreasing size.
Definition: ForwardDeclarations.h:445
HouseholderQR(EigenBase< InputType > &matrix)
Constructs a QR factorization from a given matrix.
Definition: HouseholderQR.h:135
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:82
MatrixType::RealScalar logAbsDeterminant() const
Definition: HouseholderQR.h:338
const Solve< HouseholderQR, Rhs > solve(const MatrixBase< Rhs > &b) const
Definition: Constants.h:440
constexpr Derived & derived()
Definition: EigenBase.h:49
ComputationInfo info() const
Reports whether the QR factorization was successful.
Definition: HouseholderQR.h:84
const HCoeffsType & hCoeffs() const
Definition: HouseholderQR.h:266
SolverBase()
Definition: SolverBase.h:105
Householder QR decomposition of a matrix.
Definition: ForwardDeclarations.h:427
HouseholderQR(Index rows, Index cols)
Default Constructor with memory preallocation.
Definition: HouseholderQR.h:103
void computeInPlace()
Definition: HouseholderQR.h:518
Definition: Constants.h:320
Pseudo expression representing a solving operation.
Definition: Solve.h:62
The matrix class, also used for vectors and row-vectors.
Definition: Matrix.h:186
ComputationInfo
Definition: Constants.h:438
HouseholderQR()
Default Constructor.
Definition: HouseholderQR.h:95
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:52
HouseholderQR(const EigenBase< InputType > &matrix)
Constructs a QR factorization from a given matrix.
Definition: HouseholderQR.h:119