Sophie

Sophie

distrib > Mageia > 4 > x86_64 > by-pkgid > 99cb5ede6a5329071fbeecc8218deb35 > files > 1091

eigen3-doc-3.2-3.mga4.noarch.rpm

<!DOCTYPE html PUBLIC "-//W3C//DTD XHTML 1.0 Transitional//EN" "http://www.w3.org/TR/xhtml1/DTD/xhtml1-transitional.dtd">
<html xmlns="http://www.w3.org/1999/xhtml">
<head>
<meta http-equiv="Content-Type" content="text/xhtml;charset=UTF-8"/>
<meta http-equiv="X-UA-Compatible" content="IE=9"/>
<meta name="generator" content="Doxygen 1.8.5"/>
<title>Eigen: Linear algebra and decompositions</title>
<link href="tabs.css" rel="stylesheet" type="text/css"/>
<script type="text/javascript" src="jquery.js"></script>
<script type="text/javascript" src="dynsections.js"></script>
<link href="navtree.css" rel="stylesheet" type="text/css"/>
<script type="text/javascript" src="resize.js"></script>
<script type="text/javascript" src="navtree.js"></script>
<script type="text/javascript">
  $(document).ready(initResizable);
  $(window).load(resizeHeight);
</script>
<link href="search/search.css" rel="stylesheet" type="text/css"/>
<script type="text/javascript" src="search/search.js"></script>
<script type="text/javascript">
  $(document).ready(function() { searchBox.OnSelectItem(0); });
</script>
<link href="doxygen.css"   rel="stylesheet" type="text/css" />
<link href="eigendoxy.css" rel="stylesheet" type="text/css">
<!--  -->
<script type="text/javascript" src="eigen_navtree_hacks.js"></script>
<!-- <script type="text/javascript"> -->
<!-- </script> -->
</head>
<body>
<div id="top"><!-- do not remove this div, it is closed by doxygen! -->
<!-- <a name="top"></a> -->
<div id="titlearea">
<table cellspacing="0" cellpadding="0">
 <tbody>
 <tr style="height: 56px;">
  <td id="projectlogo"><img alt="Logo" src="Eigen_Silly_Professor_64x64.png"/></td>
  <td style="padding-left: 0.5em;">
   <div id="projectname"><a href="http://eigen.tuxfamily.org">Eigen</a>
   &#160;<span id="projectnumber">3.2.0</span>
   </div>
  </td>
   <td>        <div id="MSearchBox" class="MSearchBoxInactive">
        <span class="left">
          <img id="MSearchSelect" src="search/mag_sel.png"
               onmouseover="return searchBox.OnSearchSelectShow()"
               onmouseout="return searchBox.OnSearchSelectHide()"
               alt=""/>
          <input type="text" id="MSearchField" value="Search" accesskey="S"
               onfocus="searchBox.OnSearchFieldFocus(true)" 
               onblur="searchBox.OnSearchFieldFocus(false)" 
               onkeyup="searchBox.OnSearchFieldChange(event)"/>
          </span><span class="right">
            <a id="MSearchClose" href="javascript:searchBox.CloseResultsWindow()"><img id="MSearchCloseImg" border="0" src="search/close.png" alt=""/></a>
          </span>
        </div>
</td>
 </tr>
 </tbody>
</table>
</div>
<!-- end header part -->
<!-- Generated by Doxygen 1.8.5 -->
<script type="text/javascript">
var searchBox = new SearchBox("searchBox", "search",false,'Search');
</script>
</div><!-- top -->
<div id="side-nav" class="ui-resizable side-nav-resizable">
  <div id="nav-tree">
    <div id="nav-tree-contents">
      <div id="nav-sync" class="sync"></div>
    </div>
  </div>
  <div id="splitbar" style="-moz-user-select:none;" 
       class="ui-resizable-handle">
  </div>
</div>
<script type="text/javascript">
$(document).ready(function(){initNavTree('group__TutorialLinearAlgebra.html','');});
</script>
<div id="doc-content">
<!-- window showing the filter options -->
<div id="MSearchSelectWindow"
     onmouseover="return searchBox.OnSearchSelectShow()"
     onmouseout="return searchBox.OnSearchSelectHide()"
     onkeydown="return searchBox.OnSearchSelectKey(event)">
<a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(0)"><span class="SelectionMark">&#160;</span>All</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(1)"><span class="SelectionMark">&#160;</span>Classes</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(2)"><span class="SelectionMark">&#160;</span>Namespaces</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(3)"><span class="SelectionMark">&#160;</span>Functions</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(4)"><span class="SelectionMark">&#160;</span>Variables</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(5)"><span class="SelectionMark">&#160;</span>Typedefs</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(6)"><span class="SelectionMark">&#160;</span>Enumerations</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(7)"><span class="SelectionMark">&#160;</span>Enumerator</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(8)"><span class="SelectionMark">&#160;</span>Friends</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(9)"><span class="SelectionMark">&#160;</span>Groups</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(10)"><span class="SelectionMark">&#160;</span>Pages</a></div>

