11 #ifndef EIGEN_HOUSEHOLDER_SEQUENCE_H 12 #define EIGEN_HOUSEHOLDER_SEQUENCE_H 59 template<
typename VectorsType,
typename CoeffsType,
int S
ide>
62 typedef typename VectorsType::Scalar
Scalar;
68 ColsAtCompileTime = RowsAtCompileTime,
71 MaxColsAtCompileTime = MaxRowsAtCompileTime,
78 template<
typename VectorsType,
typename CoeffsType,
int S
ide>
85 template<
typename VectorsType,
typename CoeffsType,
int S
ide>
97 template<
typename VectorsType,
typename CoeffsType>
113 typedef Matrix<
ResultScalar, MatrixType::RowsAtCompileTime, MatrixType::ColsAtCompileTime,
114 0, MatrixType::MaxRowsAtCompileTime, MatrixType::MaxColsAtCompileTime>
Type;
120 :
public EigenBase<HouseholderSequence<VectorsType,CoeffsType,Side> >
137 typename internal::conditional<NumTraits<Scalar>::IsComplex,
161 : m_vectors(v), m_coeffs(h), m_trans(false), m_length(v.diagonalSize()),
168 : m_vectors(other.m_vectors),
169 m_coeffs(other.m_coeffs),
170 m_trans(other.m_trans),
171 m_length(other.m_length),
172 m_shift(other.m_shift)
233 template<
typename DestType>
inline void evalTo(DestType& dst)
const 235 Matrix<Scalar, DestType::RowsAtCompileTime, 1,
237 evalTo(dst, workspace);
241 template<
typename Dest,
typename Workspace>
242 void evalTo(Dest& dst, Workspace& workspace)
const 244 workspace.resize(rows());
245 Index vecs = m_length;
249 dst.diagonal().setOnes();
250 dst.template triangularView<StrictlyUpper>().
setZero();
251 for(
Index k = vecs-1; k >= 0; --k)
253 Index cornerSize = rows() - k - m_shift;
255 dst.bottomRightCorner(cornerSize, cornerSize)
256 .applyHouseholderOnTheRight(essentialVector(k), m_coeffs.coeff(k), workspace.data());
258 dst.bottomRightCorner(cornerSize, cornerSize)
259 .applyHouseholderOnTheLeft(essentialVector(k), m_coeffs.coeff(k), workspace.data());
262 dst.col(k).tail(rows()-k-1).setZero();
265 for(
Index k = 0; k<cols()-vecs ; ++k)
266 dst.col(k).tail(rows()-k-1).setZero();
270 dst.setIdentity(rows(), rows());
271 for(
Index k = vecs-1; k >= 0; --k)
273 Index cornerSize = rows() - k - m_shift;
275 dst.bottomRightCorner(cornerSize, cornerSize)
276 .applyHouseholderOnTheRight(essentialVector(k), m_coeffs.coeff(k), &workspace.coeffRef(0));
278 dst.bottomRightCorner(cornerSize, cornerSize)
279 .applyHouseholderOnTheLeft(essentialVector(k), m_coeffs.coeff(k), &workspace.coeffRef(0));
288 applyThisOnTheRight(dst, workspace);
292 template<
typename Dest,
typename Workspace>
295 workspace.resize(dst.rows());
296 for(
Index k = 0; k < m_length; ++k)
298 Index actual_k = m_trans ? m_length-k-1 : k;
299 dst.rightCols(rows()-m_shift-actual_k)
300 .applyHouseholderOnTheRight(essentialVector(actual_k), m_coeffs.coeff(actual_k), workspace.data());
308 applyThisOnTheLeft(dst, workspace);
312 template<
typename Dest,
typename Workspace>
315 const Index BlockSize = 48;
317 if(m_length>=BlockSize && dst.cols()>1)
319 for(
Index i = 0; i < m_length; i+=BlockSize)
321 Index end = m_trans ? (std::min)(m_length,i+BlockSize) : m_length-i;
324 Index start = k + m_shift;
327 SubVectorsType sub_vecs1(m_vectors.const_cast_derived(), Side==
OnTheRight ? k : start,
329 Side==
OnTheRight ? bs : m_vectors.rows()-start,
330 Side==
OnTheRight ? m_vectors.cols()-start : bs);
338 workspace.resize(dst.cols());
339 for(
Index k = 0; k < m_length; ++k)
341 Index actual_k = m_trans ? k : m_length-k-1;
342 dst.bottomRows(rows()-m_shift-actual_k)
343 .applyHouseholderOnTheLeft(essentialVector(actual_k), m_coeffs.coeff(actual_k), workspace.data());
355 template<
typename OtherDerived>
360 applyThisOnTheLeft(res);
420 bool trans()
const {
return m_trans; }
437 template<
typename OtherDerived,
typename VectorsType,
typename CoeffsType,
int S
ide>
450 template<
typename VectorsType,
typename CoeffsType>
462 template<
typename VectorsType,
typename CoeffsType>
470 #endif // EIGEN_HOUSEHOLDER_SEQUENCE_H HouseholderSequence< VectorsType, CoeffsType, OnTheRight > HouseholderSequenceType
HouseholderSequence< VectorsType, CoeffsType > householderSequence(const VectorsType &v, const CoeffsType &h)
Convenience function for constructing a Householder sequence.
HouseholderSequenceShape Shape
Index length() const
Returns the length of the Householder sequence.
void apply_block_householder_on_the_left(MatrixType &mat, const VectorsType &vectors, const CoeffsType &hCoeffs, bool forward)
HouseholderSequence< typename internal::conditional< NumTraits< Scalar >::IsComplex, typename internal::remove_all< typename VectorsType::ConjugateReturnType >::type, VectorsType >::type, typename internal::conditional< NumTraits< Scalar >::IsComplex, typename internal::remove_all< typename CoeffsType::ConjugateReturnType >::type, CoeffsType >::type, Side > ConjugateReturnType
bool is_same_dense(const T1 &mat1, const T2 &mat2, typename enable_if< has_direct_access< T1 >::ret &&has_direct_access< T2 >::ret, T1 >::type *=0)
Matrix< ResultScalar, MatrixType::RowsAtCompileTime, MatrixType::ColsAtCompileTime, 0, MatrixType::MaxRowsAtCompileTime, MatrixType::MaxColsAtCompileTime > Type
HouseholderSequence & setLength(Index length)
Sets the length of the Householder sequence.
ScalarBinaryOpTraits< OtherScalarType, typename MatrixType::Scalar >::ReturnType ResultScalar
internal::hseq_side_dependent_impl< VectorsType, CoeffsType, Side >::EssentialVectorType EssentialVectorType
Index cols() const
Number of columns of transformation viewed as a matrix.
Expression of the transpose of a matrix.
HouseholderSequence(const HouseholderSequence &other)
Copy constructor.
void applyThisOnTheLeft(Dest &dst) const
VectorsType::Scalar Scalar
HouseholderSequence(const VectorsType &v, const CoeffsType &h)
Constructor.
Eigen::Index Index
The interface type of indices.
VectorsType::StorageIndex StorageIndex
HouseholderSequence transpose() const
Transpose of the Householder sequence.
Block< const VectorsType, Dynamic, 1 > EssentialVectorType
void applyThisOnTheRight(Dest &dst, Workspace &workspace) const
Sequence of Householder reflections acting on subspaces with decreasing size.
EIGEN_STRONG_INLINE EIGEN_DEVICE_FUNC half() max(const half &a, const half &b)
Transpose< Block< const VectorsType, 1, Dynamic > > EssentialVectorType
static const EssentialVectorType essentialVector(const HouseholderSequenceType &h, Index k)
HouseholderSequence< VectorsType, CoeffsType, OnTheRight > rightHouseholderSequence(const VectorsType &v, const CoeffsType &h)
Convenience function for constructing a Householder sequence.
HouseholderSequence & setTrans(bool trans)
Sets the transpose flag.
EIGEN_DEVICE_FUNC ConjugateReturnType conjugate() const
internal::matrix_type_times_scalar_type< Scalar, OtherDerived >::Type operator*(const MatrixBase< OtherDerived > &other) const
Computes the product of a Householder sequence with a matrix.
ConjugateReturnType conjugate() const
Complex conjugate of the Householder sequence.
VectorsType::Nested m_vectors
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
void applyThisOnTheLeft(Dest &dst, Workspace &workspace) const
void evalTo(Dest &dst, Workspace &workspace) const
CoeffsType::Nested m_coeffs
Index rows() const
Number of rows of transformation viewed as a matrix.
HouseholderSequence & setShift(Index shift)
Sets the shift of the Householder sequence.
static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE TensorUInt128< uint64_t, uint64_t > operator*(const TensorUInt128< HL, LL > &lhs, const TensorUInt128< HR, LR > &rhs)
void applyThisOnTheRight(Dest &dst) const
VectorsType::StorageKind StorageKind
ConjugateReturnType adjoint() const
Adjoint (conjugate transpose) of the Householder sequence.
void evalTo(DestType &dst) const
Expression of a fixed-size or dynamic-size block.
EIGEN_DEVICE_FUNC CastXpr< NewType >::Type cast() const
ConjugateReturnType inverse() const
Inverse of the Householder sequence (equals the adjoint).
const EssentialVectorType essentialVector(Index k) const
Essential part of a Householder vector.
bool trans() const
Returns the transpose flag.
Determines whether the given binary operation of two numeric types is allowed and what the scalar ret...
internal::conditional< NumTraits< Scalar >::IsComplex, const CwiseUnaryOp< internal::scalar_conjugate_op< Scalar >, const Derived >, const Derived & >::type ConjugateReturnType
Index shift() const
Returns the shift of the Householder sequence.
The matrix class, also used for vectors and row-vectors.
HouseholderSequence< VectorsType, CoeffsType, OnTheLeft > HouseholderSequenceType
Base class for all dense matrices, vectors, and expressions.
static const EssentialVectorType essentialVector(const HouseholderSequenceType &h, Index k)
internal::traits< HouseholderSequence >::Scalar Scalar