ColPivHouseholderQR_LAPACKE.h
Go to the documentation of this file.
1 /*
2  Copyright (c) 2011, Intel Corporation. All rights reserved.
3 
4  Redistribution and use in source and binary forms, with or without modification,
5  are permitted provided that the following conditions are met:
6 
7  * Redistributions of source code must retain the above copyright notice, this
8  list of conditions and the following disclaimer.
9  * Redistributions in binary form must reproduce the above copyright notice,
10  this list of conditions and the following disclaimer in the documentation
11  and/or other materials provided with the distribution.
12  * Neither the name of Intel Corporation nor the names of its contributors may
13  be used to endorse or promote products derived from this software without
14  specific prior written permission.
15 
16  THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS" AND
17  ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE IMPLIED
18  WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE ARE
19  DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT OWNER OR CONTRIBUTORS BE LIABLE FOR
20  ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES
21  (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES;
22  LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON
23  ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT
24  (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
25  SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
26 
27  ********************************************************************************
28  * Content : Eigen bindings to LAPACKe
29  * Householder QR decomposition of a matrix with column pivoting based on
30  * LAPACKE_?geqp3 function.
31  ********************************************************************************
32 */
33 
34 #ifndef EIGEN_COLPIVOTINGHOUSEHOLDERQR_LAPACKE_H
35 #define EIGEN_COLPIVOTINGHOUSEHOLDERQR_LAPACKE_H
36 
37 #include "./InternalHeaderCheck.h"
38 
39 namespace Eigen {
40 
41 #if defined(EIGEN_USE_LAPACKE)
42 
43  template<typename Scalar>
44  inline lapack_int call_geqp3(int matrix_layout, lapack_int m, lapack_int n, Scalar* a, lapack_int lda, lapack_int* jpvt, Scalar* tau);
45  template<>
46  inline lapack_int call_geqp3(int matrix_layout, lapack_int m, lapack_int n, float* a, lapack_int lda, lapack_int* jpvt, float* tau)
47  { return LAPACKE_sgeqp3(matrix_layout, m, n, a, lda, jpvt, tau); }
48  template<>
49  inline lapack_int call_geqp3(int matrix_layout, lapack_int m, lapack_int n, double* a, lapack_int lda, lapack_int* jpvt, double* tau)
50  { return LAPACKE_dgeqp3(matrix_layout, m, n, a, lda, jpvt, tau); }
51  template<>
52  inline lapack_int call_geqp3(int matrix_layout, lapack_int m, lapack_int n, lapack_complex_float* a, lapack_int lda, lapack_int* jpvt, lapack_complex_float* tau)
53  { return LAPACKE_cgeqp3(matrix_layout, m, n, a, lda, jpvt, tau); }
54  template<>
55  inline lapack_int call_geqp3(int matrix_layout, lapack_int m, lapack_int n, lapack_complex_double* a, lapack_int lda, lapack_int* jpvt, lapack_complex_double* tau)
56  { return LAPACKE_zgeqp3(matrix_layout, m, n, a, lda, jpvt, tau); }
57 
58  template <typename MatrixType>
59  struct ColPivHouseholderQR_LAPACKE_impl {
60  typedef typename MatrixType::Scalar Scalar;
61  typedef typename MatrixType::RealScalar RealScalar;
62  typedef typename internal::lapacke_helpers::translate_type_imp<Scalar>::type LapackeType;
63  static constexpr int LapackeStorage = MatrixType::IsRowMajor ? (LAPACK_ROW_MAJOR) : (LAPACK_COL_MAJOR);
64 
65  typedef typename internal::plain_diag_type<MatrixType>::type HCoeffsType;
66  typedef PermutationMatrix<Dynamic, Dynamic, lapack_int> PermutationType;
67 
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) {
71 
72  isInitialized = false;
73  hCoeffs.resize(qr.diagonalSize());
74  nonzero_pivots = 0;
75  maxpivot = RealScalar(0);
76  colsPermutation.resize(qr.cols());
77  colsPermutation.indices().setZero();
78 
79  lapack_int rows = internal::lapacke_helpers::to_lapack(qr.rows());
80  lapack_int cols = internal::lapacke_helpers::to_lapack(qr.cols());
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());
85 
86  lapack_int info = call_geqp3(LapackeStorage, rows, cols, qr_data, lda, perm_data, hCoeffs_data);
87  if (info != 0) return;
88 
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();
97  isInitialized = true;
98  };
99 
100  static void init(Index rows, Index cols, HCoeffsType& hCoeffs, PermutationType& colsPermutation,
101  bool& usePrescribedThreshold, bool& isInitialized) {
102 
103  Index diag = numext::mini(rows, cols);
104  hCoeffs.resize(diag);
105  colsPermutation.resize(cols);
106  usePrescribedThreshold = false;
107  isInitialized = false;
108  }
109  };
110 
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); } \
116 
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); } \
121 
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>) \
127 
128  typedef Matrix<float, Dynamic, Dynamic, ColMajor> MatrixXfC;
129  typedef Matrix<double, Dynamic, Dynamic, ColMajor> MatrixXdC;
130  typedef Matrix<std::complex<float>, Dynamic, Dynamic, ColMajor> MatrixXcfC;
131  typedef Matrix<std::complex<double>, Dynamic, Dynamic, ColMajor> MatrixXcdC;
132  typedef Matrix<float, Dynamic, Dynamic, RowMajor> MatrixXfR;
133  typedef Matrix<double, Dynamic, Dynamic, RowMajor> MatrixXdR;
134  typedef Matrix<std::complex<float>, Dynamic, Dynamic, RowMajor> MatrixXcfR;
135  typedef Matrix<std::complex<double>, Dynamic, Dynamic, RowMajor> MatrixXcdR;
136 
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)
145 
146 #endif
147 } // end namespace Eigen
148 
149 #endif // EIGEN_COLPIVOTINGHOUSEHOLDERQR_LAPACKE_H
Matrix3f m
int n
HouseholderQR< MatrixXf > qr(A)
Matrix< float, 1, Dynamic > MatrixType
NumTraits< Scalar >::Real RealScalar
Definition: DenseBase.h:68
internal::traits< Derived >::Scalar Scalar
Definition: DenseBase.h:61
@ ColMajor
Definition: Constants.h:321
@ RowMajor
Definition: Constants.h:323
#define LAPACK_COL_MAJOR
Definition: lapacke.h:125
lapack_int LAPACKE_sgeqp3(int matrix_order, lapack_int m, lapack_int n, float *a, lapack_int lda, lapack_int *jpvt, float *tau)
#define lapack_int
Definition: lapacke.h:52
#define LAPACK_ROW_MAJOR
Definition: lapacke.h:124
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
Definition: lapacke.h:94
#define lapack_complex_float
Definition: lapacke.h:79
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)
: InteropHeaders
Definition: Core:139
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:82
const int Dynamic
Definition: Constants.h:24