42 template<
typename _MatrixType>
class HouseholderQR
46 typedef _MatrixType MatrixType;
48 RowsAtCompileTime = MatrixType::RowsAtCompileTime,
49 ColsAtCompileTime = MatrixType::ColsAtCompileTime,
50 Options = MatrixType::Options,
51 MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
52 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
54 typedef typename MatrixType::Scalar Scalar;
55 typedef typename MatrixType::RealScalar RealScalar;
56 typedef typename MatrixType::Index Index;
57 typedef Matrix<Scalar, RowsAtCompileTime, RowsAtCompileTime, (MatrixType::Flags&RowMajorBit) ? RowMajor : ColMajor, MaxRowsAtCompileTime, MaxRowsAtCompileTime> MatrixQType;
58 typedef typename internal::plain_diag_type<MatrixType>::type HCoeffsType;
59 typedef typename internal::plain_row_type<MatrixType>::type RowVectorType;
60 typedef HouseholderSequence<MatrixType,typename internal::remove_all<typename HCoeffsType::ConjugateReturnType>::type> HouseholderSequenceType;
68 HouseholderQR() : m_qr(), m_hCoeffs(), m_temp(), m_isInitialized(false) {}
78 m_hCoeffs((std::min)(rows,cols)),
80 m_isInitialized(false) {}
95 : m_qr(matrix.rows(), matrix.cols()),
96 m_hCoeffs((std::min)(matrix.rows(),matrix.cols())),
97 m_temp(matrix.cols()),
98 m_isInitialized(false)
120 template<
typename Rhs>
121 inline const internal::solve_retval<HouseholderQR, Rhs>
124 eigen_assert(m_isInitialized &&
"HouseholderQR is not initialized.");
125 return internal::solve_retval<HouseholderQR, Rhs>(*
this, b.derived());
138 eigen_assert(m_isInitialized &&
"HouseholderQR is not initialized.");
147 eigen_assert(m_isInitialized &&
"HouseholderQR is not initialized.");
182 inline Index rows()
const {
return m_qr.rows(); }
183 inline Index cols()
const {
return m_qr.cols(); }
189 const HCoeffsType&
hCoeffs()
const {
return m_hCoeffs; }
193 static void check_template_parameters()
195 EIGEN_STATIC_ASSERT_NON_INTEGER(Scalar);
199 HCoeffsType m_hCoeffs;
200 RowVectorType m_temp;
201 bool m_isInitialized;
204 template<
typename MatrixType>
208 eigen_assert(m_isInitialized &&
"HouseholderQR is not initialized.");
209 eigen_assert(m_qr.rows() == m_qr.cols() &&
"You can't take the determinant of a non-square matrix!");
210 return abs(m_qr.diagonal().prod());
213 template<
typename MatrixType>
216 eigen_assert(m_isInitialized &&
"HouseholderQR is not initialized.");
217 eigen_assert(m_qr.rows() == m_qr.cols() &&
"You can't take the determinant of a non-square matrix!");
218 return m_qr.diagonal().cwiseAbs().array().log().sum();
224 template<
typename MatrixQR,
typename HCoeffs>
225 void householder_qr_inplace_unblocked(MatrixQR& mat, HCoeffs& hCoeffs,
typename MatrixQR::Scalar* tempData = 0)
227 typedef typename MatrixQR::Index Index;
228 typedef typename MatrixQR::Scalar Scalar;
229 typedef typename MatrixQR::RealScalar RealScalar;
230 Index rows = mat.rows();
231 Index cols = mat.cols();
232 Index size = (std::min)(rows,cols);
234 eigen_assert(hCoeffs.size() == size);
241 tempData = tempVector.data();
244 for(Index k = 0; k < size; ++k)
246 Index remainingRows = rows - k;
247 Index remainingCols = cols - k - 1;
250 mat.col(k).tail(remainingRows).makeHouseholderInPlace(hCoeffs.coeffRef(k), beta);
251 mat.coeffRef(k,k) = beta;
254 mat.bottomRightCorner(remainingRows, remainingCols)
255 .applyHouseholderOnTheLeft(mat.col(k).tail(remainingRows-1), hCoeffs.coeffRef(k), tempData+k+1);
260 template<
typename MatrixQR,
typename HCoeffs,
261 typename MatrixQRScalar =
typename MatrixQR::Scalar,
262 bool InnerStrideIsOne = (MatrixQR::InnerStrideAtCompileTime == 1 && HCoeffs::InnerStrideAtCompileTime == 1)>
263 struct householder_qr_inplace_blocked
266 static void run(MatrixQR& mat, HCoeffs& hCoeffs,
267 typename MatrixQR::Index maxBlockSize=32,
268 typename MatrixQR::Scalar* tempData = 0)
270 typedef typename MatrixQR::Index Index;
271 typedef typename MatrixQR::Scalar Scalar;
272 typedef Block<MatrixQR,Dynamic,Dynamic> BlockType;
274 Index rows = mat.rows();
275 Index cols = mat.cols();
276 Index size = (std::min)(rows, cols);
278 typedef Matrix<Scalar,Dynamic,1,ColMajor,MatrixQR::MaxColsAtCompileTime,1> TempType;
282 tempVector.resize(cols);
283 tempData = tempVector.data();
286 Index blockSize = (std::min)(maxBlockSize,size);
289 for (k = 0; k < size; k += blockSize)
291 Index bs = (std::min)(size-k,blockSize);
292 Index tcols = cols - k - bs;
293 Index brows = rows-k;
303 BlockType A11_21 = mat.block(k,k,brows,bs);
304 Block<HCoeffs,Dynamic,1> hCoeffsSegment = hCoeffs.segment(k,bs);
306 householder_qr_inplace_unblocked(A11_21, hCoeffsSegment, tempData);
310 BlockType A21_22 = mat.block(k,k+bs,brows,tcols);
311 apply_block_householder_on_the_left(A21_22,A11_21,hCoeffsSegment.adjoint());
317 template<
typename _MatrixType,
typename Rhs>
318 struct solve_retval<HouseholderQR<_MatrixType>, Rhs>
319 : solve_retval_base<HouseholderQR<_MatrixType>, Rhs>
321 EIGEN_MAKE_SOLVE_HELPERS(HouseholderQR<_MatrixType>,Rhs)
323 template<typename Dest>
void evalTo(Dest& dst)
const
325 const Index rows = dec().rows(), cols = dec().cols();
326 const Index rank = (std::min)(rows, cols);
327 eigen_assert(rhs().rows() == rows);
329 typename Rhs::PlainObject c(rhs());
333 dec().matrixQR().leftCols(rank),
334 dec().hCoeffs().head(rank)).transpose()
338 .topLeftCorner(rank, rank)
339 .template triangularView<Upper>()
340 .solveInPlace(c.topRows(rank));
342 dst.topRows(rank) = c.topRows(rank);
343 dst.bottomRows(cols-rank).setZero();
355 template<
typename MatrixType>
358 check_template_parameters();
360 Index rows = matrix.rows();
361 Index cols = matrix.cols();
362 Index size = (std::min)(rows,cols);
365 m_hCoeffs.resize(size);
369 internal::householder_qr_inplace_blocked<MatrixType, HCoeffsType>::run(m_qr, m_hCoeffs, 48, m_temp.data());
371 m_isInitialized =
true;
379 template<
typename Derived>
HouseholderQR(const MatrixType &matrix)
Constructs a QR factorization from a given matrix.
Definition: HouseholderQR.h:94
MatrixType::RealScalar logAbsDeterminant() const
Definition: HouseholderQR.h:214
HouseholderSequence< VectorsType, CoeffsType > householderSequence(const VectorsType &v, const CoeffsType &h)
Convenience function for constructing a Householder sequence.
Definition: HouseholderSequence.h:422
Sequence of Householder reflections acting on subspaces with decreasing size.
Definition: ForwardDeclarations.h:227
HouseholderSequenceType householderQ() const
Definition: HouseholderQR.h:136
const HouseholderQR< PlainObject > householderQr() const
Definition: HouseholderQR.h:381
const internal::solve_retval< HouseholderQR, Rhs > solve(const MatrixBase< Rhs > &b) const
Definition: HouseholderQR.h:122
HouseholderQR()
Default Constructor.
Definition: HouseholderQR.h:68
MatrixType::RealScalar absDeterminant() const
Definition: HouseholderQR.h:205
const HCoeffsType & hCoeffs() const
Definition: HouseholderQR.h:189
HouseholderQR & compute(const MatrixType &matrix)
Definition: HouseholderQR.h:356
Householder QR decomposition of a matrix.
Definition: ForwardDeclarations.h:221
HouseholderQR(Index rows, Index cols)
Default Constructor with memory preallocation.
Definition: HouseholderQR.h:76
void resize(Index nbRows, Index nbCols)
Definition: PlainObjectBase.h:235
const MatrixType & matrixQR() const
Definition: HouseholderQR.h:145
The matrix class, also used for vectors and row-vectors.
Definition: Matrix.h:127
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:48