<!-- iframe showing the search results (closed by default) -->
<div id="MSearchResultsWindow">
<iframe src="javascript:void(0)" frameborder="0" 
        name="MSearchResults" id="MSearchResults">
</iframe>
</div>

<div class="header">
  <div class="headertitle">
<div class="title">Linear algebra and decompositions<div class="ingroups"><a class="el" href="group__DenseLinearSolvers__chapter.html">Dense linear problems and decompositions</a></div></div>  </div>
</div><!--header-->
<div class="contents">
<p>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 <a class="el" href="group__TopicLinearAlgebraDecompositions.html">catalogue </a> of dense matrix decompositions.</p>
<h1><a class="anchor" id="TutorialLinAlgBasicSolve"></a>
Basic linear solving</h1>
<p><b>The</b> <b>problem:</b> You have a system of equations, that you have written as a single matrix equation </p>
<p class="formulaDsp">
<img class="formulaDsp" alt="\[ Ax \: = \: b \]" src="form_185.png"/>
</p>
<p> Where <em>A</em> and <em>b</em> are matrices (<em>b</em> could be a vector, as a special case). You want to find a solution <em>x</em>.</p>
<p><b>The</b> <b>solution:</b> You can choose between various decompositions, depending on what your matrix <em>A</em> looks like, 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: </p>
<table  class="example">
<tr>
<th>Example:</th><th>Output: </th></tr>
<tr>
<td><div class="fragment"><div class="line"><span class="preprocessor">#include &lt;iostream&gt;</span></div>
<div class="line"><span class="preprocessor">#include &lt;Eigen/Dense&gt;</span></div>
<div class="line"></div>
<div class="line"><span class="keyword">using namespace </span>std;</div>
<div class="line"><span class="keyword">using namespace </span>Eigen;</div>
<div class="line"></div>
<div class="line"><span class="keywordtype">int</span> main()</div>
<div class="line">{</div>
<div class="line">   <a class="code" href="classEigen_1_1Matrix.html">Matrix3f</a> A;</div>
<div class="line">   <a class="code" href="classEigen_1_1Matrix.html">Vector3f</a> b;</div>
<div class="line">   A &lt;&lt; 1,2,3,  4,5,6,  7,8,10;</div>
<div class="line">   b &lt;&lt; 3, 3, 4;</div>
<div class="line">   cout &lt;&lt; <span class="stringliteral">&quot;Here is the matrix A:\n&quot;</span> &lt;&lt; A &lt;&lt; endl;</div>
<div class="line">   cout &lt;&lt; <span class="stringliteral">&quot;Here is the vector b:\n&quot;</span> &lt;&lt; b &lt;&lt; endl;</div>
<div class="line">   <a class="code" href="classEigen_1_1Matrix.html">Vector3f</a> x = A.<a class="code" href="classEigen_1_1MatrixBase.html#a05afed751d3a7277951d1918468e0872">colPivHouseholderQr</a>().solve(b);</div>
<div class="line">   cout &lt;&lt; <span class="stringliteral">&quot;The solution is:\n&quot;</span> &lt;&lt; x &lt;&lt; endl;</div>
<div class="line">}</div>
</div><!-- fragment -->  </td><td><pre class="fragment">Here is the matrix A:
 1  2  3
 4  5  6
 7  8 10
