diff options
Diffstat (limited to 'doc/AsciiQuickReference.txt')
-rw-r--r-- | doc/AsciiQuickReference.txt | 109 |
1 files changed, 73 insertions, 36 deletions
diff --git a/doc/AsciiQuickReference.txt b/doc/AsciiQuickReference.txt index d2e973f5d..1e74e0528 100644 --- a/doc/AsciiQuickReference.txt +++ b/doc/AsciiQuickReference.txt @@ -1,8 +1,7 @@ // A simple quickref for Eigen. Add anything that's missing. // Main author: Keir Mierle -#include <Eigen/Core> -#include <Eigen/Array> +#include <Eigen/Dense> Matrix<double, 3, 3> A; // Fixed rows and cols. Same as Matrix3d. Matrix<double, 3, Dynamic> B; // Fixed rows, dynamic cols. @@ -11,13 +10,14 @@ 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. +VectorXd v; // Dynamic column vector of doubles double s; // 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 +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) // @@ -31,9 +31,19 @@ A << 1, 2, 3, // Initialize A. The elements can also be 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. -A.setRandom(); // Fill A with uniform random numbers in (-1, 1). - // Requires #include <Eigen/Array>. -A.setIdentity(); // Fill A with the identity. + +// 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::LinSpaced(size,low,high) // linspace(low,high,size)' +v.setLinSpaced(size,low,high) // v = linspace(low,high,size)' + // 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 @@ -41,20 +51,34 @@ A.setIdentity(); // Fill A with the identity. // Eigen // Matlab 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) +x.tail(n) // x(end - n + 1: end) +x.tail<n>() // x(end - n + 1: end) 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) +P.row(i) // P(i+1, :) +P.col(j) // P(:, j+1) +P.leftCols<cols>() // P(:, 1:cols) +P.leftCols(cols) // P(:, 1:cols) +P.middleCols<cols>(j) // P(:, j+1:j+cols) +P.middleCols(j, cols) // P(:, j+1:j+cols) +P.rightCols<cols>() // P(:, end-cols+1:end) +P.rightCols(cols) // P(:, end-cols+1:end) +P.topRows<rows>() // P(1:rows, :) +P.topRows(rows) // P(1:rows, :) +P.middleRows<rows>(i) // P(:, i+1:i+rows) +P.middleRows(i, rows) // P(:, i+1:i+rows) +P.bottomRows<rows>() // P(:, end-rows+1:end) +P.bottomRows(rows) // P(:, end-rows+1:end) 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.topRightCorner(rows, cols) // P(1:rows, end-cols+1:end) +P.bottomLeftCorner(rows, cols) // P(end-rows+1:end, 1:cols) +P.bottomRightCorner(rows, cols) // P(end-rows+1:end, end-cols+1:end) 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.topRightCorner<rows,cols>() // P(1:rows, end-cols+1:end) +P.bottomLeftCorner<rows,cols>() // P(end-rows+1:end, 1:cols) +P.bottomRightCorner<rows,cols>() // P(end-rows+1:end, end-cols+1:end) // Of particular note is Eigen's swap function which is highly optimized. // Eigen // Matlab @@ -67,6 +91,8 @@ R.adjoint() // R' R.transpose() // R.' or conj(R') R.diagonal() // diag(R) x.asDiagonal() // diag(x) +R.transpose().colwise().reverse(); // rot90(R) +R.conjugate() // conj(R) // All the same as Matlab, but matlab doesn't have *= style operators. // Matrix-vector. Matrix-matrix. Matrix-scalar. @@ -77,8 +103,7 @@ a *= M; R = P + Q; R = P/s; R += Q; R *= s; R -= Q; R /= s; - // Vectorized operations on each element independently - // (most require #include <Eigen/Array>) +// Vectorized operations on each element independently // Eigen // Matlab R = P.cwiseProduct(Q); // R = P .* Q R = P.array() * s.array();// R = P .* s @@ -116,16 +141,14 @@ 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); - // r = bb(dd); c = dd; s = cc -s = R.maxCoeff(&r, &c) // [aa, bb] = max(R); [cc, dd] = max(aa); - // row = bb(dd); col = dd; s = cc +s = R.minCoeff(&r, &c) // [s, i] = min(R(:)); [r, c] = ind2sub(size(R), i); +s = R.maxCoeff(&r, &c) // [s, i] = max(R(:)); [r, c] = ind2sub(size(R), i); R.sum() // sum(R(:)) -R.colwise.sum() // sum(R) -R.rowwise.sum() // sum(R, 2) or sum(R')' +R.colwise().sum() // sum(R) +R.rowwise().sum() // sum(R, 2) or sum(R')' R.prod() // prod(R(:)) -R.colwise.prod() // prod(R) -R.rowwise.prod() // prod(R, 2) or prod(R')' +R.colwise().prod() // prod(R) +R.rowwise().prod() // prod(R, 2) or prod(R')' R.trace() // trace(R) R.all() // all(R(:)) R.colwise().all() // all(R) @@ -141,21 +164,34 @@ x.squaredNorm() // dot(x, x) Note the equivalence is not true for co x.dot(y) // dot(x, y) x.cross(y) // cross(x, y) Requires #include <Eigen/Geometry> +//// Type conversion +// Eigen // Matlab +A.cast<double>(); // double(A) +A.cast<float>(); // single(A) +A.cast<int>(); // int32(A) +A.real(); // real(A) +A.imag(); // imag(A) +// if the original type equals destination type, no work is done + +// Note that for most operations Eigen requires all operands to have the same type: +MatrixXf F = MatrixXf::Zero(3,3); +A += F; // illegal in Eigen. In Matlab A = A+F is allowed +A += F.cast<double>(); // F converted to double and then added (generally, conversion happens on-the-fly) + // 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); +Vector3f::Map(array).fill(10); // create a temporary Map over array and sets entries to 10 +int data[4] = {1, 2, 3, 4}; +Matrix2i mat2x2(data); // copies data into mat2x2 +Matrix2i::Map(data) = 2*mat2x2; // overwrite elements of data with 2*mat2x2 +MatrixXi::Map(data, 2, 2) += mat2x2; // adds mat2x2 to elements of data (alternative syntax if size is not know at compile time) // Solve Ax = b. Result stored in x. Matlab: x = A \ b. -bool solved; -solved = A.ldlt().solve(b, &x)); // A sym. p.s.d. #include <Eigen/Cholesky> -solved = A.llt() .solve(b, &x)); // A sym. p.d. #include <Eigen/Cholesky> -solved = A.lu() .solve(b, &x)); // Stable and fast. #include <Eigen/LU> -solved = A.qr() .solve(b, &x)); // No pivoting. #include <Eigen/QR> -solved = A.svd() .solve(b, &x)); // Stable, slowest. #include <Eigen/SVD> +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> // .ldlt() -> .matrixL() and .matrixD() // .llt() -> .matrixL() // .lu() -> .matrixL() and .matrixU() @@ -168,3 +204,4 @@ A.eigenvalues(); // eig(A); EigenSolver<Matrix3d> eig(A); // [vec val] = eig(A) eig.eigenvalues(); // diag(val) eig.eigenvectors(); // vec +// For self-adjoint matrices use SelfAdjointEigenSolver<> |