Sophie

Sophie

distrib > Mageia > 4 > i586 > by-pkgid > 99cb5ede6a5329071fbeecc8218deb35 > files > 474

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: LLT&lt; _MatrixType, _UpLo &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_1LLT.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_1LLT-members.html">List of all members</a> &#124;
<a href="#pub-methods">Public Member Functions</a>  </div>
  <div class="headertitle">
<div class="title">LLT&lt; _MatrixType, _UpLo &gt; Class Template Reference<div class="ingroups"><a class="el" href="group__Cholesky__Module.html">Cholesky 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, int _UpLo&gt;<br/>
class Eigen::LLT&lt; _MatrixType, _UpLo &gt;</h3>

<p>Standard Cholesky decomposition (LL^T) of a matrix and associated 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 LL^T Cholesky decomposition </td></tr>
    <tr><td class="paramname">UpLo</td><td>the triangular part that will be used for the decompositon: Lower (default) or Upper. The other triangular part won't be read.</td></tr>
  </table>
  </dd>
</dl>
<p>This class performs a LL^T Cholesky decomposition of a symmetric, positive definite matrix A such that A = LL^* = U^*U, where L is lower triangular.</p>
<p>While the Cholesky decomposition is particularly useful to solve selfadjoint problems like D^*D x = b, for that purpose, we recommend the Cholesky decomposition without square root which is more stable and even faster. Nevertheless, this standard Cholesky decomposition remains useful in many other situations like generalised eigen problems with hermitian matrices.</p>
<p>Remember that Cholesky decompositions are not rank-revealing. This <a class="el" href="classEigen_1_1LLT.html" title="Standard Cholesky decomposition (LL^T) of a matrix and associated features. ">LLT</a> decomposition is only stable on positive definite matrices, use <a class="el" href="classEigen_1_1LDLT.html" title="Robust Cholesky decomposition of a matrix with pivoting. ">LDLT</a> instead for the semidefinite case. Also, do not use a Cholesky decomposition to determine whether a system of equations has a solution.</p>
<p>Example: </p>
<div class="fragment"><div class="line">MatrixXd A(3,3);</div>
<div class="line">A &lt;&lt; 4,-1,2, -1,6,0, 2,0,5;</div>
<div class="line">cout &lt;&lt; <span class="stringliteral">&quot;The matrix A is&quot;</span> &lt;&lt; endl &lt;&lt; A &lt;&lt; endl;</div>
<div class="line"></div>
<div class="line">LLT&lt;MatrixXd&gt; lltOfA(A); <span class="comment">// compute the Cholesky decomposition of A</span></div>
<div class="line">MatrixXd L = lltOfA.matrixL(); <span class="comment">// retrieve factor L  in the decomposition</span></div>
<div class="line"><span class="comment">// The previous two lines can also be written as &quot;L = A.llt().matrixL()&quot;</span></div>
<div class="line"></div>
<div class="line">cout &lt;&lt; <span class="stringliteral">&quot;The Cholesky factor L is&quot;</span> &lt;&lt; endl &lt;&lt; L &lt;&lt; endl;</div>
<div class="line">cout &lt;&lt; <span class="stringliteral">&quot;To check this, let us compute L * L.transpose()&quot;</span> &lt;&lt; endl;</div>
<div class="line">cout &lt;&lt; L * L.transpose() &lt;&lt; endl;</div>
<div class="line">cout &lt;&lt; <span class="stringliteral">&quot;This should equal the matrix A&quot;</span> &lt;&lt; endl;</div>
</div><!-- fragment --><p> Output: </p>
<pre class="fragment">The matrix A is
 4 -1  2
-1  6  0
 2  0  5
The Cholesky factor L is
    2     0     0
 -0.5   2.4     0
    1 0.209  1.99
To check this, let us compute L * L.transpose()
 4 -1  2
-1  6  0
 2  0  5
