![]() |
Eigen
3.2.7
|
Since Eigen version 3.1 and later, users can benefit from built-in Intel MKL optimizations with an installed copy of Intel MKL 10.3 (or later). Intel MKL provides highly optimized multi-threaded mathematical routines for x86-compatible architectures. Intel MKL is available on Linux, Mac and Windows for both Intel64 and IA32 architectures.
Using Intel MKL through Eigen is easy:
EIGEN_USE_MKL_ALL
macro before including any Eigen's headerWhen doing so, a number of Eigen's algorithms are silently substituted with calls to Intel MKL 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.
In addition you can coarsely select choose which parts will be substituted by defining one or multiple of the following macros:
EIGEN_USE_BLAS | Enables the use of external BLAS level 2 and 3 routines (currently works with Intel MKL only) |
EIGEN_USE_LAPACKE | Enables the use of external Lapack routines via the Intel Lapacke C interface to Lapack (currently works with Intel MKL only) |
EIGEN_USE_LAPACKE_STRICT | Same as EIGEN_USE_LAPACKE but algorithm of lower robustness are disabled. This currently concerns only JacobiSVD which otherwise would be replaced by gesvd that is less robust than Jacobi rotations. |
EIGEN_USE_MKL_VML | Enables the use of Intel VML (vector operations) |
EIGEN_USE_MKL_ALL | Defines EIGEN_USE_BLAS , EIGEN_USE_LAPACKE , and EIGEN_USE_MKL_VML |
Finally, the PARDISO sparse solver shipped with Intel MKL can be used through the PardisoLU, PardisoLLT and PardisoLDLT classes of the PardisoSupport module.
The breadth of Eigen functionality covered by Intel MKL is listed in the table below.
Functional domain | Code example | MKL routines |
---|---|---|
Matrix-matrix operations EIGEN_USE_BLAS | ?gemm
?symm/?hemm
?trmm
dsyrk/ssyrk
| |
Matrix-vector operations EIGEN_USE_BLAS | ?gemv
?symv/?hemv
?trmv
| |
LU decomposition EIGEN_USE_LAPACKE EIGEN_USE_LAPACKE_STRICT | v1 = m1.lu().solve(v2);
| ?getrf
|
Cholesky decomposition EIGEN_USE_LAPACKE EIGEN_USE_LAPACKE_STRICT | v1 = m2.selfadjointView<Upper>().llt().solve(v2);
| ?potrf
|
QR decomposition EIGEN_USE_LAPACKE EIGEN_USE_LAPACKE_STRICT | m1.householderQr();
m1.colPivHouseholderQr();
| ?geqrf
?geqp3
|
Singular value decomposition EIGEN_USE_LAPACKE | JacobiSVD<MatrixXd> svd;
svd.compute(m1, ComputeThinV);
| ?gesvd
|
Eigen-value decompositions EIGEN_USE_LAPACKE EIGEN_USE_LAPACKE_STRICT | ?gees
?gees
?syev/?heev
?syev/?heev,
?potrf
| |
Schur decomposition EIGEN_USE_LAPACKE EIGEN_USE_LAPACKE_STRICT | RealSchur<MatrixXd> schurR(m1);
ComplexSchur<MatrixXcd> schurC(m1);
| ?gees
|
Vector Math EIGEN_USE_MKL_VML | v2=v1.array().sin();
v2=v1.array().asin();
v2=v1.array().cos();
v2=v1.array().acos();
v2=v1.array().tan();
v2=v1.array().exp();
v2=v1.array().log();
v2=v1.array().sqrt();
v2=v1.array().square();
v2=v1.array().pow(1.5);
| v?Sin
v?Asin
v?Cos
v?Acos
v?Tan
v?Exp
v?Ln
v?Sqrt
v?Sqr
v?Powx
|
In the examples, m1 and m2 are dense matrices and v1 and v2 are dense vectors.