34 #ifndef EIGEN_COLPIVOTINGHOUSEHOLDERQR_LAPACKE_H
35 #define EIGEN_COLPIVOTINGHOUSEHOLDERQR_LAPACKE_H
41 #if defined(EIGEN_USE_LAPACKE)
43 template<
typename Scalar>
58 template <
typename MatrixType>
59 struct ColPivHouseholderQR_LAPACKE_impl {
62 typedef typename internal::lapacke_helpers::translate_type_imp<Scalar>::type LapackeType;
65 typedef typename internal::plain_diag_type<MatrixType>::type HCoeffsType;
66 typedef PermutationMatrix<Dynamic, Dynamic, lapack_int> PermutationType;
68 static void run(
MatrixType&
qr, HCoeffsType& hCoeffs, PermutationType& colsPermutation,
Index& nonzero_pivots,
69 RealScalar& maxpivot,
bool usePrescribedThreshold, RealScalar prescribedThreshold,
Index& det_p,
70 bool& isInitialized) {
72 isInitialized =
false;
73 hCoeffs.resize(
qr.diagonalSize());
75 maxpivot = RealScalar(0);
76 colsPermutation.resize(
qr.cols());
77 colsPermutation.indices().setZero();
81 LapackeType* qr_data = (LapackeType*)(
qr.data());
82 lapack_int lda = internal::lapacke_helpers::to_lapack(
qr.outerStride());
83 lapack_int* perm_data = colsPermutation.indices().data();
84 LapackeType* hCoeffs_data = (LapackeType*)(hCoeffs.data());
86 lapack_int info = call_geqp3(LapackeStorage,
rows,
cols, qr_data, lda, perm_data, hCoeffs_data);
87 if (info != 0)
return;
89 maxpivot =
qr.diagonal().cwiseAbs().maxCoeff();
90 hCoeffs.adjointInPlace();
91 RealScalar defaultThreshold = NumTraits<RealScalar>::epsilon() * RealScalar(
qr.diagonalSize());
92 RealScalar threshold = usePrescribedThreshold ? prescribedThreshold : defaultThreshold;
93 RealScalar premultiplied_threshold = maxpivot * threshold;
94 nonzero_pivots = (
qr.diagonal().cwiseAbs().array() > premultiplied_threshold).count();
95 colsPermutation.indices().array() -= 1;
96 det_p = colsPermutation.determinant();
100 static void init(
Index rows,
Index cols, HCoeffsType& hCoeffs, PermutationType& colsPermutation,
101 bool& usePrescribedThreshold,
bool& isInitialized) {
104 hCoeffs.resize(diag);
105 colsPermutation.resize(
cols);
106 usePrescribedThreshold =
false;
107 isInitialized =
false;
111 #define COLPIVQR_LAPACKE_COMPUTEINPLACE(EIGTYPE) \
112 template <> inline void ColPivHouseholderQR<EIGTYPE, lapack_int>::computeInPlace() { \
113 ColPivHouseholderQR_LAPACKE_impl<MatrixType>::run(m_qr, m_hCoeffs, m_colsPermutation, m_nonzero_pivots, \
114 m_maxpivot, m_usePrescribedThreshold, m_prescribedThreshold, \
115 m_det_p, m_isInitialized); } \
117 #define COLPIVQR_LAPACKE_INIT(EIGTYPE) \
118 template <> inline void ColPivHouseholderQR<EIGTYPE, lapack_int>::init(Index rows, Index cols) { \
119 ColPivHouseholderQR_LAPACKE_impl<MatrixType>::init(rows, cols, m_hCoeffs, m_colsPermutation, m_isInitialized, \
120 m_usePrescribedThreshold); } \
122 #define COLPIVQR_LAPACKE(EIGTYPE) \
123 COLPIVQR_LAPACKE_COMPUTEINPLACE(EIGTYPE) \
124 COLPIVQR_LAPACKE_INIT(EIGTYPE) \
125 COLPIVQR_LAPACKE_COMPUTEINPLACE(Ref<EIGTYPE>) \
126 COLPIVQR_LAPACKE_INIT(Ref<EIGTYPE>) \
128 typedef Matrix<float, Dynamic, Dynamic, ColMajor> MatrixXfC;
129 typedef Matrix<double, Dynamic, Dynamic, ColMajor> MatrixXdC;
132 typedef Matrix<float, Dynamic, Dynamic, RowMajor> MatrixXfR;
133 typedef Matrix<double, Dynamic, Dynamic, RowMajor> MatrixXdR;
137 COLPIVQR_LAPACKE(MatrixXfC)
138 COLPIVQR_LAPACKE(MatrixXdC)
139 COLPIVQR_LAPACKE(MatrixXcfC)
140 COLPIVQR_LAPACKE(MatrixXcdC)
141 COLPIVQR_LAPACKE(MatrixXfR)
142 COLPIVQR_LAPACKE(MatrixXdR)
143 COLPIVQR_LAPACKE(MatrixXcfR)
144 COLPIVQR_LAPACKE(MatrixXcdR)
HouseholderQR< MatrixXf > qr(A)
Matrix< float, 1, Dynamic > MatrixType
NumTraits< Scalar >::Real RealScalar
internal::traits< Derived >::Scalar Scalar
lapack_int LAPACKE_sgeqp3(int matrix_order, lapack_int m, lapack_int n, float *a, lapack_int lda, lapack_int *jpvt, float *tau)
lapack_int LAPACKE_zgeqp3(int matrix_order, lapack_int m, lapack_int n, lapack_complex_double *a, lapack_int lda, lapack_int *jpvt, lapack_complex_double *tau)
lapack_int LAPACKE_dgeqp3(int matrix_order, lapack_int m, lapack_int n, double *a, lapack_int lda, lapack_int *jpvt, double *tau)
#define lapack_complex_double
#define lapack_complex_float
lapack_int LAPACKE_cgeqp3(int matrix_order, lapack_int m, lapack_int n, lapack_complex_float *a, lapack_int lda, lapack_int *jpvt, lapack_complex_float *tau)
EIGEN_ALWAYS_INLINE T mini(const T &x, const T &y)
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.