Sophie

Sophie

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

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: FullPivLU&lt; MatrixType &gt; Class Template Reference</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('classEigen_1_1FullPivLU.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="summary">
<a href="classEigen_1_1FullPivLU-members.html">List of all members</a> &#124;
<a href="#pub-methods">Public Member Functions</a>  </div>
  <div class="headertitle">
<div class="title">FullPivLU&lt; MatrixType &gt; Class Template Reference<div class="ingroups"><a class="el" href="group__LU__Module.html">LU module</a></div></div>  </div>
</div><!--header-->
<div class="contents">
<a name="details" id="details"></a><h2 class="groupheader">Detailed Description</h2>
<div class="textblock"><h3>template&lt;typename MatrixType&gt;<br/>
class Eigen::FullPivLU&lt; MatrixType &gt;</h3>

<p>LU decomposition of a matrix with complete pivoting, and related features. </p>
<dl class="params"><dt>Parameters</dt><dd>
  <table class="params">
    <tr><td class="paramname">MatrixType</td><td>the type of the matrix of which we are computing the LU decomposition</td></tr>
  </table>
  </dd>
</dl>
<p>This class represents a LU decomposition of any matrix, with complete pivoting: the matrix A is decomposed as A = PLUQ where L is unit-lower-triangular, U is upper-triangular, and P and Q are permutation matrices. This is a rank-revealing LU decomposition. The eigenvalues (diagonal coefficients) of U are sorted in such a way that any zeros are at the end.</p>
<p>This decomposition provides the generic approach to solving systems of linear equations, computing the rank, invertibility, inverse, kernel, and determinant.</p>
<p>This LU decomposition is very stable and well tested with large matrices. However there are use cases where the SVD decomposition is inherently more stable and/or flexible. For example, when computing the kernel of a matrix, working with the SVD allows to select the smallest singular values of the matrix, something that the LU decomposition doesn't see.</p>
<p>The data of the LU decomposition can be directly accessed through the methods <a class="el" href="classEigen_1_1FullPivLU.html#ad69664a62ab4d3026566d0d4a261b187">matrixLU()</a>, <a class="el" href="classEigen_1_1FullPivLU.html#a38416f9985b9c7ad9dc3bd355479dd67">permutationP()</a>, <a class="el" href="classEigen_1_1FullPivLU.html#ae85e85c3d1488b5882a6ccd63678a4d1">permutationQ()</a>.</p>
<p>As an exemple, here is how the original matrix can be retrieved: </p>
<div class="fragment"><div class="line"><span class="keyword">typedef</span> Matrix&lt;double, 5, 3&gt; Matrix5x3;</div>
<div class="line"><span class="keyword">typedef</span> Matrix&lt;double, 5, 5&gt; Matrix5x5;</div>
<div class="line">Matrix5x3 m = Matrix5x3::Random();</div>
<div class="line">cout &lt;&lt; <span class="stringliteral">&quot;Here is the matrix m:&quot;</span> &lt;&lt; endl &lt;&lt; m &lt;&lt; endl;</div>
<div class="line"><a class="code" href="classEigen_1_1FullPivLU.html">Eigen::FullPivLU&lt;Matrix5x3&gt;</a> lu(m);</div>
<div class="line">cout &lt;&lt; <span class="stringliteral">&quot;Here is, up to permutations, its LU decomposition matrix:&quot;</span></div>
<div class="line">     &lt;&lt; endl &lt;&lt; lu.matrixLU() &lt;&lt; endl;</div>
<div class="line">cout &lt;&lt; <span class="stringliteral">&quot;Here is the L part:&quot;</span> &lt;&lt; endl;</div>
<div class="line">Matrix5x5 l = Matrix5x5::Identity();</div>
<div class="line">l.block&lt;5,3&gt;(0,0).triangularView&lt;StrictlyLower&gt;() = lu.matrixLU();</div>
<div class="line">cout &lt;&lt; l &lt;&lt; endl;</div>
<div class="line">cout &lt;&lt; <span class="stringliteral">&quot;Here is the U part:&quot;</span> &lt;&lt; endl;</div>
<div class="line">Matrix5x3 u = lu.matrixLU().triangularView&lt;<a class="code" href="group__enums.html#ggab59c1bec446b10af208f977a871d910bae70afef0d3ff7aca74e17e85ff6c9f2e">Upper</a>&gt;();</div>
<div class="line">cout &lt;&lt; u &lt;&lt; endl;</div>
<div class="line">cout &lt;&lt; <span class="stringliteral">&quot;Let us now reconstruct the original matrix m:&quot;</span> &lt;&lt; endl;</div>
<div class="line">cout &lt;&lt; lu.permutationP().inverse() * l * u * lu.permutationQ().inverse() &lt;&lt; endl;</div>
</div><!-- fragment --><p> Output: </p>
<pre class="fragment">Here is the matrix m:
   0.68  -0.605 -0.0452
 -0.211   -0.33   0.258
  0.566   0.536   -0.27
  0.597  -0.444  0.0268
  0.823   0.108   0.904
Here is, up to permutations, its LU decomposition matrix:
 0.904  0.823  0.108
-0.299  0.812  0.569
 -0.05  0.888   -1.1
0.0296  0.705  0.768
 0.285 -0.549 0.0436
Here is the L part:
     1      0      0      0      0
-0.299      1      0      0      0
 -0.05  0.888      1      0      0
0.0296  0.705  0.768      1      0
 0.285 -0.549 0.0436      0      1
Here is the U part:
0.904 0.823 0.108
    0 0.812 0.569
    0     0  -1.1
    0     0     0
    0     0     0
Let us now reconstruct the original matrix m:
   0.68  -0.605 -0.0452
 -0.211   -0.33   0.258
  0.566   0.536   -0.27
  0.597  -0.444  0.0268
  0.823   0.108   0.904
</pre><dl class="section see"><dt>See Also</dt><dd><a class="el" href="classEigen_1_1MatrixBase.html#a0844e94f8f95ae01a2cd88dbbf5cbf91">MatrixBase::fullPivLu()</a>, <a class="el" href="classEigen_1_1MatrixBase.html#ad63cea11a4bf220298dce4489a1704c7">MatrixBase::determinant()</a>, <a class="el" href="classEigen_1_1MatrixBase.html#aa2834da4c855fa35fed8c4030f79f9da">MatrixBase::inverse()</a> </dd></dl>
</div><table class="memberdecls">
<tr class="heading"><td colspan="2"><h2 class="groupheader"><a name="pub-methods"></a>
Public Member Functions</h2></td></tr>
<tr class="memitem:ac1c60e27ca149bb599662fb554a0949a"><td class="memItemLeft" align="right" valign="top"><a class="el" href="classEigen_1_1FullPivLU.html">FullPivLU</a> &amp;&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1FullPivLU.html#ac1c60e27ca149bb599662fb554a0949a">compute</a> (const MatrixType &amp;matrix)</td></tr>
<tr class="separator:ac1c60e27ca149bb599662fb554a0949a"><td class="memSeparator" colspan="2">&#160;</td></tr>
<tr class="memitem:aa47f041dae554fe1f135e2794ae914a7"><td class="memItemLeft" align="right" valign="top">internal::traits&lt; MatrixType &gt;<br class="typebreak"/>
::Scalar&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1FullPivLU.html#aa47f041dae554fe1f135e2794ae914a7">determinant</a> () const </td></tr>
<tr class="separator:aa47f041dae554fe1f135e2794ae914a7"><td class="memSeparator" colspan="2">&#160;</td></tr>
<tr class="memitem:a7c6323871c4f080fc6e2d3ad7fc607fc"><td class="memItemLeft" align="right" valign="top">Index&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1FullPivLU.html#a7c6323871c4f080fc6e2d3ad7fc607fc">dimensionOfKernel</a> () const </td></tr>
<tr class="separator:a7c6323871c4f080fc6e2d3ad7fc607fc"><td class="memSeparator" colspan="2">&#160;</td></tr>
<tr class="memitem:a10dec2fa1767ac0712af1efa732b2046"><td class="memItemLeft" align="right" valign="top">&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1FullPivLU.html#a10dec2fa1767ac0712af1efa732b2046">FullPivLU</a> ()</td></tr>
<tr class="memdesc:a10dec2fa1767ac0712af1efa732b2046"><td class="mdescLeft">&#160;</td><td class="mdescRight">Default Constructor.  <a href="#a10dec2fa1767ac0712af1efa732b2046">More...</a><br/></td></tr>
<tr class="separator:a10dec2fa1767ac0712af1efa732b2046"><td class="memSeparator" colspan="2">&#160;</td></tr>
<tr class="memitem:ab6c78ecd953ba165ef421fa67c05753a"><td class="memItemLeft" align="right" valign="top">&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1FullPivLU.html#ab6c78ecd953ba165ef421fa67c05753a">FullPivLU</a> (Index rows, Index cols)</td></tr>
<tr class="memdesc:ab6c78ecd953ba165ef421fa67c05753a"><td class="mdescLeft">&#160;</td><td class="mdescRight">Default Constructor with memory preallocation.  <a href="#ab6c78ecd953ba165ef421fa67c05753a">More...</a><br/></td></tr>
<tr class="separator:ab6c78ecd953ba165ef421fa67c05753a"><td class="memSeparator" colspan="2">&#160;</td></tr>
<tr class="memitem:ab66f95f8f6b6455a0a38d95486330808"><td class="memItemLeft" align="right" valign="top">&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1FullPivLU.html#ab66f95f8f6b6455a0a38d95486330808">FullPivLU</a> (const MatrixType &amp;matrix)</td></tr>
<tr class="separator:ab66f95f8f6b6455a0a38d95486330808"><td class="memSeparator" colspan="2">&#160;</td></tr>
<tr class="memitem:ad854843c6dc601252dc107e7b29133e9"><td class="memItemLeft" align="right" valign="top">const internal::image_retval<br class="typebreak"/>
&lt; <a class="el" href="classEigen_1_1FullPivLU.html">FullPivLU</a> &gt;&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1FullPivLU.html#ad854843c6dc601252dc107e7b29133e9">image</a> (const MatrixType &amp;originalMatrix) const </td></tr>
<tr class="separator:ad854843c6dc601252dc107e7b29133e9"><td class="memSeparator" colspan="2">&#160;</td></tr>
<tr class="memitem:a5ca20e3802e96fb14c2be37039afcae9"><td class="memItemLeft" align="right" valign="top">const internal::solve_retval<br class="typebreak"/>
&lt; <a class="el" href="classEigen_1_1FullPivLU.html">FullPivLU</a>, typename <br class="typebreak"/>
MatrixType::IdentityReturnType &gt;&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1FullPivLU.html#a5ca20e3802e96fb14c2be37039afcae9">inverse</a> () const </td></tr>
<tr class="separator:a5ca20e3802e96fb14c2be37039afcae9"><td class="memSeparator" colspan="2">&#160;</td></tr>
<tr class="memitem:a1e119085e53eca65e9ba15451c102d40"><td class="memItemLeft" align="right" valign="top">bool&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1FullPivLU.html#a1e119085e53eca65e9ba15451c102d40">isInjective</a> () const </td></tr>
<tr class="separator:a1e119085e53eca65e9ba15451c102d40"><td class="memSeparator" colspan="2">&#160;</td></tr>
<tr class="memitem:ab60c7d993c9eba31668fb8886d621094"><td class="memItemLeft" align="right" valign="top">bool&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1FullPivLU.html#ab60c7d993c9eba31668fb8886d621094">isInvertible</a> () const </td></tr>
<tr class="separator:ab60c7d993c9eba31668fb8886d621094"><td class="memSeparator" colspan="2">&#160;</td></tr>
<tr class="memitem:a5968b9ca46303b3cc7250e7b120ab7e6"><td class="memItemLeft" align="right" valign="top">bool&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1FullPivLU.html#a5968b9ca46303b3cc7250e7b120ab7e6">isSurjective</a> () const </td></tr>
<tr class="separator:a5968b9ca46303b3cc7250e7b120ab7e6"><td class="memSeparator" colspan="2">&#160;</td></tr>
<tr class="memitem:a6e8f1d2fcbd86d3dc5a8a013b6e7200a"><td class="memItemLeft" align="right" valign="top">const internal::kernel_retval<br class="typebreak"/>
&lt; <a class="el" href="classEigen_1_1FullPivLU.html">FullPivLU</a> &gt;&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1FullPivLU.html#a6e8f1d2fcbd86d3dc5a8a013b6e7200a">kernel</a> () const </td></tr>
<tr class="separator:a6e8f1d2fcbd86d3dc5a8a013b6e7200a"><td class="memSeparator" colspan="2">&#160;</td></tr>
<tr class="memitem:ad69664a62ab4d3026566d0d4a261b187"><td class="memItemLeft" align="right" valign="top">const MatrixType &amp;&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1FullPivLU.html#ad69664a62ab4d3026566d0d4a261b187">matrixLU</a> () const </td></tr>
<tr class="separator:ad69664a62ab4d3026566d0d4a261b187"><td class="memSeparator" colspan="2">&#160;</td></tr>
<tr class="memitem:a067e9d4143ce0558fc684b736128a4ed"><td class="memItemLeft" align="right" valign="top">RealScalar&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1FullPivLU.html#a067e9d4143ce0558fc684b736128a4ed">maxPivot</a> () const </td></tr>
<tr class="separator:a067e9d4143ce0558fc684b736128a4ed"><td class="memSeparator" colspan="2">&#160;</td></tr>
<tr class="memitem:a6dda9285f13dec9f49e9c17229a89988"><td class="memItemLeft" align="right" valign="top">Index&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1FullPivLU.html#a6dda9285f13dec9f49e9c17229a89988">nonzeroPivots</a> () const </td></tr>
<tr class="separator:a6dda9285f13dec9f49e9c17229a89988"><td class="memSeparator" colspan="2">&#160;</td></tr>
<tr class="memitem:a38416f9985b9c7ad9dc3bd355479dd67"><td class="memItemLeft" align="right" valign="top">const <a class="el" href="classEigen_1_1PermutationMatrix.html">PermutationPType</a> &amp;&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1FullPivLU.html#a38416f9985b9c7ad9dc3bd355479dd67">permutationP</a> () const </td></tr>
<tr class="separator:a38416f9985b9c7ad9dc3bd355479dd67"><td class="memSeparator" colspan="2">&#160;</td></tr>
<tr class="memitem:ae85e85c3d1488b5882a6ccd63678a4d1"><td class="memItemLeft" align="right" valign="top">const <a class="el" href="classEigen_1_1PermutationMatrix.html">PermutationQType</a> &amp;&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1FullPivLU.html#ae85e85c3d1488b5882a6ccd63678a4d1">permutationQ</a> () const </td></tr>
<tr class="separator:ae85e85c3d1488b5882a6ccd63678a4d1"><td class="memSeparator" colspan="2">&#160;</td></tr>
<tr class="memitem:a363d1c09d77f09d6ea2d2789776e7be3"><td class="memItemLeft" align="right" valign="top">Index&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1FullPivLU.html#a363d1c09d77f09d6ea2d2789776e7be3">rank</a> () const </td></tr>
<tr class="separator:a363d1c09d77f09d6ea2d2789776e7be3"><td class="memSeparator" colspan="2">&#160;</td></tr>
<tr class="memitem:ac36925ac693435a090efee1cb5d6d16a"><td class="memItemLeft" align="right" valign="top">MatrixType&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1FullPivLU.html#ac36925ac693435a090efee1cb5d6d16a">reconstructedMatrix</a> () const </td></tr>
<tr class="separator:ac36925ac693435a090efee1cb5d6d16a"><td class="memSeparator" colspan="2">&#160;</td></tr>
<tr class="memitem:ae47399e3cf7943075ce18ef89fe17f21"><td class="memItemLeft" align="right" valign="top"><a class="el" href="classEigen_1_1FullPivLU.html">FullPivLU</a> &amp;&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1FullPivLU.html#ae47399e3cf7943075ce18ef89fe17f21">setThreshold</a> (const RealScalar &amp;<a class="el" href="classEigen_1_1FullPivLU.html#aa5a87faaa92a3081045d1f934e292ef0">threshold</a>)</td></tr>
<tr class="separator:ae47399e3cf7943075ce18ef89fe17f21"><td class="memSeparator" colspan="2">&#160;</td></tr>
<tr class="memitem:a02b74b9b823e060825443fd82cd25fd4"><td class="memItemLeft" align="right" valign="top"><a class="el" href="classEigen_1_1FullPivLU.html">FullPivLU</a> &amp;&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1FullPivLU.html#a02b74b9b823e060825443fd82cd25fd4">setThreshold</a> (Default_t)</td></tr>
<tr class="separator:a02b74b9b823e060825443fd82cd25fd4"><td class="memSeparator" colspan="2">&#160;</td></tr>
<tr class="memitem:ac1a642d728c059c7625863f126e2e718"><td class="memTemplParams" colspan="2">template&lt;typename Rhs &gt; </td></tr>
<tr class="memitem:ac1a642d728c059c7625863f126e2e718"><td class="memTemplItemLeft" align="right" valign="top">const internal::solve_retval<br class="typebreak"/>
&lt; <a class="el" href="classEigen_1_1FullPivLU.html">FullPivLU</a>, Rhs &gt;&#160;</td><td class="memTemplItemRight" valign="bottom"><a class="el" href="classEigen_1_1FullPivLU.html#ac1a642d728c059c7625863f126e2e718">solve</a> (const <a class="el" href="classEigen_1_1MatrixBase.html">MatrixBase</a>&lt; Rhs &gt; &amp;b) const </td></tr>
<tr class="separator:ac1a642d728c059c7625863f126e2e718"><td class="memSeparator" colspan="2">&#160;</td></tr>
<tr class="memitem:aa5a87faaa92a3081045d1f934e292ef0"><td class="memItemLeft" align="right" valign="top">RealScalar&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1FullPivLU.html#aa5a87faaa92a3081045d1f934e292ef0">threshold</a> () const </td></tr>
<tr class="separator:aa5a87faaa92a3081045d1f934e292ef0"><td class="memSeparator" colspan="2">&#160;</td></tr>
</table>
<h2 class="groupheader">Constructor &amp; Destructor Documentation</h2>
<a class="anchor" id="a10dec2fa1767ac0712af1efa732b2046"></a>
<div class="memitem">
<div class="memproto">
      <table class="memname">
        <tr>
          <td class="memname"><a class="el" href="classEigen_1_1FullPivLU.html">FullPivLU</a> </td>
          <td>(</td>
          <td class="paramname"></td><td>)</td>
          <td></td>
        </tr>
      </table>
</div><div class="memdoc">

<p>Default Constructor. </p>
<p>The default constructor is useful in cases in which the user intends to perform decompositions via LU::compute(const MatrixType&amp;). </p>

</div>
</div>
<a class="anchor" id="ab6c78ecd953ba165ef421fa67c05753a"></a>
<div class="memitem">
<div class="memproto">
      <table class="memname">
        <tr>
          <td class="memname"><a class="el" href="classEigen_1_1FullPivLU.html">FullPivLU</a> </td>
          <td>(</td>
          <td class="paramtype">Index&#160;</td>
          <td class="paramname"><em>rows</em>, </td>
        </tr>
        <tr>
          <td class="paramkey"></td>
          <td></td>
          <td class="paramtype">Index&#160;</td>
          <td class="paramname"><em>cols</em>&#160;</td>
        </tr>
        <tr>
          <td></td>
          <td>)</td>
          <td></td><td></td>
        </tr>
      </table>
</div><div class="memdoc">

<p>Default Constructor with memory preallocation. </p>
<p>Like the default constructor but with preallocation of the internal data according to the specified problem <em>size</em>. </p>
<dl class="section see"><dt>See Also</dt><dd><a class="el" href="classEigen_1_1FullPivLU.html#a10dec2fa1767ac0712af1efa732b2046" title="Default Constructor. ">FullPivLU()</a> </dd></dl>

</div>
</div>
<a class="anchor" id="ab66f95f8f6b6455a0a38d95486330808"></a>
<div class="memitem">
<div class="memproto">
      <table class="memname">
        <tr>
          <td class="memname"><a class="el" href="classEigen_1_1FullPivLU.html">FullPivLU</a> </td>
          <td>(</td>
          <td class="paramtype">const MatrixType &amp;&#160;</td>
          <td class="paramname"><em>matrix</em></td><td>)</td>
          <td></td>
        </tr>
      </table>
</div><div class="memdoc">
<p>Constructor.</p>
<dl class="params"><dt>Parameters</dt><dd>
  <table class="params">
    <tr><td class="paramname">matrix</td><td>the matrix of which to compute the LU decomposition. It is required to be nonzero. </td></tr>
  </table>
  </dd>
</dl>

<p>References <a class="el" href="classEigen_1_1FullPivLU.html#ac1c60e27ca149bb599662fb554a0949a">FullPivLU&lt; MatrixType &gt;::compute()</a>.</p>

</div>
</div>
<h2 class="groupheader">Member Function Documentation</h2>
<a class="anchor" id="ac1c60e27ca149bb599662fb554a0949a"></a>
<div class="memitem">
<div class="memproto">
      <table class="memname">
        <tr>
          <td class="memname"><a class="el" href="classEigen_1_1FullPivLU.html">FullPivLU</a>&lt; MatrixType &gt; &amp; compute </td>
          <td>(</td>
          <td class="paramtype">const MatrixType &amp;&#160;</td>
          <td class="paramname"><em>matrix</em></td><td>)</td>
          <td></td>
        </tr>
      </table>
</div><div class="memdoc">
<p>Computes the LU decomposition of the given matrix.</p>
<dl class="params"><dt>Parameters</dt><dd>
  <table class="params">
    <tr><td class="paramname">matrix</td><td>the matrix of which to compute the LU decomposition. It is required to be nonzero.</td></tr>
  </table>
  </dd>
</dl>
<dl class="section return"><dt>Returns</dt><dd>a reference to *this </dd></dl>

<p>Referenced by <a class="el" href="classEigen_1_1FullPivLU.html#ab66f95f8f6b6455a0a38d95486330808">FullPivLU&lt; MatrixType &gt;::FullPivLU()</a>.</p>

</div>
</div>
<a class="anchor" id="aa47f041dae554fe1f135e2794ae914a7"></a>
<div class="memitem">
<div class="memproto">
      <table class="memname">
        <tr>
          <td class="memname">internal::traits&lt; MatrixType &gt;::Scalar determinant </td>
          <td>(</td>
          <td class="paramname"></td><td>)</td>
          <td> const</td>
        </tr>
      </table>
</div><div class="memdoc">
<dl class="section return"><dt>Returns</dt><dd>the determinant of the matrix of which *this is the LU decomposition. It has only linear complexity (that is, O(n) where n is the dimension of the square matrix) as the LU decomposition has already been computed.</dd></dl>
<dl class="section note"><dt>Note</dt><dd>This is only for square matrices.</dd>
<dd>
For fixed-size matrices of size up to 4, <a class="el" href="classEigen_1_1MatrixBase.html#ad63cea11a4bf220298dce4489a1704c7">MatrixBase::determinant()</a> offers optimized paths.</dd></dl>
<dl class="section warning"><dt>Warning</dt><dd>a determinant can be very big or small, so for matrices of large enough dimension, there is a risk of overflow/underflow.</dd></dl>
<dl class="section see"><dt>See Also</dt><dd><a class="el" href="classEigen_1_1MatrixBase.html#ad63cea11a4bf220298dce4489a1704c7">MatrixBase::determinant()</a> </dd></dl>

</div>
</div>
<a class="anchor" id="a7c6323871c4f080fc6e2d3ad7fc607fc"></a>
<div class="memitem">
<div class="memproto">
<table class="mlabels">
  <tr>
  <td class="mlabels-left">
      <table class="memname">
        <tr>
          <td class="memname">Index dimensionOfKernel </td>
          <td>(</td>
          <td class="paramname"></td><td>)</td>
          <td> const</td>
        </tr>
      </table>
  </td>
  <td class="mlabels-right">
<span class="mlabels"><span class="mlabel">inline</span></span>  </td>
  </tr>
</table>
</div><div class="memdoc">
<dl class="section return"><dt>Returns</dt><dd>the dimension of the kernel of the matrix of which *this is the LU decomposition.</dd></dl>
<dl class="section note"><dt>Note</dt><dd>This method has to determine which pivots should be considered nonzero. For that, it uses the threshold value that you can control by calling <a class="el" href="classEigen_1_1FullPivLU.html#ae47399e3cf7943075ce18ef89fe17f21">setThreshold(const RealScalar&amp;)</a>. </dd></dl>

<p>References <a class="el" href="classEigen_1_1FullPivLU.html#a363d1c09d77f09d6ea2d2789776e7be3">FullPivLU&lt; MatrixType &gt;::rank()</a>.</p>

</div>
</div>
<a class="anchor" id="ad854843c6dc601252dc107e7b29133e9"></a>
<div class="memitem">
<div class="memproto">
<table class="mlabels">
  <tr>
  <td class="mlabels-left">
      <table class="memname">
        <tr>
          <td class="memname">const internal::image_retval&lt;<a class="el" href="classEigen_1_1FullPivLU.html">FullPivLU</a>&gt; image </td>
          <td>(</td>
          <td class="paramtype">const MatrixType &amp;&#160;</td>
          <td class="paramname"><em>originalMatrix</em></td><td>)</td>
          <td> const</td>
        </tr>
      </table>
  </td>
  <td class="mlabels-right">
<span class="mlabels"><span class="mlabel">inline</span></span>  </td>
  </tr>
</table>
</div><div class="memdoc">
<dl class="section return"><dt>Returns</dt><dd>the image of the matrix, also called its column-space. The columns of the returned matrix will form a basis of the kernel.</dd></dl>
<dl class="params"><dt>Parameters</dt><dd>
  <table class="params">
    <tr><td class="paramname">originalMatrix</td><td>the original matrix, of which *this is the LU decomposition. The reason why it is needed to pass it here, is that this allows a large optimization, as otherwise this method would need to reconstruct it from the LU decomposition.</td></tr>
  </table>
  </dd>
</dl>
<dl class="section note"><dt>Note</dt><dd>If the image has dimension zero, then the returned matrix is a column-vector filled with zeros.</dd>
<dd>
This method has to determine which pivots should be considered nonzero. For that, it uses the threshold value that you can control by calling <a class="el" href="classEigen_1_1FullPivLU.html#ae47399e3cf7943075ce18ef89fe17f21">setThreshold(const RealScalar&amp;)</a>.</dd></dl>
<p>Example: </p>
<div class="fragment"><div class="line"><a class="code" href="group__matrixtypedefs.html#ga45a14b423c687c3e2e8325f148e27bf3">Matrix3d</a> m;</div>
<div class="line">m &lt;&lt; 1,1,0,</div>
<div class="line">     1,3,2,</div>
<div class="line">     0,1,1;</div>
<div class="line">cout &lt;&lt; <span class="stringliteral">&quot;Here is the matrix m:&quot;</span> &lt;&lt; endl &lt;&lt; m &lt;&lt; endl;</div>
<div class="line">cout &lt;&lt; <span class="stringliteral">&quot;Notice that the middle column is the sum of the two others, so the &quot;</span></div>
<div class="line">     &lt;&lt; <span class="stringliteral">&quot;columns are linearly dependent.&quot;</span> &lt;&lt; endl;</div>
<div class="line">cout &lt;&lt; <span class="stringliteral">&quot;Here is a matrix whose columns have the same span but are linearly independent:&quot;</span></div>
<div class="line">     &lt;&lt; endl &lt;&lt; m.fullPivLu().image(m) &lt;&lt; endl;</div>
</div><!-- fragment --><p> Output: </p>
<pre class="fragment">Here is the matrix m:
1 1 0
1 3 2
0 1 1
Notice that the middle column is the sum of the two others, so the columns are linearly dependent.
Here is a matrix whose columns have the same span but are linearly independent:
1 1
3 1
1 0
</pre><dl class="section see"><dt>See Also</dt><dd><a class="el" href="classEigen_1_1FullPivLU.html#a6e8f1d2fcbd86d3dc5a8a013b6e7200a">kernel()</a> </dd></dl>

</div>
</div>
<a class="anchor" id="a5ca20e3802e96fb14c2be37039afcae9"></a>
<div class="memitem">
<div class="memproto">
<table class="mlabels">
  <tr>
  <td class="mlabels-left">
      <table class="memname">
        <tr>
          <td class="memname">const internal::solve_retval&lt;<a class="el" href="classEigen_1_1FullPivLU.html">FullPivLU</a>,typename MatrixType::IdentityReturnType&gt; inverse </td>
          <td>(</td>
          <td class="paramname"></td><td>)</td>
          <td> const</td>
        </tr>
      </table>
  </td>
  <td class="mlabels-right">
<span class="mlabels"><span class="mlabel">inline</span></span>  </td>
  </tr>
</table>
</div><div class="memdoc">
<dl class="section return"><dt>Returns</dt><dd>the inverse of the matrix of which *this is the LU decomposition.</dd></dl>
<dl class="section note"><dt>Note</dt><dd>If this matrix is not invertible, the returned matrix has undefined coefficients. Use <a class="el" href="classEigen_1_1FullPivLU.html#ab60c7d993c9eba31668fb8886d621094">isInvertible()</a> to first determine whether this matrix is invertible.</dd></dl>
<dl class="section see"><dt>See Also</dt><dd><a class="el" href="classEigen_1_1MatrixBase.html#aa2834da4c855fa35fed8c4030f79f9da">MatrixBase::inverse()</a> </dd></dl>

</div>
</div>
<a class="anchor" id="a1e119085e53eca65e9ba15451c102d40"></a>
<div class="memitem">
<div class="memproto">
<table class="mlabels">
  <tr>
  <td class="mlabels-left">
      <table class="memname">
        <tr>
          <td class="memname">bool isInjective </td>
          <td>(</td>
          <td class="paramname"></td><td>)</td>
          <td> const</td>
        </tr>
      </table>
  </td>
  <td class="mlabels-right">
<span class="mlabels"><span class="mlabel">inline</span></span>  </td>
  </tr>
</table>
</div><div class="memdoc">
<dl class="section return"><dt>Returns</dt><dd>true if the matrix of which *this is the LU decomposition represents an injective linear map, i.e. has trivial kernel; false otherwise.</dd></dl>
<dl class="section note"><dt>Note</dt><dd>This method has to determine which pivots should be considered nonzero. For that, it uses the threshold value that you can control by calling <a class="el" href="classEigen_1_1FullPivLU.html#ae47399e3cf7943075ce18ef89fe17f21">setThreshold(const RealScalar&amp;)</a>. </dd></dl>

<p>References <a class="el" href="classEigen_1_1FullPivLU.html#a363d1c09d77f09d6ea2d2789776e7be3">FullPivLU&lt; MatrixType &gt;::rank()</a>.</p>

<p>Referenced by <a class="el" href="classEigen_1_1FullPivLU.html#ab60c7d993c9eba31668fb8886d621094">FullPivLU&lt; MatrixType &gt;::isInvertible()</a>.</p>

</div>
</div>
<a class="anchor" id="ab60c7d993c9eba31668fb8886d621094"></a>
<div class="memitem">
<div class="memproto">
<table class="mlabels">
  <tr>
  <td class="mlabels-left">
      <table class="memname">
        <tr>
          <td class="memname">bool isInvertible </td>
          <td>(</td>
          <td class="paramname"></td><td>)</td>
          <td> const</td>
        </tr>
      </table>
  </td>
  <td class="mlabels-right">
<span class="mlabels"><span class="mlabel">inline</span></span>  </td>
  </tr>
</table>
</div><div class="memdoc">
<dl class="section return"><dt>Returns</dt><dd>true if the matrix of which *this is the LU decomposition is invertible.</dd></dl>
<dl class="section note"><dt>Note</dt><dd>This method has to determine which pivots should be considered nonzero. For that, it uses the threshold value that you can control by calling <a class="el" href="classEigen_1_1FullPivLU.html#ae47399e3cf7943075ce18ef89fe17f21">setThreshold(const RealScalar&amp;)</a>. </dd></dl>

<p>References <a class="el" href="classEigen_1_1FullPivLU.html#a1e119085e53eca65e9ba15451c102d40">FullPivLU&lt; MatrixType &gt;::isInjective()</a>.</p>

</div>
</div>
<a class="anchor" id="a5968b9ca46303b3cc7250e7b120ab7e6"></a>
<div class="memitem">
<div class="memproto">
<table class="mlabels">
  <tr>
  <td class="mlabels-left">
      <table class="memname">
        <tr>
          <td class="memname">bool isSurjective </td>
          <td>(</td>
          <td class="paramname"></td><td>)</td>
          <td> const</td>
        </tr>
      </table>
  </td>
  <td class="mlabels-right">
<span class="mlabels"><span class="mlabel">inline</span></span>  </td>
  </tr>
</table>
</div><div class="memdoc">
<dl class="section return"><dt>Returns</dt><dd>true if the matrix of which *this is the LU decomposition represents a surjective linear map; false otherwise.</dd></dl>
<dl class="section note"><dt>Note</dt><dd>This method has to determine which pivots should be considered nonzero. For that, it uses the threshold value that you can control by calling <a class="el" href="classEigen_1_1FullPivLU.html#ae47399e3cf7943075ce18ef89fe17f21">setThreshold(const RealScalar&amp;)</a>. </dd></dl>

<p>References <a class="el" href="classEigen_1_1FullPivLU.html#a363d1c09d77f09d6ea2d2789776e7be3">FullPivLU&lt; MatrixType &gt;::rank()</a>.</p>

</div>
</div>
<a class="anchor" id="a6e8f1d2fcbd86d3dc5a8a013b6e7200a"></a>
<div class="memitem">
<div class="memproto">
<table class="mlabels">
  <tr>
  <td class="mlabels-left">
      <table class="memname">
        <tr>
          <td class="memname">const internal::kernel_retval&lt;<a class="el" href="classEigen_1_1FullPivLU.html">FullPivLU</a>&gt; kernel </td>
          <td>(</td>
          <td class="paramname"></td><td>)</td>
          <td> const</td>
        </tr>
      </table>
  </td>
  <td class="mlabels-right">
<span class="mlabels"><span class="mlabel">inline</span></span>  </td>
  </tr>
</table>
</div><div class="memdoc">
<dl class="section return"><dt>Returns</dt><dd>the kernel of the matrix, also called its null-space. The columns of the returned matrix will form a basis of the kernel.</dd></dl>
<dl class="section note"><dt>Note</dt><dd>If the kernel has dimension zero, then the returned matrix is a column-vector filled with zeros.</dd>
<dd>
This method has to determine which pivots should be considered nonzero. For that, it uses the threshold value that you can control by calling <a class="el" href="classEigen_1_1FullPivLU.html#ae47399e3cf7943075ce18ef89fe17f21">setThreshold(const RealScalar&amp;)</a>.</dd></dl>
<p>Example: </p>
<div class="fragment"><div class="line"><a class="code" href="group__matrixtypedefs.html#gabab09c32e96cfa9829a88400627af162">MatrixXf</a> m = <a class="code" href="classEigen_1_1DenseBase.html#a8e759dafdd9ecc446d397b7f5435f60a">MatrixXf::Random</a>(3,5);</div>
<div class="line">cout &lt;&lt; <span class="stringliteral">&quot;Here is the matrix m:&quot;</span> &lt;&lt; endl &lt;&lt; m &lt;&lt; endl;</div>
<div class="line"><a class="code" href="group__matrixtypedefs.html#gabab09c32e96cfa9829a88400627af162">MatrixXf</a> ker = m.fullPivLu().kernel();</div>
<div class="line">cout &lt;&lt; <span class="stringliteral">&quot;Here is a matrix whose columns form a basis of the kernel of m:&quot;</span></div>
<div class="line">     &lt;&lt; endl &lt;&lt; ker &lt;&lt; endl;</div>
<div class="line">cout &lt;&lt; <span class="stringliteral">&quot;By definition of the kernel, m*ker is zero:&quot;</span></div>
<div class="line">     &lt;&lt; endl &lt;&lt; m*ker &lt;&lt; endl;</div>
</div><!-- fragment --><p> Output: </p>
<pre class="fragment">Here is the matrix m:
   0.68   0.597   -0.33   0.108   -0.27
 -0.211   0.823   0.536 -0.0452  0.0268
  0.566  -0.605  -0.444   0.258   0.904
Here is a matrix whose columns form a basis of the kernel of m:
-0.219  0.763
0.00335 -0.447
     0      1
     1      0
-0.145 -0.285
By definition of the kernel, m*ker is zero:
-1.12e-08  1.49e-08
 -1.4e-09 -4.05e-08
 1.49e-08 -2.98e-08
</pre><dl class="section see"><dt>See Also</dt><dd><a class="el" href="classEigen_1_1FullPivLU.html#ad854843c6dc601252dc107e7b29133e9">image()</a> </dd></dl>

</div>
</div>
<a class="anchor" id="ad69664a62ab4d3026566d0d4a261b187"></a>
<div class="memitem">
<div class="memproto">
<table class="mlabels">
  <tr>
  <td class="mlabels-left">
      <table class="memname">
        <tr>
          <td class="memname">const MatrixType&amp; matrixLU </td>
          <td>(</td>
          <td class="paramname"></td><td>)</td>
          <td> const</td>
        </tr>
      </table>
  </td>
  <td class="mlabels-right">
<span class="mlabels"><span class="mlabel">inline</span></span>  </td>
  </tr>
</table>
</div><div class="memdoc">
<dl class="section return"><dt>Returns</dt><dd>the LU decomposition matrix: the upper-triangular part is U, the unit-lower-triangular part is L (at least for square matrices; in the non-square case, special care is needed, see the documentation of class <a class="el" href="classEigen_1_1FullPivLU.html" title="LU decomposition of a matrix with complete pivoting, and related features. ">FullPivLU</a>).</dd></dl>
<dl class="section see"><dt>See Also</dt><dd>matrixL(), matrixU() </dd></dl>

</div>
</div>
<a class="anchor" id="a067e9d4143ce0558fc684b736128a4ed"></a>
<div class="memitem">
<div class="memproto">
<table class="mlabels">
  <tr>
  <td class="mlabels-left">
      <table class="memname">
        <tr>
          <td class="memname">RealScalar maxPivot </td>
          <td>(</td>
          <td class="paramname"></td><td>)</td>
          <td> const</td>
        </tr>
      </table>
  </td>
  <td class="mlabels-right">
<span class="mlabels"><span class="mlabel">inline</span></span>  </td>
  </tr>
</table>
</div><div class="memdoc">
<dl class="section return"><dt>Returns</dt><dd>the absolute value of the biggest pivot, i.e. the biggest diagonal coefficient of U. </dd></dl>

</div>
</div>
<a class="anchor" id="a6dda9285f13dec9f49e9c17229a89988"></a>
<div class="memitem">
<div class="memproto">
<table class="mlabels">
  <tr>
  <td class="mlabels-left">
      <table class="memname">
        <tr>
          <td class="memname">Index nonzeroPivots </td>
          <td>(</td>
          <td class="paramname"></td><td>)</td>
          <td> const</td>
        </tr>
      </table>
  </td>
  <td class="mlabels-right">
<span class="mlabels"><span class="mlabel">inline</span></span>  </td>
  </tr>
</table>
</div><div class="memdoc">
<dl class="section return"><dt>Returns</dt><dd>the number of nonzero pivots in the LU decomposition. Here nonzero is meant in the exact sense, not in a fuzzy sense. So that notion isn't really intrinsically interesting, but it is still useful when implementing algorithms.</dd></dl>
<dl class="section see"><dt>See Also</dt><dd><a class="el" href="classEigen_1_1FullPivLU.html#a363d1c09d77f09d6ea2d2789776e7be3">rank()</a> </dd></dl>

</div>
</div>
<a class="anchor" id="a38416f9985b9c7ad9dc3bd355479dd67"></a>
<div class="memitem">
<div class="memproto">
<table class="mlabels">
  <tr>
  <td class="mlabels-left">
      <table class="memname">
        <tr>
          <td class="memname">const <a class="el" href="classEigen_1_1PermutationMatrix.html">PermutationPType</a>&amp; permutationP </td>
          <td>(</td>
          <td class="paramname"></td><td>)</td>
          <td> const</td>
        </tr>
      </table>
  </td>
  <td class="mlabels-right">
<span class="mlabels"><span class="mlabel">inline</span></span>  </td>
  </tr>
</table>
</div><div class="memdoc">
<dl class="section return"><dt>Returns</dt><dd>the permutation matrix P</dd></dl>
<dl class="section see"><dt>See Also</dt><dd><a class="el" href="classEigen_1_1FullPivLU.html#ae85e85c3d1488b5882a6ccd63678a4d1">permutationQ()</a> </dd></dl>

</div>
</div>
<a class="anchor" id="ae85e85c3d1488b5882a6ccd63678a4d1"></a>
<div class="memitem">
<div class="memproto">
<table class="mlabels">
  <tr>
  <td class="mlabels-left">
      <table class="memname">
        <tr>
          <td class="memname">const <a class="el" href="classEigen_1_1PermutationMatrix.html">PermutationQType</a>&amp; permutationQ </td>
          <td>(</td>
          <td class="paramname"></td><td>)</td>
          <td> const</td>
        </tr>
      </table>
  </td>
  <td class="mlabels-right">
<span class="mlabels"><span class="mlabel">inline</span></span>  </td>
  </tr>
</table>
</div><div class="memdoc">
<dl class="section return"><dt>Returns</dt><dd>the permutation matrix Q</dd></dl>
<dl class="section see"><dt>See Also</dt><dd><a class="el" href="classEigen_1_1FullPivLU.html#a38416f9985b9c7ad9dc3bd355479dd67">permutationP()</a> </dd></dl>

</div>
</div>
<a class="anchor" id="a363d1c09d77f09d6ea2d2789776e7be3"></a>
<div class="memitem">
<div class="memproto">
<table class="mlabels">
  <tr>
  <td class="mlabels-left">
      <table class="memname">
        <tr>
          <td class="memname">Index rank </td>
          <td>(</td>
          <td class="paramname"></td><td>)</td>
          <td> const</td>
        </tr>
      </table>
  </td>
  <td class="mlabels-right">
<span class="mlabels"><span class="mlabel">inline</span></span>  </td>
  </tr>
</table>
</div><div class="memdoc">
<dl class="section return"><dt>Returns</dt><dd>the rank of the matrix of which *this is the LU decomposition.</dd></dl>
<dl class="section note"><dt>Note</dt><dd>This method has to determine which pivots should be considered nonzero. For that, it uses the threshold value that you can control by calling <a class="el" href="classEigen_1_1FullPivLU.html#ae47399e3cf7943075ce18ef89fe17f21">setThreshold(const RealScalar&amp;)</a>. </dd></dl>

<p>References <a class="el" href="classEigen_1_1FullPivLU.html#aa5a87faaa92a3081045d1f934e292ef0">FullPivLU&lt; MatrixType &gt;::threshold()</a>.</p>

<p>Referenced by <a class="el" href="classEigen_1_1FullPivLU.html#a7c6323871c4f080fc6e2d3ad7fc607fc">FullPivLU&lt; MatrixType &gt;::dimensionOfKernel()</a>, <a class="el" href="classEigen_1_1FullPivLU.html#a1e119085e53eca65e9ba15451c102d40">FullPivLU&lt; MatrixType &gt;::isInjective()</a>, and <a class="el" href="classEigen_1_1FullPivLU.html#a5968b9ca46303b3cc7250e7b120ab7e6">FullPivLU&lt; MatrixType &gt;::isSurjective()</a>.</p>

</div>
</div>
<a class="anchor" id="ac36925ac693435a090efee1cb5d6d16a"></a>
<div class="memitem">
<div class="memproto">
      <table class="memname">
        <tr>
          <td class="memname">MatrixType reconstructedMatrix </td>
          <td>(</td>
          <td class="paramname"></td><td>)</td>
          <td> const</td>
        </tr>
      </table>
</div><div class="memdoc">
<dl class="section return"><dt>Returns</dt><dd>the matrix represented by the decomposition, i.e., it returns the product: P^{-1} L U Q^{-1}. This function is provided for debug purpose. </dd></dl>

</div>
</div>
<a class="anchor" id="ae47399e3cf7943075ce18ef89fe17f21"></a>
<div class="memitem">
<div class="memproto">
<table class="mlabels">
  <tr>
  <td class="mlabels-left">
      <table class="memname">
        <tr>
          <td class="memname"><a class="el" href="classEigen_1_1FullPivLU.html">FullPivLU</a>&amp; setThreshold </td>
          <td>(</td>
          <td class="paramtype">const RealScalar &amp;&#160;</td>
          <td class="paramname"><em>threshold</em></td><td>)</td>
          <td></td>
        </tr>
      </table>
  </td>
  <td class="mlabels-right">
<span class="mlabels"><span class="mlabel">inline</span></span>  </td>
  </tr>
</table>
</div><div class="memdoc">
<p>Allows to prescribe a threshold to be used by certain methods, such as <a class="el" href="classEigen_1_1FullPivLU.html#a363d1c09d77f09d6ea2d2789776e7be3">rank()</a>, who need to determine when pivots are to be considered nonzero. This is not used for the LU decomposition itself.</p>
<p>When it needs to get the threshold value, <a class="el" href="namespaceEigen.html" title="Namespace containing all symbols from the Eigen library. ">Eigen</a> calls <a class="el" href="classEigen_1_1FullPivLU.html#aa5a87faaa92a3081045d1f934e292ef0">threshold()</a>. By default, this uses a formula to automatically determine a reasonable threshold. Once you have called the present method <a class="el" href="classEigen_1_1FullPivLU.html#ae47399e3cf7943075ce18ef89fe17f21">setThreshold(const RealScalar&amp;)</a>, your value is used instead.</p>
<dl class="params"><dt>Parameters</dt><dd>
  <table class="params">
    <tr><td class="paramname">threshold</td><td>The new value to use as the threshold.</td></tr>
  </table>
  </dd>
</dl>
<p>A pivot will be considered nonzero if its absolute value is strictly greater than <img class="formulaInl" alt="$ \vert pivot \vert \leqslant threshold \times \vert maxpivot \vert $" src="form_156.png"/> where maxpivot is the biggest pivot.</p>
<p>If you want to come back to the default behavior, call <a class="el" href="classEigen_1_1FullPivLU.html#a02b74b9b823e060825443fd82cd25fd4">setThreshold(Default_t)</a> </p>

<p>References <a class="el" href="classEigen_1_1FullPivLU.html#aa5a87faaa92a3081045d1f934e292ef0">FullPivLU&lt; MatrixType &gt;::threshold()</a>.</p>

</div>
</div>
<a class="anchor" id="a02b74b9b823e060825443fd82cd25fd4"></a>
<div class="memitem">
<div class="memproto">
<table class="mlabels">
  <tr>
  <td class="mlabels-left">
      <table class="memname">
        <tr>
          <td class="memname"><a class="el" href="classEigen_1_1FullPivLU.html">FullPivLU</a>&amp; setThreshold </td>
          <td>(</td>
          <td class="paramtype">Default_t&#160;</td>
          <td class="paramname"></td><td>)</td>
          <td></td>
        </tr>
      </table>
  </td>
  <td class="mlabels-right">
<span class="mlabels"><span class="mlabel">inline</span></span>  </td>
  </tr>
</table>
</div><div class="memdoc">
<p>Allows to come back to the default behavior, letting <a class="el" href="namespaceEigen.html" title="Namespace containing all symbols from the Eigen library. ">Eigen</a> use its default formula for determining the threshold.</p>
<p>You should pass the special object Eigen::Default as parameter here. </p>
<div class="fragment"><div class="line">lu.setThreshold(Eigen::Default); </div>
</div><!-- fragment --><p>See the documentation of <a class="el" href="classEigen_1_1FullPivLU.html#ae47399e3cf7943075ce18ef89fe17f21">setThreshold(const RealScalar&amp;)</a>. </p>

</div>
</div>
<a class="anchor" id="ac1a642d728c059c7625863f126e2e718"></a>
<div class="memitem">
<div class="memproto">
<table class="mlabels">
  <tr>
  <td class="mlabels-left">
      <table class="memname">
        <tr>
          <td class="memname">const internal::solve_retval&lt;<a class="el" href="classEigen_1_1FullPivLU.html">FullPivLU</a>, Rhs&gt; solve </td>
          <td>(</td>
          <td class="paramtype">const <a class="el" href="classEigen_1_1MatrixBase.html">MatrixBase</a>&lt; Rhs &gt; &amp;&#160;</td>
          <td class="paramname"><em>b</em></td><td>)</td>
          <td> const</td>
        </tr>
      </table>
  </td>
  <td class="mlabels-right">
<span class="mlabels"><span class="mlabel">inline</span></span>  </td>
  </tr>
</table>
</div><div class="memdoc">
<dl class="section return"><dt>Returns</dt><dd>a solution x to the equation Ax=b, where A is the matrix of which *this is the LU decomposition.</dd></dl>
<dl class="params"><dt>Parameters</dt><dd>
  <table class="params">
    <tr><td class="paramname">b</td><td>the right-hand-side of the equation to solve. Can be a vector or a matrix, the only requirement in order for the equation to make sense is that b.rows()==A.rows(), where A is the matrix of which *this is the LU decomposition.</td></tr>
  </table>
  </dd>
</dl>
<dl class="section return"><dt>Returns</dt><dd>a solution.</dd></dl>
<p>This method just tries to find as good a solution as possible. If you want to check whether a solution exists or if it is accurate, just call this function to get a result and then compute the error of this result, or use <a class="el" href="classEigen_1_1DenseBase.html#a158c2184951e6e415c2e9b98db8e8966">MatrixBase::isApprox()</a> directly, for instance like this:</p>
<div class="fragment"><div class="line"><span class="keywordtype">bool</span> a_solution_exists = (A*result).isApprox(b, precision); </div>
</div><!-- fragment --><p> This method avoids dividing by zero, so that the non-existence of a solution doesn't by itself mean that you'll get <code>inf</code> or <code>nan</code> values.</p>
<p>If there exists more than one solution, this method will arbitrarily choose one. If you need a complete analysis of the space of solutions, take the one solution obtained by this method and add to it elements of the kernel, as determined by <a class="el" href="classEigen_1_1FullPivLU.html#a6e8f1d2fcbd86d3dc5a8a013b6e7200a">kernel()</a>.</p>
<p>Example: </p>
<div class="fragment"><div class="line">Matrix&lt;float,2,3&gt; m = <a class="code" href="classEigen_1_1DenseBase.html#a8e759dafdd9ecc446d397b7f5435f60a">Matrix&lt;float,2,3&gt;::Random</a>();</div>
<div class="line"><a class="code" href="group__matrixtypedefs.html#ga535a919504bb3bc463b8995c196c1eed">Matrix2f</a> y = <a class="code" href="classEigen_1_1DenseBase.html#a8e759dafdd9ecc446d397b7f5435f60a">Matrix2f::Random</a>();</div>
<div class="line">cout &lt;&lt; <span class="stringliteral">&quot;Here is the matrix m:&quot;</span> &lt;&lt; endl &lt;&lt; m &lt;&lt; endl;</div>
<div class="line">cout &lt;&lt; <span class="stringliteral">&quot;Here is the matrix y:&quot;</span> &lt;&lt; endl &lt;&lt; y &lt;&lt; endl;</div>
<div class="line">Matrix&lt;float,3,2&gt; x = m.fullPivLu().solve(y);</div>
<div class="line"><span class="keywordflow">if</span>((m*x).isApprox(y))</div>
<div class="line">{</div>
<div class="line">  cout &lt;&lt; <span class="stringliteral">&quot;Here is a solution x to the equation mx=y:&quot;</span> &lt;&lt; endl &lt;&lt; x &lt;&lt; endl;</div>
<div class="line">}</div>
<div class="line"><span class="keywordflow">else</span></div>
<div class="line">  cout &lt;&lt; <span class="stringliteral">&quot;The equation mx=y does not have any solution.&quot;</span> &lt;&lt; endl;</div>
</div><!-- fragment --><p> Output: </p>
<pre class="fragment">Here is the matrix m:
  0.68  0.566  0.823
-0.211  0.597 -0.605
Here is the matrix y:
 -0.33 -0.444
 0.536  0.108
Here is a solution x to the equation mx=y:
     0      0
 0.291 -0.216
  -0.6 -0.391
</pre><dl class="section see"><dt>See Also</dt><dd>TriangularView::solve(), <a class="el" href="classEigen_1_1FullPivLU.html#a6e8f1d2fcbd86d3dc5a8a013b6e7200a">kernel()</a>, <a class="el" href="classEigen_1_1FullPivLU.html#a5ca20e3802e96fb14c2be37039afcae9">inverse()</a> </dd></dl>

</div>
</div>
<a class="anchor" id="aa5a87faaa92a3081045d1f934e292ef0"></a>
<div class="memitem">
<div class="memproto">
<table class="mlabels">
  <tr>
  <td class="mlabels-left">
      <table class="memname">
        <tr>
          <td class="memname">RealScalar threshold </td>
          <td>(</td>
          <td class="paramname"></td><td>)</td>
          <td> const</td>
        </tr>
      </table>
  </td>
  <td class="mlabels-right">
<span class="mlabels"><span class="mlabel">inline</span></span>  </td>
  </tr>
</table>
</div><div class="memdoc">
<p>Returns the threshold that will be used by certain methods such as <a class="el" href="classEigen_1_1FullPivLU.html#a363d1c09d77f09d6ea2d2789776e7be3">rank()</a>.</p>
<p>See the documentation of <a class="el" href="classEigen_1_1FullPivLU.html#ae47399e3cf7943075ce18ef89fe17f21">setThreshold(const RealScalar&amp;)</a>. </p>

<p>Referenced by <a class="el" href="classEigen_1_1FullPivLU.html#a363d1c09d77f09d6ea2d2789776e7be3">FullPivLU&lt; MatrixType &gt;::rank()</a>, and <a class="el" href="classEigen_1_1FullPivLU.html#ae47399e3cf7943075ce18ef89fe17f21">FullPivLU&lt; MatrixType &gt;::setThreshold()</a>.</p>

</div>
</div>
<hr/>The documentation for this class was generated from the following files:<ul>
<li><a class="el" href="ForwardDeclarations_8h_source.html">ForwardDeclarations.h</a></li>
<li><a class="el" href="FullPivLU_8h_source.html">FullPivLU.h</a></li>
</ul>
</div><!-- contents -->
</div><!-- doc-content -->
<!-- start footer part -->
<div id="nav-path" class="navpath"><!-- id is needed for treeview function! -->
  <ul>
    <li class="navelem"><a class="el" href="namespaceEigen.html">Eigen</a></li><li class="navelem"><a class="el" href="classEigen_1_1FullPivLU.html">FullPivLU</a></li>
    <li class="footer">Generated on Mon Oct 28 2013 11:04:29 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>