Here is the vector b:
3
3
4
The solution is:
-2
1
1
</pre>   </td></tr>
</table>
<p>In this example, the colPivHouseholderQr() method returns an object of class <a class="el" href="classEigen_1_1ColPivHouseholderQR.html" title="Householder rank-revealing QR decomposition of a matrix with column-pivoting. ">ColPivHouseholderQR</a>. Since here the matrix is of type Matrix3f, this line could have been replaced by: </p>
<div class="fragment"><div class="line">ColPivHouseholderQR&lt;Matrix3f&gt; dec(A);</div>
<div class="line"><a class="code" href="group__matrixtypedefs.html#ga5ac9fb0df0c4858477890cce1f998096">Vector3f</a> x = dec.solve(b);</div>
</div><!-- fragment --><p>Here, <a class="el" href="classEigen_1_1ColPivHouseholderQR.html" title="Householder rank-revealing QR decomposition of a matrix with column-pivoting. ">ColPivHouseholderQR</a> 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, depending on your matrix and the trade-off you want to make:</p>
<table  class="manual">
<tr>
<th>Decomposition </th><th>Method </th><th>Requirements on the matrix </th><th>Speed </th><th>Accuracy  </th></tr>
<tr>
<td><a class="el" href="classEigen_1_1PartialPivLU.html" title="LU decomposition of a matrix with partial pivoting, and related features. ">PartialPivLU</a> </td><td>partialPivLu() </td><td>Invertible </td><td>++ </td><td>+  </td></tr>
<tr class="alt">
<td><a class="el" href="classEigen_1_1FullPivLU.html" title="LU decomposition of a matrix with complete pivoting, and related features. ">FullPivLU</a> </td><td>fullPivLu() </td><td>None </td><td>- </td><td>+++  </td></tr>
<tr>
<td><a class="el" href="classEigen_1_1HouseholderQR.html" title="Householder QR decomposition of a matrix. ">HouseholderQR</a> </td><td>householderQr() </td><td>None </td><td>++ </td><td>+  </td></tr>
<tr class="alt">
<td><a class="el" href="classEigen_1_1ColPivHouseholderQR.html" title="Householder rank-revealing QR decomposition of a matrix with column-pivoting. ">ColPivHouseholderQR</a> </td><td>colPivHouseholderQr() </td><td>None </td><td>+ </td><td>++  </td></tr>
<tr>
<td><a class="el" href="classEigen_1_1FullPivHouseholderQR.html" title="Householder rank-revealing QR decomposition of a matrix with full pivoting. ">FullPivHouseholderQR</a> </td><td>fullPivHouseholderQr() </td><td>None </td><td>- </td><td>+++  </td></tr>
<tr class="alt">
<td><a class="el" href="classEigen_1_1LLT.html" title="Standard Cholesky decomposition (LL^T) of a matrix and associated features. ">LLT</a> </td><td>llt() </td><td>Positive definite </td><td>+++ </td><td>+  </td></tr>
<tr>
<td><a class="el" href="classEigen_1_1LDLT.html" title="Robust Cholesky decomposition of a matrix with pivoting. ">LDLT</a> </td><td>ldlt() </td><td>Positive or negative semidefinite </td><td>+++ </td><td>++  </td></tr>
</table>
<p>All of these decompositions offer a solve() method that works as in the above example.</p>
<p>For example, if your matrix is positive definite, the above table says that a very good choice is then the <a class="el" href="classEigen_1_1LDLT.html" title="Robust Cholesky decomposition of a matrix with pivoting. ">LDLT</a> decomposition. Here's an example, also demonstrating that using a general matrix (not a vector) as right hand side is possible.</p>
<table  class="example">
<tr>
<th>Example:</th><th>Output: </th></tr>
<tr>
<td><div class="fragment"><div class="line"><span class="preprocessor">#include &lt;iostream&gt;</span></div>
<div class="line"><span class="preprocessor">#include &lt;Eigen/Dense&gt;</span></div>
<div class="line"></div>
<div class="line"><span class="keyword">using namespace </span>std;</div>
<div class="line"><span class="keyword">using namespace </span>Eigen;</div>
<div class="line"></div>
<div class="line"><span class="keywordtype">int</span> main()</div>
<div class="line">{</div>
<div class="line">   <a class="code" href="classEigen_1_1Matrix.html">Matrix2f</a> A, b;</div>
<div class="line">   A &lt;&lt; 2, -1, -1, 3;</div>
<div class="line">   b &lt;&lt; 1, 2, 3, 1;</div>
<div class="line">   cout &lt;&lt; <span class="stringliteral">&quot;Here is the matrix A:\n&quot;</span> &lt;&lt; A &lt;&lt; endl;</div>
<div class="line">   cout &lt;&lt; <span class="stringliteral">&quot;Here is the right hand side b:\n&quot;</span> &lt;&lt; b &lt;&lt; endl;</div>
<div class="line">   <a class="code" href="classEigen_1_1Matrix.html">Matrix2f</a> x = A.<a class="code" href="classEigen_1_1MatrixBase.html#a197a04cda6b4606ec2416fd3f950371f">ldlt</a>().solve(b);</div>
<div class="line">   cout &lt;&lt; <span class="stringliteral">&quot;The solution is:\n&quot;</span> &lt;&lt; x &lt;&lt; endl;</div>
<div class="line">}</div>
</div><!-- fragment -->  </td><td><pre class="fragment">Here is the matrix A:
 2 -1
