1namespace Eigen {
2
3/** \page Eigen2ToEigen3 Porting from Eigen2 to Eigen3
4
5<div class="bigwarning">Eigen2 support is deprecated in Eigen 3.2.x and it will be removed in Eigen 3.3.</div>
6
7This page lists the most important API changes between Eigen2 and Eigen3,
8and gives tips to help porting your application from Eigen2 to Eigen3.
9
10\eigenAutoToc
11
12\section CompatibilitySupport Eigen2 compatibility support
13
14In order to ease the switch from Eigen2 to Eigen3, Eigen3 features \subpage Eigen2SupportModes "Eigen2 support modes".
15
16The quick way to enable this is to define the \c EIGEN2_SUPPORT preprocessor token \b before including any Eigen header (typically it should be set in your project options).
17
18A more powerful, \em staged migration path is also provided, which may be useful to migrate larger projects from Eigen2 to Eigen3. This is explained in the \ref Eigen2SupportModes "Eigen 2 support modes" page. 
19
20\section Using The USING_PART_OF_NAMESPACE_EIGEN macro
21
22The USING_PART_OF_NAMESPACE_EIGEN macro has been removed. In Eigen 3, just do:
23\code
24using namespace Eigen;
25\endcode
26
27\section ComplexDot Dot products over complex numbers
28
29This is the single trickiest change between Eigen 2 and Eigen 3. It only affects code using \c std::complex numbers as scalar type.
30
31Eigen 2's dot product was linear in the first variable. Eigen 3's dot product is linear in the second variable. In other words, the Eigen 2 code \code x.dot(y) \endcode is equivalent to the Eigen 3 code \code y.dot(x) \endcode In yet other words, dot products are complex-conjugated in Eigen 3 compared to Eigen 2. The switch to the new convention was commanded by common usage, especially with the notation \f$ x^Ty \f$ for dot products of column-vectors.
32
33\section VectorBlocks Vector blocks
34
35<table class="manual">
36<tr><th>Eigen 2</th><th>Eigen 3</th></th>
37<tr><td>\code
38vector.start(length)
39vector.start<length>()
40vector.end(length)
41vector.end<length>()
42\endcode</td><td>\code
43vector.head(length)
44vector.head<length>()
45vector.tail(length)
46vector.tail<length>()
47\endcode</td></tr>
48</table>
49
50
51\section Corners Matrix Corners
52
53<table class="manual">
54<tr><th>Eigen 2</th><th>Eigen 3</th></th>
55<tr><td>\code
56matrix.corner(TopLeft,r,c)
57matrix.corner(TopRight,r,c)
58matrix.corner(BottomLeft,r,c)
59matrix.corner(BottomRight,r,c)
60matrix.corner<r,c>(TopLeft)
61matrix.corner<r,c>(TopRight)
62matrix.corner<r,c>(BottomLeft)
63matrix.corner<r,c>(BottomRight)
64\endcode</td><td>\code
65matrix.topLeftCorner(r,c)
66matrix.topRightCorner(r,c)
67matrix.bottomLeftCorner(r,c)
68matrix.bottomRightCorner(r,c)
69matrix.topLeftCorner<r,c>()
70matrix.topRightCorner<r,c>()
71matrix.bottomLeftCorner<r,c>()
72matrix.bottomRightCorner<r,c>()
73\endcode</td>
74</tr>
75</table>
76
77Notice that Eigen3 also provides these new convenience methods: topRows(), bottomRows(), leftCols(), rightCols(). See in class DenseBase.
78
79\section CoefficientWiseOperations Coefficient wise operations
80
81In Eigen2, coefficient wise operations which have no proper mathematical definition (as a coefficient wise product)
82were achieved using the .cwise() prefix, e.g.:
83\code a.cwise() * b \endcode
84In Eigen3 this .cwise() prefix has been superseded by a new kind of matrix type called
85Array for which all operations are performed coefficient wise. You can easily view a matrix as an array and vice versa using
86the MatrixBase::array() and ArrayBase::matrix() functions respectively. Here is an example:
87\code
88Vector4f a, b, c;
89c = a.array() * b.array();
90\endcode
91Note that the .array() function is not at all a synonym of the deprecated .cwise() prefix.
92While the .cwise() prefix changed the behavior of the following operator, the array() function performs
93a permanent conversion to the array world. Therefore, for binary operations such as the coefficient wise product,
94both sides must be converted to an \em array as in the above example. On the other hand, when you
95concatenate multiple coefficient wise operations you only have to do the conversion once, e.g.:
96\code
97Vector4f a, b, c;
98c = a.array().abs().pow(3) * b.array().abs().sin();
99\endcode
100With Eigen2 you would have written:
101\code
102c = (a.cwise().abs().cwise().pow(3)).cwise() * (b.cwise().abs().cwise().sin());
103\endcode
104
105\section PartAndExtract Triangular and self-adjoint matrices
106
107In 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:
108
109<table class="manual">
110<tr><th>Eigen 2</th><th>Eigen 3</th></tr>
111<tr><td>\code
112A.part<UpperTriangular>();
113A.part<StrictlyLowerTriangular>(); \endcode</td>
114<td>\code
115A.triangularView<Upper>()
116A.triangularView<StrictlyLower>()\endcode</td></tr>
117<tr><td>\code
118A.extract<UpperTriangular>();
119A.extract<StrictlyLowerTriangular>();\endcode</td>
120<td>\code
121A.triangularView<Upper>()
122A.triangularView<StrictlyLower>()\endcode</td></tr>
123<tr><td>\code
124A.marked<UpperTriangular>();
125A.marked<StrictlyLowerTriangular>();\endcode</td>
126<td>\code
127A.triangularView<Upper>()
128A.triangularView<StrictlyLower>()\endcode</td></tr>
129<tr><td colspan="2"></td></tr>
130<tr><td>\code
131A.part<SelfAdfjoint|UpperTriangular>();
132A.extract<SelfAdfjoint|LowerTriangular>();\endcode</td>
133<td>\code
134A.selfadjointView<Upper>()
135A.selfadjointView<Lower>()\endcode</td></tr>
136<tr><td colspan="2"></td></tr>
137<tr><td>\code
138UpperTriangular
139LowerTriangular
140UnitUpperTriangular
141UnitLowerTriangular
142StrictlyUpperTriangular
143StrictlyLowerTriangular
144\endcode</td><td>\code
145Upper
146Lower
147UnitUpper
148UnitLower
149StrictlyUpper
150StrictlyLower
151\endcode</td>
152</tr>
153</table>
154
155\sa class TriangularView, class SelfAdjointView
156
157\section TriangularSolveInPlace Triangular in-place solving
158
159<table class="manual">
160<tr><th>Eigen 2</th><th>Eigen 3</th></tr>
161<tr><td>\code A.triangularSolveInPlace<XxxTriangular>(Y);\endcode</td><td>\code A.triangularView<Xxx>().solveInPlace(Y);\endcode</td></tr>
162</table>
163
164
165\section Decompositions Matrix decompositions
166
167Some of Eigen 2's matrix decompositions have been renamed in Eigen 3, while some others have been removed and are replaced by other decompositions in Eigen 3.
168
169<table class="manual">
170  <tr>
171    <th>Eigen 2</th>
172    <th>Eigen 3</th>
173    <th>Notes</th>
174  </tr>
175  <tr>
176    <td>LU</td>
177    <td>FullPivLU</td>
178    <td class="alt">See also the new PartialPivLU, it's much faster</td>
179  </tr>
180  <tr>
181    <td>QR</td>
182    <td>HouseholderQR</td>
183    <td class="alt">See also the new ColPivHouseholderQR, it's more reliable</td>
184  </tr>
185  <tr>
186    <td>SVD</td>
187    <td>JacobiSVD</td>
188    <td class="alt">We currently don't have a bidiagonalizing SVD; of course this is planned.</td>
189  </tr>
190  <tr>
191    <td>EigenSolver and friends</td>
192    <td>\code #include<Eigen/Eigenvalues> \endcode </td>
193    <td class="alt">Moved to separate module</td>
194  </tr>
195</table>
196
197\section LinearSolvers Linear solvers
198
199<table class="manual">
200<tr><th>Eigen 2</th><th>Eigen 3</th><th>Notes</th></tr>
201<tr><td>\code A.lu();\endcode</td>
202<td>\code A.fullPivLu();\endcode</td>
203<td class="alt">Now A.lu() returns a PartialPivLU</td></tr>
204<tr><td>\code A.lu().solve(B,&X);\endcode</td>
205<td>\code X = A.lu().solve(B);
206 X = A.fullPivLu().solve(B);\endcode</td>
207<td class="alt">The returned by value is fully optimized</td></tr>
208<tr><td>\code A.llt().solve(B,&X);\endcode</td>
209<td>\code X = A.llt().solve(B);
210 X = A.selfadjointView<Lower>.llt().solve(B);
211 X = A.selfadjointView<Upper>.llt().solve(B);\endcode</td>
212<td class="alt">The returned by value is fully optimized and \n
213the selfadjointView API allows you to select the \n
214triangular part to work on (default is lower part)</td></tr>
215<tr><td>\code A.llt().solveInPlace(B);\endcode</td>
216<td>\code B = A.llt().solve(B);
217 B = A.selfadjointView<Lower>.llt().solve(B);
218 B = A.selfadjointView<Upper>.llt().solve(B);\endcode</td>
219<td class="alt">In place solving</td></tr>
220<tr><td>\code A.ldlt().solve(B,&X);\endcode</td>
221<td>\code X = A.ldlt().solve(B);
222 X = A.selfadjointView<Lower>.ldlt().solve(B);
223 X = A.selfadjointView<Upper>.ldlt().solve(B);\endcode</td>
224<td class="alt">The returned by value is fully optimized and \n
225the selfadjointView API allows you to select the \n
226triangular part to work on</td></tr>
227</table>
228
229\section GeometryModule Changes in the Geometry module
230
231The Geometry module is the one that changed the most. If you rely heavily on it, it's probably a good idea to use the \ref Eigen2SupportModes "Eigen 2 support modes" to perform your migration.
232
233\section Transform The Transform class
234
235In Eigen 2, the Transform class didn't really know whether it was a projective or affine transformation. In Eigen 3, it takes a new \a Mode template parameter, which indicates whether it's \a Projective or \a Affine transform. There is no default value.
236
237The Transform3f (etc) typedefs are no more. In Eigen 3, the Transform typedefs explicitly refer to the \a Projective and \a Affine modes:
238
239<table class="manual">
240<tr><th>Eigen 2</th><th>Eigen 3</th><th>Notes</th></tr>
241<tr>
242  <td> Transform3f </td>
243  <td> Affine3f or Projective3f </td>
244  <td> Of course 3f is just an example here </td>
245</tr>
246</table>
247
248
249\section LazyVsNoalias Lazy evaluation and noalias
250
251In Eigen all operations are performed in a lazy fashion except the matrix products which are always evaluated into a temporary by default.
252In Eigen2, lazy evaluation could be enforced by tagging a product using the .lazy() function. However, in complex expressions it was not
253easy to determine where to put the lazy() function. In Eigen3, the lazy() feature has been superseded by the MatrixBase::noalias() function
254which can be used on the left hand side of an assignment when no aliasing can occur. Here is an example:
255\code
256MatrixXf a, b, c;
257...
258c.noalias() += 2 * a.transpose() * b;
259\endcode
260However, the noalias mechanism does not cover all the features of the old .lazy(). Indeed, in some extremely rare cases,
261it 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,
262just like any other expressions. To this end you can use the MatrixBase::lazyProduct() function, however we strongly discourage you to
263use it unless you are sure of what you are doing, i.e., you have rigourosly measured a speed improvement.
264
265\section AlignMacros Alignment-related macros
266
267The EIGEN_ALIGN_128 macro has been renamed to EIGEN_ALIGN16. Don't be surprised, it's just that we switched to counting in bytes ;-)
268
269The EIGEN_DONT_ALIGN option still exists in Eigen 3, but it has a new cousin: EIGEN_DONT_ALIGN_STATICALLY. It allows to get rid of all static alignment issues while keeping alignment of dynamic-size heap-allocated arrays, thus keeping vectorization for dynamic-size objects.
270
271\section AlignedMap Aligned Map objects
272
273A common issue with Eigen 2 was that when mapping an array with Map, there was no way to tell Eigen that your array was aligned. There was a ForceAligned option but it didn't mean that; it was just confusing and has been removed.
274
275New in Eigen3 is the #Aligned option. See the documentation of class Map. Use it like this:
276\code
277Map<Vector4f, Aligned> myMappedVector(some_aligned_array);
278\endcode
279There also are related convenience static methods, which actually are the preferred way as they take care of such things as constness:
280\code
281result = Vector4f::MapAligned(some_aligned_array);
282\endcode
283
284\section StdContainers STL Containers
285
286In Eigen2, <tt>#include<Eigen/StdVector></tt> tweaked std::vector to automatically align elements. The problem was that that was quite invasive. In Eigen3, we only override standard behavior if you use Eigen::aligned_allocator<T> as your allocator type. So for example, if you use std::vector<Matrix4f>, you need to do the following change (note that aligned_allocator is under namespace Eigen):
287
288<table class="manual">
289<tr><th>Eigen 2</th><th>Eigen 3</th></tr>
290<tr>
291  <td> \code std::vector<Matrix4f> \endcode </td>
292  <td> \code std::vector<Matrix4f, aligned_allocator<Matrix4f> > \endcode </td>
293</tr>
294</table>
295
296\section eiPrefix Internal ei_ prefix
297
298In Eigen2, global internal functions and structures were prefixed by \c ei_. In Eigen3, they all have been moved into the more explicit \c internal namespace. So, e.g., \c ei_sqrt(x) now becomes \c internal::sqrt(x). Of course it is not recommended to rely on Eigen's internal features.
299
300
301
302*/
303
304}
305