| 
									
										
										
										
											2010-06-29 22:02:33 +08:00
										 |  |  | namespace Eigen { | 
					
						
							|  |  |  | 
 | 
					
						
							| 
									
										
											  
											
												Big changes in Eigen documentation:
- Organize the documentation into "chapters".
  - Each chapter include many documentation pages, reference pages organized as modules, and a quick reference page.
  - The "Chapters" tree is created using the defgroup/ingroup mechanism, even for the documentation pages (i.e., .dox files for which I added an \eigenManualPage macro that we can switch between \page or \defgroup ).
  - Add a "General topics" entry for all pages that do not fit well in the previous "chapters".
  - The highlevel struture is managed by a new eigendoxy_layout.xml file.
- remove the "index" and quite useless pages (namespace list, class hierarchy, member list, file list, etc.)
- add the javascript search-engine.
- add the "treeview" panel.
- remove \tableofcontents (replace them by a custom \eigenAutoToc macro to be able to easily re-enable if needed).
- add javascript to automatically generate a TOC from the h1/h2 tags of the current page, and put the TOC in the left side panel.
- overload various javascript function generated by doxygen to:
  - remove the root of the treeview
  - remove links to section/subsection from the treeview
  - automatically expand the "Chapters" section
  - automatically expand the current section
  - adjust the height of the treeview to take into account the TOC
- always use the default .css file, eigendoxy.css now only includes our modifications
- use Doxyfile to specify our logo
- remove cross references to unsupported modules (temporarily)
											
										 
											2013-01-05 23:37:11 +08:00
										 |  |  | /** \eigenManualPage TutorialLinearAlgebra Linear algebra and decompositions | 
					
						
							| 
									
										
										
										
											2010-06-29 22:02:33 +08:00
										 |  |  | 
 | 
					
						
							| 
									
										
										
										
											2013-01-07 06:48:59 +08:00
										 |  |  | This page explains how to solve linear systems, compute various decompositions such as LU, | 
					
						
							|  |  |  | QR, %SVD, eigendecompositions... After reading this page, don't miss our | 
					
						
							|  |  |  | \link TopicLinearAlgebraDecompositions catalogue \endlink of dense matrix decompositions. | 
					
						
							| 
									
										
										
										
											2010-06-29 22:02:33 +08:00
										 |  |  | 
 | 
					
						
							| 
									
										
											  
											
												Big changes in Eigen documentation:
- Organize the documentation into "chapters".
  - Each chapter include many documentation pages, reference pages organized as modules, and a quick reference page.
  - The "Chapters" tree is created using the defgroup/ingroup mechanism, even for the documentation pages (i.e., .dox files for which I added an \eigenManualPage macro that we can switch between \page or \defgroup ).
  - Add a "General topics" entry for all pages that do not fit well in the previous "chapters".
  - The highlevel struture is managed by a new eigendoxy_layout.xml file.
- remove the "index" and quite useless pages (namespace list, class hierarchy, member list, file list, etc.)
- add the javascript search-engine.
- add the "treeview" panel.
- remove \tableofcontents (replace them by a custom \eigenAutoToc macro to be able to easily re-enable if needed).
- add javascript to automatically generate a TOC from the h1/h2 tags of the current page, and put the TOC in the left side panel.
- overload various javascript function generated by doxygen to:
  - remove the root of the treeview
  - remove links to section/subsection from the treeview
  - automatically expand the "Chapters" section
  - automatically expand the current section
  - adjust the height of the treeview to take into account the TOC
- always use the default .css file, eigendoxy.css now only includes our modifications
- use Doxyfile to specify our logo
- remove cross references to unsupported modules (temporarily)
											
										 
											2013-01-05 23:37:11 +08:00
										 |  |  | \eigenAutoToc | 
					
						
							| 
									
										
										
										
											2010-07-23 04:52:04 +08:00
										 |  |  | 
 | 
					
						
							| 
									
										
										
										
											2010-06-30 22:11:55 +08:00
										 |  |  | \section TutorialLinAlgBasicSolve Basic linear solving | 
					
						
							| 
									
										
										
										
											2010-06-29 22:02:33 +08:00
										 |  |  | 
 | 
					
						
							|  |  |  | \b The \b problem: You have a system of equations, that you have written as a single matrix equation | 
					
						
							|  |  |  |     \f[ Ax \: = \: b \f] | 
					
						
							|  |  |  | Where \a A and \a b are matrices (\a b could be a vector, as a special case). You want to find a solution \a x. | 
					
						
							|  |  |  | 
 | 
					
						
							| 
									
										
										
										
											2021-08-17 05:56:18 +08:00
										 |  |  | \b The \b solution: You can choose between various decompositions, depending on the properties of your matrix \a A, | 
					
						
							| 
									
										
										
										
											2010-06-29 22:02:33 +08:00
										 |  |  | and depending on whether you favor speed or accuracy. However, let's start with an example that works in all cases, | 
					
						
							|  |  |  | and is a good compromise: | 
					
						
							| 
									
										
										
										
											2010-10-19 17:40:49 +08:00
										 |  |  | <table class="example"> | 
					
						
							|  |  |  | <tr><th>Example:</th><th>Output:</th></tr> | 
					
						
							| 
									
										
										
										
											2010-06-29 22:02:33 +08:00
										 |  |  | <tr> | 
					
						
							|  |  |  |   <td>\include TutorialLinAlgExSolveColPivHouseholderQR.cpp </td> | 
					
						
							| 
									
										
										
										
											2010-10-19 17:40:49 +08:00
										 |  |  |   <td>\verbinclude TutorialLinAlgExSolveColPivHouseholderQR.out </td> | 
					
						
							| 
									
										
										
										
											2010-06-29 22:02:33 +08:00
										 |  |  | </tr> | 
					
						
							|  |  |  | </table> | 
					
						
							|  |  |  | 
 | 
					
						
							| 
									
										
										
										
											2010-06-30 22:11:55 +08:00
										 |  |  | In this example, the colPivHouseholderQr() method returns an object of class ColPivHouseholderQR. Since here the | 
					
						
							|  |  |  | matrix is of type Matrix3f, this line could have been replaced by: | 
					
						
							| 
									
										
										
										
											2010-06-29 22:02:33 +08:00
										 |  |  | \code | 
					
						
							| 
									
										
										
										
											2010-06-30 22:11:55 +08:00
										 |  |  | ColPivHouseholderQR<Matrix3f> dec(A); | 
					
						
							| 
									
										
										
										
											2010-06-29 22:02:33 +08:00
										 |  |  | Vector3f x = dec.solve(b); | 
					
						
							|  |  |  | \endcode | 
					
						
							|  |  |  | 
 | 
					
						
							|  |  |  | Here, ColPivHouseholderQR is a QR decomposition with column pivoting. It's a good compromise for this tutorial, as it | 
					
						
							|  |  |  | works for all matrices while being quite fast. Here is a table of some other decompositions that you can choose from, | 
					
						
							| 
									
										
										
										
											2021-08-17 05:56:18 +08:00
										 |  |  | depending on your matrix, the problem you are trying to solve, and the trade-off you want to make: | 
					
						
							| 
									
										
										
										
											2010-06-29 22:02:33 +08:00
										 |  |  | 
 | 
					
						
							| 
									
										
										
										
											2010-10-19 17:40:49 +08:00
										 |  |  | <table class="manual"> | 
					
						
							| 
									
										
										
										
											2010-06-29 22:02:33 +08:00
										 |  |  |     <tr> | 
					
						
							| 
									
										
										
										
											2010-10-19 17:40:49 +08:00
										 |  |  |         <th>Decomposition</th> | 
					
						
							|  |  |  |         <th>Method</th> | 
					
						
							| 
									
										
										
										
											2014-06-17 15:37:07 +08:00
										 |  |  |         <th>Requirements<br/>on the matrix</th> | 
					
						
							|  |  |  |         <th>Speed<br/> (small-to-medium)</th> | 
					
						
							|  |  |  |         <th>Speed<br/> (large)</th> | 
					
						
							| 
									
										
										
										
											2010-10-19 17:40:49 +08:00
										 |  |  |         <th>Accuracy</th> | 
					
						
							| 
									
										
										
										
											2010-06-29 22:02:33 +08:00
										 |  |  |     </tr> | 
					
						
							|  |  |  |     <tr> | 
					
						
							|  |  |  |         <td>PartialPivLU</td> | 
					
						
							|  |  |  |         <td>partialPivLu()</td> | 
					
						
							|  |  |  |         <td>Invertible</td> | 
					
						
							|  |  |  |         <td>++</td> | 
					
						
							| 
									
										
										
										
											2014-06-17 15:37:07 +08:00
										 |  |  |         <td>++</td> | 
					
						
							| 
									
										
										
										
											2010-06-29 22:02:33 +08:00
										 |  |  |         <td>+</td> | 
					
						
							|  |  |  |     </tr> | 
					
						
							| 
									
										
										
										
											2010-10-19 17:40:49 +08:00
										 |  |  |     <tr class="alt"> | 
					
						
							| 
									
										
										
										
											2010-06-29 22:02:33 +08:00
										 |  |  |         <td>FullPivLU</td> | 
					
						
							|  |  |  |         <td>fullPivLu()</td> | 
					
						
							|  |  |  |         <td>None</td> | 
					
						
							|  |  |  |         <td>-</td> | 
					
						
							| 
									
										
										
										
											2014-06-17 15:37:07 +08:00
										 |  |  |         <td>- -</td> | 
					
						
							| 
									
										
										
										
											2010-06-29 22:02:33 +08:00
										 |  |  |         <td>+++</td> | 
					
						
							|  |  |  |     </tr> | 
					
						
							|  |  |  |     <tr> | 
					
						
							|  |  |  |         <td>HouseholderQR</td> | 
					
						
							|  |  |  |         <td>householderQr()</td> | 
					
						
							|  |  |  |         <td>None</td> | 
					
						
							|  |  |  |         <td>++</td> | 
					
						
							| 
									
										
										
										
											2014-06-17 15:37:07 +08:00
										 |  |  |         <td>++</td> | 
					
						
							| 
									
										
										
										
											2010-06-29 22:02:33 +08:00
										 |  |  |         <td>+</td> | 
					
						
							|  |  |  |     </tr> | 
					
						
							| 
									
										
										
										
											2010-10-19 17:40:49 +08:00
										 |  |  |     <tr class="alt"> | 
					
						
							| 
									
										
										
										
											2010-06-29 22:02:33 +08:00
										 |  |  |         <td>ColPivHouseholderQR</td> | 
					
						
							|  |  |  |         <td>colPivHouseholderQr()</td> | 
					
						
							|  |  |  |         <td>None</td> | 
					
						
							| 
									
										
										
										
											2018-04-11 16:46:11 +08:00
										 |  |  |         <td>+</td> | 
					
						
							| 
									
										
										
										
											2014-06-17 15:37:07 +08:00
										 |  |  |         <td>-</td> | 
					
						
							|  |  |  |         <td>+++</td> | 
					
						
							| 
									
										
										
										
											2010-06-29 22:02:33 +08:00
										 |  |  |     </tr> | 
					
						
							|  |  |  |     <tr> | 
					
						
							|  |  |  |         <td>FullPivHouseholderQR</td> | 
					
						
							|  |  |  |         <td>fullPivHouseholderQr()</td> | 
					
						
							|  |  |  |         <td>None</td> | 
					
						
							|  |  |  |         <td>-</td> | 
					
						
							| 
									
										
										
										
											2014-06-17 15:37:07 +08:00
										 |  |  |         <td>- -</td> | 
					
						
							| 
									
										
										
										
											2010-06-29 22:02:33 +08:00
										 |  |  |         <td>+++</td> | 
					
						
							|  |  |  |     </tr> | 
					
						
							| 
									
										
										
										
											2018-04-11 16:46:11 +08:00
										 |  |  |     <tr class="alt"> | 
					
						
							|  |  |  |         <td>CompleteOrthogonalDecomposition</td> | 
					
						
							|  |  |  |         <td>completeOrthogonalDecomposition()</td> | 
					
						
							|  |  |  |         <td>None</td> | 
					
						
							|  |  |  |         <td>+</td> | 
					
						
							|  |  |  |         <td>-</td> | 
					
						
							|  |  |  |         <td>+++</td> | 
					
						
							|  |  |  |     </tr> | 
					
						
							| 
									
										
										
										
											2010-10-19 17:40:49 +08:00
										 |  |  |     <tr class="alt"> | 
					
						
							| 
									
										
										
										
											2010-06-29 22:02:33 +08:00
										 |  |  |         <td>LLT</td> | 
					
						
							|  |  |  |         <td>llt()</td> | 
					
						
							|  |  |  |         <td>Positive definite</td> | 
					
						
							|  |  |  |         <td>+++</td> | 
					
						
							| 
									
										
										
										
											2014-06-17 15:37:07 +08:00
										 |  |  |         <td>+++</td> | 
					
						
							| 
									
										
										
										
											2010-06-29 22:02:33 +08:00
										 |  |  |         <td>+</td> | 
					
						
							|  |  |  |     </tr> | 
					
						
							|  |  |  |     <tr> | 
					
						
							|  |  |  |         <td>LDLT</td> | 
					
						
							|  |  |  |         <td>ldlt()</td> | 
					
						
							| 
									
										
										
										
											2014-06-17 15:37:07 +08:00
										 |  |  |         <td>Positive or negative<br/> semidefinite</td> | 
					
						
							| 
									
										
										
										
											2010-06-29 22:02:33 +08:00
										 |  |  |         <td>+++</td> | 
					
						
							| 
									
										
										
										
											2014-06-17 15:37:07 +08:00
										 |  |  |         <td>+</td> | 
					
						
							| 
									
										
										
										
											2010-06-29 22:02:33 +08:00
										 |  |  |         <td>++</td> | 
					
						
							|  |  |  |     </tr> | 
					
						
							| 
									
										
										
										
											2018-04-11 16:46:11 +08:00
										 |  |  |     <tr class="alt"> | 
					
						
							|  |  |  |         <td>BDCSVD</td> | 
					
						
							|  |  |  |         <td>bdcSvd()</td> | 
					
						
							|  |  |  |         <td>None</td> | 
					
						
							|  |  |  |         <td>-</td> | 
					
						
							|  |  |  |         <td>-</td> | 
					
						
							|  |  |  |         <td>+++</td> | 
					
						
							|  |  |  |     </tr> | 
					
						
							| 
									
										
										
										
											2014-06-17 15:37:07 +08:00
										 |  |  |     <tr class="alt"> | 
					
						
							|  |  |  |         <td>JacobiSVD</td> | 
					
						
							|  |  |  |         <td>jacobiSvd()</td> | 
					
						
							|  |  |  |         <td>None</td> | 
					
						
							| 
									
										
										
										
											2018-04-11 16:46:11 +08:00
										 |  |  |         <td>-</td> | 
					
						
							| 
									
										
										
										
											2014-06-17 15:37:07 +08:00
										 |  |  |         <td>- - -</td> | 
					
						
							|  |  |  |         <td>+++</td> | 
					
						
							|  |  |  |     </tr> | 
					
						
							| 
									
										
										
										
											2010-06-29 22:02:33 +08:00
										 |  |  | </table> | 
					
						
							| 
									
										
										
										
											2018-04-11 16:46:11 +08:00
										 |  |  | To get an overview of the true relative speed of the different decompositions, check this \link DenseDecompositionBenchmark benchmark \endlink. | 
					
						
							| 
									
										
										
										
											2010-06-29 22:02:33 +08:00
										 |  |  | 
 | 
					
						
							| 
									
										
										
										
											2021-08-17 05:56:18 +08:00
										 |  |  | All of these decompositions offer a solve() method that works as in the above example.  | 
					
						
							| 
									
										
										
										
											2010-06-29 22:02:33 +08:00
										 |  |  | 
 | 
					
						
							| 
									
										
										
										
											2021-08-17 05:56:18 +08:00
										 |  |  | If you know more about the properties of your matrix, you can use the above table to select the best method. | 
					
						
							|  |  |  | For example, a good choice for solving linear systems with a non-symmetric matrix of full rank is PartialPivLU. | 
					
						
							|  |  |  | If you know that your matrix is also symmetric and positive definite, the above table says that | 
					
						
							|  |  |  | a very good choice is the LLT or LDLT decomposition. Here's an example, also demonstrating that using a general | 
					
						
							|  |  |  | matrix (not a vector) as right hand side is possible: | 
					
						
							| 
									
										
										
										
											2010-06-30 22:11:55 +08:00
										 |  |  | 
 | 
					
						
							| 
									
										
										
										
											2010-10-19 17:40:49 +08:00
										 |  |  | <table class="example"> | 
					
						
							|  |  |  | <tr><th>Example:</th><th>Output:</th></tr> | 
					
						
							| 
									
										
										
										
											2010-06-30 22:11:55 +08:00
										 |  |  | <tr> | 
					
						
							|  |  |  |   <td>\include TutorialLinAlgExSolveLDLT.cpp </td> | 
					
						
							| 
									
										
										
										
											2010-10-19 17:40:49 +08:00
										 |  |  |   <td>\verbinclude TutorialLinAlgExSolveLDLT.out </td> | 
					
						
							| 
									
										
										
										
											2010-06-30 22:11:55 +08:00
										 |  |  | </tr> | 
					
						
							|  |  |  | </table> | 
					
						
							|  |  |  | 
 | 
					
						
							|  |  |  | For a \ref TopicLinearAlgebraDecompositions "much more complete table" comparing all decompositions supported by Eigen (notice that Eigen | 
					
						
							| 
									
										
										
										
											2010-06-29 22:02:33 +08:00
										 |  |  | supports many other decompositions), see our special page on | 
					
						
							|  |  |  | \ref TopicLinearAlgebraDecompositions "this topic". | 
					
						
							|  |  |  | 
 | 
					
						
							| 
									
										
										
										
											2021-08-17 05:56:18 +08:00
										 |  |  | 
 | 
					
						
							|  |  |  | \section TutorialLinAlgLeastsquares Least squares solving | 
					
						
							|  |  |  | 
 | 
					
						
							|  |  |  | The most general and accurate method to solve under- or over-determined linear systems | 
					
						
							|  |  |  | in the least squares sense, is the SVD decomposition. Eigen provides two implementations. | 
					
						
							|  |  |  | The recommended one is the BDCSVD class, which scales well for large problems | 
					
						
							|  |  |  | and automatically falls back to the JacobiSVD class for smaller problems. | 
					
						
							|  |  |  | For both classes, their solve() method solved the linear system in the least-squares | 
					
						
							|  |  |  | sense.  | 
					
						
							|  |  |  | 
 | 
					
						
							|  |  |  | Here is an example: | 
					
						
							|  |  |  | <table class="example"> | 
					
						
							|  |  |  | <tr><th>Example:</th><th>Output:</th></tr> | 
					
						
							|  |  |  | <tr> | 
					
						
							|  |  |  |   <td>\include TutorialLinAlgSVDSolve.cpp </td> | 
					
						
							|  |  |  |   <td>\verbinclude TutorialLinAlgSVDSolve.out </td> | 
					
						
							|  |  |  | </tr> | 
					
						
							|  |  |  | </table> | 
					
						
							|  |  |  | 
 | 
					
						
							|  |  |  | An alternative to the SVD, which is usually faster and about as accurate, is CompleteOrthogonalDecomposition.  | 
					
						
							|  |  |  | 
 | 
					
						
							|  |  |  | Again, if you know more about the problem, the table above contains methods that are potentially faster. | 
					
						
							|  |  |  | If your matrix is full rank, HouseHolderQR is the method of choice. If your matrix is full rank and well conditioned, | 
					
						
							|  |  |  | using the Cholesky decomposition (LLT) on the matrix of the normal equations can be faster still. | 
					
						
							|  |  |  | Our page on \link LeastSquares least squares solving \endlink has more details. | 
					
						
							|  |  |  | 
 | 
					
						
							|  |  |  | 
 | 
					
						
							|  |  |  | \section TutorialLinAlgSolutionExists Checking if a matrix is singular | 
					
						
							| 
									
										
										
										
											2010-06-30 22:11:55 +08:00
										 |  |  | 
 | 
					
						
							|  |  |  | Only you know what error margin you want to allow for a solution to be considered valid. | 
					
						
							|  |  |  | So Eigen lets you do this computation for yourself, if you want to, as in this example: | 
					
						
							|  |  |  | 
 | 
					
						
							| 
									
										
										
										
											2010-10-19 17:40:49 +08:00
										 |  |  | <table class="example"> | 
					
						
							|  |  |  | <tr><th>Example:</th><th>Output:</th></tr> | 
					
						
							| 
									
										
										
										
											2010-06-30 22:11:55 +08:00
										 |  |  | <tr> | 
					
						
							|  |  |  |   <td>\include TutorialLinAlgExComputeSolveError.cpp </td> | 
					
						
							| 
									
										
										
										
											2010-10-19 17:40:49 +08:00
										 |  |  |   <td>\verbinclude TutorialLinAlgExComputeSolveError.out </td> | 
					
						
							| 
									
										
										
										
											2010-06-30 22:11:55 +08:00
										 |  |  | </tr> | 
					
						
							|  |  |  | </table> | 
					
						
							|  |  |  | 
 | 
					
						
							|  |  |  | \section TutorialLinAlgEigensolving Computing eigenvalues and eigenvectors | 
					
						
							|  |  |  | 
 | 
					
						
							|  |  |  | You need an eigendecomposition here, see available such decompositions on \ref TopicLinearAlgebraDecompositions "this page". | 
					
						
							|  |  |  | Make sure to check if your matrix is self-adjoint, as is often the case in these problems. Here's an example using | 
					
						
							|  |  |  | SelfAdjointEigenSolver, it could easily be adapted to general matrices using EigenSolver or ComplexEigenSolver. | 
					
						
							|  |  |  | 
 | 
					
						
							| 
									
										
										
										
											2011-11-08 01:14:06 +08:00
										 |  |  | The computation of eigenvalues and eigenvectors does not necessarily converge, but such failure to converge is | 
					
						
							|  |  |  | very rare. The call to info() is to check for this possibility. | 
					
						
							|  |  |  | 
 | 
					
						
							| 
									
										
										
										
											2010-10-19 17:40:49 +08:00
										 |  |  | <table class="example"> | 
					
						
							|  |  |  | <tr><th>Example:</th><th>Output:</th></tr> | 
					
						
							| 
									
										
										
										
											2010-06-30 22:11:55 +08:00
										 |  |  | <tr> | 
					
						
							|  |  |  |   <td>\include TutorialLinAlgSelfAdjointEigenSolver.cpp </td> | 
					
						
							| 
									
										
										
										
											2010-10-19 17:40:49 +08:00
										 |  |  |   <td>\verbinclude TutorialLinAlgSelfAdjointEigenSolver.out </td> | 
					
						
							| 
									
										
										
										
											2010-06-30 22:11:55 +08:00
										 |  |  | </tr> | 
					
						
							|  |  |  | </table> | 
					
						
							|  |  |  | 
 | 
					
						
							| 
									
										
										
										
											2010-07-23 04:52:04 +08:00
										 |  |  | \section TutorialLinAlgInverse Computing inverse and determinant | 
					
						
							| 
									
										
										
										
											2010-06-30 22:11:55 +08:00
										 |  |  | 
 | 
					
						
							|  |  |  | First of all, make sure that you really want this. While inverse and determinant are fundamental mathematical concepts, | 
					
						
							| 
									
										
										
										
											2021-08-17 05:56:18 +08:00
										 |  |  | in \em numerical linear algebra they are not as useful as in pure mathematics. Inverse computations are often | 
					
						
							| 
									
										
										
										
											2010-06-30 22:11:55 +08:00
										 |  |  | advantageously replaced by solve() operations, and the determinant is often \em not a good way of checking if a matrix | 
					
						
							|  |  |  | is invertible. | 
					
						
							|  |  |  | 
 | 
					
						
							| 
									
										
										
										
											2021-08-17 05:56:18 +08:00
										 |  |  | However, for \em very \em small matrices, the above may not be true, and inverse and determinant can be very useful. | 
					
						
							| 
									
										
										
										
											2010-06-30 22:11:55 +08:00
										 |  |  | 
 | 
					
						
							|  |  |  | While certain decompositions, such as PartialPivLU and FullPivLU, offer inverse() and determinant() methods, you can also | 
					
						
							|  |  |  | call inverse() and determinant() directly on a matrix. If your matrix is of a very small fixed size (at most 4x4) this | 
					
						
							|  |  |  | allows Eigen to avoid performing a LU decomposition, and instead use formulas that are more efficient on such small matrices. | 
					
						
							|  |  |  | 
 | 
					
						
							|  |  |  | Here is an example: | 
					
						
							| 
									
										
										
										
											2010-10-19 17:40:49 +08:00
										 |  |  | <table class="example"> | 
					
						
							|  |  |  | <tr><th>Example:</th><th>Output:</th></tr> | 
					
						
							| 
									
										
										
										
											2010-06-30 22:11:55 +08:00
										 |  |  | <tr> | 
					
						
							|  |  |  |   <td>\include TutorialLinAlgInverseDeterminant.cpp </td> | 
					
						
							| 
									
										
										
										
											2010-10-19 17:40:49 +08:00
										 |  |  |   <td>\verbinclude TutorialLinAlgInverseDeterminant.out </td> | 
					
						
							| 
									
										
										
										
											2010-06-30 22:11:55 +08:00
										 |  |  | </tr> | 
					
						
							|  |  |  | </table> | 
					
						
							|  |  |  | 
 | 
					
						
							|  |  |  | \section TutorialLinAlgSeparateComputation Separating the computation from the construction | 
					
						
							|  |  |  | 
 | 
					
						
							|  |  |  | In the above examples, the decomposition was computed at the same time that the decomposition object was constructed. | 
					
						
							|  |  |  | There are however situations where you might want to separate these two things, for example if you don't know, | 
					
						
							|  |  |  | at the time of the construction, the matrix that you will want to decompose; or if you want to reuse an existing | 
					
						
							|  |  |  | decomposition object. | 
					
						
							|  |  |  | 
 | 
					
						
							|  |  |  | What makes this possible is that: | 
					
						
							|  |  |  | \li all decompositions have a default constructor, | 
					
						
							|  |  |  | \li all decompositions have a compute(matrix) method that does the computation, and that may be called again | 
					
						
							|  |  |  |     on an already-computed decomposition, reinitializing it. | 
					
						
							|  |  |  | 
 | 
					
						
							|  |  |  | For example: | 
					
						
							|  |  |  | 
 | 
					
						
							| 
									
										
										
										
											2010-10-19 17:40:49 +08:00
										 |  |  | <table class="example"> | 
					
						
							|  |  |  | <tr><th>Example:</th><th>Output:</th></tr> | 
					
						
							| 
									
										
										
										
											2010-06-30 22:11:55 +08:00
										 |  |  | <tr> | 
					
						
							|  |  |  |   <td>\include TutorialLinAlgComputeTwice.cpp </td> | 
					
						
							| 
									
										
										
										
											2010-10-19 17:40:49 +08:00
										 |  |  |   <td>\verbinclude TutorialLinAlgComputeTwice.out </td> | 
					
						
							| 
									
										
										
										
											2010-06-30 22:11:55 +08:00
										 |  |  | </tr> | 
					
						
							|  |  |  | </table> | 
					
						
							|  |  |  | 
 | 
					
						
							|  |  |  | Finally, you can tell the decomposition constructor to preallocate storage for decomposing matrices of a given size, | 
					
						
							|  |  |  | so that when you subsequently decompose such matrices, no dynamic memory allocation is performed (of course, if you | 
					
						
							|  |  |  | are using fixed-size matrices, no dynamic memory allocation happens at all). This is done by just | 
					
						
							|  |  |  | passing the size to the decomposition constructor, as in this example: | 
					
						
							|  |  |  | \code | 
					
						
							|  |  |  | HouseholderQR<MatrixXf> qr(50,50); | 
					
						
							|  |  |  | MatrixXf A = MatrixXf::Random(50,50); | 
					
						
							|  |  |  | qr.compute(A); // no dynamic memory allocation | 
					
						
							|  |  |  | \endcode | 
					
						
							|  |  |  | 
 | 
					
						
							|  |  |  | \section TutorialLinAlgRankRevealing Rank-revealing decompositions | 
					
						
							|  |  |  | 
 | 
					
						
							|  |  |  | Certain decompositions are rank-revealing, i.e. are able to compute the rank of a matrix. These are typically | 
					
						
							|  |  |  | also the decompositions that behave best in the face of a non-full-rank matrix (which in the square case means a | 
					
						
							|  |  |  | singular matrix). On \ref TopicLinearAlgebraDecompositions "this table" you can see for all our decompositions | 
					
						
							|  |  |  | whether they are rank-revealing or not. | 
					
						
							|  |  |  | 
 | 
					
						
							|  |  |  | Rank-revealing decompositions offer at least a rank() method. They can also offer convenience methods such as isInvertible(), | 
					
						
							|  |  |  | and some are also providing methods to compute the kernel (null-space) and image (column-space) of the matrix, as is the | 
					
						
							|  |  |  | case with FullPivLU: | 
					
						
							|  |  |  | 
 | 
					
						
							| 
									
										
										
										
											2010-10-19 17:40:49 +08:00
										 |  |  | <table class="example"> | 
					
						
							|  |  |  | <tr><th>Example:</th><th>Output:</th></tr> | 
					
						
							| 
									
										
										
										
											2010-06-30 22:11:55 +08:00
										 |  |  | <tr> | 
					
						
							|  |  |  |   <td>\include TutorialLinAlgRankRevealing.cpp </td> | 
					
						
							| 
									
										
										
										
											2010-10-19 17:40:49 +08:00
										 |  |  |   <td>\verbinclude TutorialLinAlgRankRevealing.out </td> | 
					
						
							| 
									
										
										
										
											2010-06-30 22:11:55 +08:00
										 |  |  | </tr> | 
					
						
							|  |  |  | </table> | 
					
						
							|  |  |  | 
 | 
					
						
							|  |  |  | Of course, any rank computation depends on the choice of an arbitrary threshold, since practically no | 
					
						
							|  |  |  | floating-point matrix is \em exactly rank-deficient. Eigen picks a sensible default threshold, which depends | 
					
						
							|  |  |  | on the decomposition but is typically the diagonal size times machine epsilon. While this is the best default we | 
					
						
							|  |  |  | could pick, only you know what is the right threshold for your application. You can set this by calling setThreshold() | 
					
						
							| 
									
										
										
										
											2010-07-01 07:27:30 +08:00
										 |  |  | on your decomposition object before calling rank() or any other method that needs to use such a threshold. | 
					
						
							|  |  |  | The decomposition itself, i.e. the compute() method, is independent of the threshold. You don't need to recompute the | 
					
						
							|  |  |  | decomposition after you've changed the threshold. | 
					
						
							| 
									
										
										
										
											2010-06-30 22:11:55 +08:00
										 |  |  | 
 | 
					
						
							| 
									
										
										
										
											2010-10-19 17:40:49 +08:00
										 |  |  | <table class="example"> | 
					
						
							|  |  |  | <tr><th>Example:</th><th>Output:</th></tr> | 
					
						
							| 
									
										
										
										
											2010-06-30 22:11:55 +08:00
										 |  |  | <tr> | 
					
						
							|  |  |  |   <td>\include TutorialLinAlgSetThreshold.cpp </td> | 
					
						
							| 
									
										
										
										
											2010-10-19 17:40:49 +08:00
										 |  |  |   <td>\verbinclude TutorialLinAlgSetThreshold.out </td> | 
					
						
							| 
									
										
										
										
											2010-06-30 22:11:55 +08:00
										 |  |  | </tr> | 
					
						
							|  |  |  | </table> | 
					
						
							|  |  |  | 
 | 
					
						
							| 
									
										
										
										
											2010-06-29 22:02:33 +08:00
										 |  |  | */ | 
					
						
							|  |  |  | 
 | 
					
						
							|  |  |  | } |