BDCSVD_LAPACKE.h
Go to the documentation of this file.
1 // This file is part of Eigen, a lightweight C++ template library
2 // for linear algebra.
3 //
4 // Copyright (C) 2022 Melven Roehrig-Zoellner <Melven.Roehrig-Zoellner@DLR.de>
5 // Copyright (c) 2011, Intel Corporation. All rights reserved.
6 //
7 // This file is based on the JacobiSVD_LAPACKE.h originally from Intel -
8 // see license notice below:
9 /*
10  Redistribution and use in source and binary forms, with or without modification,
11  are permitted provided that the following conditions are met:
12 
13  * Redistributions of source code must retain the above copyright notice, this
14  list of conditions and the following disclaimer.
15  * Redistributions in binary form must reproduce the above copyright notice,
16  this list of conditions and the following disclaimer in the documentation
17  and/or other materials provided with the distribution.
18  * Neither the name of Intel Corporation nor the names of its contributors may
19  be used to endorse or promote products derived from this software without
20  specific prior written permission.
21 
22  THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS" AND
23  ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE IMPLIED
24  WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE ARE
25  DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT OWNER OR CONTRIBUTORS BE LIABLE FOR
26  ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES
27  (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES;
28  LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON
29  ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT
30  (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
31  SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
32 
33  ********************************************************************************
34  * Content : Eigen bindings to LAPACKe
35  * Singular Value Decomposition - SVD (divide and conquer variant)
36  ********************************************************************************
37 */
38 #ifndef EIGEN_BDCSVD_LAPACKE_H
39 #define EIGEN_BDCSVD_LAPACKE_H
40 
41 namespace Eigen {
42 
43 namespace internal {
44 
45 namespace lapacke_helpers {
46 
49 // defining a derived class to allow access to protected members
50 template <typename MatrixType_, int Options>
51 class BDCSVD_LAPACKE : public BDCSVD<MatrixType_, Options> {
52  typedef BDCSVD<MatrixType_, Options> SVD;
53  typedef typename SVD::MatrixType MatrixType;
54  typedef typename SVD::Scalar Scalar;
55  typedef typename SVD::RealScalar RealScalar;
56 
57 public:
58  // construct this by moving from a parent object
59  BDCSVD_LAPACKE(SVD&& svd) : SVD(std::move(svd)) {}
60 
61  void compute_impl_lapacke(const MatrixType& matrix, unsigned int computationOptions) {
62 
63  SVD::allocate(matrix.rows(), matrix.cols(), computationOptions);
64 
66 
67  // prepare arguments to ?gesdd
68  const lapack_int matrix_order = lapack_storage_of(matrix);
69  const char jobz = (SVD::m_computeFullU || SVD::m_computeFullV) ? 'A' : (SVD::m_computeThinU || SVD::m_computeThinV) ? 'S' : 'N';
70  const lapack_int u_cols = (jobz == 'A') ? to_lapack(SVD::m_rows) : (jobz == 'S') ? to_lapack(SVD::m_diagSize) : 1;
71  const lapack_int vt_rows = (jobz == 'A') ? to_lapack(SVD::m_cols) : (jobz == 'S') ? to_lapack(SVD::m_diagSize) : 1;
72  lapack_int ldu, ldvt;
73  Scalar *u, *vt, dummy;
74  MatrixType localU;
76  ldu = to_lapack(SVD::m_matrixU.outerStride());
77  u = SVD::m_matrixU.data();
78  } else if (SVD::computeV()) {
79  localU.resize(SVD::m_rows, u_cols);
80  ldu = to_lapack(localU.outerStride());
81  u = localU.data();
82  } else { ldu=1; u=&dummy; }
83  MatrixType localV;
84  if (SVD::computeU() || SVD::computeV()) {
85  localV.resize(vt_rows, SVD::m_cols);
86  ldvt = to_lapack(localV.outerStride());
87  vt = localV.data();
88  } else { ldvt=1; vt=&dummy; }
89  MatrixType temp; temp = matrix;
90 
91  // actual call to ?gesdd
92  lapack_int info = gesdd( matrix_order, jobz, to_lapack(SVD::m_rows), to_lapack(SVD::m_cols),
93  to_lapack(temp.data()), to_lapack(temp.outerStride()), (RealScalar*)SVD::m_singularValues.data(),
94  to_lapack(u), ldu, to_lapack(vt), ldvt);
95 
96  // Check the result of the LAPACK call
97  if (info < 0 || !SVD::m_singularValues.allFinite()) {
98  // this includes info == -4 => NaN entry in A
100  } else if (info > 0 ) {
102  } else {
105  SVD::m_matrixU = localU.leftCols(SVD::m_matrixU.cols());
106  }
107  if (SVD::computeV()) {
108  SVD::m_matrixV = localV.adjoint().leftCols(SVD::m_matrixV.cols());
109  }
110  }
111  SVD::m_isInitialized = true;
112  }
113 };
114 
115 template<typename MatrixType_, int Options>
116 BDCSVD<MatrixType_, Options>& BDCSVD_wrapper(BDCSVD<MatrixType_, Options>& svd, const MatrixType_& matrix, int computationOptions)
117 {
118  // we need to move to the wrapper type and back
119  BDCSVD_LAPACKE<MatrixType_, Options> tmpSvd(std::move(svd));
120  tmpSvd.compute_impl_lapacke(matrix, computationOptions);
121  svd = std::move(tmpSvd);
122  return svd;
123 }
124 
125 } // end namespace lapacke_helpers
126 
127 } // end namespace internal
128 
129 #define EIGEN_LAPACKE_SDD(EIGTYPE, EIGCOLROW, OPTIONS) \
130 template<> inline \
131 BDCSVD<Matrix<EIGTYPE, Dynamic, Dynamic, EIGCOLROW, Dynamic, Dynamic>, OPTIONS>& \
132 BDCSVD<Matrix<EIGTYPE, Dynamic, Dynamic, EIGCOLROW, Dynamic, Dynamic>, OPTIONS>::compute_impl(const Matrix<EIGTYPE, Dynamic, Dynamic, EIGCOLROW, Dynamic, Dynamic>& matrix, unsigned int computationOptions) {\
133  return internal::lapacke_helpers::BDCSVD_wrapper(*this, matrix, computationOptions); \
134 }
135 
136 #define EIGEN_LAPACK_SDD_OPTIONS(OPTIONS) \
137  EIGEN_LAPACKE_SDD(double, ColMajor, OPTIONS) \
138  EIGEN_LAPACKE_SDD(float, ColMajor, OPTIONS) \
139  EIGEN_LAPACKE_SDD(dcomplex, ColMajor, OPTIONS) \
140  EIGEN_LAPACKE_SDD(scomplex, ColMajor, OPTIONS) \
141 \
142  EIGEN_LAPACKE_SDD(double, RowMajor, OPTIONS) \
143  EIGEN_LAPACKE_SDD(float, RowMajor, OPTIONS) \
144  EIGEN_LAPACKE_SDD(dcomplex, RowMajor, OPTIONS) \
145  EIGEN_LAPACKE_SDD(scomplex, RowMajor, OPTIONS)
146 
156 
157 #undef EIGEN_LAPACK_SDD_OPTIONS
158 
159 #undef EIGEN_LAPACKE_SDD
160 
161 } // end namespace Eigen
162 
163 #endif // EIGEN_BDCSVD_LAPACKE_H
#define EIGEN_LAPACK_SDD_OPTIONS(OPTIONS)
cout<< "Here is the matrix m:"<< endl<< m<< endl;JacobiSVD< MatrixXf, ComputeThinU|ComputeThinV > svd(m)
Matrix< float, 1, Dynamic > MatrixType
void allocate(Index rows, Index cols, unsigned int computationOptions)
Definition: BDCSVD.h:287
const AdjointReturnType adjoint() const
Definition: Transpose.h:223
const Scalar * data() const
ComputationInfo info() const
Reports whether previous computation was successful.
Definition: SVDBase.h:310
@ InvalidInput
Definition: Constants.h:453
@ Success
Definition: Constants.h:446
@ NoConvergence
Definition: Constants.h:450
@ ComputeFullV
Definition: Constants.h:399
@ ComputeThinV
Definition: Constants.h:401
@ ComputeFullU
Definition: Constants.h:395
@ ComputeThinU
Definition: Constants.h:397
#define lapack_int
Definition: lapacke.h:52
: InteropHeaders
Definition: Core:139
Definition: BFloat16.h:222