34 #ifndef EIGEN_COLPIVOTINGHOUSEHOLDERQR_LAPACKE_H 35 #define EIGEN_COLPIVOTINGHOUSEHOLDERQR_LAPACKE_H 38 #include "./InternalHeaderCheck.h" 42 #if defined(EIGEN_USE_LAPACKE) 44 template <
typename Scalar>
45 inline lapack_int call_geqp3(
int matrix_layout, lapack_int m, lapack_int n, Scalar* a, lapack_int lda, lapack_int* jpvt,
48 inline lapack_int call_geqp3(
int matrix_layout, lapack_int m, lapack_int n,
float* a, lapack_int lda, lapack_int* jpvt,
50 return LAPACKE_sgeqp3(matrix_layout, m, n, a, lda, jpvt, tau);
53 inline lapack_int call_geqp3(
int matrix_layout, lapack_int m, lapack_int n,
double* a, lapack_int lda, lapack_int* jpvt,
55 return LAPACKE_dgeqp3(matrix_layout, m, n, a, lda, jpvt, tau);
58 inline lapack_int call_geqp3(
int matrix_layout, lapack_int m, lapack_int n, lapack_complex_float* a, lapack_int lda,
59 lapack_int* jpvt, lapack_complex_float* tau) {
60 return LAPACKE_cgeqp3(matrix_layout, m, n, a, lda, jpvt, tau);
63 inline lapack_int call_geqp3(
int matrix_layout, lapack_int m, lapack_int n, lapack_complex_double* a, lapack_int lda,
64 lapack_int* jpvt, lapack_complex_double* tau) {
65 return LAPACKE_zgeqp3(matrix_layout, m, n, a, lda, jpvt, tau);
68 template <
typename MatrixType>
69 struct ColPivHouseholderQR_LAPACKE_impl {
70 typedef typename MatrixType::Scalar Scalar;
71 typedef typename MatrixType::RealScalar RealScalar;
72 typedef typename internal::lapacke_helpers::translate_type_imp<Scalar>::type LapackeType;
73 static constexpr
int LapackeStorage = MatrixType::IsRowMajor ? (LAPACK_ROW_MAJOR) : (LAPACK_COL_MAJOR);
75 typedef typename internal::plain_diag_type<MatrixType>::type HCoeffsType;
76 typedef PermutationMatrix<Dynamic, Dynamic, lapack_int> PermutationType;
78 static void run(MatrixType& qr, HCoeffsType& hCoeffs, PermutationType& colsPermutation,
Index& nonzero_pivots,
79 RealScalar& maxpivot,
bool usePrescribedThreshold, RealScalar prescribedThreshold,
Index& det_p,
80 bool& isInitialized) {
81 isInitialized =
false;
82 hCoeffs.resize(qr.diagonalSize());
84 maxpivot = RealScalar(0);
85 colsPermutation.resize(qr.cols());
86 colsPermutation.indices().setZero();
88 lapack_int rows = internal::lapacke_helpers::to_lapack(qr.rows());
89 lapack_int cols = internal::lapacke_helpers::to_lapack(qr.cols());
90 LapackeType* qr_data = (LapackeType*)(qr.data());
91 lapack_int lda = internal::lapacke_helpers::to_lapack(qr.outerStride());
92 lapack_int* perm_data = colsPermutation.indices().data();
93 LapackeType* hCoeffs_data = (LapackeType*)(hCoeffs.data());
95 lapack_int info = call_geqp3(LapackeStorage, rows, cols, qr_data, lda, perm_data, hCoeffs_data);
96 if (info != 0)
return;
98 maxpivot = qr.diagonal().cwiseAbs().maxCoeff();
99 hCoeffs.adjointInPlace();
100 RealScalar defaultThreshold = NumTraits<RealScalar>::epsilon() * RealScalar(qr.diagonalSize());
101 RealScalar threshold = usePrescribedThreshold ? prescribedThreshold : defaultThreshold;
102 RealScalar premultiplied_threshold = maxpivot * threshold;
103 nonzero_pivots = (qr.diagonal().cwiseAbs().array() > premultiplied_threshold).count();
104 colsPermutation.indices().array() -= 1;
105 det_p = colsPermutation.determinant();
106 isInitialized =
true;
109 static void init(
Index rows,
Index cols, HCoeffsType& hCoeffs, PermutationType& colsPermutation,
110 bool& usePrescribedThreshold,
bool& isInitialized) {
111 Index diag = numext::mini(rows, cols);
112 hCoeffs.resize(diag);
113 colsPermutation.resize(cols);
114 usePrescribedThreshold =
false;
115 isInitialized =
false;
119 #define COLPIVQR_LAPACKE_COMPUTEINPLACE(EIGTYPE) \ 121 inline void ColPivHouseholderQR<EIGTYPE, lapack_int>::computeInPlace() { \ 122 ColPivHouseholderQR_LAPACKE_impl<MatrixType>::run(m_qr, m_hCoeffs, m_colsPermutation, m_nonzero_pivots, \ 123 m_maxpivot, m_usePrescribedThreshold, m_prescribedThreshold, \ 124 m_det_p, m_isInitialized); \ 127 #define COLPIVQR_LAPACKE_INIT(EIGTYPE) \ 129 inline void ColPivHouseholderQR<EIGTYPE, lapack_int>::init(Index rows, Index cols) { \ 130 ColPivHouseholderQR_LAPACKE_impl<MatrixType>::init(rows, cols, m_hCoeffs, m_colsPermutation, m_isInitialized, \ 131 m_usePrescribedThreshold); \ 134 #define COLPIVQR_LAPACKE(EIGTYPE) \ 135 COLPIVQR_LAPACKE_COMPUTEINPLACE(EIGTYPE) \ 136 COLPIVQR_LAPACKE_INIT(EIGTYPE) \ 137 COLPIVQR_LAPACKE_COMPUTEINPLACE(Ref<EIGTYPE>) \ 138 COLPIVQR_LAPACKE_INIT(Ref<EIGTYPE>) 140 typedef Matrix<float, Dynamic, Dynamic, ColMajor> MatrixXfC;
141 typedef Matrix<double, Dynamic, Dynamic, ColMajor> MatrixXdC;
144 typedef Matrix<float, Dynamic, Dynamic, RowMajor> MatrixXfR;
145 typedef Matrix<double, Dynamic, Dynamic, RowMajor> MatrixXdR;
149 COLPIVQR_LAPACKE(MatrixXfC)
150 COLPIVQR_LAPACKE(MatrixXdC)
151 COLPIVQR_LAPACKE(MatrixXcfC)
152 COLPIVQR_LAPACKE(MatrixXcdC)
153 COLPIVQR_LAPACKE(MatrixXfR)
154 COLPIVQR_LAPACKE(MatrixXdR)
155 COLPIVQR_LAPACKE(MatrixXcfR)
156 COLPIVQR_LAPACKE(MatrixXcdR)
161 #endif // EIGEN_COLPIVOTINGHOUSEHOLDERQR_LAPACKE_H Definition: Constants.h:318
Namespace containing all symbols from the Eigen library.
Definition: B01_Experimental.dox:1
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:82
Definition: Constants.h:320
const int Dynamic
Definition: Constants.h:25