QuickReference.dox 29 KB

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