2009-02-20 00:46:35 +08:00
|
|
|
// A simple quickref for Eigen. Add anything that's missing.
|
|
|
|
// Main author: Keir Mierle
|
|
|
|
|
2013-03-12 04:20:12 +08:00
|
|
|
#include <Eigen/Dense>
|
2009-02-20 00:46:35 +08:00
|
|
|
|
|
|
|
Matrix<double, 3, 3> A; // Fixed rows and cols. Same as Matrix3d.
|
|
|
|
Matrix<double, 3, Dynamic> B; // Fixed rows, dynamic cols.
|
|
|
|
Matrix<double, Dynamic, Dynamic> C; // Full dynamic. Same as MatrixXd.
|
|
|
|
Matrix<double, 3, 3, RowMajor> E; // Row major; default is column-major.
|
|
|
|
Matrix3f P, Q, R; // 3x3 float matrix.
|
|
|
|
Vector3f x, y, z; // 3x1 float matrix.
|
|
|
|
RowVector3f a, b, c; // 1x3 float matrix.
|
2013-03-12 04:20:12 +08:00
|
|
|
VectorXd v; // Dynamic column vector of doubles
|
2009-02-20 00:46:35 +08:00
|
|
|
double s;
|
2009-02-18 18:27:18 +08:00
|
|
|
|
2009-02-19 00:35:16 +08:00
|
|
|
// Basic usage
|
|
|
|
// Eigen // Matlab // comments
|
|
|
|
x.size() // length(x) // vector size
|
|
|
|
C.rows() // size(C)(1) // number of rows
|
|
|
|
C.cols() // size(C)(2) // number of columns
|
|
|
|
x(i) // x(i+1) // Matlab is 1-based
|
|
|
|
C(i,j) // C(i+1,j+1) //
|
|
|
|
|
2009-02-20 00:46:35 +08:00
|
|
|
A.resize(4, 4); // Runtime error if assertions are on.
|
|
|
|
B.resize(4, 9); // Runtime error if assertions are on.
|
|
|
|
A.resize(3, 3); // Ok; size didn't change.
|
|
|
|
B.resize(3, 9); // Ok; only dynamic cols changed.
|
|
|
|
|
|
|
|
A << 1, 2, 3, // Initialize A. The elements can also be
|
|
|
|
4, 5, 6, // matrices, which are stacked along cols
|
|
|
|
7, 8, 9; // and then the rows are stacked.
|
|
|
|
B << A, A, A; // B is three horizontally stacked A's.
|
|
|
|
A.fill(10); // Fill A with all 10's.
|
2013-03-12 04:20:12 +08:00
|
|
|
|
|
|
|
// Eigen // Matlab
|
|
|
|
MatrixXd::Identity(rows,cols) // eye(rows,cols)
|
|
|
|
C.setIdentity(rows,cols) // C = eye(rows,cols)
|
|
|
|
MatrixXd::Zero(rows,cols) // zeros(rows,cols)
|
|
|
|
C.setZero(rows,cols) // C = ones(rows,cols)
|
|
|
|
MatrixXd::Ones(rows,cols) // ones(rows,cols)
|
|
|
|
C.setOnes(rows,cols) // C = ones(rows,cols)
|
|
|
|
MatrixXd::Random(rows,cols) // rand(rows,cols)*2-1 // MatrixXd::Random returns uniform random numbers in (-1, 1).
|
|
|
|
C.setRandom(rows,cols) // C = rand(rows,cols)*2-1
|
|
|
|
VectorXd::LinSpace(size,low,high) // linspace(low,high,size)'
|
|
|
|
v.setLinSpace(size,low,high) // v = linspace(low,high,size)'
|
|
|
|
|
2009-02-18 18:27:18 +08:00
|
|
|
|
2009-02-20 00:46:35 +08:00
|
|
|
// Matrix slicing and blocks. All expressions listed here are read/write.
|
|
|
|
// Templated size versions are faster. Note that Matlab is 1-based (a size N
|
|
|
|
// vector is x(1)...x(N)).
|
|
|
|
// Eigen // Matlab
|
2010-01-06 20:56:04 +08:00
|
|
|
x.head(n) // x(1:n)
|
|
|
|
x.head<n>() // x(1:n)
|
|
|
|
x.tail(n) // N = rows(x); x(N - n: N)
|
|
|
|
x.tail<n>() // N = rows(x); x(N - n: N)
|
2009-02-20 00:46:35 +08:00
|
|
|
x.segment(i, n) // x(i+1 : i+n)
|
|
|
|
x.segment<n>(i) // x(i+1 : i+n)
|
|
|
|
P.block(i, j, rows, cols) // P(i+1 : i+rows, j+1 : j+cols)
|
|
|
|
P.block<rows, cols>(i, j) // P(i+1 : i+rows, j+1 : j+cols)
|
2010-04-23 02:11:18 +08:00
|
|
|
P.topLeftCorner(rows, cols) // P(1:rows, 1:cols)
|
|
|
|
P.topRightCorner(rows, cols) // [m n]=size(P); P(1:rows, n-cols+1:n)
|
|
|
|
P.bottomLeftCorner(rows, cols) // [m n]=size(P); P(m-rows+1:m, 1:cols)
|
|
|
|
P.bottomRightCorner(rows, cols) // [m n]=size(P); P(m-rows+1:m, n-cols+1:n)
|
|
|
|
P.topLeftCorner<rows,cols>() // P(1:rows, 1:cols)
|
|
|
|
P.topRightCorner<rows,cols>() // [m n]=size(P); P(1:rows, n-cols+1:n)
|
|
|
|
P.bottomLeftCorner<rows,cols>() // [m n]=size(P); P(m-rows+1:m, 1:cols)
|
|
|
|
P.bottomRightCorner<rows,cols>() // [m n]=size(P); P(m-rows+1:m, n-cols+1:n)
|
2009-02-18 18:27:18 +08:00
|
|
|
|
2009-02-20 00:46:35 +08:00
|
|
|
// Of particular note is Eigen's swap function which is highly optimized.
|
|
|
|
// Eigen // Matlab
|
|
|
|
R.row(i) = P.col(j); // R(i, :) = P(:, i)
|
|
|
|
R.col(j1).swap(mat1.col(j2)); // R(:, [j1 j2]) = R(:, [j2, j1])
|
2009-02-18 18:27:18 +08:00
|
|
|
|
2009-02-20 00:46:35 +08:00
|
|
|
// Views, transpose, etc; all read-write except for .adjoint().
|
|
|
|
// Eigen // Matlab
|
|
|
|
R.adjoint() // R'
|
|
|
|
R.transpose() // R.' or conj(R')
|
|
|
|
R.diagonal() // diag(R)
|
|
|
|
x.asDiagonal() // diag(x)
|
2009-02-18 18:27:18 +08:00
|
|
|
|
2009-02-20 00:46:35 +08:00
|
|
|
// All the same as Matlab, but matlab doesn't have *= style operators.
|
|
|
|
// Matrix-vector. Matrix-matrix. Matrix-scalar.
|
|
|
|
y = M*x; R = P*Q; R = P*s;
|
|
|
|
a = b*M; R = P - Q; R = s*P;
|
|
|
|
a *= M; R = P + Q; R = P/s;
|
|
|
|
R *= Q; R = s*P;
|
|
|
|
R += Q; R *= s;
|
|
|
|
R -= Q; R /= s;
|
2009-02-18 18:27:18 +08:00
|
|
|
|
2013-03-12 04:20:12 +08:00
|
|
|
// Vectorized operations on each element independently
|
2009-02-20 00:46:35 +08:00
|
|
|
// Eigen // Matlab
|
2010-01-06 20:56:04 +08:00
|
|
|
R = P.cwiseProduct(Q); // R = P .* Q
|
|
|
|
R = P.array() * s.array();// R = P .* s
|
|
|
|
R = P.cwiseQuotient(Q); // R = P ./ Q
|
|
|
|
R = P.array() / Q.array();// R = P ./ Q
|
|
|
|
R = P.array() + s.array();// R = P + s
|
|
|
|
R = P.array() - s.array();// R = P - s
|
|
|
|
R.array() += s; // R = R + s
|
|
|
|
R.array() -= s; // R = R - s
|
|
|
|
R.array() < Q.array(); // R < Q
|
|
|
|
R.array() <= Q.array(); // R <= Q
|
|
|
|
R.cwiseInverse(); // 1 ./ P
|
|
|
|
R.array().inverse(); // 1 ./ P
|
|
|
|
R.array().sin() // sin(P)
|
|
|
|
R.array().cos() // cos(P)
|
|
|
|
R.array().pow(s) // P .^ s
|
|
|
|
R.array().square() // P .^ 2
|
|
|
|
R.array().cube() // P .^ 3
|
|
|
|
R.cwiseSqrt() // sqrt(P)
|
|
|
|
R.array().sqrt() // sqrt(P)
|
|
|
|
R.array().exp() // exp(P)
|
|
|
|
R.array().log() // log(P)
|
|
|
|
R.cwiseMax(P) // max(R, P)
|
|
|
|
R.array().max(P.array()) // max(R, P)
|
|
|
|
R.cwiseMin(P) // min(R, P)
|
|
|
|
R.array().min(P.array()) // min(R, P)
|
|
|
|
R.cwiseAbs() // abs(P)
|
|
|
|
R.array().abs() // abs(P)
|
|
|
|
R.cwiseAbs2() // abs(P.^2)
|
|
|
|
R.array().abs2() // abs(P.^2)
|
|
|
|
(R.array() < s).select(P,Q); // (R < s ? P : Q)
|
2009-02-18 18:27:18 +08:00
|
|
|
|
2009-02-20 00:46:35 +08:00
|
|
|
// Reductions.
|
|
|
|
int r, c;
|
|
|
|
// Eigen // Matlab
|
|
|
|
R.minCoeff() // min(R(:))
|
|
|
|
R.maxCoeff() // max(R(:))
|
|
|
|
s = R.minCoeff(&r, &c) // [aa, bb] = min(R); [cc, dd] = min(aa);
|
2009-02-18 18:27:18 +08:00
|
|
|
// r = bb(dd); c = dd; s = cc
|
2009-02-20 00:46:35 +08:00
|
|
|
s = R.maxCoeff(&r, &c) // [aa, bb] = max(R); [cc, dd] = max(aa);
|
2009-02-18 18:27:18 +08:00
|
|
|
// row = bb(dd); col = dd; s = cc
|
2009-02-20 00:46:35 +08:00
|
|
|
R.sum() // sum(R(:))
|
2012-10-03 15:48:33 +08:00
|
|
|
R.colwise().sum() // sum(R)
|
|
|
|
R.rowwise().sum() // sum(R, 2) or sum(R')'
|
2009-10-14 15:40:22 +08:00
|
|
|
R.prod() // prod(R(:))
|
2012-10-03 15:48:33 +08:00
|
|
|
R.colwise().prod() // prod(R)
|
|
|
|
R.rowwise().prod() // prod(R, 2) or prod(R')'
|
2009-02-20 00:46:35 +08:00
|
|
|
R.trace() // trace(R)
|
|
|
|
R.all() // all(R(:))
|
|
|
|
R.colwise().all() // all(R)
|
|
|
|
R.rowwise().all() // all(R, 2)
|
|
|
|
R.any() // any(R(:))
|
|
|
|
R.colwise().any() // any(R)
|
|
|
|
R.rowwise().any() // any(R, 2)
|
2009-02-18 18:27:18 +08:00
|
|
|
|
2009-02-20 00:46:35 +08:00
|
|
|
// Dot products, norms, etc.
|
|
|
|
// Eigen // Matlab
|
|
|
|
x.norm() // norm(x). Note that norm(R) doesn't work in Eigen.
|
|
|
|
x.squaredNorm() // dot(x, x) Note the equivalence is not true for complex
|
|
|
|
x.dot(y) // dot(x, y)
|
|
|
|
x.cross(y) // cross(x, y) Requires #include <Eigen/Geometry>
|
2009-02-18 18:27:18 +08:00
|
|
|
|
2009-02-20 00:46:35 +08:00
|
|
|
// Eigen can map existing memory into Eigen matrices.
|
|
|
|
float array[3];
|
|
|
|
Map<Vector3f>(array, 3).fill(10);
|
|
|
|
int data[4] = 1, 2, 3, 4;
|
|
|
|
Matrix2i mat2x2(data);
|
|
|
|
MatrixXi mat2x2 = Map<Matrix2i>(data);
|
|
|
|
MatrixXi mat2x2 = Map<MatrixXi>(data, 2, 2);
|
2009-02-18 18:27:18 +08:00
|
|
|
|
2009-02-20 00:46:35 +08:00
|
|
|
// Solve Ax = b. Result stored in x. Matlab: x = A \ b.
|
2013-03-12 04:20:12 +08:00
|
|
|
x = A.ldlt().solve(b)); // A sym. p.s.d. #include <Eigen/Cholesky>
|
|
|
|
x = A.llt() .solve(b)); // A sym. p.d. #include <Eigen/Cholesky>
|
|
|
|
x = A.lu() .solve(b)); // Stable and fast. #include <Eigen/LU>
|
|
|
|
x = A.qr() .solve(b)); // No pivoting. #include <Eigen/QR>
|
|
|
|
x = A.svd() .solve(b)); // Stable, slowest. #include <Eigen/SVD>
|
2009-02-20 00:46:35 +08:00
|
|
|
// .ldlt() -> .matrixL() and .matrixD()
|
|
|
|
// .llt() -> .matrixL()
|
|
|
|
// .lu() -> .matrixL() and .matrixU()
|
|
|
|
// .qr() -> .matrixQ() and .matrixR()
|
|
|
|
// .svd() -> .matrixU(), .singularValues(), and .matrixV()
|
2009-02-18 18:27:18 +08:00
|
|
|
|
2009-02-20 00:46:35 +08:00
|
|
|
// Eigenvalue problems
|
|
|
|
// Eigen // Matlab
|
|
|
|
A.eigenvalues(); // eig(A);
|
|
|
|
EigenSolver<Matrix3d> eig(A); // [vec val] = eig(A)
|
|
|
|
eig.eigenvalues(); // diag(val)
|
|
|
|
eig.eigenvectors(); // vec
|
2013-03-12 04:20:12 +08:00
|
|
|
// For self-adjoint matrices use SelfAdjointEigenSolver<>
|