mirror of
https://gitlab.com/libeigen/eigen.git
synced 2024-12-15 07:10:37 +08:00
134 lines
6.5 KiB
Plaintext
134 lines
6.5 KiB
Plaintext
/*
|
|
Copyright (c) 2011, Intel Corporation. All rights reserved.
|
|
Copyright (C) 2011-2016 Gael Guennebaud <gael.guennebaud@inria.fr>
|
|
|
|
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 <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.
|
|
|
|
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 <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):
|
|
|
|
\note For Mac users, in order to use the lapack version shipped with the Accelerate framework, you also need the lapacke library.
|
|
Using <a href="https://www.macports.org/">MacPorts</a>, 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
|
|
|
|
<table class="manual">
|
|
<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>
|
|
<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>
|
|
<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>
|
|
</table>
|
|
|
|
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<float>, and \c complex<double>.
|
|
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.
|
|
<table class="manual">
|
|
<tr><th>Functional domain</th><th>Code example</th><th>BLAS/LAPACK routines</th></tr>
|
|
<tr><td>Matrix-matrix operations \n \c EIGEN_USE_BLAS </td><td>\code
|
|
m1*m2.transpose();
|
|
m1.selfadjointView<Lower>()*m2;
|
|
m1*m2.triangularView<Upper>();
|
|
m1.selfadjointView<Lower>().rankUpdate(m2,1.0);
|
|
\endcode</td><td>\code
|
|
?gemm
|
|
?symm/?hemm
|
|
?trmm
|
|
dsyrk/ssyrk
|
|
\endcode</td></tr>
|
|
<tr class="alt"><td>Matrix-vector operations \n \c EIGEN_USE_BLAS </td><td>\code
|
|
m1.adjoint()*b;
|
|
m1.selfadjointView<Lower>()*b;
|
|
m1.triangularView<Upper>()*b;
|
|
\endcode</td><td>\code
|
|
?gemv
|
|
?symv/?hemv
|
|
?trmv
|
|
\endcode</td></tr>
|
|
<tr><td>LU decomposition \n \c EIGEN_USE_LAPACKE \n \c EIGEN_USE_LAPACKE_STRICT </td><td>\code
|
|
v1 = m1.lu().solve(v2);
|
|
\endcode</td><td>\code
|
|
?getrf
|
|
\endcode</td></tr>
|
|
<tr class="alt"><td>Cholesky decomposition \n \c EIGEN_USE_LAPACKE \n \c EIGEN_USE_LAPACKE_STRICT </td><td>\code
|
|
v1 = m2.selfadjointView<Upper>().llt().solve(v2);
|
|
\endcode</td><td>\code
|
|
?potrf
|
|
\endcode</td></tr>
|
|
<tr><td>QR decomposition \n \c EIGEN_USE_LAPACKE \n \c EIGEN_USE_LAPACKE_STRICT </td><td>\code
|
|
m1.householderQr();
|
|
m1.colPivHouseholderQr();
|
|
\endcode</td><td>\code
|
|
?geqrf
|
|
?geqp3
|
|
\endcode</td></tr>
|
|
<tr class="alt"><td>Singular value decomposition \n \c EIGEN_USE_LAPACKE </td><td>\code
|
|
JacobiSVD<MatrixXd> svd;
|
|
svd.compute(m1, ComputeThinV);
|
|
\endcode</td><td>\code
|
|
?gesvd
|
|
\endcode</td></tr>
|
|
<tr><td>Eigen-value decompositions \n \c EIGEN_USE_LAPACKE \n \c EIGEN_USE_LAPACKE_STRICT </td><td>\code
|
|
EigenSolver<MatrixXd> es(m1);
|
|
ComplexEigenSolver<MatrixXcd> ces(m1);
|
|
SelfAdjointEigenSolver<MatrixXd> saes(m1+m1.transpose());
|
|
GeneralizedSelfAdjointEigenSolver<MatrixXd>
|
|
gsaes(m1+m1.transpose(),m2+m2.transpose());
|
|
\endcode</td><td>\code
|
|
?gees
|
|
?gees
|
|
?syev/?heev
|
|
?syev/?heev,
|
|
?potrf
|
|
\endcode</td></tr>
|
|
<tr class="alt"><td>Schur decomposition \n \c EIGEN_USE_LAPACKE \n \c EIGEN_USE_LAPACKE_STRICT </td><td>\code
|
|
RealSchur<MatrixXd> schurR(m1);
|
|
ComplexSchur<MatrixXcd> schurC(m1);
|
|
\endcode</td><td>\code
|
|
?gees
|
|
\endcode</td></tr>
|
|
</table>
|
|
In the examples, m1 and m2 are dense matrices and v1 and v2 are dense vectors.
|
|
|
|
*/
|
|
|
|
}
|