Brian Silverman | 72890c2 | 2015-09-19 14:37:37 -0400 | [diff] [blame] | 1 | // A simple quickref for Eigen. Add anything that's missing. |
| 2 | // Main author: Keir Mierle |
| 3 | |
| 4 | #include <Eigen/Dense> |
| 5 | |
| 6 | Matrix<double, 3, 3> A; // Fixed rows and cols. Same as Matrix3d. |
| 7 | Matrix<double, 3, Dynamic> B; // Fixed rows, dynamic cols. |
| 8 | Matrix<double, Dynamic, Dynamic> C; // Full dynamic. Same as MatrixXd. |
| 9 | Matrix<double, 3, 3, RowMajor> E; // Row major; default is column-major. |
| 10 | Matrix3f P, Q, R; // 3x3 float matrix. |
| 11 | Vector3f x, y, z; // 3x1 float matrix. |
| 12 | RowVector3f a, b, c; // 1x3 float matrix. |
| 13 | VectorXd v; // Dynamic column vector of doubles |
| 14 | double s; |
| 15 | |
| 16 | // Basic usage |
| 17 | // Eigen // Matlab // comments |
| 18 | x.size() // length(x) // vector size |
| 19 | C.rows() // size(C,1) // number of rows |
| 20 | C.cols() // size(C,2) // number of columns |
| 21 | x(i) // x(i+1) // Matlab is 1-based |
| 22 | C(i,j) // C(i+1,j+1) // |
| 23 | |
| 24 | A.resize(4, 4); // Runtime error if assertions are on. |
| 25 | B.resize(4, 9); // Runtime error if assertions are on. |
| 26 | A.resize(3, 3); // Ok; size didn't change. |
| 27 | B.resize(3, 9); // Ok; only dynamic cols changed. |
| 28 | |
| 29 | A << 1, 2, 3, // Initialize A. The elements can also be |
| 30 | 4, 5, 6, // matrices, which are stacked along cols |
| 31 | 7, 8, 9; // and then the rows are stacked. |
| 32 | B << A, A, A; // B is three horizontally stacked A's. |
| 33 | A.fill(10); // Fill A with all 10's. |
| 34 | |
Austin Schuh | 189376f | 2018-12-20 22:11:15 +1100 | [diff] [blame] | 35 | // Eigen // Matlab |
| 36 | MatrixXd::Identity(rows,cols) // eye(rows,cols) |
| 37 | C.setIdentity(rows,cols) // C = eye(rows,cols) |
| 38 | MatrixXd::Zero(rows,cols) // zeros(rows,cols) |
| 39 | C.setZero(rows,cols) // C = zeros(rows,cols) |
| 40 | MatrixXd::Ones(rows,cols) // ones(rows,cols) |
| 41 | C.setOnes(rows,cols) // C = ones(rows,cols) |
| 42 | MatrixXd::Random(rows,cols) // rand(rows,cols)*2-1 // MatrixXd::Random returns uniform random numbers in (-1, 1). |
| 43 | C.setRandom(rows,cols) // C = rand(rows,cols)*2-1 |
| 44 | VectorXd::LinSpaced(size,low,high) // linspace(low,high,size)' |
| 45 | v.setLinSpaced(size,low,high) // v = linspace(low,high,size)' |
| 46 | VectorXi::LinSpaced(((hi-low)/step)+1, // low:step:hi |
| 47 | low,low+step*(size-1)) // |
Brian Silverman | 72890c2 | 2015-09-19 14:37:37 -0400 | [diff] [blame] | 48 | |
| 49 | |
| 50 | // Matrix slicing and blocks. All expressions listed here are read/write. |
| 51 | // Templated size versions are faster. Note that Matlab is 1-based (a size N |
| 52 | // vector is x(1)...x(N)). |
| 53 | // Eigen // Matlab |
| 54 | x.head(n) // x(1:n) |
| 55 | x.head<n>() // x(1:n) |
| 56 | x.tail(n) // x(end - n + 1: end) |
| 57 | x.tail<n>() // x(end - n + 1: end) |
| 58 | x.segment(i, n) // x(i+1 : i+n) |
| 59 | x.segment<n>(i) // x(i+1 : i+n) |
| 60 | P.block(i, j, rows, cols) // P(i+1 : i+rows, j+1 : j+cols) |
| 61 | P.block<rows, cols>(i, j) // P(i+1 : i+rows, j+1 : j+cols) |
| 62 | P.row(i) // P(i+1, :) |
| 63 | P.col(j) // P(:, j+1) |
| 64 | P.leftCols<cols>() // P(:, 1:cols) |
| 65 | P.leftCols(cols) // P(:, 1:cols) |
| 66 | P.middleCols<cols>(j) // P(:, j+1:j+cols) |
| 67 | P.middleCols(j, cols) // P(:, j+1:j+cols) |
| 68 | P.rightCols<cols>() // P(:, end-cols+1:end) |
| 69 | P.rightCols(cols) // P(:, end-cols+1:end) |
| 70 | P.topRows<rows>() // P(1:rows, :) |
| 71 | P.topRows(rows) // P(1:rows, :) |
| 72 | P.middleRows<rows>(i) // P(i+1:i+rows, :) |
| 73 | P.middleRows(i, rows) // P(i+1:i+rows, :) |
| 74 | P.bottomRows<rows>() // P(end-rows+1:end, :) |
| 75 | P.bottomRows(rows) // P(end-rows+1:end, :) |
| 76 | P.topLeftCorner(rows, cols) // P(1:rows, 1:cols) |
| 77 | P.topRightCorner(rows, cols) // P(1:rows, end-cols+1:end) |
| 78 | P.bottomLeftCorner(rows, cols) // P(end-rows+1:end, 1:cols) |
| 79 | P.bottomRightCorner(rows, cols) // P(end-rows+1:end, end-cols+1:end) |
| 80 | P.topLeftCorner<rows,cols>() // P(1:rows, 1:cols) |
| 81 | P.topRightCorner<rows,cols>() // P(1:rows, end-cols+1:end) |
| 82 | P.bottomLeftCorner<rows,cols>() // P(end-rows+1:end, 1:cols) |
| 83 | P.bottomRightCorner<rows,cols>() // P(end-rows+1:end, end-cols+1:end) |
| 84 | |
| 85 | // Of particular note is Eigen's swap function which is highly optimized. |
| 86 | // Eigen // Matlab |
Austin Schuh | 189376f | 2018-12-20 22:11:15 +1100 | [diff] [blame] | 87 | R.row(i) = P.col(j); // R(i, :) = P(:, j) |
Brian Silverman | 72890c2 | 2015-09-19 14:37:37 -0400 | [diff] [blame] | 88 | R.col(j1).swap(mat1.col(j2)); // R(:, [j1 j2]) = R(:, [j2, j1]) |
| 89 | |
Austin Schuh | 189376f | 2018-12-20 22:11:15 +1100 | [diff] [blame] | 90 | // Views, transpose, etc; |
Brian Silverman | 72890c2 | 2015-09-19 14:37:37 -0400 | [diff] [blame] | 91 | // Eigen // Matlab |
| 92 | R.adjoint() // R' |
Austin Schuh | 189376f | 2018-12-20 22:11:15 +1100 | [diff] [blame] | 93 | R.transpose() // R.' or conj(R') // Read-write |
| 94 | R.diagonal() // diag(R) // Read-write |
Brian Silverman | 72890c2 | 2015-09-19 14:37:37 -0400 | [diff] [blame] | 95 | x.asDiagonal() // diag(x) |
Austin Schuh | 189376f | 2018-12-20 22:11:15 +1100 | [diff] [blame] | 96 | R.transpose().colwise().reverse() // rot90(R) // Read-write |
| 97 | R.rowwise().reverse() // fliplr(R) |
| 98 | R.colwise().reverse() // flipud(R) |
| 99 | R.replicate(i,j) // repmat(P,i,j) |
| 100 | |
Brian Silverman | 72890c2 | 2015-09-19 14:37:37 -0400 | [diff] [blame] | 101 | |
| 102 | // All the same as Matlab, but matlab doesn't have *= style operators. |
| 103 | // Matrix-vector. Matrix-matrix. Matrix-scalar. |
| 104 | y = M*x; R = P*Q; R = P*s; |
| 105 | a = b*M; R = P - Q; R = s*P; |
| 106 | a *= M; R = P + Q; R = P/s; |
| 107 | R *= Q; R = s*P; |
| 108 | R += Q; R *= s; |
| 109 | R -= Q; R /= s; |
| 110 | |
| 111 | // Vectorized operations on each element independently |
Austin Schuh | 189376f | 2018-12-20 22:11:15 +1100 | [diff] [blame] | 112 | // Eigen // Matlab |
| 113 | R = P.cwiseProduct(Q); // R = P .* Q |
| 114 | R = P.array() * s.array(); // R = P .* s |
| 115 | R = P.cwiseQuotient(Q); // R = P ./ Q |
| 116 | R = P.array() / Q.array(); // R = P ./ Q |
| 117 | R = P.array() + s.array(); // R = P + s |
| 118 | R = P.array() - s.array(); // R = P - s |
| 119 | R.array() += s; // R = R + s |
| 120 | R.array() -= s; // R = R - s |
| 121 | R.array() < Q.array(); // R < Q |
| 122 | R.array() <= Q.array(); // R <= Q |
| 123 | R.cwiseInverse(); // 1 ./ P |
| 124 | R.array().inverse(); // 1 ./ P |
| 125 | R.array().sin() // sin(P) |
| 126 | R.array().cos() // cos(P) |
| 127 | R.array().pow(s) // P .^ s |
| 128 | R.array().square() // P .^ 2 |
| 129 | R.array().cube() // P .^ 3 |
| 130 | R.cwiseSqrt() // sqrt(P) |
| 131 | R.array().sqrt() // sqrt(P) |
| 132 | R.array().exp() // exp(P) |
| 133 | R.array().log() // log(P) |
| 134 | R.cwiseMax(P) // max(R, P) |
| 135 | R.array().max(P.array()) // max(R, P) |
| 136 | R.cwiseMin(P) // min(R, P) |
| 137 | R.array().min(P.array()) // min(R, P) |
| 138 | R.cwiseAbs() // abs(P) |
| 139 | R.array().abs() // abs(P) |
| 140 | R.cwiseAbs2() // abs(P.^2) |
| 141 | R.array().abs2() // abs(P.^2) |
| 142 | (R.array() < s).select(P,Q ); // (R < s ? P : Q) |
| 143 | R = (Q.array()==0).select(P,R) // R(Q==0) = P(Q==0) |
| 144 | R = P.unaryExpr(ptr_fun(func)) // R = arrayfun(func, P) // with: scalar func(const scalar &x); |
| 145 | |
Brian Silverman | 72890c2 | 2015-09-19 14:37:37 -0400 | [diff] [blame] | 146 | |
| 147 | // Reductions. |
| 148 | int r, c; |
| 149 | // Eigen // Matlab |
| 150 | R.minCoeff() // min(R(:)) |
| 151 | R.maxCoeff() // max(R(:)) |
| 152 | s = R.minCoeff(&r, &c) // [s, i] = min(R(:)); [r, c] = ind2sub(size(R), i); |
| 153 | s = R.maxCoeff(&r, &c) // [s, i] = max(R(:)); [r, c] = ind2sub(size(R), i); |
| 154 | R.sum() // sum(R(:)) |
| 155 | R.colwise().sum() // sum(R) |
| 156 | R.rowwise().sum() // sum(R, 2) or sum(R')' |
| 157 | R.prod() // prod(R(:)) |
| 158 | R.colwise().prod() // prod(R) |
| 159 | R.rowwise().prod() // prod(R, 2) or prod(R')' |
| 160 | R.trace() // trace(R) |
| 161 | R.all() // all(R(:)) |
| 162 | R.colwise().all() // all(R) |
| 163 | R.rowwise().all() // all(R, 2) |
| 164 | R.any() // any(R(:)) |
| 165 | R.colwise().any() // any(R) |
| 166 | R.rowwise().any() // any(R, 2) |
| 167 | |
| 168 | // Dot products, norms, etc. |
| 169 | // Eigen // Matlab |
| 170 | x.norm() // norm(x). Note that norm(R) doesn't work in Eigen. |
| 171 | x.squaredNorm() // dot(x, x) Note the equivalence is not true for complex |
| 172 | x.dot(y) // dot(x, y) |
| 173 | x.cross(y) // cross(x, y) Requires #include <Eigen/Geometry> |
| 174 | |
| 175 | //// Type conversion |
Austin Schuh | 189376f | 2018-12-20 22:11:15 +1100 | [diff] [blame] | 176 | // Eigen // Matlab |
| 177 | A.cast<double>(); // double(A) |
| 178 | A.cast<float>(); // single(A) |
| 179 | A.cast<int>(); // int32(A) |
| 180 | A.real(); // real(A) |
| 181 | A.imag(); // imag(A) |
Brian Silverman | 72890c2 | 2015-09-19 14:37:37 -0400 | [diff] [blame] | 182 | // if the original type equals destination type, no work is done |
| 183 | |
| 184 | // Note that for most operations Eigen requires all operands to have the same type: |
| 185 | MatrixXf F = MatrixXf::Zero(3,3); |
| 186 | A += F; // illegal in Eigen. In Matlab A = A+F is allowed |
| 187 | A += F.cast<double>(); // F converted to double and then added (generally, conversion happens on-the-fly) |
| 188 | |
| 189 | // Eigen can map existing memory into Eigen matrices. |
| 190 | float array[3]; |
| 191 | Vector3f::Map(array).fill(10); // create a temporary Map over array and sets entries to 10 |
| 192 | int data[4] = {1, 2, 3, 4}; |
| 193 | Matrix2i mat2x2(data); // copies data into mat2x2 |
| 194 | Matrix2i::Map(data) = 2*mat2x2; // overwrite elements of data with 2*mat2x2 |
| 195 | MatrixXi::Map(data, 2, 2) += mat2x2; // adds mat2x2 to elements of data (alternative syntax if size is not know at compile time) |
| 196 | |
| 197 | // Solve Ax = b. Result stored in x. Matlab: x = A \ b. |
| 198 | x = A.ldlt().solve(b)); // A sym. p.s.d. #include <Eigen/Cholesky> |
| 199 | x = A.llt() .solve(b)); // A sym. p.d. #include <Eigen/Cholesky> |
| 200 | x = A.lu() .solve(b)); // Stable and fast. #include <Eigen/LU> |
| 201 | x = A.qr() .solve(b)); // No pivoting. #include <Eigen/QR> |
| 202 | x = A.svd() .solve(b)); // Stable, slowest. #include <Eigen/SVD> |
| 203 | // .ldlt() -> .matrixL() and .matrixD() |
| 204 | // .llt() -> .matrixL() |
| 205 | // .lu() -> .matrixL() and .matrixU() |
| 206 | // .qr() -> .matrixQ() and .matrixR() |
| 207 | // .svd() -> .matrixU(), .singularValues(), and .matrixV() |
| 208 | |
| 209 | // Eigenvalue problems |
| 210 | // Eigen // Matlab |
| 211 | A.eigenvalues(); // eig(A); |
| 212 | EigenSolver<Matrix3d> eig(A); // [vec val] = eig(A) |
| 213 | eig.eigenvalues(); // diag(val) |
| 214 | eig.eigenvectors(); // vec |
| 215 | // For self-adjoint matrices use SelfAdjointEigenSolver<> |