-1  3
Here is the right hand side b:
1 2
3 1
The solution is:
1.2 1.4
1.4 0.8
</pre>   </td></tr>
</table>
<p>For a <a class="el" href="group__TopicLinearAlgebraDecompositions.html">much more complete table</a> comparing all decompositions supported by <a class="el" href="namespaceEigen.html" title="Namespace containing all symbols from the Eigen library. ">Eigen</a> (notice that <a class="el" href="namespaceEigen.html" title="Namespace containing all symbols from the Eigen library. ">Eigen</a> supports many other decompositions), see our special page on <a class="el" href="group__TopicLinearAlgebraDecompositions.html">this topic</a>.</p>
<h1><a class="anchor" id="TutorialLinAlgSolutionExists"></a>
Checking if a solution really exists</h1>
<p>Only you know what error margin you want to allow for a solution to be considered valid. So <a class="el" href="namespaceEigen.html" title="Namespace containing all symbols from the Eigen library. ">Eigen</a> lets you do this computation for yourself, if you want to, as in this example:</p>
<table  class="example">
<tr>
<th>Example:</th><th>Output: </th></tr>
<tr>
<td><div class="fragment"><div class="line"><span class="preprocessor">#include &lt;iostream&gt;</span></div>
<div class="line"><span class="preprocessor">#include &lt;Eigen/Dense&gt;</span></div>
<div class="line"></div>
<div class="line"><span class="keyword">using namespace </span>std;</div>
<div class="line"><span class="keyword">using namespace </span>Eigen;</div>
<div class="line"></div>
<div class="line"><span class="keywordtype">int</span> main()</div>
<div class="line">{</div>
<div class="line">   <a class="code" href="classEigen_1_1Matrix.html">MatrixXd</a> A = <a class="code" href="classEigen_1_1DenseBase.html#a8e759dafdd9ecc446d397b7f5435f60a">MatrixXd::Random</a>(100,100);</div>
<div class="line">   <a class="code" href="classEigen_1_1Matrix.html">MatrixXd</a> b = <a class="code" href="classEigen_1_1DenseBase.html#a8e759dafdd9ecc446d397b7f5435f60a">MatrixXd::Random</a>(100,50);</div>
<div class="line">   <a class="code" href="classEigen_1_1Matrix.html">MatrixXd</a> x = A.<a class="code" href="classEigen_1_1MatrixBase.html#a0844e94f8f95ae01a2cd88dbbf5cbf91">fullPivLu</a>().solve(b);</div>
<div class="line">   <span class="keywordtype">double</span> relative_error = (A*x - b).norm() / b.<a class="code" href="classEigen_1_1MatrixBase.html#a0be1b433c65ce9d92c81a4718daf54e5">norm</a>(); <span class="comment">// norm() is L2 norm</span></div>
<div class="line">   cout &lt;&lt; <span class="stringliteral">&quot;The relative error is:\n&quot;</span> &lt;&lt; relative_error &lt;&lt; endl;</div>
<div class="line">}</div>
</div><!-- fragment -->  </td><td><pre class="fragment">The relative error is:
2.63205e-14
</pre>   </td></tr>
</table>
<h1><a class="anchor" id="TutorialLinAlgEigensolving"></a>
Computing eigenvalues and eigenvectors</h1>
<p>You need an eigendecomposition here, see available such decompositions on <a class="el" href="group__TopicLinearAlgebraDecompositions.html">this page</a>. Make sure to check if your matrix is self-adjoint, as is often the case in these problems. Here's an example using <a class="el" href="classEigen_1_1SelfAdjointEigenSolver.html" title="Computes eigenvalues and eigenvectors of selfadjoint matrices. ">SelfAdjointEigenSolver</a>, it could easily be adapted to general matrices using <a class="el" href="classEigen_1_1EigenSolver.html" title="Computes eigenvalues and eigenvectors of general matrices. ">EigenSolver</a> or <a class="el" href="classEigen_1_1ComplexEigenSolver.html" title="Computes eigenvalues and eigenvectors of general complex matrices. ">ComplexEigenSolver</a>.</p>
<p>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.</p>
<table  class="example">
<tr>
<th>Example:</th><th>Output: </th></tr>
<tr>
<td><div class="fragment"><div class="line"><span class="preprocessor">#include &lt;iostream&gt;</span></div>
<div class="line"><span class="preprocessor">#include &lt;Eigen/Dense&gt;</span></div>
<div class="line"></div>
<div class="line"><span class="keyword">using namespace </span>std;</div>
<div class="line"><span class="keyword">using namespace </span>Eigen;</div>
<div class="line"></div>
<div class="line"><span class="keywordtype">int</span> main()</div>
<div class="line">{</div>
<div class="line">   <a class="code" href="classEigen_1_1Matrix.html">Matrix2f</a> A;</div>
<div class="line">   A &lt;&lt; 1, 2, 2, 3;</div>
<div class="line">   cout &lt;&lt; <span class="stringliteral">&quot;Here is the matrix A:\n&quot;</span> &lt;&lt; A &lt;&lt; endl;</div>
<div class="line">   <a class="code" href="classEigen_1_1SelfAdjointEigenSolver.html">SelfAdjointEigenSolver&lt;Matrix2f&gt;</a> eigensolver(A);</div>
<div class="line">   <span class="keywordflow">if</span> (eigensolver.info() != <a class="code" href="group__enums.html#gga51bc1ac16f26ebe51eae1abb77bd037bafdfbdf3247bd36a1f17270d5cec74c9c">Success</a>) abort();</div>
<div class="line">   cout &lt;&lt; <span class="stringliteral">&quot;The eigenvalues of A are:\n&quot;</span> &lt;&lt; eigensolver.eigenvalues() &lt;&lt; endl;</div>
<div class="line">   cout &lt;&lt; <span class="stringliteral">&quot;Here&#39;s a matrix whose columns are eigenvectors of A \n&quot;</span></div>
<div class="line">        &lt;&lt; <span class="stringliteral">&quot;corresponding to these eigenvalues:\n&quot;</span></div>
<div class="line">        &lt;&lt; eigensolver.eigenvectors() &lt;&lt; endl;</div>
<div class="line">}</div>
</div><!-- fragment -->  </td><td><pre class="fragment">Here is the matrix A:
1 2
2 3
The eigenvalues of A are:
-0.236
4.24
Here's a matrix whose columns are eigenvectors of A 
corresponding to these eigenvalues:
-0.851 -0.526
 0.526 -0.851
