You can not select more than 25 topics Topics must start with a letter or number, can include dashes ('-') and can be up to 35 characters long.

167 lines
7.5 KiB

  1. /*
  2. Copyright (c) 2011, Intel Corporation. All rights reserved.
  3. Copyright (C) 2011 Gael Guennebaud <gael.guennebaud@inria.fr>
  4. Redistribution and use in source and binary forms, with or without modification,
  5. are permitted provided that the following conditions are met:
  6. * Redistributions of source code must retain the above copyright notice, this
  7. list of conditions and the following disclaimer.
  8. * Redistributions in binary form must reproduce the above copyright notice,
  9. this list of conditions and the following disclaimer in the documentation
  10. and/or other materials provided with the distribution.
  11. * Neither the name of Intel Corporation nor the names of its contributors may
  12. be used to endorse or promote products derived from this software without
  13. specific prior written permission.
  14. THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS" AND
  15. ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE IMPLIED
  16. WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE ARE
  17. DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT OWNER OR CONTRIBUTORS BE LIABLE FOR
  18. ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES
  19. (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES;
  20. LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON
  21. ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT
  22. (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
  23. SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
  24. ********************************************************************************
  25. * Content : Documentation on the use of Intel MKL through Eigen
  26. ********************************************************************************
  27. */
  28. namespace Eigen {
  29. /** \page TopicUsingIntelMKL Using Intel® Math Kernel Library from Eigen
  30. \section TopicUsingIntelMKL_Intro Eigen and Intel® Math Kernel Library (Intel® MKL)
  31. 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).
  32. <a href="http://eigen.tuxfamily.org/Counter/redirect_to_mkl.php"> Intel MKL </a> provides highly optimized multi-threaded mathematical routines for x86-compatible architectures.
  33. Intel MKL is available on Linux, Mac and Windows for both Intel64 and IA32 architectures.
  34. \warning Be aware that Intel® MKL is a proprietary software. It is the responsibility of the users to buy MKL licenses for their products. Moreover, the license of the user product has to allow linking to proprietary software that excludes any unmodified versions of the GPL. As a consequence, this also means that Eigen has to be used through the LGPL3+ license.
  35. Using Intel MKL through Eigen is easy:
  36. -# define the \c EIGEN_USE_MKL_ALL macro before including any Eigen's header
  37. -# link your program to MKL libraries (see the <a href="http://software.intel.com/en-us/articles/intel-mkl-link-line-advisor/">MKL linking advisor</a>)
  38. -# on a 64bits system, you must use the LP64 interface (not the ILP64 one)
  39. When doing so, a number of Eigen's algorithms are silently substituted with calls to Intel MKL routines.
  40. These substitutions apply only for \b Dynamic \b or \b large enough objects with one of the following four standard scalar types: \c float, \c double, \c complex<float>, and \c complex<double>.
  41. Operations on other scalar types or mixing reals and complexes will continue to use the built-in algorithms.
  42. In addition you can coarsely select choose which parts will be substituted by defining one or multiple of the following macros:
  43. <table class="manual">
  44. <tr><td>\c EIGEN_USE_BLAS </td><td>Enables the use of external BLAS level 2 and 3 routines (currently works with Intel MKL only)</td></tr>
  45. <tr class="alt"><td>\c EIGEN_USE_LAPACKE </td><td>Enables the use of external Lapack routines via the <a href="http://www.netlib.org/lapack/lapacke.html">Intel Lapacke</a> C interface to Lapack (currently works with Intel MKL only)</td></tr>
  46. <tr><td>\c EIGEN_USE_LAPACKE_STRICT </td><td>Same as \c EIGEN_USE_LAPACKE but algorithm of lower robustness are disabled. This currently concerns only JacobiSVD which otherwise would be replaced by \c gesvd that is less robust than Jacobi rotations.</td></tr>
  47. <tr class="alt"><td>\c EIGEN_USE_MKL_VML </td><td>Enables the use of Intel VML (vector operations)</td></tr>
  48. <tr><td>\c EIGEN_USE_MKL_ALL </td><td>Defines \c EIGEN_USE_BLAS, \c EIGEN_USE_LAPACKE, and \c EIGEN_USE_MKL_VML </td></tr>
  49. </table>
  50. Finally, the PARDISO sparse solver shipped with Intel MKL can be used through the \ref PardisoLU, \ref PardisoLLT and \ref PardisoLDLT classes of the \ref PARDISOSupport_Module.
  51. \section TopicUsingIntelMKL_SupportedFeatures List of supported features
  52. The breadth of Eigen functionality covered by Intel MKL is listed in the table below.
  53. <table class="manual">
  54. <tr><th>Functional domain</th><th>Code example</th><th>MKL routines</th></tr>
  55. <tr><td>Matrix-matrix operations \n \c EIGEN_USE_BLAS </td><td>\code
  56. m1*m2.transpose();
  57. m1.selfadjointView<Lower>()*m2;
  58. m1*m2.triangularView<Upper>();
  59. m1.selfadjointView<Lower>().rankUpdate(m2,1.0);
  60. \endcode</td><td>\code
  61. ?gemm
  62. ?symm/?hemm
  63. ?trmm
  64. dsyrk/ssyrk
  65. \endcode</td></tr>
  66. <tr class="alt"><td>Matrix-vector operations \n \c EIGEN_USE_BLAS </td><td>\code
  67. m1.adjoint()*b;
  68. m1.selfadjointView<Lower>()*b;
  69. m1.triangularView<Upper>()*b;
  70. \endcode</td><td>\code
  71. ?gemv
  72. ?symv/?hemv
  73. ?trmv
  74. \endcode</td></tr>
  75. <tr><td>LU decomposition \n \c EIGEN_USE_LAPACKE \n \c EIGEN_USE_LAPACKE_STRICT </td><td>\code
  76. v1 = m1.lu().solve(v2);
  77. \endcode</td><td>\code
  78. ?getrf
  79. \endcode</td></tr>
  80. <tr class="alt"><td>Cholesky decomposition \n \c EIGEN_USE_LAPACKE \n \c EIGEN_USE_LAPACKE_STRICT </td><td>\code
  81. v1 = m2.selfadjointView<Upper>().llt().solve(v2);
  82. \endcode</td><td>\code
  83. ?potrf
  84. \endcode</td></tr>
  85. <tr><td>QR decomposition \n \c EIGEN_USE_LAPACKE \n \c EIGEN_USE_LAPACKE_STRICT </td><td>\code
  86. m1.householderQr();
  87. m1.colPivHouseholderQr();
  88. \endcode</td><td>\code
  89. ?geqrf
  90. ?geqp3
  91. \endcode</td></tr>
  92. <tr class="alt"><td>Singular value decomposition \n \c EIGEN_USE_LAPACKE </td><td>\code
  93. JacobiSVD<MatrixXd> svd;
  94. svd.compute(m1, ComputeThinV);
  95. \endcode</td><td>\code
  96. ?gesvd
  97. \endcode</td></tr>
  98. <tr><td>Eigen-value decompositions \n \c EIGEN_USE_LAPACKE \n \c EIGEN_USE_LAPACKE_STRICT </td><td>\code
  99. EigenSolver<MatrixXd> es(m1);
  100. ComplexEigenSolver<MatrixXcd> ces(m1);
  101. SelfAdjointEigenSolver<MatrixXd> saes(m1+m1.transpose());
  102. GeneralizedSelfAdjointEigenSolver<MatrixXd>
  103. gsaes(m1+m1.transpose(),m2+m2.transpose());
  104. \endcode</td><td>\code
  105. ?gees
  106. ?gees
  107. ?syev/?heev
  108. ?syev/?heev,
  109. ?potrf
  110. \endcode</td></tr>
  111. <tr class="alt"><td>Schur decomposition \n \c EIGEN_USE_LAPACKE \n \c EIGEN_USE_LAPACKE_STRICT </td><td>\code
  112. RealSchur<MatrixXd> schurR(m1);
  113. ComplexSchur<MatrixXcd> schurC(m1);
  114. \endcode</td><td>\code
  115. ?gees
  116. \endcode</td></tr>
  117. <tr><td>Vector Math \n \c EIGEN_USE_MKL_VML </td><td>\code
  118. v2=v1.array().sin();
  119. v2=v1.array().asin();
  120. v2=v1.array().cos();
  121. v2=v1.array().acos();
  122. v2=v1.array().tan();
  123. v2=v1.array().exp();
  124. v2=v1.array().log();
  125. v2=v1.array().sqrt();
  126. v2=v1.array().square();
  127. v2=v1.array().pow(1.5);
  128. \endcode</td><td>\code
  129. v?Sin
  130. v?Asin
  131. v?Cos
  132. v?Acos
  133. v?Tan
  134. v?Exp
  135. v?Ln
  136. v?Sqrt
  137. v?Sqr
  138. v?Powx
  139. \endcode</td></tr>
  140. </table>
  141. In the examples, m1 and m2 are dense matrices and v1 and v2 are dense vectors.
  142. \section TopicUsingIntelMKL_Links Links
  143. - Intel MKL can be purchased and downloaded <a href="http://eigen.tuxfamily.org/Counter/redirect_to_mkl.php">here</a>.
  144. - Intel MKL is also bundled with <a href="http://software.intel.com/en-us/articles/intel-composer-xe/">Intel Composer XE</a>.
  145. */
  146. }