UsingBlasLapackBackends.dox 6.5 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133
  1. /*
  2. Copyright (c) 2011, Intel Corporation. All rights reserved.
  3. Copyright (C) 2011-2016 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 BLAS/LAPACK libraries through Eigen
  26. ********************************************************************************
  27. */
  28. namespace Eigen {
  29. /** \page TopicUsingBlasLapack Using BLAS/LAPACK from %Eigen
  30. 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.
  31. For instance, one can use <a href="http://eigen.tuxfamily.org/Counter/redirect_to_mkl.php">Intel® MKL</a>, Apple's Accelerate framework on OSX, <a href="http://www.openblas.net/">OpenBLAS</a>, <a href="http://www.netlib.org/lapack">Netlib LAPACK</a>, etc.
  32. Do not miss this \link TopicUsingIntelMKL page \endlink for further discussions on the specific use of Intel® MKL (also includes VML, PARDISO, etc.)
  33. In order to use an external BLAS and/or LAPACK library, you must link you own application to the respective libraries and their dependencies.
  34. For LAPACK, you must also link to the standard <a href="http://www.netlib.org/lapack/lapacke.html">Lapacke</a> 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 (\b before including any %Eigen's header):
  35. \note For Mac users, in order to use the lapack version shipped with the Accelerate framework, you also need the lapacke library.
  36. Using <a href="https://www.macports.org/">MacPorts</a>, this is as easy as:
  37. \code
  38. sudo port install lapack
  39. \endcode
  40. and then use the following link flags: \c -framework \c Accelerate \c /opt/local/lib/lapack/liblapacke.dylib
  41. <table class="manual">
  42. <tr><td>\c EIGEN_USE_BLAS </td><td>Enables the use of external BLAS level 2 and 3 routines (compatible with any F77 BLAS interface)</td></tr>
  43. <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">Lapacke</a> C interface to Lapack (compatible with any F77 LAPACK interface)</td></tr>
  44. <tr><td>\c EIGEN_USE_LAPACKE_STRICT </td><td>Same as \c EIGEN_USE_LAPACKE but algorithms of lower numerical robustness are disabled. \n This currently concerns only JacobiSVD which otherwise would be replaced by \c gesvd that is less robust than Jacobi rotations.</td></tr>
  45. </table>
  46. When doing so, a number of %Eigen's algorithms are silently substituted with calls to BLAS or LAPACK routines.
  47. 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>.
  48. Operations on other scalar types or mixing reals and complexes will continue to use the built-in algorithms.
  49. The breadth of %Eigen functionality that can be substituted is listed in the table below.
  50. <table class="manual">
  51. <tr><th>Functional domain</th><th>Code example</th><th>BLAS/LAPACK routines</th></tr>
  52. <tr><td>Matrix-matrix operations \n \c EIGEN_USE_BLAS </td><td>\code
  53. m1*m2.transpose();
  54. m1.selfadjointView<Lower>()*m2;
  55. m1*m2.triangularView<Upper>();
  56. m1.selfadjointView<Lower>().rankUpdate(m2,1.0);
  57. \endcode</td><td>\code
  58. ?gemm
  59. ?symm/?hemm
  60. ?trmm
  61. dsyrk/ssyrk
  62. \endcode</td></tr>
  63. <tr class="alt"><td>Matrix-vector operations \n \c EIGEN_USE_BLAS </td><td>\code
  64. m1.adjoint()*b;
  65. m1.selfadjointView<Lower>()*b;
  66. m1.triangularView<Upper>()*b;
  67. \endcode</td><td>\code
  68. ?gemv
  69. ?symv/?hemv
  70. ?trmv
  71. \endcode</td></tr>
  72. <tr><td>LU decomposition \n \c EIGEN_USE_LAPACKE \n \c EIGEN_USE_LAPACKE_STRICT </td><td>\code
  73. v1 = m1.lu().solve(v2);
  74. \endcode</td><td>\code
  75. ?getrf
  76. \endcode</td></tr>
  77. <tr class="alt"><td>Cholesky decomposition \n \c EIGEN_USE_LAPACKE \n \c EIGEN_USE_LAPACKE_STRICT </td><td>\code
  78. v1 = m2.selfadjointView<Upper>().llt().solve(v2);
  79. \endcode</td><td>\code
  80. ?potrf
  81. \endcode</td></tr>
  82. <tr><td>QR decomposition \n \c EIGEN_USE_LAPACKE \n \c EIGEN_USE_LAPACKE_STRICT </td><td>\code
  83. m1.householderQr();
  84. m1.colPivHouseholderQr();
  85. \endcode</td><td>\code
  86. ?geqrf
  87. ?geqp3
  88. \endcode</td></tr>
  89. <tr class="alt"><td>Singular value decomposition \n \c EIGEN_USE_LAPACKE </td><td>\code
  90. JacobiSVD<MatrixXd> svd;
  91. svd.compute(m1, ComputeThinV);
  92. \endcode</td><td>\code
  93. ?gesvd
  94. \endcode</td></tr>
  95. <tr><td>Eigen-value decompositions \n \c EIGEN_USE_LAPACKE \n \c EIGEN_USE_LAPACKE_STRICT </td><td>\code
  96. EigenSolver<MatrixXd> es(m1);
  97. ComplexEigenSolver<MatrixXcd> ces(m1);
  98. SelfAdjointEigenSolver<MatrixXd> saes(m1+m1.transpose());
  99. GeneralizedSelfAdjointEigenSolver<MatrixXd>
  100. gsaes(m1+m1.transpose(),m2+m2.transpose());
  101. \endcode</td><td>\code
  102. ?gees
  103. ?gees
  104. ?syev/?heev
  105. ?syev/?heev,
  106. ?potrf
  107. \endcode</td></tr>
  108. <tr class="alt"><td>Schur decomposition \n \c EIGEN_USE_LAPACKE \n \c EIGEN_USE_LAPACKE_STRICT </td><td>\code
  109. RealSchur<MatrixXd> schurR(m1);
  110. ComplexSchur<MatrixXcd> schurC(m1);
  111. \endcode</td><td>\code
  112. ?gees
  113. \endcode</td></tr>
  114. </table>
  115. In the examples, m1 and m2 are dense matrices and v1 and v2 are dense vectors.
  116. */
  117. }