Home | History | Annotate | Download | only in doc
      1 /*
      2  Copyright (c) 2011, Intel Corporation. All rights reserved.
      3  Copyright (C) 2011-2016 Gael Guennebaud <gael.guennebaud (at) inria.fr>
      4 
      5  Redistribution and use in source and binary forms, with or without modification,
      6  are permitted provided that the following conditions are met:
      7 
      8  * Redistributions of source code must retain the above copyright notice, this
      9    list of conditions and the following disclaimer.
     10  * Redistributions in binary form must reproduce the above copyright notice,
     11    this list of conditions and the following disclaimer in the documentation
     12    and/or other materials provided with the distribution.
     13  * Neither the name of Intel Corporation nor the names of its contributors may
     14    be used to endorse or promote products derived from this software without
     15    specific prior written permission.
     16 
     17  THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS" AND
     18  ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE IMPLIED
     19  WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE ARE
     20  DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT OWNER OR CONTRIBUTORS BE LIABLE FOR
     21  ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES
     22  (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES;
     23  LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON
     24  ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT
     25  (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
     26  SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
     27 
     28  ********************************************************************************
     29  *   Content : Documentation on the use of BLAS/LAPACK libraries through Eigen
     30  ********************************************************************************
     31 */
     32 
     33 namespace Eigen {
     34 
     35 /** \page TopicUsingBlasLapack Using BLAS/LAPACK from %Eigen
     36 
     37 
     38 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.
     39 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.
     40 
     41 Do not miss this \link TopicUsingIntelMKL page \endlink for further discussions on the specific use of Intel MKL (also includes VML, PARDISO, etc.)
     42 
     43 In order to use an external BLAS and/or LAPACK library, you must link you own application to the respective libraries and their dependencies.
     44 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):
     45 
     46 \note For Mac users, in order to use the lapack version shipped with the Accelerate framework, you also need the lapacke library.
     47 Using <a href="https://www.macports.org/">MacPorts</a>, this is as easy as:
     48 \code
     49 sudo port install lapack
     50 \endcode
     51 and then use the following link flags: \c -framework \c Accelerate \c /opt/local/lib/lapack/liblapacke.dylib
     52 
     53 <table class="manual">
     54 <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>
     55 <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>
     56 <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>
     57 </table>
     58 
     59 When doing so, a number of %Eigen's algorithms are silently substituted with calls to BLAS or LAPACK routines.
     60 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>.
     61 Operations on other scalar types or mixing reals and complexes will continue to use the built-in algorithms.
     62 
     63 The breadth of %Eigen functionality that can be substituted is listed in the table below.
     64 <table class="manual">
     65 <tr><th>Functional domain</th><th>Code example</th><th>BLAS/LAPACK routines</th></tr>
     66 <tr><td>Matrix-matrix operations \n \c EIGEN_USE_BLAS </td><td>\code
     67 m1*m2.transpose();
     68 m1.selfadjointView<Lower>()*m2;
     69 m1*m2.triangularView<Upper>();
     70 m1.selfadjointView<Lower>().rankUpdate(m2,1.0);
     71 \endcode</td><td>\code
     72 ?gemm
     73 ?symm/?hemm
     74 ?trmm
     75 dsyrk/ssyrk
     76 \endcode</td></tr>
     77 <tr class="alt"><td>Matrix-vector operations \n \c EIGEN_USE_BLAS </td><td>\code
     78 m1.adjoint()*b;
     79 m1.selfadjointView<Lower>()*b;
     80 m1.triangularView<Upper>()*b;
     81 \endcode</td><td>\code
     82 ?gemv
     83 ?symv/?hemv
     84 ?trmv
     85 \endcode</td></tr>
     86 <tr><td>LU decomposition \n \c EIGEN_USE_LAPACKE \n \c EIGEN_USE_LAPACKE_STRICT </td><td>\code
     87 v1 = m1.lu().solve(v2);
     88 \endcode</td><td>\code
     89 ?getrf
     90 \endcode</td></tr>
     91 <tr class="alt"><td>Cholesky decomposition \n \c EIGEN_USE_LAPACKE \n \c EIGEN_USE_LAPACKE_STRICT </td><td>\code
     92 v1 = m2.selfadjointView<Upper>().llt().solve(v2);
     93 \endcode</td><td>\code
     94 ?potrf
     95 \endcode</td></tr>
     96 <tr><td>QR decomposition \n \c EIGEN_USE_LAPACKE \n \c EIGEN_USE_LAPACKE_STRICT </td><td>\code
     97 m1.householderQr();
     98 m1.colPivHouseholderQr();
     99 \endcode</td><td>\code
    100 ?geqrf
    101 ?geqp3
    102 \endcode</td></tr>
    103 <tr class="alt"><td>Singular value decomposition \n \c EIGEN_USE_LAPACKE </td><td>\code
    104 JacobiSVD<MatrixXd> svd;
    105 svd.compute(m1, ComputeThinV);
    106 \endcode</td><td>\code
    107 ?gesvd
    108 \endcode</td></tr>
    109 <tr><td>Eigen-value decompositions \n \c EIGEN_USE_LAPACKE \n \c EIGEN_USE_LAPACKE_STRICT </td><td>\code
    110 EigenSolver<MatrixXd> es(m1);
    111 ComplexEigenSolver<MatrixXcd> ces(m1);
    112 SelfAdjointEigenSolver<MatrixXd> saes(m1+m1.transpose());
    113 GeneralizedSelfAdjointEigenSolver<MatrixXd>
    114     gsaes(m1+m1.transpose(),m2+m2.transpose());
    115 \endcode</td><td>\code
    116 ?gees
    117 ?gees
    118 ?syev/?heev
    119 ?syev/?heev,
    120 ?potrf
    121 \endcode</td></tr>
    122 <tr class="alt"><td>Schur decomposition \n \c EIGEN_USE_LAPACKE \n \c EIGEN_USE_LAPACKE_STRICT </td><td>\code
    123 RealSchur<MatrixXd> schurR(m1);
    124 ComplexSchur<MatrixXcd> schurC(m1);
    125 \endcode</td><td>\code
    126 ?gees
    127 \endcode</td></tr>
    128 </table>
    129 In the examples, m1 and m2 are dense matrices and v1 and v2 are dense vectors.
    130 
    131 */
    132 
    133 }
    134