</pre>   </td></tr>
</table>
<h1><a class="anchor" id="TutorialLinAlgInverse"></a>
Computing inverse and determinant</h1>
<p>First of all, make sure that you really want this. While inverse and determinant are fundamental mathematical concepts, in <em>numerical</em> linear algebra they are not as popular as in pure mathematics. Inverse computations are often advantageously replaced by solve() operations, and the determinant is often <em>not</em> a good way of checking if a matrix is invertible.</p>
<p>However, for <em>very</em> <em>small</em> matrices, the above is not true, and inverse and determinant can be very useful.</p>
<p>While certain decompositions, such as <a class="el" href="classEigen_1_1PartialPivLU.html" title="LU decomposition of a matrix with partial pivoting, and related features. ">PartialPivLU</a> and <a class="el" href="classEigen_1_1FullPivLU.html" title="LU decomposition of a matrix with complete pivoting, and related features. ">FullPivLU</a>, 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 <a class="el" href="namespaceEigen.html" title="Namespace containing all symbols from the Eigen library. ">Eigen</a> to avoid performing a LU decomposition, and instead use formulas that are more efficient on such small matrices.</p>
<p>Here is an example: </p>
<table  class="example">
<tr>
<th>Example:</th><th>Output: </th></tr>
<tr>
<td><div class="fragment"><div class="line"><span class="preprocessor">#include &lt;iostream&gt;</span></div>
<div class="line"><span class="preprocessor">#include &lt;Eigen/Dense&gt;</span></div>
<div class="line"></div>
<div class="line"><span class="keyword">using namespace </span>std;</div>
<div class="line"><span class="keyword">using namespace </span>Eigen;</div>
<div class="line"></div>
<div class="line"><span class="keywordtype">int</span> main()</div>
<div class="line">{</div>
<div class="line">   <a class="code" href="classEigen_1_1Matrix.html">Matrix3f</a> A;</div>
<div class="line">   A &lt;&lt; 1, 2, 1,</div>
<div class="line">        2, 1, 0,</div>
<div class="line">        -1, 1, 2;</div>
<div class="line">   cout &lt;&lt; <span class="stringliteral">&quot;Here is the matrix A:\n&quot;</span> &lt;&lt; A &lt;&lt; endl;</div>
<div class="line">   cout &lt;&lt; <span class="stringliteral">&quot;The determinant of A is &quot;</span> &lt;&lt; A.<a class="code" href="classEigen_1_1MatrixBase.html#ad63cea11a4bf220298dce4489a1704c7">determinant</a>() &lt;&lt; endl;</div>
<div class="line">   cout &lt;&lt; <span class="stringliteral">&quot;The inverse of A is:\n&quot;</span> &lt;&lt; A.<a class="code" href="classEigen_1_1MatrixBase.html#aa2834da4c855fa35fed8c4030f79f9da">inverse</a>() &lt;&lt; endl;</div>
<div class="line">}</div>
</div><!-- fragment -->  </td><td><pre class="fragment">Here is the matrix A:
1 2 1
2 1 0
-1 1 2
The determinant of A is -3
The inverse of A is:
-0.667      1  0.333
  1.33     -1 -0.667
    -1      1      1
