You can not select more than 25 topics Topics must start with a letter or number, can include dashes ('-') and can be up to 35 characters long.

737 lines
26 KiB

  1. namespace Eigen {
  2. /** \page QuickRefPage Quick reference guide
  3. \b Table \b of \b contents
  4. - \ref QuickRef_Headers
  5. - \ref QuickRef_Types
  6. - \ref QuickRef_Map
  7. - \ref QuickRef_ArithmeticOperators
  8. - \ref QuickRef_Coeffwise
  9. - \ref QuickRef_Reductions
  10. - \ref QuickRef_Blocks
  11. - \ref QuickRef_Misc
  12. - \ref QuickRef_DiagTriSymm
  13. \n
  14. <hr>
  15. <a href="#" class="top">top</a>
  16. \section QuickRef_Headers Modules and Header files
  17. The Eigen library is divided in a Core module and several additional modules. Each module has a corresponding header file which has to be included in order to use the module. The \c %Dense and \c Eigen header files are provided to conveniently gain access to several modules at once.
  18. <table class="manual">
  19. <tr><th>Module</th><th>Header file</th><th>Contents</th></tr>
  20. <tr><td>\link Core_Module Core \endlink</td><td>\code#include <Eigen/Core>\endcode</td><td>Matrix and Array classes, basic linear algebra (including triangular and selfadjoint products), array manipulation</td></tr>
  21. <tr class="alt"><td>\link Geometry_Module Geometry \endlink</td><td>\code#include <Eigen/Geometry>\endcode</td><td>Transform, Translation, Scaling, Rotation2D and 3D rotations (Quaternion, AngleAxis)</td></tr>
  22. <tr><td>\link LU_Module LU \endlink</td><td>\code#include <Eigen/LU>\endcode</td><td>Inverse, determinant, LU decompositions with solver (FullPivLU, PartialPivLU)</td></tr>
  23. <tr><td>\link Cholesky_Module Cholesky \endlink</td><td>\code#include <Eigen/Cholesky>\endcode</td><td>LLT and LDLT Cholesky factorization with solver</td></tr>
  24. <tr class="alt"><td>\link Householder_Module Householder \endlink</td><td>\code#include <Eigen/Householder>\endcode</td><td>Householder transformations; this module is used by several linear algebra modules</td></tr>
  25. <tr><td>\link SVD_Module SVD \endlink</td><td>\code#include <Eigen/SVD>\endcode</td><td>SVD decomposition with least-squares solver (JacobiSVD)</td></tr>
  26. <tr class="alt"><td>\link QR_Module QR \endlink</td><td>\code#include <Eigen/QR>\endcode</td><td>QR decomposition with solver (HouseholderQR, ColPivHouseholderQR, FullPivHouseholderQR)</td></tr>
  27. <tr><td>\link Eigenvalues_Module Eigenvalues \endlink</td><td>\code#include <Eigen/Eigenvalues>\endcode</td><td>Eigenvalue, eigenvector decompositions (EigenSolver, SelfAdjointEigenSolver, ComplexEigenSolver)</td></tr>
  28. <tr class="alt"><td>\link Sparse_Module Sparse \endlink</td><td>\code#include <Eigen/Sparse>\endcode</td><td>%Sparse matrix storage and related basic linear algebra (SparseMatrix, DynamicSparseMatrix, SparseVector)</td></tr>
  29. <tr><td></td><td>\code#include <Eigen/Dense>\endcode</td><td>Includes Core, Geometry, LU, Cholesky, SVD, QR, and Eigenvalues header files</td></tr>
  30. <tr class="alt"><td></td><td>\code#include <Eigen/Eigen>\endcode</td><td>Includes %Dense and %Sparse header files (the whole Eigen library)</td></tr>
  31. </table>
  32. <a href="#" class="top">top</a>
  33. \section QuickRef_Types Array, matrix and vector types
  34. \b Recall: Eigen provides two kinds of dense objects: mathematical matrices and vectors which are both represented by the template class Matrix, and general 1D and 2D arrays represented by the template class Array:
  35. \code
  36. typedef Matrix<Scalar, RowsAtCompileTime, ColsAtCompileTime, Options> MyMatrixType;
  37. typedef Array<Scalar, RowsAtCompileTime, ColsAtCompileTime, Options> MyArrayType;
  38. \endcode
  39. \li \c Scalar is the scalar type of the coefficients (e.g., \c float, \c double, \c bool, \c int, etc.).
  40. \li \c RowsAtCompileTime and \c ColsAtCompileTime are the number of rows and columns of the matrix as known at compile-time or \c Dynamic.
  41. \li \c Options can be \c ColMajor or \c RowMajor, default is \c ColMajor. (see class Matrix for more options)
  42. All combinations are allowed: you can have a matrix with a fixed number of rows and a dynamic number of columns, etc. The following are all valid:
  43. \code
  44. Matrix<double, 6, Dynamic> // Dynamic number of columns (heap allocation)
  45. Matrix<double, Dynamic, 2> // Dynamic number of rows (heap allocation)
  46. Matrix<double, Dynamic, Dynamic, RowMajor> // Fully dynamic, row major (heap allocation)
  47. Matrix<double, 13, 3> // Fully fixed (static allocation)
  48. \endcode
  49. In most cases, you can simply use one of the convenience typedefs for \ref matrixtypedefs "matrices" and \ref arraytypedefs "arrays". Some examples:
  50. <table class="example">
  51. <tr><th>Matrices</th><th>Arrays</th></tr>
  52. <tr><td>\code
  53. Matrix<float,Dynamic,Dynamic> <=> MatrixXf
  54. Matrix<double,Dynamic,1> <=> VectorXd
  55. Matrix<int,1,Dynamic> <=> RowVectorXi
  56. Matrix<float,3,3> <=> Matrix3f
  57. Matrix<float,4,1> <=> Vector4f
  58. \endcode</td><td>\code
  59. Array<float,Dynamic,Dynamic> <=> ArrayXXf
  60. Array<double,Dynamic,1> <=> ArrayXd
  61. Array<int,1,Dynamic> <=> RowArrayXi
  62. Array<float,3,3> <=> Array33f
  63. Array<float,4,1> <=> Array4f
  64. \endcode</td></tr>
  65. </table>
  66. Conversion between the matrix and array worlds:
  67. \code
  68. Array44f a1, a1;
  69. Matrix4f m1, m2;
  70. m1 = a1 * a2; // coeffwise product, implicit conversion from array to matrix.
  71. a1 = m1 * m2; // matrix product, implicit conversion from matrix to array.
  72. a2 = a1 + m1.array(); // mixing array and matrix is forbidden
  73. m2 = a1.matrix() + m1; // and explicit conversion is required.
  74. ArrayWrapper<Matrix4f> m1a(m1); // m1a is an alias for m1.array(), they share the same coefficients
  75. MatrixWrapper<Array44f> a1m(a1);
  76. \endcode
  77. In the rest of this document we will use the following symbols to emphasize the features which are specifics to a given kind of object:
  78. \li <a name="matrixonly"><a/>\matrixworld linear algebra matrix and vector only
  79. \li <a name="arrayonly"><a/>\arrayworld array objects only
  80. \subsection QuickRef_Basics Basic matrix manipulation
  81. <table class="manual">
  82. <tr><th></th><th>1D objects</th><th>2D objects</th><th>Notes</th></tr>
  83. <tr><td>Constructors</td>
  84. <td>\code
  85. Vector4d v4;
  86. Vector2f v1(x, y);
  87. Array3i v2(x, y, z);
  88. Vector4d v3(x, y, z, w);
  89. VectorXf v5; // empty object
  90. ArrayXf v6(size);
  91. \endcode</td><td>\code
  92. Matrix4f m1;
  93. MatrixXf m5; // empty object
  94. MatrixXf m6(nb_rows, nb_columns);
  95. \endcode</td><td class="note">
  96. By default, the coefficients \n are left uninitialized</td></tr>
  97. <tr class="alt"><td>Comma initializer</td>
  98. <td>\code
  99. Vector3f v1; v1 << x, y, z;
  100. ArrayXf v2(4); v2 << 1, 2, 3, 4;
  101. \endcode</td><td>\code
  102. Matrix3f m1; m1 << 1, 2, 3,
  103. 4, 5, 6,
  104. 7, 8, 9;
  105. \endcode</td><td></td></tr>
  106. <tr><td>Comma initializer (bis)</td>
  107. <td colspan="2">
  108. \include Tutorial_commainit_02.cpp
  109. </td>
  110. <td>
  111. output:
  112. \verbinclude Tutorial_commainit_02.out
  113. </td>
  114. </tr>
  115. <tr class="alt"><td>Runtime info</td>
  116. <td>\code
  117. vector.size();
  118. vector.innerStride();
  119. vector.data();
  120. \endcode</td><td>\code
  121. matrix.rows(); matrix.cols();
  122. matrix.innerSize(); matrix.outerSize();
  123. matrix.innerStride(); matrix.outerStride();
  124. matrix.data();
  125. \endcode</td><td class="note">Inner/Outer* are storage order dependent</td></tr>
  126. <tr><td>Compile-time info</td>
  127. <td colspan="2">\code
  128. ObjectType::Scalar ObjectType::RowsAtCompileTime
  129. ObjectType::RealScalar ObjectType::ColsAtCompileTime
  130. ObjectType::Index ObjectType::SizeAtCompileTime
  131. \endcode</td><td></td></tr>
  132. <tr class="alt"><td>Resizing</td>
  133. <td>\code
  134. vector.resize(size);
  135. vector.resizeLike(other_vector);
  136. vector.conservativeResize(size);
  137. \endcode</td><td>\code
  138. matrix.resize(nb_rows, nb_cols);
  139. matrix.resize(Eigen::NoChange, nb_cols);
  140. matrix.resize(nb_rows, Eigen::NoChange);
  141. matrix.resizeLike(other_matrix);
  142. matrix.conservativeResize(nb_rows, nb_cols);
  143. \endcode</td><td class="note">no-op if the new sizes match,<br/>otherwise data are lost<br/><br/>resizing with data preservation</td></tr>
  144. <tr><td>Coeff access with \n range checking</td>
  145. <td>\code
  146. vector(i) vector.x()
  147. vector[i] vector.y()
  148. vector.z()
  149. vector.w()
  150. \endcode</td><td>\code
  151. matrix(i,j)
  152. \endcode</td><td class="note">Range checking is disabled if \n NDEBUG or EIGEN_NO_DEBUG is defined</td></tr>
  153. <tr class="alt"><td>Coeff access without \n range checking</td>
  154. <td>\code
  155. vector.coeff(i)
  156. vector.coeffRef(i)
  157. \endcode</td><td>\code
  158. matrix.coeff(i,j)
  159. matrix.coeffRef(i,j)
  160. \endcode</td><td></td></tr>
  161. <tr><td>Assignment/copy</td>
  162. <td colspan="2">\code
  163. object = expression;
  164. object_of_float = expression_of_double.cast<float>();
  165. \endcode</td><td class="note">the destination is automatically resized (if possible)</td></tr>
  166. </table>
  167. \subsection QuickRef_PredefMat Predefined Matrices
  168. <table class="manual">
  169. <tr>
  170. <th>Fixed-size matrix or vector</th>
  171. <th>Dynamic-size matrix</th>
  172. <th>Dynamic-size vector</th>
  173. </tr>
  174. <tr style="border-bottom-style: none;">
  175. <td>
  176. \code
  177. typedef {Matrix3f|Array33f} FixedXD;
  178. FixedXD x;
  179. x = FixedXD::Zero();
  180. x = FixedXD::Ones();
  181. x = FixedXD::Constant(value);
  182. x = FixedXD::Random();
  183. x = FixedXD::LinSpaced(size, low, high);
  184. x.setZero();
  185. x.setOnes();
  186. x.setConstant(value);
  187. x.setRandom();
  188. x.setLinSpaced(size, low, high);
  189. \endcode
  190. </td>
  191. <td>
  192. \code
  193. typedef {MatrixXf|ArrayXXf} Dynamic2D;
  194. Dynamic2D x;
  195. x = Dynamic2D::Zero(rows, cols);
  196. x = Dynamic2D::Ones(rows, cols);
  197. x = Dynamic2D::Constant(rows, cols, value);
  198. x = Dynamic2D::Random(rows, cols);
  199. N/A
  200. x.setZero(rows, cols);
  201. x.setOnes(rows, cols);
  202. x.setConstant(rows, cols, value);
  203. x.setRandom(rows, cols);
  204. N/A
  205. \endcode
  206. </td>
  207. <td>
  208. \code
  209. typedef {VectorXf|ArrayXf} Dynamic1D;
  210. Dynamic1D x;
  211. x = Dynamic1D::Zero(size);
  212. x = Dynamic1D::Ones(size);
  213. x = Dynamic1D::Constant(size, value);
  214. x = Dynamic1D::Random(size);
  215. x = Dynamic1D::LinSpaced(size, low, high);
  216. x.setZero(size);
  217. x.setOnes(size);
  218. x.setConstant(size, value);
  219. x.setRandom(size);
  220. x.setLinSpaced(size, low, high);
  221. \endcode
  222. </td>
  223. </tr>
  224. <tr><td colspan="3">Identity and \link MatrixBase::Unit basis vectors \endlink \matrixworld</td></tr>
  225. <tr style="border-bottom-style: none;">
  226. <td>
  227. \code
  228. x = FixedXD::Identity();
  229. x.setIdentity();
  230. Vector3f::UnitX() // 1 0 0
  231. Vector3f::UnitY() // 0 1 0
  232. Vector3f::UnitZ() // 0 0 1
  233. \endcode
  234. </td>
  235. <td>
  236. \code
  237. x = Dynamic2D::Identity(rows, cols);
  238. x.setIdentity(rows, cols);
  239. N/A
  240. \endcode
  241. </td>
  242. <td>\code
  243. N/A
  244. VectorXf::Unit(size,i)
  245. VectorXf::Unit(4,1) == Vector4f(0,1,0,0)
  246. == Vector4f::UnitY()
  247. \endcode
  248. </td>
  249. </tr>
  250. </table>
  251. \subsection QuickRef_Map Mapping external arrays
  252. <table class="manual">
  253. <tr>
  254. <td>Contiguous \n memory</td>
  255. <td>\code
  256. float data[] = {1,2,3,4};
  257. Map<Vector3f> v1(data); // uses v1 as a Vector3f object
  258. Map<ArrayXf> v2(data,3); // uses v2 as a ArrayXf object
  259. Map<Array22f> m1(data); // uses m1 as a Array22f object
  260. Map<MatrixXf> m2(data,2,2); // uses m2 as a MatrixXf object
  261. \endcode</td>
  262. </tr>
  263. <tr>
  264. <td>Typical usage \n of strides</td>
  265. <td>\code
  266. float data[] = {1,2,3,4,5,6,7,8,9};
  267. Map<VectorXf,0,InnerStride<2> > v1(data,3); // = [1,3,5]
  268. Map<VectorXf,0,InnerStride<> > v2(data,3,InnerStride<>(3)); // = [1,4,7]
  269. Map<MatrixXf,0,OuterStride<3> > m2(data,2,3); // both lines |1,4,7|
  270. Map<MatrixXf,0,OuterStride<> > m1(data,2,3,OuterStride<>(3)); // are equal to: |2,5,8|
  271. \endcode</td>
  272. </tr>
  273. </table>
  274. <a href="#" class="top">top</a>
  275. \section QuickRef_ArithmeticOperators Arithmetic Operators
  276. <table class="manual">
  277. <tr><td>
  278. add \n subtract</td><td>\code
  279. mat3 = mat1 + mat2; mat3 += mat1;
  280. mat3 = mat1 - mat2; mat3 -= mat1;\endcode
  281. </td></tr>
  282. <tr class="alt"><td>
  283. scalar product</td><td>\code
  284. mat3 = mat1 * s1; mat3 *= s1; mat3 = s1 * mat1;
  285. mat3 = mat1 / s1; mat3 /= s1;\endcode
  286. </td></tr>
  287. <tr><td>
  288. matrix/vector \n products \matrixworld</td><td>\code
  289. col2 = mat1 * col1;
  290. row2 = row1 * mat1; row1 *= mat1;
  291. mat3 = mat1 * mat2; mat3 *= mat1; \endcode
  292. </td></tr>
  293. <tr class="alt"><td>
  294. transposition \n adjoint \matrixworld</td><td>\code
  295. mat1 = mat2.transpose(); mat1.transposeInPlace();
  296. mat1 = mat2.adjoint(); mat1.adjointInPlace();
  297. \endcode
  298. </td></tr>
  299. <tr><td>
  300. \link MatrixBase::dot() dot \endlink product \n inner product \matrixworld</td><td>\code
  301. scalar = vec1.dot(vec2);
  302. scalar = col1.adjoint() * col2;
  303. scalar = (col1.adjoint() * col2).value();\endcode
  304. </td></tr>
  305. <tr class="alt"><td>
  306. outer product \matrixworld</td><td>\code
  307. mat = col1 * col2.transpose();\endcode
  308. </td></tr>
  309. <tr><td>
  310. \link MatrixBase::norm() norm \endlink \n \link MatrixBase::normalized() normalization \endlink \matrixworld</td><td>\code
  311. scalar = vec1.norm(); scalar = vec1.squaredNorm()
  312. vec2 = vec1.normalized(); vec1.normalize(); // inplace \endcode
  313. </td></tr>
  314. <tr class="alt"><td>
  315. \link MatrixBase::cross() cross product \endlink \matrixworld</td><td>\code
  316. #include <Eigen/Geometry>
  317. vec3 = vec1.cross(vec2);\endcode</td></tr>
  318. </table>
  319. <a href="#" class="top">top</a>
  320. \section QuickRef_Coeffwise Coefficient-wise \& Array operators
  321. Coefficient-wise operators for matrices and vectors:
  322. <table class="manual">
  323. <tr><th>Matrix API \matrixworld</th><th>Via Array conversions</th></tr>
  324. <tr><td>\code
  325. mat1.cwiseMin(mat2)
  326. mat1.cwiseMax(mat2)
  327. mat1.cwiseAbs2()
  328. mat1.cwiseAbs()
  329. mat1.cwiseSqrt()
  330. mat1.cwiseProduct(mat2)
  331. mat1.cwiseQuotient(mat2)\endcode
  332. </td><td>\code
  333. mat1.array().min(mat2.array())
  334. mat1.array().max(mat2.array())
  335. mat1.array().abs2()
  336. mat1.array().abs()
  337. mat1.array().sqrt()
  338. mat1.array() * mat2.array()
  339. mat1.array() / mat2.array()
  340. \endcode</td></tr>
  341. </table>
  342. It is also very simple to apply any user defined function \c foo using DenseBase::unaryExpr together with std::ptr_fun:
  343. \code mat1.unaryExpr(std::ptr_fun(foo))\endcode
  344. Array operators:\arrayworld
  345. <table class="manual">
  346. <tr><td>Arithmetic operators</td><td>\code
  347. array1 * array2 array1 / array2 array1 *= array2 array1 /= array2
  348. array1 + scalar array1 - scalar array1 += scalar array1 -= scalar
  349. \endcode</td></tr>
  350. <tr><td>Comparisons</td><td>\code
  351. array1 < array2 array1 > array2 array1 < scalar array1 > scalar
  352. array1 <= array2 array1 >= array2 array1 <= scalar array1 >= scalar
  353. array1 == array2 array1 != array2 array1 == scalar array1 != scalar
  354. \endcode</td></tr>
  355. <tr><td>Trigo, power, and \n misc functions \n and the STL variants</td><td>\code
  356. array1.min(array2)
  357. array1.max(array2)
  358. array1.abs2()
  359. array1.abs() std::abs(array1)
  360. array1.sqrt() std::sqrt(array1)
  361. array1.log() std::log(array1)
  362. array1.exp() std::exp(array1)
  363. array1.pow(exponent) std::pow(array1,exponent)
  364. array1.square()
  365. array1.cube()
  366. array1.inverse()
  367. array1.sin() std::sin(array1)
  368. array1.cos() std::cos(array1)
  369. array1.tan() std::tan(array1)
  370. array1.asin() std::asin(array1)
  371. array1.acos() std::acos(array1)
  372. \endcode
  373. </td></tr>
  374. </table>
  375. <a href="#" class="top">top</a>
  376. \section QuickRef_Reductions Reductions
  377. Eigen provides several reduction methods such as:
  378. \link DenseBase::minCoeff() minCoeff() \endlink, \link DenseBase::maxCoeff() maxCoeff() \endlink,
  379. \link DenseBase::sum() sum() \endlink, \link DenseBase::prod() prod() \endlink,
  380. \link MatrixBase::trace() trace() \endlink \matrixworld,
  381. \link MatrixBase::norm() norm() \endlink \matrixworld, \link MatrixBase::squaredNorm() squaredNorm() \endlink \matrixworld,
  382. \link DenseBase::all() all() \endlink, and \link DenseBase::any() any() \endlink.
  383. All reduction operations can be done matrix-wise,
  384. \link DenseBase::colwise() column-wise \endlink or
  385. \link DenseBase::rowwise() row-wise \endlink. Usage example:
  386. <table class="manual">
  387. <tr><td rowspan="3" style="border-right-style:dashed;vertical-align:middle">\code
  388. 5 3 1
  389. mat = 2 7 8
  390. 9 4 6 \endcode
  391. </td> <td>\code mat.minCoeff(); \endcode</td><td>\code 1 \endcode</td></tr>
  392. <tr class="alt"><td>\code mat.colwise().minCoeff(); \endcode</td><td>\code 2 3 1 \endcode</td></tr>
  393. <tr style="vertical-align:middle"><td>\code mat.rowwise().minCoeff(); \endcode</td><td>\code
  394. 1
  395. 2
  396. 4
  397. \endcode</td></tr>
  398. </table>
  399. Special versions of \link DenseBase::minCoeff(Index*,Index*) minCoeff \endlink and \link DenseBase::maxCoeff(Index*,Index*) maxCoeff \endlink:
  400. \code
  401. int i, j;
  402. s = vector.minCoeff(&i); // s == vector[i]
  403. s = matrix.maxCoeff(&i, &j); // s == matrix(i,j)
  404. \endcode
  405. Typical use cases of all() and any():
  406. \code
  407. if((array1 > 0).all()) ... // if all coefficients of array1 are greater than 0 ...
  408. if((array1 < array2).any()) ... // if there exist a pair i,j such that array1(i,j) < array2(i,j) ...
  409. \endcode
  410. <a href="#" class="top">top</a>\section QuickRef_Blocks Sub-matrices
  411. Read-write access to a \link DenseBase::col(Index) column \endlink
  412. or a \link DenseBase::row(Index) row \endlink of a matrix (or array):
  413. \code
  414. mat1.row(i) = mat2.col(j);
  415. mat1.col(j1).swap(mat1.col(j2));
  416. \endcode
  417. Read-write access to sub-vectors:
  418. <table class="manual">
  419. <tr>
  420. <th>Default versions</th>
  421. <th>Optimized versions when the size \n is known at compile time</th></tr>
  422. <th></th>
  423. <tr><td>\code vec1.head(n)\endcode</td><td>\code vec1.head<n>()\endcode</td><td>the first \c n coeffs </td></tr>
  424. <tr><td>\code vec1.tail(n)\endcode</td><td>\code vec1.tail<n>()\endcode</td><td>the last \c n coeffs </td></tr>
  425. <tr><td>\code vec1.segment(pos,n)\endcode</td><td>\code vec1.segment<n>(pos)\endcode</td>
  426. <td>the \c n coeffs in the \n range [\c pos : \c pos + \c n - 1]</td></tr>
  427. <tr class="alt"><td colspan="3">
  428. Read-write access to sub-matrices:</td></tr>
  429. <tr>
  430. <td>\code mat1.block(i,j,rows,cols)\endcode
  431. \link DenseBase::block(Index,Index,Index,Index) (more) \endlink</td>
  432. <td>\code mat1.block<rows,cols>(i,j)\endcode
  433. \link DenseBase::block(Index,Index) (more) \endlink</td>
  434. <td>the \c rows x \c cols sub-matrix \n starting from position (\c i,\c j)</td></tr>
  435. <tr><td>\code
  436. mat1.topLeftCorner(rows,cols)
  437. mat1.topRightCorner(rows,cols)
  438. mat1.bottomLeftCorner(rows,cols)
  439. mat1.bottomRightCorner(rows,cols)\endcode
  440. <td>\code
  441. mat1.topLeftCorner<rows,cols>()
  442. mat1.topRightCorner<rows,cols>()
  443. mat1.bottomLeftCorner<rows,cols>()
  444. mat1.bottomRightCorner<rows,cols>()\endcode
  445. <td>the \c rows x \c cols sub-matrix \n taken in one of the four corners</td></tr>
  446. <tr><td>\code
  447. mat1.topRows(rows)
  448. mat1.bottomRows(rows)
  449. mat1.leftCols(cols)
  450. mat1.rightCols(cols)\endcode
  451. <td>\code
  452. mat1.topRows<rows>()
  453. mat1.bottomRows<rows>()
  454. mat1.leftCols<cols>()
  455. mat1.rightCols<cols>()\endcode
  456. <td>specialized versions of block() \n when the block fit two corners</td></tr>
  457. </table>
  458. <a href="#" class="top">top</a>\section QuickRef_Misc Miscellaneous operations
  459. \subsection QuickRef_Reverse Reverse
  460. Vectors, rows, and/or columns of a matrix can be reversed (see DenseBase::reverse(), DenseBase::reverseInPlace(), VectorwiseOp::reverse()).
  461. \code
  462. vec.reverse() mat.colwise().reverse() mat.rowwise().reverse()
  463. vec.reverseInPlace()
  464. \endcode
  465. \subsection QuickRef_Replicate Replicate
  466. Vectors, matrices, rows, and/or columns can be replicated in any direction (see DenseBase::replicate(), VectorwiseOp::replicate())
  467. \code
  468. vec.replicate(times) vec.replicate<Times>
  469. mat.replicate(vertical_times, horizontal_times) mat.replicate<VerticalTimes, HorizontalTimes>()
  470. mat.colwise().replicate(vertical_times, horizontal_times) mat.colwise().replicate<VerticalTimes, HorizontalTimes>()
  471. mat.rowwise().replicate(vertical_times, horizontal_times) mat.rowwise().replicate<VerticalTimes, HorizontalTimes>()
  472. \endcode
  473. <a href="#" class="top">top</a>\section QuickRef_DiagTriSymm Diagonal, Triangular, and Self-adjoint matrices
  474. (matrix world \matrixworld)
  475. \subsection QuickRef_Diagonal Diagonal matrices
  476. <table class="example">
  477. <tr><th>Operation</th><th>Code</th></tr>
  478. <tr><td>
  479. view a vector \link MatrixBase::asDiagonal() as a diagonal matrix \endlink \n </td><td>\code
  480. mat1 = vec1.asDiagonal();\endcode
  481. </td></tr>
  482. <tr><td>
  483. Declare a diagonal matrix</td><td>\code
  484. DiagonalMatrix<Scalar,SizeAtCompileTime> diag1(size);
  485. diag1.diagonal() = vector;\endcode
  486. </td></tr>
  487. <tr><td>Access the \link MatrixBase::diagonal() diagonal \endlink and \link MatrixBase::diagonal(Index) super/sub diagonals \endlink of a matrix as a vector (read/write)</td>
  488. <td>\code
  489. vec1 = mat1.diagonal(); mat1.diagonal() = vec1; // main diagonal
  490. vec1 = mat1.diagonal(+n); mat1.diagonal(+n) = vec1; // n-th super diagonal
  491. vec1 = mat1.diagonal(-n); mat1.diagonal(-n) = vec1; // n-th sub diagonal
  492. vec1 = mat1.diagonal<1>(); mat1.diagonal<1>() = vec1; // first super diagonal
  493. vec1 = mat1.diagonal<-2>(); mat1.diagonal<-2>() = vec1; // second sub diagonal
  494. \endcode</td>
  495. </tr>
  496. <tr><td>Optimized products and inverse</td>
  497. <td>\code
  498. mat3 = scalar * diag1 * mat1;
  499. mat3 += scalar * mat1 * vec1.asDiagonal();
  500. mat3 = vec1.asDiagonal().inverse() * mat1
  501. mat3 = mat1 * diag1.inverse()
  502. \endcode</td>
  503. </tr>
  504. </table>
  505. \subsection QuickRef_TriangularView Triangular views
  506. TriangularView gives a view on a triangular part of a dense matrix and allows to perform optimized operations on it. The opposite triangular part is never referenced and can be used to store other information.
  507. \note The .triangularView() template member function requires the \c template keyword if it is used on an
  508. object of a type that depends on a template parameter; see \ref TopicTemplateKeyword for details.
  509. <table class="example">
  510. <tr><th>Operation</th><th>Code</th></tr>
  511. <tr><td>
  512. Reference to a triangular with optional \n
  513. unit or null diagonal (read/write):
  514. </td><td>\code
  515. m.triangularView<Xxx>()
  516. \endcode \n
  517. \c Xxx = ::Upper, ::Lower, ::StrictlyUpper, ::StrictlyLower, ::UnitUpper, ::UnitLower
  518. </td></tr>
  519. <tr><td>
  520. Writing to a specific triangular part:\n (only the referenced triangular part is evaluated)
  521. </td><td>\code
  522. m1.triangularView<Eigen::Lower>() = m2 + m3 \endcode
  523. </td></tr>
  524. <tr><td>
  525. Conversion to a dense matrix setting the opposite triangular part to zero:
  526. </td><td>\code
  527. m2 = m1.triangularView<Eigen::UnitUpper>()\endcode
  528. </td></tr>
  529. <tr><td>
  530. Products:
  531. </td><td>\code
  532. m3 += s1 * m1.adjoint().triangularView<Eigen::UnitUpper>() * m2
  533. m3 -= s1 * m2.conjugate() * m1.adjoint().triangularView<Eigen::Lower>() \endcode
  534. </td></tr>
  535. <tr><td>
  536. Solving linear equations:\n
  537. \f$ M_2 := L_1^{-1} M_2 \f$ \n
  538. \f$ M_3 := {L_1^*}^{-1} M_3 \f$ \n
  539. \f$ M_4 := M_4 U_1^{-1} \f$
  540. </td><td>\n \code
  541. L1.triangularView<Eigen::UnitLower>().solveInPlace(M2)
  542. L1.triangularView<Eigen::Lower>().adjoint().solveInPlace(M3)
  543. U1.triangularView<Eigen::Upper>().solveInPlace<OnTheRight>(M4)\endcode
  544. </td></tr>
  545. </table>
  546. \subsection QuickRef_SelfadjointMatrix Symmetric/selfadjoint views
  547. Just as for triangular matrix, you can reference any triangular part of a square matrix to see it as a selfadjoint
  548. matrix and perform special and optimized operations. Again the opposite triangular part is never referenced and can be
  549. used to store other information.
  550. \note The .selfadjointView() template member function requires the \c template keyword if it is used on an
  551. object of a type that depends on a template parameter; see \ref TopicTemplateKeyword for details.
  552. <table class="example">
  553. <tr><th>Operation</th><th>Code</th></tr>
  554. <tr><td>
  555. Conversion to a dense matrix:
  556. </td><td>\code
  557. m2 = m.selfadjointView<Eigen::Lower>();\endcode
  558. </td></tr>
  559. <tr><td>
  560. Product with another general matrix or vector:
  561. </td><td>\code
  562. m3 = s1 * m1.conjugate().selfadjointView<Eigen::Upper>() * m3;
  563. m3 -= s1 * m3.adjoint() * m1.selfadjointView<Eigen::Lower>();\endcode
  564. </td></tr>
  565. <tr><td>
  566. Rank 1 and rank K update: \n
  567. \f$ upper(M_1) \mathrel{{+}{=}} s_1 M_2 M_2^* \f$ \n
  568. \f$ lower(M_1) \mathbin{{-}{=}} M_2^* M_2 \f$
  569. </td><td>\n \code
  570. M1.selfadjointView<Eigen::Upper>().rankUpdate(M2,s1);
  571. M1.selfadjointView<Eigen::Lower>().rankUpdate(M2.adjoint(),-1); \endcode
  572. </td></tr>
  573. <tr><td>
  574. Rank 2 update: (\f$ M \mathrel{{+}{=}} s u v^* + s v u^* \f$)
  575. </td><td>\code
  576. M.selfadjointView<Eigen::Upper>().rankUpdate(u,v,s);
  577. \endcode
  578. </td></tr>
  579. <tr><td>
  580. Solving linear equations:\n(\f$ M_2 := M_1^{-1} M_2 \f$)
  581. </td><td>\code
  582. // via a standard Cholesky factorization
  583. m2 = m1.selfadjointView<Eigen::Upper>().llt().solve(m2);
  584. // via a Cholesky factorization with pivoting
  585. m2 = m1.selfadjointView<Eigen::Lower>().ldlt().solve(m2);
  586. \endcode
  587. </td></tr>
  588. </table>
  589. */
  590. /*
  591. <table class="tutorial_code">
  592. <tr><td>
  593. \link MatrixBase::asDiagonal() make a diagonal matrix \endlink \n from a vector </td><td>\code
  594. mat1 = vec1.asDiagonal();\endcode
  595. </td></tr>
  596. <tr><td>
  597. Declare a diagonal matrix</td><td>\code
  598. DiagonalMatrix<Scalar,SizeAtCompileTime> diag1(size);
  599. diag1.diagonal() = vector;\endcode
  600. </td></tr>
  601. <tr><td>Access \link MatrixBase::diagonal() the diagonal and super/sub diagonals of a matrix \endlink as a vector (read/write)</td>
  602. <td>\code
  603. vec1 = mat1.diagonal(); mat1.diagonal() = vec1; // main diagonal
  604. vec1 = mat1.diagonal(+n); mat1.diagonal(+n) = vec1; // n-th super diagonal
  605. vec1 = mat1.diagonal(-n); mat1.diagonal(-n) = vec1; // n-th sub diagonal
  606. vec1 = mat1.diagonal<1>(); mat1.diagonal<1>() = vec1; // first super diagonal
  607. vec1 = mat1.diagonal<-2>(); mat1.diagonal<-2>() = vec1; // second sub diagonal
  608. \endcode</td>
  609. </tr>
  610. <tr><td>View on a triangular part of a matrix (read/write)</td>
  611. <td>\code
  612. mat2 = mat1.triangularView<Xxx>();
  613. // Xxx = Upper, Lower, StrictlyUpper, StrictlyLower, UnitUpper, UnitLower
  614. mat1.triangularView<Upper>() = mat2 + mat3; // only the upper part is evaluated and referenced
  615. \endcode</td></tr>
  616. <tr><td>View a triangular part as a symmetric/self-adjoint matrix (read/write)</td>
  617. <td>\code
  618. mat2 = mat1.selfadjointView<Xxx>(); // Xxx = Upper or Lower
  619. mat1.selfadjointView<Upper>() = mat2 + mat2.adjoint(); // evaluated and write to the upper triangular part only
  620. \endcode</td></tr>
  621. </table>
  622. Optimized products:
  623. \code
  624. mat3 += scalar * vec1.asDiagonal() * mat1
  625. mat3 += scalar * mat1 * vec1.asDiagonal()
  626. mat3.noalias() += scalar * mat1.triangularView<Xxx>() * mat2
  627. mat3.noalias() += scalar * mat2 * mat1.triangularView<Xxx>()
  628. mat3.noalias() += scalar * mat1.selfadjointView<Upper or Lower>() * mat2
  629. mat3.noalias() += scalar * mat2 * mat1.selfadjointView<Upper or Lower>()
  630. mat1.selfadjointView<Upper or Lower>().rankUpdate(mat2);
  631. mat1.selfadjointView<Upper or Lower>().rankUpdate(mat2.adjoint(), scalar);
  632. \endcode
  633. Inverse products: (all are optimized)
  634. \code
  635. mat3 = vec1.asDiagonal().inverse() * mat1
  636. mat3 = mat1 * diag1.inverse()
  637. mat1.triangularView<Xxx>().solveInPlace(mat2)
  638. mat1.triangularView<Xxx>().solveInPlace<OnTheRight>(mat2)
  639. mat2 = mat1.selfadjointView<Upper or Lower>().llt().solve(mat2)
  640. \endcode
  641. */
  642. }