1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
|
namespace Eigen {
/** \page TutorialCore Tutorial 1/3 - Core features
\ingroup Tutorial
<div class="eimainmenu">\ref index "Overview"
| \b Core \b features
| \ref TutorialGeometry "Geometry"
| \ref TutorialAdvancedLinearAlgebra "Advanced linear algebra"
</div>
\b Table \b of \b contents
- \ref TutorialCoreSimpleExampleFixedSize
- \ref TutorialCoreSimpleExampleDynamicSize
- \ref TutorialCoreMatrixTypes
- \ref TutorialCoreMatrixInitialization
- \ref TutorialCoreBasicLinearAlgebra
- \ref TutorialCoreReductions
- \ref TutorialCoreSubMatrix
- \ref TutorialCoreMatrixTransformations
- \ref TutorialCoreTriangularMatrix
- \ref TutorialCorePerformance
\n
<hr>
\section TutorialCoreSimpleExampleFixedSize Simple example with fixed-size matrices and vectors
By fixed-size, we mean that the number of rows and columns are known at compile-time. In this case, Eigen avoids dynamic memory allocation and unroll loops. This is useful for very small sizes (typically up to 4x4).
<table class="tutorial_code"><tr><td>
\include Tutorial_simple_example_fixed_size.cpp
</td>
<td>
output:
\include Tutorial_simple_example_fixed_size.out
</td></tr></table>
<a href="#" class="top">top</a>\section TutorialCoreSimpleExampleDynamicSize Simple example with dynamic-size matrices and vectors
Dynamic-size means that the number of rows and columns are not known at compile-time. In this case, they are stored as runtime variables and the arrays are dynamically allocated.
<table class="tutorial_code"><tr><td>
\include Tutorial_simple_example_dynamic_size.cpp
</td>
<td>
output:
\include Tutorial_simple_example_dynamic_size.out
</td></tr></table>
<a href="#" class="top">top</a>\section TutorialCoreMatrixTypes Matrix and vector types
In Eigen, all kinds of dense matrices and vectors are represented by the template class Matrix. In most cases you can simply use one of the \ref matrixtypedefs "several convenient typedefs".
The template class Matrix takes a number of template parameters, but for now it is enough to understand the 3 first ones (and the others can then be left unspecified):
\code Matrix<Scalar, RowsAtCompileTime, ColsAtCompileTime> \endcode
\li \c Scalar is the scalar type, i.e. the type of the coefficients. That is, if you want a vector of floats, choose \c float here.
\li \c RowsAtCompileTime and \c ColsAtCompileTime are the number of rows and columns of the matrix as known at compile-time.
For example, \c Vector3d is a typedef for \code Matrix<double, 3, 1> \endcode
What if the matrix has dynamic-size i.e. the number of rows or cols isn't known at compile-time? Then use the special value Eigen::Dynamic. For example, \c VectorXd is a typedef for \code Matrix<double, Dynamic, 1> \endcode
<a href="#" class="top">top</a>\section TutorialCoreMatrixInitialization Matrix and vector creation and initialization
\subsection TutorialPredefMat PredefinedMatrix
Eigen offers several methods to create or set matrices with coefficients equals to either a constant value, the identity matrix or even random values:
<table class="tutorial_code">
<tr>
<td>Fixed-size matrix or vector</td>
<td>Dynamic-size matrix</td>
<td>Dynamic-size vector</td>
</tr>
<tr>
<td>
\code
Matrix3f x;
x = Matrix3f::Zero();
x = Matrix3f::Ones();
x = Matrix3f::Constant(value);
x = Matrix3f::Identity();
x = Matrix3f::Random();
x.setZero();
x.setOnes();
x.setIdentity();
x.setConstant(value);
x.setRandom();
\endcode
</td>
<td>
\code
MatrixXf x;
x = MatrixXf::Zero(rows, cols);
x = MatrixXf::Ones(rows, cols);
x = MatrixXf::Constant(rows, cols, value);
x = MatrixXf::Identity(rows, cols);
x = MatrixXf::Random(rows, cols);
x.setZero(rows, cols);
x.setOnes(rows, cols);
x.setConstant(rows, cols, value);
x.setIdentity(rows, cols);
x.setRandom(rows, cols);
\endcode
</td>
<td>
\code
VectorXf x;
x = VectorXf::Zero(size);
x = VectorXf::Ones(size);
x = VectorXf::Constant(size, value);
x = VectorXf::Identity(size);
x = VectorXf::Random(size);
x.setZero(size);
x.setOnes(size);
x.setConstant(size, value);
x.setIdentity(size);
x.setRandom(size);
\endcode
</td>
</tr>
<tr><td colspan="3">Basis vectors \link MatrixBase::Unit [details]\endlink</td></tr>
<tr><td>\code
Vector3f::UnixX() // 1 0 0
Vector3f::UnixY() // 0 1 0
Vector3f::UnixZ() // 0 0 1
\endcode</td><td></td><td>\code
VectorXf::Unit(size,i)
VectorXf::Unit(4,1) == Vector4f(0,1,0,0)
== Vector4f::UnitY()
\endcode
</table>
Here is an usage example:
<table class="tutorial_code"><tr><td>
\code
cout << MatrixXf::Constant(2, 3, sqrt(2)) << endl;
RowVector3i v;
v.setConstant(6);
cout << "v = " << v << endl;
\endcode
</td>
<td>
output:
\code
1.41 1.41 1.41
1.41 1.41 1.41
v = 6 6 6
\endcode
</td></tr></table>
\subsection TutorialMap Map
Any memory buffer can be mapped as an Eigen's expression:
<table class="tutorial_code"><tr><td>
\code
std::vector<float> stlarray(10);
Map<VectorXf>(&stlarray[0], stlarray.size()).setOnes();
int data[4] = 1, 2, 3, 4;
Matrix2i mat2x2(data);
MatrixXi mat2x2 = Map<Matrix2i>(data);
MatrixXi mat2x2 = Map<MatrixXi>(data,2,2);
\endcode
</td></tr></table>
\subsection TutorialCommaInit CommaInitializer
Eigen also offer a comma initializer syntax which allows you to set all the coefficients of a matrix to specific values:
<table class="tutorial_code"><tr><td>
\include Tutorial_commainit_01.cpp
</td>
<td>
output:
\verbinclude Tutorial_commainit_01.out
</td></tr></table>
Feel the above example boring ? Look at the following example where the matrix is set per block:
<table class="tutorial_code"><tr><td>
\include Tutorial_commainit_02.cpp
</td>
<td>
output:
\verbinclude Tutorial_commainit_02.out
</td></tr></table>
<span class="note">\b Side \b note: here .finished() is used to get the actual matrix object once the comma initialization
of our temporary submatrix is done. Note that despite the appearant complexity of such an expression
Eigen's comma initializer usually yields to very optimized code without any overhead.</span>
<a href="#" class="top">top</a>\section TutorialCoreBasicLinearAlgebra Basic Linear Algebra
In short all mathematically well defined operators can be used right away as in the following example:
\code
mat4 -= mat1*1.5 + mat2 * mat3/4;
\endcode
which includes two matrix scalar products ("mat1*1.5" and "mat3/4"), a matrix-matrix product ("mat2 * mat3/4"),
a matrix addition ("+") and substraction with assignment ("-=").
<table class="tutorial_code">
<tr><td>
matrix/vector product</td><td>\code
col2 = mat1 * col1;
row2 = row1 * mat1; row1 *= mat1;
mat3 = mat1 * mat2; mat3 *= mat1; \endcode
</td></tr>
<tr><td>
add/subtract</td><td>\code
mat3 = mat1 + mat2; mat3 += mat1;
mat3 = mat1 - mat2; mat3 -= mat1;\endcode
</td></tr>
<tr><td>
scalar product</td><td>\code
mat3 = mat1 * s1; mat3 = s1 * mat1; mat3 *= s1;
mat3 = mat1 / s1; mat3 /= s1;\endcode
</td></tr>
<tr><td>
\link MatrixBase::dot() dot product \endlink (inner product)</td><td>\code
scalar = vec1.dot(vec2);\endcode
</td></tr>
<tr><td>
outer product</td><td>\code
mat = vec1 * vec2.transpose();\endcode
</td></tr>
<tr><td>
\link MatrixBase::cross() cross product \endlink</td><td>\code
#include <Eigen/Geometry>
vec3 = vec1.cross(vec2);\endcode</td></tr>
</table>
In Eigen only mathematically well defined operators can be used right away,
but don't worry, thanks to the \link Cwise .cwise() \endlink operator prefix,
Eigen's matrices also provide a very powerful numerical container supporting
most common coefficient wise operators:
<table class="noborder">
<tr><td>
<table class="tutorial_code" style="margin-right:10pt">
<tr><td>Coefficient wise product</td>
<td>\code mat3 = mat1.cwise() * mat2; \endcode
</td></tr>
<tr><td>
Add a scalar to all coefficients</td><td>\code
mat3 = mat1.cwise() + scalar;
mat3.cwise() += scalar;
mat3.cwise() -= scalar;
\endcode
</td></tr>
<tr><td>
Coefficient wise division</td><td>\code
mat3 = mat1.cwise() / mat2; \endcode
</td></tr>
<tr><td>
Coefficient wise reciprocal</td><td>\code
mat3 = mat1.cwise().inverse(); \endcode
</td></tr>
<tr><td>
Coefficient wise comparisons \n
(support all operators)</td><td>\code
mat3 = mat1.cwise() < mat2;
mat3 = mat1.cwise() <= mat2;
mat3 = mat1.cwise() > mat2;
etc.
\endcode
</td></tr></table>
</td>
<td><table class="tutorial_code">
<tr><td>
Trigo:\n sin, cos, tan</td><td>\code
mat3 = mat1.cwise().sin();
etc.
\endcode
</td></tr>
<tr><td>
Power:\n pow, square, cube,\n sqrt, exp, log</td><td>\code
mat3 = mat1.cwise().square();
mat3 = mat1.cwise().pow(5);
mat3 = mat1.cwise().log();
etc.
\endcode
</td></tr>
<tr><td>
min, max, absolute value</td><td>\code
mat3 = mat1.cwise().min(mat2);
mat3 = mat1.cwise().max(mat2);
mat3 = mat1.cwise().abs(mat2);
mat3 = mat1.cwise().abs2(mat2);
\endcode</td></tr>
</table>
</td></tr></table>
<span class="note">\b Side \b note: If you feel the \c .cwise() syntax is too verbose for your taste and don't bother to have non mathematical operator directly available feel free to extend MatrixBase as described \ref ExtendingMatrixBase "here".</span>
<a href="#" class="top">top</a>\section TutorialCoreReductions Reductions
Eigen provides several several reduction methods such as:
\link MatrixBase::minCoeff() minCoeff() \endlink, \link MatrixBase::maxCoeff() maxCoeff() \endlink,
\link MatrixBase::sum() sum() \endlink, \link MatrixBase::trace() trace() \endlink,
\link MatrixBase::norm() norm() \endlink, \link MatrixBase::norm2() norm2() \endlink,
\link MatrixBase::all() all() \endlink,and \link MatrixBase::any() any() \endlink.
All reduction operations can be done matrix-wise,
\link MatrixBase::colwise() column-wise \endlink or
\link MatrixBase::rowwise() row-wise \endlink. Usage example:
<table class="tutorial_code">
<tr><td rowspan="3" style="border-right-style:dashed">\code
5 3 1
mat = 2 7 8
9 4 6 \endcode
</td> <td>\code mat.minCoeff(); \endcode</td><td>\code 1 \endcode</td></tr>
<tr><td>\code mat.colwise().minCoeff(); \endcode</td><td>\code 2 3 1 \endcode</td></tr>
<tr><td>\code mat.rowwise().minCoeff(); \endcode</td><td>\code
1
2
4
\endcode</td></tr>
</table>
<span class="note">\b Side \b note: The all() and any() functions are especially useful in combinaison with coeff-wise comparison operators (\ref CwiseAll "example").</span>
<a href="#" class="top">top</a>\section TutorialCoreSubMatrix Sub matrices
Read-write access to a \link MatrixBase::col(int) column \endlink
or a \link MatrixBase::row(int) row \endlink of a matrix:
\code
mat1.row(i) = mat2.col(j);
mat1.col(j1).swap(mat1.col(j2));
\endcode
Read-write access to sub-vectors:
<table class="tutorial_code">
<tr>
<td>Default versions</td>
<td>Optimized versions when the size is known at compile time</td></tr>
<td></td>
<tr><td>\code vec1.start(n)\endcode</td><td>\code vec1.start<n>()\endcode</td><td>the first \c n coeffs </td></tr>
<tr><td>\code vec1.end(n)\endcode</td><td>\code vec1.end<n>()\endcode</td><td>the last \c n coeffs </td></tr>
<tr><td>\code vec1.block(pos,n)\endcode</td><td>\code vec1.block<n>(pos)\endcode</td>
<td>the \c size coeffs in the range [\c pos : \c pos + \c n [</td></tr>
</table>
Read-write access to sub-matrices:
<table class="tutorial_code">
<tr><td>Default versions</td>
<td>Optimized versions when the size is known at compile time</td><td></td></tr>
<tr>
<td>\code mat1.block(i,j,rows,cols)\endcode
\link MatrixBase::block(int,int,int,int) (more) \endlink</td>
<td>\code mat1.block<rows,cols>(i,j)\endcode
\link MatrixBase::block(int,int) (more) \endlink</td>
<td>the \c rows x \c cols sub-matrix starting from position (\c i,\c j) </td>
</tr>
<tr>
<td>\code
mat1.corner(TopLeft,rows,cols)
mat1.corner(TopRight,rows,cols)
mat1.corner(BottomLeft,rows,cols)
mat1.corner(BottomRight,rows,cols)\endcode
\link MatrixBase::corner(CornerType,int,int) (more) \endlink</td>
<td>\code
mat1.corner<rows,cols>(TopLeft)
mat1.corner<rows,cols>(TopRight)
mat1.corner<rows,cols>(BottomLeft)
mat1.corner<rows,cols>(BottomRight)\endcode
\link MatrixBase::corner(CornerType) (more) \endlink</td>
<td>the \c rows x \c cols sub-matrix \n taken in one of the four corners</td></tr>
<tr>
<td>\code
vec1 = mat1.diagonal();
mat1.diagonal() = vec1;
\endcode
\link MatrixBase::diagonal() (more) \endlink</td></td>
<td></td>
</table>
<a href="#" class="top">top</a>\section TutorialCoreMatrixTransformations Matrix transformations
<table class="tutorial_code">
<tr><td>
\link MatrixBase::transpose() transposition \endlink (read-write)</td><td>\code
mat3 = mat1.transpose() * mat2;
mat3.transpose() = mat1 * mat2.transpose();
\endcode
</td></tr>
<tr><td>
\link MatrixBase::adjoint() adjoint \endlink (read only)\n</td><td>\code
mat3 = mat1.adjoint() * mat2;
mat3 = mat1.conjugate().transpose() * mat2;
\endcode
</td></tr>
<tr><td>
\link MatrixBase::asDiagonal() make a diagonal matrix \endlink from a vector \n
\b Note: this product is automatically optimized !</td><td>\code
mat3 = mat1 * vec2.asDiagonal();\endcode
</td></tr>
<tr><td>
\link MatrixBase::minor() minor \endlink (read-write)</td><td>\code
mat4x4.minor(i,j) = mat3x3;
mat3x3 = mat4x4.minor(i,j);\endcode
</td></tr>
</table>
<a href="#" class="top">top</a>\section TutorialCoreTriangularMatrix Dealing with triangular matrices
todo
<a href="#" class="top">top</a>\section TutorialCorePerformance Notes on performances
<table class="tutorial_code">
<tr><td>\code
m4 = m4 * m4;\endcode</td><td>
auto-evaluates so no aliasing problem (performance penalty is low)</td></tr>
<tr><td>\code
Matrix4f other = (m4 * m4).lazy();\endcode</td><td>
forces lazy evaluation</td></tr>
<tr><td>\code
m4 = m4 + m4;\endcode</td><td>
here Eigen goes for lazy evaluation, as with most expressions</td></tr>
<tr><td>\code
m4 = -m4 + m4 + 5 * m4;\endcode</td><td>
same here, Eigen chooses lazy evaluation for all that.</td></tr>
<tr><td>\code
m4 = m4 * (m4 + m4);\endcode</td><td>
here Eigen chooses to first evaluate m4 + m4 into a temporary.
indeed, here it is an optimization to cache this intermediate result.</td></tr>
<tr><td>\code
m3 = m3 * m4.block<3,3>(1,1);\endcode</td><td>
here Eigen chooses \b not to evaluate block() into a temporary
because accessing coefficients of that block expression is not more costly than accessing
coefficients of a plain matrix.</td></tr>
<tr><td>\code
m4 = m4 * m4.transpose();\endcode</td><td>
same here, lazy evaluation of the transpose.</td></tr>
<tr><td>\code
m4 = m4 * m4.transpose().eval();\endcode</td><td>
forces immediate evaluation of the transpose</td></tr>
</table>
*/
/** \page TutorialGeometry Tutorial 2/3 - Geometry
\ingroup Tutorial
<div class="eimainmenu">\ref index "Overview"
| \ref TutorialCore "Core features"
| \b Geometry
| \ref TutorialAdvancedLinearAlgebra "Advanced linear algebra"
</div>
In this tutorial chapter we will shortly introduce the many possibilities offered by the \ref GeometryModule "geometry module",
namely 2D and 3D rotations and affine transformations.
\b Table \b of \b contents
- \ref TutorialGeoRotations
- \ref TutorialGeoTransformation
<a href="#" class="top">top</a>\section TutorialGeoRotations 2D and 3D Rotations
\subsection TutorialGeoRotationTypes Rotation types
<table class="tutorial_code">
<tr><td>Rotation type</td><td>Typical initialization code</td><td>Recommended usage</td></tr>
<tr><td>2D rotation from an angle</td><td>\code
Rotation2D<float> rot2(angle_in_radian);\endcode</td><td></td></tr>
<tr><td>2D rotation matrix</td><td>\code
Matrix2f rotmat2 = Rotation2Df(angle_in_radian);\endcode</td><td></td></tr>
<tr><td>3D rotation as an angle + axis</td><td>\code
AngleAxis<float> aa(angle_in_radian, Vector3f(ax,ay,az));\endcode</td><td></td></tr>
<tr><td>3D rotation as a quaternion</td><td>\code
Quaternion<float> q = AngleAxis<float>(angle_in_radian, axis);\endcode</td><td></td></tr>
<tr><td>3D rotation matrix</td><td>\code
Matrix3f rotmat3 = AngleAxis<float>(angle_in_radian, axis);\endcode</td><td></td></tr>
</table>
To transform more than a single vector the prefered representations are rotation matrices,
for other usage Rotation2D and Quaternion are the representations of choice as they are
more compact, fast and stable. AngleAxis are only useful to create other rotation objects.
\subsection TutorialGeoCommonRotationAPI Common API across rotation types
To some extent, Eigen's \ref Geometry_Module "geometry module" allows you to write
generic algorithms working on both 2D and 3D rotations of any of the five above types.
The following operation are supported:
<table class="tutorial_code">
<tr><td>Convertion from and to any types (of same space dimension)</td><td>\code
RotType2 a = RotType1();\endcode</td></tr>
<tr><td>Concatenation of two rotations</td><td>\code
rot3 = rot1 * rot2;\endcode</td></tr>
<tr><td>Apply the rotation to a vector</td><td>\code
vec2 = rot1 * vec1;\endcode</td></tr>
<tr><td>Get the inverse rotation \n (not always the most effient choice)</td><td>\code
rot2 = rot1.inverse();\endcode</td></tr>
<tr><td>Spherical interpolation \n (Rotation2D and Quaternion only)</td><td>\code
rot3 = rot1.slerp(alpha,rot2);\endcode</td></tr>
</table>
\subsection TutorialGeoEulerAngles Euler angles
<table class="tutorial_code">
<tr><td style="max-width:30em;">
Euler angles might be convenient to create rotation object.
Since there exist 24 differents convensions, they are one
the ahand pretty confusing to use. This example shows how
to create a rotation matrix according to the 2-1-2 convention.</td><td>\code
Matrix3f m;
m = AngleAxisf(angle1, Vector3f::UnitZ())
* AngleAxisf(angle2, Vector3f::UnitY())
* AngleAxisf(angle3, Vector3f::UnitZ());
\endcode</td></tr>
</table>
<a href="#" class="top">top</a>\section TutorialGeoTransformation Affine transformations
In Eigen we have chosen to not distinghish between points and vectors such that all points are
actually represented by displacement vector from the origine (pt \~ pt-0). With that in mind,
real points and vector distinguish when the rotation is applied.
<table class="tutorial_code">
<tr><td></td><td>\b 3D </td><td>\b 2D </td></tr>
<tr><td>Creation \n <span class="note">rot2D can also be an angle in radian</span></td><td>\code
Transform3f t;
t.fromPositionOrientationScale(
pos,any_3D_rotation,Vector3f(sx,sy,sz)); \endcode</td><td>\code
Transform2f t;
t.fromPositionOrientationScale(
pos,any_2D_rotation,Vector2f(sx,sy)); \endcode</td></tr>
<tr><td>Apply the transformation to a \b point </td><td>\code
Vector3f p1, p2;
p2 = t * p1;\endcode</td><td>\code
Vector2f p1, p2;
p2 = t * p1;\endcode</td></tr>
<tr><td>Apply the transformation to a \b vector </td><td>\code
Vector3f v1, v2;
v2 = t.linear() * v1;\endcode</td><td>\code
Vector2f v1, v2;
v2 = t.linear() * v1;\endcode</td></tr>
<tr><td>Concatenate two transformations</td><td>\code
t3 = t1 * t2;\endcode</td><td>\code
t3 = t1 * t2;\endcode</td></tr>
<tr><td>OpenGL compatibility</td><td>\code
glLoadMatrixf(t.data());\endcode</td><td>\code
Transform3f aux(Transform3f::Identity);
aux.linear().corner<2,2>(TopLeft) = t.linear();
aux.translation().start<2>() = t.translation();
glLoadMatrixf(aux.data());\endcode</td></tr>
<tr><td colspan="3">\b Component \b accessors</td></tr>
<tr><td>translation part</td><td>\code
t.translation() = vec3;
vec3 = t.translation();
\endcode</td><td>\code
t.translation() = vec2;
vec2 = t.translation();
\endcode</td></tr>
<tr><td>linear part</td><td>\code
t.linear() = mat3x3;
mat3x3 = t.linear();
\endcode</td><td>\code
t.linear() = mat2x2;
mat2x2 = t.linear();
\endcode</td></tr>
<tr><td colspan="3">\b Editing \b shortcuts</td></tr>
<tr><td>Applies a translation</td><td>\code
t.translate(Vector3f(tx, ty, tz));
t.pretranslate(Vector3f(tx, ty, tz));
\endcode</td><td>\code
t.translate(Vector2f(tx, ty));
t.pretranslate(Vector2f(tx, ty));
\endcode</td></tr>
<tr><td>Applies a rotation \n <span class="note">rot2D can also be an angle in radian</span></td><td>\code
t.rotate(rot3D);
t.prerotate(rot3D);
\endcode</td><td>\code
t.rotate(rot2D);
t.prerotate(rot2D);
\endcode</td></tr>
<tr><td>Applies a scaling</td><td>\code
t.scale(Vector3f(sx, sy, sz));
t.scale(Vector3f::Constant(s));
t.prescale(Vector3f(sx, sy, sz));
\endcode</td><td>\code
t.scale(Vector2f(tx, ty));
t.scale(Vector2f::Constant(s));
t.prescale(Vector2f(tx, ty));
\endcode</td></tr>
<tr><td>Applies a shear transformation \n(2D only)</td><td></td><td>\code
t.shear(sx,sy);
t.preshear(sx,sy);
\endcode</td></tr>
</table>
*/
/** \page TutorialAdvancedLinearAlgebra Tutorial 3/3 - Advanced linear algebra
\ingroup Tutorial
<div class="eimainmenu">\ref index "Overview"
| \ref TutorialCore "Core features"
| \ref TutorialGeometry "Geometry"
| \b Advanced \b linear \b algebra
</div>
\b Table \b of \b contents
- \ref TutorialAdvLinearSolvers
- \ref TutorialAdvLU
- \ref TutorialAdvCholesky
- \ref TutorialAdvQR
- \ref TutorialAdvEigenProblems
\section TutorialAdvLinearSolvers Solving linear problems
todo
<a href="#" class="top">top</a>\section TutorialAdvLU LU
todo
<a href="#" class="top">top</a>\section TutorialAdvCholesky Cholesky
todo
<a href="#" class="top">top</a>\section TutorialAdvQR QR
todo
<a href="#" class="top">top</a>\section TutorialAdvEigenProblems Eigen value problems
todo
*/
}
|