</pre>   </td></tr>
</table>
<h1><a class="anchor" id="TutorialLinAlgLeastsquares"></a>
Least squares solving</h1>
<p>The best way to do least squares solving is with a SVD decomposition. <a class="el" href="namespaceEigen.html" title="Namespace containing all symbols from the Eigen library. ">Eigen</a> provides one as the <a class="el" href="classEigen_1_1JacobiSVD.html" title="Two-sided Jacobi SVD decomposition of a rectangular matrix. ">JacobiSVD</a> class, and its solve() is doing least-squares solving.</p>
<p>Here is an example: </p>
<table  class="example">
<tr>
<th>Example:</th><th>Output: </th></tr>
<tr>
<td><div class="fragment"><div class="line"><span class="preprocessor">#include &lt;iostream&gt;</span></div>
<div class="line"><span class="preprocessor">#include &lt;Eigen/Dense&gt;</span></div>
<div class="line"></div>
<div class="line"><span class="keyword">using namespace </span>std;</div>
<div class="line"><span class="keyword">using namespace </span>Eigen;</div>
<div class="line"></div>
<div class="line"><span class="keywordtype">int</span> main()</div>
<div class="line">{</div>
<div class="line">   <a class="code" href="classEigen_1_1Matrix.html">MatrixXf</a> A = <a class="code" href="classEigen_1_1DenseBase.html#a8e759dafdd9ecc446d397b7f5435f60a">MatrixXf::Random</a>(3, 2);</div>
<div class="line">   cout &lt;&lt; <span class="stringliteral">&quot;Here is the matrix A:\n&quot;</span> &lt;&lt; A &lt;&lt; endl;</div>
<div class="line">   VectorXf b = <a class="code" href="classEigen_1_1DenseBase.html#a8e759dafdd9ecc446d397b7f5435f60a">VectorXf::Random</a>(3);</div>
<div class="line">   cout &lt;&lt; <span class="stringliteral">&quot;Here is the right hand side b:\n&quot;</span> &lt;&lt; b &lt;&lt; endl;</div>
<div class="line">   cout &lt;&lt; <span class="stringliteral">&quot;The least-squares solution is:\n&quot;</span></div>
<div class="line">        &lt;&lt; A.<a class="code" href="classEigen_1_1MatrixBase.html#a27a6a2e9d8839662a29cc8294398ff97">jacobiSvd</a>(<a class="code" href="group__enums.html#gga2d78499b99ddc29b9494f7ea33864d52aa1954d61a33cbb0e9960ca88085ec487">ComputeThinU</a> | <a class="code" href="group__enums.html#gga2d78499b99ddc29b9494f7ea33864d52aa0bec1072bc1ecbbaaa436f96fe02e78">ComputeThinV</a>).solve(b) &lt;&lt; endl;</div>
<div class="line">}</div>
</div><!-- fragment -->  </td><td><pre class="fragment">Here is the matrix A:
  0.68  0.597
-0.211  0.823
 0.566 -0.605
