mirror of
https://gitlab.com/libeigen/eigen.git
synced 2024-12-21 07:19:46 +08:00
212 lines
7.3 KiB
Plaintext
212 lines
7.3 KiB
Plaintext
namespace Eigen {
|
|
|
|
/** \page Eigen2ToEigen3 Porting from Eigen2 to Eigen3
|
|
|
|
The goals of this page is to enumerate the API changes between Eigen2 and Eigen3,
|
|
and to help porting an application from Eigen2 to Eigen3.
|
|
|
|
\b Table \b of \b contents
|
|
- \ref CompatibilitySupport
|
|
- \ref VectorBlocks
|
|
- \ref CoefficientWiseOperations
|
|
- \ref PartAndExtract
|
|
- \ref TriangularSolveInPlace
|
|
- \ref Using
|
|
- \ref Corners
|
|
- \ref LazyVsNoalias
|
|
|
|
\section CompatibilitySupport Eigen2 compatibility support
|
|
|
|
In order to ease the switch from Eigen2 to Eigen3, Eigen3 features a compatibility mode which can be enabled by defining the EIGEN2_SUPPORT preprocessor token \b before including any Eigen header (typically it should be set in your project options).
|
|
|
|
\section VectorBlocks Vector blocks
|
|
|
|
<table>
|
|
<tr><td>Eigen 2</td><td>Eigen 3</td></tr>
|
|
<tr><td>\code
|
|
vector.start(length)
|
|
vector.start<length>()
|
|
vector.end(length)
|
|
vector.end<length>()
|
|
\endcode</td><td>\code
|
|
vector.head(length)
|
|
vector.head<length>()
|
|
vector.tail(length)
|
|
vector.tail<length>()
|
|
\endcode</td></tr>
|
|
</table>
|
|
|
|
|
|
\section Corners Matrix Corners
|
|
|
|
<table>
|
|
<tr><td>Eigen 2</td><td>Eigen 3</td></tr>
|
|
<tr><td>\code
|
|
matrix.corner(TopLeft,r,c)
|
|
matrix.corner(TopRight,r,c)
|
|
matrix.corner(BottomLeft,r,c)
|
|
matrix.corner(BottomRight,r,c)
|
|
matrix.corner<r,c>(TopLeft)
|
|
matrix.corner<r,c>(TopRight)
|
|
matrix.corner<r,c>(BottomLeft)
|
|
matrix.corner<r,c>(BottomRight)
|
|
\endcode</td><td>\code
|
|
matrix.topLeftCorner(r,c)
|
|
matrix.topRightCorner(r,c)
|
|
matrix.bottomLeftCorner(r,c)
|
|
matrix.bottomRightCorner(r,c)
|
|
matrix.topLeftCorner<r,c>()
|
|
matrix.topRightCorner<r,c>()
|
|
matrix.bottomLeftCorner<r,c>()
|
|
matrix.bottomRightCorner<r,c>()
|
|
\endcode</td>
|
|
</tr>
|
|
</table>
|
|
|
|
Notice that Eigen3 also provides these new convenience methods: topRows(), bottomRows(), leftCols(), rightCols(). See in class DenseBase.
|
|
|
|
\section CoefficientWiseOperations Coefficient wise operations
|
|
|
|
In Eigen2, coefficient wise operations which have no proper mathematical definition (as a coefficient wise product)
|
|
were achieved using the .cwise() prefix, e.g.:
|
|
\code a.cwise() * b \endcode
|
|
In Eigen3 this .cwise() prefix has been superseded by a new kind of matrix type called
|
|
Array for which all operations are performed coefficient wise. You can easily view a matrix as an array and vice versa using
|
|
the MatrixBase::array() and ArrayBase::matrix() functions respectively. Here is an example:
|
|
\code
|
|
Vector4f a, b, c;
|
|
c = a.array() * b.array();
|
|
\endcode
|
|
Note that the .array() function is not at all a synonym of the deprecated .cwise() prefix.
|
|
While the .cwise() prefix changed the behavior of the following operator, the array() function performs
|
|
a permanent conversion to the array world. Therefore, for binary operations such as the coefficient wise product,
|
|
both sides must be converted to an \em array as in the above example. On the other hand, when you
|
|
concatenate multiple coefficient wise operations you only have to do the conversion once, e.g.:
|
|
\code
|
|
Vector4f a, b, c;
|
|
c = a.array().abs().pow(3) * b.array().abs().sin();
|
|
\endcode
|
|
With Eigen2 you would have written:
|
|
\code
|
|
c = (a.cwise().abs().cwise().pow(3)).cwise() * (b.cwise().abs().cwise().sin());
|
|
\endcode
|
|
|
|
\section PartAndExtract Triangular and self-adjoint matrices
|
|
|
|
In Eigen 2 you had to play with the part, extract, and marked functions to deal with triangular and selfadjoint matrices. In Eigen 3, all these functions have been removed in favor of the concept of \em views:
|
|
|
|
<table>
|
|
<tr><td>Eigen 2</td><td>Eigen 3</td></tr>
|
|
<tr><td>\code
|
|
A.part<UpperTriangular>();
|
|
A.part<StrictlyLowerTriangular>(); \endcode</td>
|
|
<td>\code
|
|
A.triangularView<Upper>()
|
|
A.triangularView<StrictlyLower>()\endcode</td></tr>
|
|
<tr><td>\code
|
|
A.extract<UpperTriangular>();
|
|
A.extract<StrictlyLowerTriangular>();\endcode</td>
|
|
<td>\code
|
|
A.triangularView<Upper>()
|
|
A.triangularView<StrictlyLower>()\endcode</td></tr>
|
|
<tr><td>\code
|
|
A.marked<UpperTriangular>();
|
|
A.marked<StrictlyLowerTriangular>();\endcode</td>
|
|
<td>\code
|
|
A.triangularView<Upper>()
|
|
A.triangularView<StrictlyLower>()\endcode</td></tr>
|
|
<tr><td colspan="2"></td></tr>
|
|
<tr><td>\code
|
|
A.part<SelfAdfjoint|UpperTriangular>();
|
|
A.extract<SelfAdfjoint|LowerTriangular>();\endcode</td>
|
|
<td>\code
|
|
A.selfadjointView<Upper>()
|
|
A.selfadjointView<Lower>()\endcode</td></tr>
|
|
<tr><td colspan="2"></td></tr>
|
|
<tr><td>\code
|
|
UpperTriangular
|
|
LowerTriangular
|
|
UnitUpperTriangular
|
|
UnitLowerTriangular
|
|
StrictlyUpperTriangular
|
|
StrictlyLowerTriangular
|
|
\endcode</td><td>\code
|
|
Upper
|
|
Lower
|
|
UnitUpper
|
|
UnitLower
|
|
StrictlyUpper
|
|
StrictlyLower
|
|
\endcode</td>
|
|
</tr>
|
|
</table>
|
|
|
|
\sa class TriangularView, class SelfAdjointView
|
|
|
|
\section TriangularSolveInPlace Triangular in-place solving
|
|
|
|
<table>
|
|
<tr><td>Eigen 2</td><td>Eigen 3</td></tr>
|
|
<tr><td>\code A.triangularSolveInPlace<XxxTriangular>(Y);\endcode</td><td>\code A.triangularView<Xxx>().solveInPlace(Y);\endcode</td></tr>
|
|
</table>
|
|
|
|
\section LinearSolvers Linear solvers
|
|
|
|
<table>
|
|
<tr><td>Eigen 2</td><td>Eigen 3</td><td>Notes</td></tr>
|
|
<tr><td>\code A.lu();\endcode</td>
|
|
<td>\code A.fullPivLu();\endcode</td>
|
|
<td>Now A.lu() returns a PartialPivLU</td></tr>
|
|
<tr><td>\code A.lu().solve(B,&X);\endcode</td>
|
|
<td>\code X = A.lu().solve(B);
|
|
X = A.fullPivLu().solve(B);\endcode</td>
|
|
<td>The returned by value is fully optimized</td></tr>
|
|
<tr><td>\code A.llt().solve(B,&X);\endcode</td>
|
|
<td>\code X = A.llt().solve(B);
|
|
X = A.selfadjointView<Lower>.llt().solve(B);
|
|
X = A.selfadjointView<Upper>.llt().solve(B);\endcode</td>
|
|
<td>The returned by value is fully optimized and \n
|
|
the selfadjointView API allows you to select the \n
|
|
triangular part to work on (default is lower part)</td></tr>
|
|
<tr><td>\code A.llt().solveInPlace(B);\endcode</td>
|
|
<td>\code B = A.llt().solve(B);
|
|
B = A.selfadjointView<Lower>.llt().solve(B);
|
|
B = A.selfadjointView<Upper>.llt().solve(B);\endcode</td>
|
|
<td>In place solving</td></tr>
|
|
<tr><td>\code A.ldlt().solve(B,&X);\endcode</td>
|
|
<td>\code X = A.ldlt().solve(B);
|
|
X = A.selfadjointView<Lower>.ldlt().solve(B);
|
|
X = A.selfadjointView<Upper>.ldlt().solve(B);\endcode</td>
|
|
<td>The returned by value is fully optimized and \n
|
|
the selfadjointView API allows you to select the \n
|
|
triangular part to work on</td></tr>
|
|
</table>
|
|
|
|
\section Using The USING_PART_OF_NAMESPACE_EIGEN macro
|
|
|
|
The USING_PART_OF_NAMESPACE_EIGEN macro has been removed. In Eigen 3, just do:
|
|
\code
|
|
using namespace Eigen;
|
|
\endcode
|
|
|
|
|
|
\section LazyVsNoalias Lazy evaluation and noalias
|
|
|
|
In Eigen all operations are performed in a lazy fashion except the matrix products which are always evaluated into a temporary by default.
|
|
In Eigen2, lazy evaluation could be enforced by tagging a product using the .lazy() function. However, in complex expressions it was not
|
|
easy to determine where to put the lazy() function. In Eigen3, the lazy() feature has been superseded by the MatrixBase::noalias() function
|
|
which can be used on the left hand side of an assignment when no aliasing can occur. Here is an example:
|
|
\code
|
|
MatrixXf a, b, c;
|
|
...
|
|
c.noalias() += 2 * a.transpose() * b;
|
|
\endcode
|
|
However, the noalias mechanism does not cover all the features of the old .lazy(). Indeed, in some extremely rare cases,
|
|
it might be useful to explicit request for a lay product, i.e., for a product which will be evaluated one coefficient at once, on request,
|
|
just like any other expressions. To this end you can use the MatrixBase::lazyProduct() function, however we strongly discourage you to
|
|
use it unless you are sure of what you are doing, i.e., you have rigourosly measured a speed improvement.
|
|
|
|
*/
|
|
|
|
}
|