/* Copyright (c) 2011, Intel Corporation. All rights reserved. Copyright (C) 2011-2016 Gael Guennebaud Redistribution and use in source and binary forms, with or without modification, are permitted provided that the following conditions are met: * Redistributions of source code must retain the above copyright notice, this list of conditions and the following disclaimer. * Redistributions in binary form must reproduce the above copyright notice, this list of conditions and the following disclaimer in the documentation and/or other materials provided with the distribution. * Neither the name of Intel Corporation nor the names of its contributors may be used to endorse or promote products derived from this software without specific prior written permission. THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT OWNER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE. ******************************************************************************** * Content : Documentation on the use of BLAS/LAPACK libraries through Eigen ******************************************************************************** */ namespace Eigen { /** \page TopicUsingBlasLapack Using BLAS/LAPACK from %Eigen 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 \link TopicUsingIntelMKL page \endlink 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 (\b before including any %Eigen's header): \note For Mac users, in order to use the lapack version shipped with the Accelerate framework, you also need the lapacke library. Using MacPorts, this is as easy as: \code sudo port install lapack \endcode and then use the following link flags: \c -framework \c Accelerate \c /opt/local/lib/lapack/liblapacke.dylib
\c EIGEN_USE_BLAS Enables the use of external BLAS level 2 and 3 routines (compatible with any F77 BLAS interface)
\c EIGEN_USE_LAPACKE Enables the use of external Lapack routines via the Lapacke C interface to Lapack (compatible with any F77 LAPACK interface)
\c EIGEN_USE_LAPACKE_STRICT 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.
When doing so, a number of %Eigen's algorithms are silently substituted with calls to BLAS or LAPACK routines. 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, and \c complex. 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 domainCode exampleBLAS/LAPACK routines
Matrix-matrix operations \n \c EIGEN_USE_BLAS \code m1*m2.transpose(); m1.selfadjointView()*m2; m1*m2.triangularView(); m1.selfadjointView().rankUpdate(m2,1.0); \endcode\code ?gemm ?symm/?hemm ?trmm dsyrk/ssyrk \endcode
Matrix-vector operations \n \c EIGEN_USE_BLAS \code m1.adjoint()*b; m1.selfadjointView()*b; m1.triangularView()*b; \endcode\code ?gemv ?symv/?hemv ?trmv \endcode
LU decomposition \n \c EIGEN_USE_LAPACKE \n \c EIGEN_USE_LAPACKE_STRICT \code v1 = m1.lu().solve(v2); \endcode\code ?getrf \endcode
Cholesky decomposition \n \c EIGEN_USE_LAPACKE \n \c EIGEN_USE_LAPACKE_STRICT \code v1 = m2.selfadjointView().llt().solve(v2); \endcode\code ?potrf \endcode
QR decomposition \n \c EIGEN_USE_LAPACKE \n \c EIGEN_USE_LAPACKE_STRICT \code m1.householderQr(); m1.colPivHouseholderQr(); \endcode\code ?geqrf ?geqp3 \endcode
Singular value decomposition \n \c EIGEN_USE_LAPACKE \code JacobiSVD svd; svd.compute(m1, ComputeThinV); \endcode\code ?gesvd \endcode
Eigen-value decompositions \n \c EIGEN_USE_LAPACKE \n \c EIGEN_USE_LAPACKE_STRICT \code EigenSolver es(m1); ComplexEigenSolver ces(m1); SelfAdjointEigenSolver saes(m1+m1.transpose()); GeneralizedSelfAdjointEigenSolver gsaes(m1+m1.transpose(),m2+m2.transpose()); \endcode\code ?gees ?gees ?syev/?heev ?syev/?heev, ?potrf \endcode
Schur decomposition \n \c EIGEN_USE_LAPACKE \n \c EIGEN_USE_LAPACKE_STRICT \code RealSchur schurR(m1); ComplexSchur schurC(m1); \endcode\code ?gees \endcode
In the examples, m1 and m2 are dense matrices and v1 and v2 are dense vectors. */ }