This should equal the matrix A
</pre><dl class="section see"><dt>See Also</dt><dd><a class="el" href="classEigen_1_1MatrixBase.html#a72ca76a1df449fd21a21d15ac1d1042a">MatrixBase::llt()</a>, class <a class="el" href="classEigen_1_1LDLT.html" title="Robust Cholesky decomposition of a matrix with pivoting. ">LDLT</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:ac97b2f0b61df78a30544670b3711e3db"><td class="memItemLeft" align="right" valign="top"><a class="el" href="classEigen_1_1LLT.html">LLT</a> &amp;&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1LLT.html#ac97b2f0b61df78a30544670b3711e3db">compute</a> (const MatrixType &amp;matrix)</td></tr>
<tr class="separator:ac97b2f0b61df78a30544670b3711e3db"><td class="memSeparator" colspan="2">&#160;</td></tr>
<tr class="memitem:a0c06d5c2034ebb329c54235369643ad2"><td class="memItemLeft" align="right" valign="top"><a class="el" href="group__enums.html#ga51bc1ac16f26ebe51eae1abb77bd037b">ComputationInfo</a>&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1LLT.html#a0c06d5c2034ebb329c54235369643ad2">info</a> () const </td></tr>
<tr class="memdesc:a0c06d5c2034ebb329c54235369643ad2"><td class="mdescLeft">&#160;</td><td class="mdescRight">Reports whether previous computation was successful.  <a href="#a0c06d5c2034ebb329c54235369643ad2">More...</a><br/></td></tr>
<tr class="separator:a0c06d5c2034ebb329c54235369643ad2"><td class="memSeparator" colspan="2">&#160;</td></tr>
<tr class="memitem:a5bb9ddb29f0e0710acfe830794571f25"><td class="memItemLeft" align="right" valign="top">&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1LLT.html#a5bb9ddb29f0e0710acfe830794571f25">LLT</a> ()</td></tr>
<tr class="memdesc:a5bb9ddb29f0e0710acfe830794571f25"><td class="mdescLeft">&#160;</td><td class="mdescRight">Default Constructor.  <a href="#a5bb9ddb29f0e0710acfe830794571f25">More...</a><br/></td></tr>
<tr class="separator:a5bb9ddb29f0e0710acfe830794571f25"><td class="memSeparator" colspan="2">&#160;</td></tr>
<tr class="memitem:ad14ba9fddbad27dbf73bafd1128f4f93"><td class="memItemLeft" align="right" valign="top">&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1LLT.html#ad14ba9fddbad27dbf73bafd1128f4f93">LLT</a> (Index size)</td></tr>
<tr class="memdesc:ad14ba9fddbad27dbf73bafd1128f4f93"><td class="mdescLeft">&#160;</td><td class="mdescRight">Default Constructor with memory preallocation.  <a href="#ad14ba9fddbad27dbf73bafd1128f4f93">More...</a><br/></td></tr>
<tr class="separator:ad14ba9fddbad27dbf73bafd1128f4f93"><td class="memSeparator" colspan="2">&#160;</td></tr>
<tr class="memitem:a625f69b684b4434db1cf3cc434e86fe6"><td class="memItemLeft" align="right" valign="top">Traits::MatrixL&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1LLT.html#a625f69b684b4434db1cf3cc434e86fe6">matrixL</a> () const </td></tr>
<tr class="separator:a625f69b684b4434db1cf3cc434e86fe6"><td class="memSeparator" colspan="2">&#160;</td></tr>
<tr class="memitem:a0b0ede52ac49b086994a447a321b1bdb"><td class="memItemLeft" align="right" valign="top">const MatrixType &amp;&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1LLT.html#a0b0ede52ac49b086994a447a321b1bdb">matrixLLT</a> () const </td></tr>
<tr class="separator:a0b0ede52ac49b086994a447a321b1bdb"><td class="memSeparator" colspan="2">&#160;</td></tr>
<tr class="memitem:adc32fdba9f5b478afb2a96d53c6eacbb"><td class="memItemLeft" align="right" valign="top">Traits::MatrixU&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1LLT.html#adc32fdba9f5b478afb2a96d53c6eacbb">matrixU</a> () const </td></tr>
<tr class="separator:adc32fdba9f5b478afb2a96d53c6eacbb"><td class="memSeparator" colspan="2">&#160;</td></tr>
<tr class="memitem:a67b49e986da81c3c90fe1793681e70ad"><td class="memTemplParams" colspan="2">template&lt;typename VectorType &gt; </td></tr>
<tr class="memitem:a67b49e986da81c3c90fe1793681e70ad"><td class="memTemplItemLeft" align="right" valign="top"><a class="el" href="classEigen_1_1LLT.html">LLT</a>&lt; _MatrixType, _UpLo &gt;&#160;</td><td class="memTemplItemRight" valign="bottom"><a class="el" href="classEigen_1_1LLT.html#a67b49e986da81c3c90fe1793681e70ad">rankUpdate</a> (const VectorType &amp;v, const RealScalar &amp;sigma)</td></tr>
<tr class="separator:a67b49e986da81c3c90fe1793681e70ad"><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_1LLT.html#ac36925ac693435a090efee1cb5d6d16a">reconstructedMatrix</a> () const </td></tr>
<tr class="separator:ac36925ac693435a090efee1cb5d6d16a"><td class="memSeparator" colspan="2">&#160;</td></tr>
<tr class="memitem:a301e28d205fb0308cb4f4b04718e1685"><td class="memTemplParams" colspan="2">template&lt;typename Rhs &gt; </td></tr>
<tr class="memitem:a301e28d205fb0308cb4f4b04718e1685"><td class="memTemplItemLeft" align="right" valign="top">const internal::solve_retval<br class="typebreak"/>
&lt; <a class="el" href="classEigen_1_1LLT.html">LLT</a>, Rhs &gt;&#160;</td><td class="memTemplItemRight" valign="bottom"><a class="el" href="classEigen_1_1LLT.html#a301e28d205fb0308cb4f4b04718e1685">solve</a> (const <a class="el" href="classEigen_1_1MatrixBase.html">MatrixBase</a>&lt; Rhs &gt; &amp;b) const </td></tr>
<tr class="separator:a301e28d205fb0308cb4f4b04718e1685"><td class="memSeparator" colspan="2">&#160;</td></tr>
</table>
<h2 class="groupheader">Constructor &amp; Destructor Documentation</h2>
<a class="anchor" id="a5bb9ddb29f0e0710acfe830794571f25"></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_1LLT.html">LLT</a> </td>
          <td>(</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>Default Constructor. </p>
<p>The default constructor is useful in cases in which the user intends to perform decompositions via <a class="el" href="classEigen_1_1LLT.html#ac97b2f0b61df78a30544670b3711e3db">LLT::compute(const MatrixType&amp;)</a>. </p>

</div>
</div>
<a class="anchor" id="ad14ba9fddbad27dbf73bafd1128f4f93"></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_1LLT.html">LLT</a> </td>
          <td>(</td>
          <td class="paramtype">Index&#160;</td>
          <td class="paramname"><em>size</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>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_1LLT.html#a5bb9ddb29f0e0710acfe830794571f25" title="Default Constructor. ">LLT()</a> </dd></dl>

</div>
</div>
<h2 class="groupheader">Member Function Documentation</h2>
<a class="anchor" id="ac97b2f0b61df78a30544670b3711e3db"></a>
<div class="memitem">
<div class="memproto">
      <table class="memname">
        <tr>
          <td class="memname"><a class="el" href="classEigen_1_1LLT.html">LLT</a>&lt; MatrixType, _UpLo &gt; &amp; compute </td>
          <td>(</td>
          <td class="paramtype">const MatrixType &amp;&#160;</td>
          <td class="paramname"><em>a</em></td><td>)</td>
          <td></td>
        </tr>
      </table>
</div><div class="memdoc">
<p>Computes / recomputes the Cholesky decomposition A = LL^* = U^*U of <em>matrix</em> </p>
<dl class="section return"><dt>Returns</dt><dd>a reference to *this</dd></dl>
<p>Example: </p>
<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 --><p> Output: </p>
<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> 
<p>References <a class="el" href="group__enums.html#gga51bc1ac16f26ebe51eae1abb77bd037ba710fff14e8fc77846d4b75d8f4cc2d5c">Eigen::NumericalIssue</a>, and <a class="el" href="group__enums.html#gga51bc1ac16f26ebe51eae1abb77bd037bafdfbdf3247bd36a1f17270d5cec74c9c">Eigen::Success</a>.</p>

</div>
</div>
<a class="anchor" id="a0c06d5c2034ebb329c54235369643ad2"></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="group__enums.html#ga51bc1ac16f26ebe51eae1abb77bd037b">ComputationInfo</a> info </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>Reports whether previous computation was successful. </p>
<dl class="section return"><dt>Returns</dt><dd><code>Success</code> if computation was succesful, <code>NumericalIssue</code> if the matrix.appears to be negative. </dd></dl>

</div>
</div>
<a class="anchor" id="a625f69b684b4434db1cf3cc434e86fe6"></a>
<div class="memitem">
<div class="memproto">
<table class="mlabels">
  <tr>
  <td class="mlabels-left">
      <table class="memname">
        <tr>
          <td class="memname">Traits::MatrixL matrixL </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>a view of the lower triangular matrix L </dd></dl>

<p>Referenced by <a class="el" href="classEigen_1_1GeneralizedSelfAdjointEigenSolver.html#aaa204ef15aaefac270c0376269083ed6">GeneralizedSelfAdjointEigenSolver&lt; _MatrixType &gt;::compute()</a>.</p>

</div>
</div>
<a class="anchor" id="a0b0ede52ac49b086994a447a321b1bdb"></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; matrixLLT </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 <a class="el" href="classEigen_1_1LLT.html" title="Standard Cholesky decomposition (LL^T) of a matrix and associated features. ">LLT</a> decomposition matrix</dd></dl>
<p>TODO: document the storage layout </p>

</div>
</div>
<a class="anchor" id="adc32fdba9f5b478afb2a96d53c6eacbb"></a>
<div class="memitem">
<div class="memproto">
<table class="mlabels">
  <tr>
  <td class="mlabels-left">
      <table class="memname">
        <tr>
          <td class="memname">Traits::MatrixU matrixU </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>a view of the upper triangular matrix U </dd></dl>

<p>Referenced by <a class="el" href="classEigen_1_1GeneralizedSelfAdjointEigenSolver.html#aaa204ef15aaefac270c0376269083ed6">GeneralizedSelfAdjointEigenSolver&lt; _MatrixType &gt;::compute()</a>.</p>

</div>
</div>
<a class="anchor" id="a67b49e986da81c3c90fe1793681e70ad"></a>
<div class="memitem">
<div class="memproto">
      <table class="memname">
        <tr>
          <td class="memname"><a class="el" href="classEigen_1_1LLT.html">LLT</a>&lt;_MatrixType,_UpLo&gt; rankUpdate </td>
          <td>(</td>
          <td class="paramtype">const VectorType &amp;&#160;</td>
          <td class="paramname"><em>v</em>, </td>
        </tr>
        <tr>
          <td class="paramkey"></td>
          <td></td>
          <td class="paramtype">const RealScalar &amp;&#160;</td>
          <td class="paramname"><em>sigma</em>&#160;</td>
        </tr>
        <tr>
          <td></td>
          <td>)</td>
          <td></td><td></td>
        </tr>
      </table>
</div><div class="memdoc">
<p>Performs a rank one update (or dowdate) of the current decomposition. If A = LL^* before the rank one update, then after it we have LL^* = A + sigma * v v^* where <em>v</em> must be a vector of same dimension. </p>

<p>References <a class="el" href="group__enums.html#gga51bc1ac16f26ebe51eae1abb77bd037ba710fff14e8fc77846d4b75d8f4cc2d5c">Eigen::NumericalIssue</a>, and <a class="el" href="group__enums.html#gga51bc1ac16f26ebe51eae1abb77bd037bafdfbdf3247bd36a1f17270d5cec74c9c">Eigen::Success</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: L L^*. This function is provided for debug purpose. </dd></dl>

</div>
</div>
<a class="anchor" id="a301e28d205fb0308cb4f4b04718e1685"></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_1LLT.html">LLT</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>the solution x of <img class="formulaInl" alt="$ A x = b $" src="form_3.png"/> using the current decomposition of A.</dd></dl>
<p>Since this <a class="el" href="classEigen_1_1LLT.html" title="Standard Cholesky decomposition (LL^T) of a matrix and associated features. ">LLT</a> class assumes anyway that the matrix A is invertible, the solution theoretically exists and is unique regardless of b.</p>
<p>Example: </p>
<div class="fragment"><div class="line"><span class="keyword">typedef</span> Matrix&lt;float,Dynamic,2&gt; DataMatrix;</div>
<div class="line"><span class="comment">// let&#39;s generate some samples on the 3D plane of equation z = 2x+3y (with some noise)</span></div>
<div class="line">DataMatrix samples = DataMatrix::Random(12,2);</div>
<div class="line">VectorXf elevations = 2*samples.col(0) + 3*samples.col(1) + <a class="code" href="classEigen_1_1DenseBase.html#a8e759dafdd9ecc446d397b7f5435f60a">VectorXf::Random</a>(12)*0.1;</div>
<div class="line"><span class="comment">// and let&#39;s solve samples * [x y]^T = elevations in least square sense:</span></div>
<div class="line">Matrix&lt;float,2,1&gt; xy</div>
<div class="line"> = (samples.adjoint() * samples).llt().<a class="code" href="classEigen_1_1LLT.html#a301e28d205fb0308cb4f4b04718e1685">solve</a>((samples.adjoint()*elevations));</div>
<div class="line">cout &lt;&lt; xy &lt;&lt; endl;</div>
</div><!-- fragment --><p> Output: </p>
<pre class="fragment">2.02
2.97
</pre><dl class="section see"><dt>See Also</dt><dd>solveInPlace(), <a class="el" href="classEigen_1_1MatrixBase.html#a72ca76a1df449fd21a21d15ac1d1042a">MatrixBase::llt()</a> </dd></dl>

</div>
</div>
<hr/>The documentation for this class was generated from the following file:<ul>
<li><a class="el" href="LLT_8h_source.html">LLT.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_1LLT.html">LLT</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>