|
TR-mbed 1.0
|
Since Eigen version 3.3 and later, any F77 compatible BLAS or LAPACK libraries can be used as backends for dense matrix products and dense matrix decompositions. For instance, one can use IntelĀ® MKL, Apple's Accelerate framework on OSX, OpenBLAS, Netlib LAPACK, etc.
Do not miss this page for further discussions on the specific use of IntelĀ® MKL (also includes VML, PARDISO, etc.)
In order to use an external BLAS and/or LAPACK library, you must link you own application to the respective libraries and their dependencies. For LAPACK, you must also link to the standard Lapacke library, which is used as a convenient think layer between Eigen's C++ code and LAPACK F77 interface. Then you must activate their usage by defining one or multiple of the following macros (before including any Eigen's header):
-framework Accelerate /opt/local/lib/lapack/liblapacke.dylibEIGEN_USE_BLAS | Enables the use of external BLAS level 2 and 3 routines (compatible with any F77 BLAS interface) |
EIGEN_USE_LAPACKE | Enables the use of external Lapack routines via the Lapacke C interface to Lapack (compatible with any F77 LAPACK interface) |
EIGEN_USE_LAPACKE_STRICT | Same as EIGEN_USE_LAPACKE but algorithms of lower numerical robustness are disabled. This currently concerns only JacobiSVD which otherwise would be replaced by gesvd that is less robust than Jacobi rotations. |
When doing so, a number of Eigen's algorithms are silently substituted with calls to BLAS or LAPACK routines. These substitutions apply only for Dynamic or large enough objects with one of the following four standard scalar types: float, double, complex<float>, and complex<double>. Operations on other scalar types or mixing reals and complexes will continue to use the built-in algorithms.
The breadth of Eigen functionality that can be substituted is listed in the table below.
| Functional domain | Code example | BLAS/LAPACK routines |
|---|---|---|
Matrix-matrix operations EIGEN_USE_BLAS | ?gemm
?trmm
int BLASFUNC() ssyrk(const char *, const char *, const int *, const int *, const float *, const float *, const int *, const float *, float *, const int *) int BLASFUNC() dsyrk(const char *, const char *, const int *, const int *, const double *, const double *, const int *, const double *, double *, const int *) int EIGEN_BLAS_FUNC() symm(const char *side, const char *uplo, const int *m, const int *n, const RealScalar *palpha, const RealScalar *pa, const int *lda, const RealScalar *pb, const int *ldb, const RealScalar *pbeta, RealScalar *pc, const int *ldc) Definition level3_impl.h:287 int EIGEN_BLAS_FUNC() hemm(const char *side, const char *uplo, const int *m, const int *n, const RealScalar *palpha, const RealScalar *pa, const int *lda, const RealScalar *pb, const int *ldb, const RealScalar *pbeta, RealScalar *pc, const int *ldc) Definition level3_impl.h:505 int EIGEN_BLAS_FUNC() trmm(const char *side, const char *uplo, const char *opa, const char *diag, const int *m, const int *n, const RealScalar *palpha, const RealScalar *pa, const int *lda, RealScalar *pb, const int *ldb) Definition level3_impl.h:183 | |
Matrix-vector operations EIGEN_USE_BLAS | ?gemv
?trmv
int EIGEN_BLAS_FUNC() hemv(const char *uplo, const int *n, const RealScalar *palpha, const RealScalar *pa, const int *lda, const RealScalar *px, const int *incx, const RealScalar *pbeta, RealScalar *py, const int *incy) Definition level2_cplx_impl.h:19 int EIGEN_BLAS_FUNC() symv(const char *uplo, const int *n, const RealScalar *palpha, const RealScalar *pa, const int *lda, const RealScalar *px, const int *incx, const RealScalar *pbeta, RealScalar *py, const int *incy) Definition level2_real_impl.h:13 | |
LU decomposition EIGEN_USE_LAPACKE EIGEN_USE_LAPACKE_STRICT | ?getrf
| |
Cholesky decomposition EIGEN_USE_LAPACKE EIGEN_USE_LAPACKE_STRICT | ?potrf
| |
QR decomposition EIGEN_USE_LAPACKE EIGEN_USE_LAPACKE_STRICT | ?geqrf
?geqp3
| |
Singular value decomposition EIGEN_USE_LAPACKE | JacobiSVD<MatrixXd> svd;
cout<< "Here is the matrix m:"<< endl<< m<< endl;JacobiSVD< MatrixXf > svd(m, ComputeThinU|ComputeThinV) | ?gesvd
|
Eigen-value decompositions EIGEN_USE_LAPACKE EIGEN_USE_LAPACKE_STRICT | GeneralizedSelfAdjointEigenSolver<MatrixXd>
cout<< "Here is a random 4x4 matrix, A:"<< endl<< A<< endl<< endl;ComplexEigenSolver< MatrixXcf > ces Definition ComplexEigenSolver_compute.cpp:4 | ?gees
?gees
?syev/?heev
?syev/?heev,
?potrf
|
Schur decomposition EIGEN_USE_LAPACKE EIGEN_USE_LAPACKE_STRICT | ?gees
|
In the examples, m1 and m2 are dense matrices and v1 and v2 are dense vectors.