Here is the right hand side b:
-0.33
0.536
-0.444
The least-squares solution is:
-0.67
0.314
</pre>   </td></tr>
</table>
<p>Another way, potentially faster but less reliable, is to use a <a class="el" href="classEigen_1_1LDLT.html" title="Robust Cholesky decomposition of a matrix with pivoting. ">LDLT</a> decomposition of the normal matrix. In any case, just read any reference text on least squares, and it will be very easy for you to implement any linear least squares computation on top of <a class="el" href="namespaceEigen.html" title="Namespace containing all symbols from the Eigen library. ">Eigen</a>.</p>
<h1><a class="anchor" id="TutorialLinAlgSeparateComputation"></a>
Separating the computation from the construction</h1>
<p>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.</p>
<p>What makes this possible is that: </p>
<ul>
<li>all decompositions have a default constructor, </li>
<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.</li>
</ul>
<p>For example:</p>
<table  class="example">
<tr>
<th>Example:</th><th>Output: </th></tr>
<tr>
<td><div class="fragment"><div class="line"><span class="preprocessor">#include &lt;iostream&gt;</span></div>
<div class="line"><span class="preprocessor">#include &lt;Eigen/Dense&gt;</span></div>
<div class="line"></div>
<div class="line"><span class="keyword">using namespace </span>std;</div>
<div class="line"><span class="keyword">using namespace </span>Eigen;</div>
<div class="line"></div>
<div class="line"><span class="keywordtype">int</span> main()</div>
<div class="line">{</div>
<div class="line">   <a class="code" href="classEigen_1_1Matrix.html">Matrix2f</a> A, b;</div>
<div class="line">   <a class="code" href="classEigen_1_1LLT.html">LLT&lt;Matrix2f&gt;</a> llt;</div>
<div class="line">   A &lt;&lt; 2, -1, -1, 3;</div>
<div class="line">   b &lt;&lt; 1, 2, 3, 1;</div>
<div class="line">   cout &lt;&lt; <span class="stringliteral">&quot;Here is the matrix A:\n&quot;</span> &lt;&lt; A &lt;&lt; endl;</div>
<div class="line">   cout &lt;&lt; <span class="stringliteral">&quot;Here is the right hand side b:\n&quot;</span> &lt;&lt; b &lt;&lt; endl;</div>
<div class="line">   cout &lt;&lt; <span class="stringliteral">&quot;Computing LLT decomposition...&quot;</span> &lt;&lt; endl;</div>
<div class="line">   llt.<a class="code" href="classEigen_1_1LLT.html#ac97b2f0b61df78a30544670b3711e3db">compute</a>(A);</div>
<div class="line">   cout &lt;&lt; <span class="stringliteral">&quot;The solution is:\n&quot;</span> &lt;&lt; llt.<a class="code" href="classEigen_1_1LLT.html#a301e28d205fb0308cb4f4b04718e1685">solve</a>(b) &lt;&lt; endl;</div>
<div class="line">   A(1,1)++;</div>
<div class="line">   cout &lt;&lt; <span class="stringliteral">&quot;The matrix A is now:\n&quot;</span> &lt;&lt; A &lt;&lt; endl;</div>
<div class="line">   cout &lt;&lt; <span class="stringliteral">&quot;Computing LLT decomposition...&quot;</span> &lt;&lt; endl;</div>
<div class="line">   llt.<a class="code" href="classEigen_1_1LLT.html#ac97b2f0b61df78a30544670b3711e3db">compute</a>(A);</div>
<div class="line">   cout &lt;&lt; <span class="stringliteral">&quot;The solution is now:\n&quot;</span> &lt;&lt; llt.<a class="code" href="classEigen_1_1LLT.html#a301e28d205fb0308cb4f4b04718e1685">solve</a>(b) &lt;&lt; endl;</div>
<div class="line">}</div>
</div><!-- fragment -->  </td><td><pre class="fragment">Here is the matrix A:
 2 -1
-1  3
Here is the right hand side b:
1 2
3 1
Computing LLT decomposition...
The solution is:
1.2 1.4
1.4 0.8
The matrix A is now:
 2 -1
-1  4
Computing LLT decomposition...
The solution is now:
    1  1.29
    1 0.571
</pre>   </td></tr>
</table>
<p>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: </p>
<div class="fragment"><div class="line">HouseholderQR&lt;MatrixXf&gt; qr(50,50);</div>
<div class="line"><a class="code" href="group__matrixtypedefs.html#gabab09c32e96cfa9829a88400627af162">MatrixXf</a> A = <a class="code" href="classEigen_1_1DenseBase.html#a8e759dafdd9ecc446d397b7f5435f60a">MatrixXf::Random</a>(50,50);</div>
<div class="line">qr.compute(A); <span class="comment">// no dynamic memory allocation</span></div>
</div><!-- fragment --><h1><a class="anchor" id="TutorialLinAlgRankRevealing"></a>
Rank-revealing decompositions</h1>
<p>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 <a class="el" href="group__TopicLinearAlgebraDecompositions.html">this table</a> you can see for all our decompositions whether they are rank-revealing or not.</p>
<p>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 <a class="el" href="classEigen_1_1FullPivLU.html" title="LU decomposition of a matrix with complete pivoting, and related features. ">FullPivLU</a>:</p>
<table  class="example">
<tr>
<th>Example:</th><th>Output: </th></tr>
<tr>
<td><div class="fragment"><div class="line"><span class="preprocessor">#include &lt;iostream&gt;</span></div>
<div class="line"><span class="preprocessor">#include &lt;Eigen/Dense&gt;</span></div>
<div class="line"></div>
<div class="line"><span class="keyword">using namespace </span>std;</div>
<div class="line"><span class="keyword">using namespace </span>Eigen;</div>
<div class="line"></div>
<div class="line"><span class="keywordtype">int</span> main()</div>
<div class="line">{</div>
<div class="line">   <a class="code" href="classEigen_1_1Matrix.html">Matrix3f</a> A;</div>
<div class="line">   A &lt;&lt; 1, 2, 5,</div>
<div class="line">        2, 1, 4,</div>
<div class="line">        3, 0, 3;</div>
<div class="line">   cout &lt;&lt; <span class="stringliteral">&quot;Here is the matrix A:\n&quot;</span> &lt;&lt; A &lt;&lt; endl;</div>
<div class="line">   <a class="code" href="classEigen_1_1FullPivLU.html">FullPivLU&lt;Matrix3f&gt;</a> lu_decomp(A);</div>
<div class="line">   cout &lt;&lt; <span class="stringliteral">&quot;The rank of A is &quot;</span> &lt;&lt; lu_decomp.rank() &lt;&lt; endl;</div>
<div class="line">   cout &lt;&lt; <span class="stringliteral">&quot;Here is a matrix whose columns form a basis of the null-space of A:\n&quot;</span></div>
<div class="line">        &lt;&lt; lu_decomp.kernel() &lt;&lt; endl;</div>
<div class="line">   cout &lt;&lt; <span class="stringliteral">&quot;Here is a matrix whose columns form a basis of the column-space of A:\n&quot;</span></div>
<div class="line">        &lt;&lt; lu_decomp.image(A) &lt;&lt; endl; <span class="comment">// yes, have to pass the original A</span></div>
<div class="line">}</div>
</div><!-- fragment -->  </td><td><pre class="fragment">Here is the matrix A:
1 2 5
2 1 4
3 0 3
The rank of A is 2
Here is a matrix whose columns form a basis of the null-space of A:
0.5
1
-0.5
Here is a matrix whose columns form a basis of the column-space of A:
5 1
4 2
3 3
</pre>   </td></tr>
</table>
<p>Of course, any rank computation depends on the choice of an arbitrary threshold, since practically no floating-point matrix is <em>exactly</em> rank-deficient. <a class="el" href="namespaceEigen.html" title="Namespace containing all symbols from the Eigen library. ">Eigen</a> 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() 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.</p>
<table  class="example">
<tr>
<th>Example:</th><th>Output: </th></tr>
<tr>
<td><div class="fragment"><div class="line"><span class="preprocessor">#include &lt;iostream&gt;</span></div>
<div class="line"><span class="preprocessor">#include &lt;Eigen/Dense&gt;</span></div>
<div class="line"></div>
<div class="line"><span class="keyword">using namespace </span>std;</div>
<div class="line"><span class="keyword">using namespace </span>Eigen;</div>
<div class="line"></div>
<div class="line"><span class="keywordtype">int</span> main()</div>
<div class="line">{</div>
<div class="line">   <a class="code" href="classEigen_1_1Matrix.html">Matrix2d</a> A;</div>
<div class="line">   A &lt;&lt; 2, 1,</div>
<div class="line">        2, 0.9999999999;</div>
<div class="line">   <a class="code" href="classEigen_1_1FullPivLU.html">FullPivLU&lt;Matrix2d&gt;</a> lu(A);</div>
<div class="line">   cout &lt;&lt; <span class="stringliteral">&quot;By default, the rank of A is found to be &quot;</span> &lt;&lt; lu.rank() &lt;&lt; endl;</div>
<div class="line">   lu.setThreshold(1e-5);</div>
<div class="line">   cout &lt;&lt; <span class="stringliteral">&quot;With threshold 1e-5, the rank of A is found to be &quot;</span> &lt;&lt; lu.rank() &lt;&lt; endl;</div>
<div class="line">}</div>
</div><!-- fragment -->  </td><td><pre class="fragment">By default, the rank of A is found to be 2
With threshold 1e-5, the rank of A is found to be 1
</pre>   </td></tr>
</table>
</div><!-- contents -->
</div><!-- doc-content -->
<!-- start footer part -->
<div id="nav-path" class="navpath"><!-- id is needed for treeview function! -->
  <ul>
    <li class="footer">Generated on Mon Oct 28 2013 11:04:27 for Eigen by
    <a href="http://www.doxygen.org/index.html">
    <img class="footer" src="doxygen.png" alt="doxygen"/></a> 1.8.5 </li>
  </ul>
</div>
<!-- Piwik --> 
<!--
<script type="text/javascript">
var pkBaseURL = (("https:" == document.location.protocol) ? "https://stats.sylphide-consulting.com/piwik/" : "http://stats.sylphide-consulting.com/piwik/");
document.write(unescape("%3Cscript src='" + pkBaseURL + "piwik.js' type='text/javascript'%3E%3C/script%3E"));
</script><script type="text/javascript">
try {
var piwikTracker = Piwik.getTracker(pkBaseURL + "piwik.php", 20);
piwikTracker.trackPageView();
piwikTracker.enableLinkTracking();
} catch( err ) {}
</script><noscript><p><img src="http://stats.sylphide-consulting.com/piwik/piwik.php?idsite=20" style="border:0" alt="" /></p></noscript>
-->
<!-- End Piwik Tracking Code -->
</body>
</html>