Sophie

Sophie

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

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: Tridiagonalization&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_1Tridiagonalization.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_1Tridiagonalization-members.html">List of all members</a> &#124;
<a href="#pub-types">Public Types</a> &#124;
<a href="#pub-methods">Public Member Functions</a>  </div>
  <div class="headertitle">
<div class="title">Tridiagonalization&lt; _MatrixType &gt; Class Template Reference<div class="ingroups"><a class="el" href="group__Eigenvalues__Module.html">Eigenvalues 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::Tridiagonalization&lt; _MatrixType &gt;</h3>

<p>Tridiagonal decomposition of a selfadjoint matrix. </p>
<p>This is defined in the Eigenvalues module.</p>
<div class="fragment"><div class="line"><span class="preprocessor">#include &lt;Eigen/Eigenvalues&gt;</span> </div>
</div><!-- fragment --><dl class="tparams"><dt>Template Parameters</dt><dd>
  <table class="tparams">
    <tr><td class="paramname">_MatrixType</td><td>the type of the matrix of which we are computing the tridiagonal decomposition; this is expected to be an instantiation of the <a class="el" href="classEigen_1_1Matrix.html" title="The matrix class, also used for vectors and row-vectors. ">Matrix</a> class template.</td></tr>
  </table>
  </dd>
</dl>
<p>This class performs a tridiagonal decomposition of a selfadjoint matrix <img class="formulaInl" alt="$ A $" src="form_1.png"/> such that: <img class="formulaInl" alt="$ A = Q T Q^* $" src="form_98.png"/> where <img class="formulaInl" alt="$ Q $" src="form_73.png"/> is unitary and <img class="formulaInl" alt="$ T $" src="form_99.png"/> a real symmetric tridiagonal matrix.</p>
<p>A tridiagonal matrix is a matrix which has nonzero elements only on the main diagonal and the first diagonal below and above it. The Hessenberg decomposition of a selfadjoint matrix is in fact a tridiagonal decomposition. This class is used in <a class="el" href="classEigen_1_1SelfAdjointEigenSolver.html" title="Computes eigenvalues and eigenvectors of selfadjoint matrices. ">SelfAdjointEigenSolver</a> to compute the eigenvalues and eigenvectors of a selfadjoint matrix.</p>
<p>Call the function <a class="el" href="classEigen_1_1Tridiagonalization.html#aa69e607a4aab4fb6321ca6acbf074fc2" title="Computes tridiagonal decomposition of given matrix. ">compute()</a> to compute the tridiagonal decomposition of a given matrix. Alternatively, you can use the <a class="el" href="classEigen_1_1Tridiagonalization.html#aa9f9722d2cef9425e2c0da3553dfbac7" title="Constructor; computes tridiagonal decomposition of given matrix. ">Tridiagonalization(const MatrixType&amp;)</a> constructor which computes the tridiagonal Schur decomposition at construction time. Once the decomposition is computed, you can use the <a class="el" href="classEigen_1_1Tridiagonalization.html#ad13845d7490115664924b3dc208ec369" title="Returns the unitary matrix Q in the decomposition. ">matrixQ()</a> and <a class="el" href="classEigen_1_1Tridiagonalization.html#aceb0f16a166f4c236a1b536b7424d292" title="Returns an expression of the tridiagonal matrix T in the decomposition. ">matrixT()</a> functions to retrieve the matrices Q and T in the decomposition.</p>
<p>The documentation of <a class="el" href="classEigen_1_1Tridiagonalization.html#aa9f9722d2cef9425e2c0da3553dfbac7" title="Constructor; computes tridiagonal decomposition of given matrix. ">Tridiagonalization(const MatrixType&amp;)</a> contains an example of the typical use of this class.</p>
<dl class="section see"><dt>See Also</dt><dd>class <a class="el" href="classEigen_1_1HessenbergDecomposition.html" title="Reduces a square matrix to Hessenberg form by an orthogonal similarity transformation. ">HessenbergDecomposition</a>, class <a class="el" href="classEigen_1_1SelfAdjointEigenSolver.html" title="Computes eigenvalues and eigenvectors of selfadjoint matrices. ">SelfAdjointEigenSolver</a> </dd></dl>
</div><table class="memberdecls">
<tr class="heading"><td colspan="2"><h2 class="groupheader"><a name="pub-types"></a>
Public Types</h2></td></tr>
<tr class="memitem:aa96bdbc1b19c647e3372c31301ea4999"><td class="memItemLeft" align="right" valign="top"><a class="anchor" id="aa96bdbc1b19c647e3372c31301ea4999"></a>
typedef <a class="el" href="classEigen_1_1HouseholderSequence.html">HouseholderSequence</a><br class="typebreak"/>
&lt; <a class="el" href="classEigen_1_1Tridiagonalization.html#aeb6c0eb89cc982629305f6c7e0791caf">MatrixType</a>, typename <br class="typebreak"/>
internal::remove_all&lt; typename <br class="typebreak"/>
CoeffVectorType::ConjugateReturnType &gt;<br class="typebreak"/>
::type &gt;&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1Tridiagonalization.html#aa96bdbc1b19c647e3372c31301ea4999">HouseholderSequenceType</a></td></tr>
<tr class="memdesc:aa96bdbc1b19c647e3372c31301ea4999"><td class="mdescLeft">&#160;</td><td class="mdescRight">Return type of <a class="el" href="classEigen_1_1Tridiagonalization.html#ad13845d7490115664924b3dc208ec369" title="Returns the unitary matrix Q in the decomposition. ">matrixQ()</a> <br/></td></tr>
<tr class="separator:aa96bdbc1b19c647e3372c31301ea4999"><td class="memSeparator" colspan="2">&#160;</td></tr>
<tr class="memitem:aeb6c0eb89cc982629305f6c7e0791caf"><td class="memItemLeft" align="right" valign="top"><a class="anchor" id="aeb6c0eb89cc982629305f6c7e0791caf"></a>
typedef _MatrixType&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1Tridiagonalization.html#aeb6c0eb89cc982629305f6c7e0791caf">MatrixType</a></td></tr>
<tr class="memdesc:aeb6c0eb89cc982629305f6c7e0791caf"><td class="mdescLeft">&#160;</td><td class="mdescRight">Synonym for the template parameter <code>_MatrixType</code>. <br/></td></tr>
<tr class="separator:aeb6c0eb89cc982629305f6c7e0791caf"><td class="memSeparator" colspan="2">&#160;</td></tr>
</table><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:aa69e607a4aab4fb6321ca6acbf074fc2"><td class="memItemLeft" align="right" valign="top"><a class="el" href="classEigen_1_1Tridiagonalization.html">Tridiagonalization</a> &amp;&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1Tridiagonalization.html#aa69e607a4aab4fb6321ca6acbf074fc2">compute</a> (const <a class="el" href="classEigen_1_1Tridiagonalization.html#aeb6c0eb89cc982629305f6c7e0791caf">MatrixType</a> &amp;matrix)</td></tr>
<tr class="memdesc:aa69e607a4aab4fb6321ca6acbf074fc2"><td class="mdescLeft">&#160;</td><td class="mdescRight">Computes tridiagonal decomposition of given matrix.  <a href="#aa69e607a4aab4fb6321ca6acbf074fc2">More...</a><br/></td></tr>
<tr class="separator:aa69e607a4aab4fb6321ca6acbf074fc2"><td class="memSeparator" colspan="2">&#160;</td></tr>
<tr class="memitem:ac109eefddd733d8e82841da5bb2dd8d3"><td class="memItemLeft" align="right" valign="top">DiagonalReturnType&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1Tridiagonalization.html#ac109eefddd733d8e82841da5bb2dd8d3">diagonal</a> () const </td></tr>
<tr class="memdesc:ac109eefddd733d8e82841da5bb2dd8d3"><td class="mdescLeft">&#160;</td><td class="mdescRight">Returns the diagonal of the tridiagonal matrix T in the decomposition.  <a href="#ac109eefddd733d8e82841da5bb2dd8d3">More...</a><br/></td></tr>
<tr class="separator:ac109eefddd733d8e82841da5bb2dd8d3"><td class="memSeparator" colspan="2">&#160;</td></tr>
<tr class="memitem:a2ab889c75460c178d941ee24e371b206"><td class="memItemLeft" align="right" valign="top"><a class="el" href="classEigen_1_1Matrix.html">CoeffVectorType</a>&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1Tridiagonalization.html#a2ab889c75460c178d941ee24e371b206">householderCoefficients</a> () const </td></tr>
<tr class="memdesc:a2ab889c75460c178d941ee24e371b206"><td class="mdescLeft">&#160;</td><td class="mdescRight">Returns the Householder coefficients.  <a href="#a2ab889c75460c178d941ee24e371b206">More...</a><br/></td></tr>
<tr class="separator:a2ab889c75460c178d941ee24e371b206"><td class="memSeparator" colspan="2">&#160;</td></tr>
<tr class="memitem:ad13845d7490115664924b3dc208ec369"><td class="memItemLeft" align="right" valign="top"><a class="el" href="classEigen_1_1Tridiagonalization.html#aa96bdbc1b19c647e3372c31301ea4999">HouseholderSequenceType</a>&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1Tridiagonalization.html#ad13845d7490115664924b3dc208ec369">matrixQ</a> () const </td></tr>
<tr class="memdesc:ad13845d7490115664924b3dc208ec369"><td class="mdescLeft">&#160;</td><td class="mdescRight">Returns the unitary matrix Q in the decomposition.  <a href="#ad13845d7490115664924b3dc208ec369">More...</a><br/></td></tr>
<tr class="separator:ad13845d7490115664924b3dc208ec369"><td class="memSeparator" colspan="2">&#160;</td></tr>
<tr class="memitem:aceb0f16a166f4c236a1b536b7424d292"><td class="memItemLeft" align="right" valign="top">MatrixTReturnType&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1Tridiagonalization.html#aceb0f16a166f4c236a1b536b7424d292">matrixT</a> () const </td></tr>
<tr class="memdesc:aceb0f16a166f4c236a1b536b7424d292"><td class="mdescLeft">&#160;</td><td class="mdescRight">Returns an expression of the tridiagonal matrix T in the decomposition.  <a href="#aceb0f16a166f4c236a1b536b7424d292">More...</a><br/></td></tr>
<tr class="separator:aceb0f16a166f4c236a1b536b7424d292"><td class="memSeparator" colspan="2">&#160;</td></tr>
<tr class="memitem:a66adece364b64b26b3771662de70f2df"><td class="memItemLeft" align="right" valign="top">const <a class="el" href="classEigen_1_1Tridiagonalization.html#aeb6c0eb89cc982629305f6c7e0791caf">MatrixType</a> &amp;&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1Tridiagonalization.html#a66adece364b64b26b3771662de70f2df">packedMatrix</a> () const </td></tr>
<tr class="memdesc:a66adece364b64b26b3771662de70f2df"><td class="mdescLeft">&#160;</td><td class="mdescRight">Returns the internal representation of the decomposition.  <a href="#a66adece364b64b26b3771662de70f2df">More...</a><br/></td></tr>
<tr class="separator:a66adece364b64b26b3771662de70f2df"><td class="memSeparator" colspan="2">&#160;</td></tr>
<tr class="memitem:a8fa49216273ab7579b7bea06debb1e51"><td class="memItemLeft" align="right" valign="top">SubDiagonalReturnType&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1Tridiagonalization.html#a8fa49216273ab7579b7bea06debb1e51">subDiagonal</a> () const </td></tr>
<tr class="memdesc:a8fa49216273ab7579b7bea06debb1e51"><td class="mdescLeft">&#160;</td><td class="mdescRight">Returns the subdiagonal of the tridiagonal matrix T in the decomposition.  <a href="#a8fa49216273ab7579b7bea06debb1e51">More...</a><br/></td></tr>
<tr class="separator:a8fa49216273ab7579b7bea06debb1e51"><td class="memSeparator" colspan="2">&#160;</td></tr>
<tr class="memitem:a0698ae78b0ab6f239c475b73b9c6bbee"><td class="memItemLeft" align="right" valign="top">&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1Tridiagonalization.html#a0698ae78b0ab6f239c475b73b9c6bbee">Tridiagonalization</a> (Index size=Size==<a class="el" href="namespaceEigen.html#adc9da5be31bdce40c25a92c27999c0e3">Dynamic</a>?2:Size)</td></tr>
<tr class="memdesc:a0698ae78b0ab6f239c475b73b9c6bbee"><td class="mdescLeft">&#160;</td><td class="mdescRight">Default constructor.  <a href="#a0698ae78b0ab6f239c475b73b9c6bbee">More...</a><br/></td></tr>
<tr class="separator:a0698ae78b0ab6f239c475b73b9c6bbee"><td class="memSeparator" colspan="2">&#160;</td></tr>
<tr class="memitem:aa9f9722d2cef9425e2c0da3553dfbac7"><td class="memItemLeft" align="right" valign="top">&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1Tridiagonalization.html#aa9f9722d2cef9425e2c0da3553dfbac7">Tridiagonalization</a> (const <a class="el" href="classEigen_1_1Tridiagonalization.html#aeb6c0eb89cc982629305f6c7e0791caf">MatrixType</a> &amp;matrix)</td></tr>
<tr class="memdesc:aa9f9722d2cef9425e2c0da3553dfbac7"><td class="mdescLeft">&#160;</td><td class="mdescRight">Constructor; computes tridiagonal decomposition of given matrix.  <a href="#aa9f9722d2cef9425e2c0da3553dfbac7">More...</a><br/></td></tr>
<tr class="separator:aa9f9722d2cef9425e2c0da3553dfbac7"><td class="memSeparator" colspan="2">&#160;</td></tr>
</table>
<h2 class="groupheader">Constructor &amp; Destructor Documentation</h2>
<a class="anchor" id="a0698ae78b0ab6f239c475b73b9c6bbee"></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_1Tridiagonalization.html">Tridiagonalization</a> </td>
          <td>(</td>
          <td class="paramtype">Index&#160;</td>
          <td class="paramname"><em>size</em> = <code>Size==<a class="el" href="namespaceEigen.html#adc9da5be31bdce40c25a92c27999c0e3">Dynamic</a>&#160;?&#160;2&#160;:&#160;Size</code></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>
<dl class="params"><dt>Parameters</dt><dd>
  <table class="params">
    <tr><td class="paramdir">[in]</td><td class="paramname">size</td><td>Positive integer, size of the matrix whose tridiagonal decomposition will be computed.</td></tr>
  </table>
  </dd>
</dl>
<p>The default constructor is useful in cases in which the user intends to perform decompositions via <a class="el" href="classEigen_1_1Tridiagonalization.html#aa69e607a4aab4fb6321ca6acbf074fc2" title="Computes tridiagonal decomposition of given matrix. ">compute()</a>. The <code>size</code> parameter is only used as a hint. It is not an error to give a wrong <code>size</code>, but it may impair performance.</p>
<dl class="section see"><dt>See Also</dt><dd><a class="el" href="classEigen_1_1Tridiagonalization.html#aa69e607a4aab4fb6321ca6acbf074fc2" title="Computes tridiagonal decomposition of given matrix. ">compute()</a> for an example. </dd></dl>

</div>
</div>
<a class="anchor" id="aa9f9722d2cef9425e2c0da3553dfbac7"></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_1Tridiagonalization.html">Tridiagonalization</a> </td>
          <td>(</td>
          <td class="paramtype">const <a class="el" href="classEigen_1_1Tridiagonalization.html#aeb6c0eb89cc982629305f6c7e0791caf">MatrixType</a> &amp;&#160;</td>
          <td class="paramname"><em>matrix</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>Constructor; computes tridiagonal decomposition of given matrix. </p>
<dl class="params"><dt>Parameters</dt><dd>
  <table class="params">
    <tr><td class="paramdir">[in]</td><td class="paramname">matrix</td><td>Selfadjoint matrix whose tridiagonal decomposition is to be computed.</td></tr>
  </table>
  </dd>
</dl>
<p>This constructor calls <a class="el" href="classEigen_1_1Tridiagonalization.html#aa69e607a4aab4fb6321ca6acbf074fc2" title="Computes tridiagonal decomposition of given matrix. ">compute()</a> to compute the tridiagonal decomposition.</p>
<p>Example: </p>
<div class="fragment"><div class="line">MatrixXd X = <a class="code" href="classEigen_1_1DenseBase.html#a8e759dafdd9ecc446d397b7f5435f60a">MatrixXd::Random</a>(5,5);</div>
<div class="line">MatrixXd A = X + X.transpose();</div>
<div class="line">cout &lt;&lt; <span class="stringliteral">&quot;Here is a random symmetric 5x5 matrix:&quot;</span> &lt;&lt; endl &lt;&lt; A &lt;&lt; endl &lt;&lt; endl;</div>
<div class="line">Tridiagonalization&lt;MatrixXd&gt; triOfA(A);</div>
<div class="line">MatrixXd Q = triOfA.matrixQ();</div>
<div class="line">cout &lt;&lt; <span class="stringliteral">&quot;The orthogonal matrix Q is:&quot;</span> &lt;&lt; endl &lt;&lt; Q &lt;&lt; endl;</div>
<div class="line">MatrixXd T = triOfA.matrixT();</div>
<div class="line">cout &lt;&lt; <span class="stringliteral">&quot;The tridiagonal matrix T is:&quot;</span> &lt;&lt; endl &lt;&lt; T &lt;&lt; endl &lt;&lt; endl;</div>
<div class="line">cout &lt;&lt; <span class="stringliteral">&quot;Q * T * Q^T = &quot;</span> &lt;&lt; endl &lt;&lt; Q * T * Q.transpose() &lt;&lt; endl;</div>
</div><!-- fragment --><p> Output: </p>
<pre class="fragment">Here is a random symmetric 5x5 matrix:
  1.36 -0.816  0.521   1.43 -0.144
-0.816 -0.659  0.794 -0.173 -0.406
 0.521  0.794 -0.541  0.461  0.179
  1.43 -0.173  0.461  -1.43  0.822
-0.144 -0.406  0.179  0.822  -1.37

The orthogonal matrix Q is:
       1        0        0        0        0
       0   -0.471    0.127   -0.671   -0.558
       0    0.301   -0.195    0.437   -0.825
       0    0.825   0.0459   -0.563 -0.00872
       0  -0.0832   -0.971   -0.202   0.0922
The tridiagonal matrix T is:
  1.36   1.73      0      0      0
  1.73   -1.2 -0.966      0      0
     0 -0.966  -1.28  0.214      0
     0      0  0.214  -1.69  0.345
     0      0      0  0.345  0.164

Q * T * Q^T = 
  1.36 -0.816  0.521   1.43 -0.144
-0.816 -0.659  0.794 -0.173 -0.406
 0.521  0.794 -0.541  0.461  0.179
  1.43 -0.173  0.461  -1.43  0.822
-0.144 -0.406  0.179  0.822  -1.37
</pre> 
</div>
</div>
<h2 class="groupheader">Member Function Documentation</h2>
<a class="anchor" id="aa69e607a4aab4fb6321ca6acbf074fc2"></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_1Tridiagonalization.html">Tridiagonalization</a>&amp; compute </td>
          <td>(</td>
          <td class="paramtype">const <a class="el" href="classEigen_1_1Tridiagonalization.html#aeb6c0eb89cc982629305f6c7e0791caf">MatrixType</a> &amp;&#160;</td>
          <td class="paramname"><em>matrix</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>Computes tridiagonal decomposition of given matrix. </p>
<dl class="params"><dt>Parameters</dt><dd>
  <table class="params">
    <tr><td class="paramdir">[in]</td><td class="paramname">matrix</td><td>Selfadjoint matrix whose tridiagonal decomposition is to be computed. </td></tr>
  </table>
  </dd>
</dl>
<dl class="section return"><dt>Returns</dt><dd>Reference to <code>*this</code> </dd></dl>
<p>The tridiagonal decomposition is computed by bringing the columns of the matrix successively in the required form using Householder reflections. The cost is <img class="formulaInl" alt="$ 4n^3/3 $" src="form_95.png"/> flops, where <img class="formulaInl" alt="$ n $" src="form_45.png"/> denotes the size of the given matrix.</p>
<p>This method reuses of the allocated data in the <a class="el" href="classEigen_1_1Tridiagonalization.html" title="Tridiagonal decomposition of a selfadjoint matrix. ">Tridiagonalization</a> object, if the size of the matrix does not change.</p>
<p>Example: </p>
<div class="fragment"><div class="line">Tridiagonalization&lt;MatrixXf&gt; tri;</div>
<div class="line"><a class="code" href="group__matrixtypedefs.html#gabab09c32e96cfa9829a88400627af162">MatrixXf</a> X = <a class="code" href="classEigen_1_1DenseBase.html#a8e759dafdd9ecc446d397b7f5435f60a">MatrixXf::Random</a>(4,4);</div>
<div class="line"><a class="code" href="group__matrixtypedefs.html#gabab09c32e96cfa9829a88400627af162">MatrixXf</a> A = X + X.transpose();</div>
<div class="line">tri.compute(A);</div>
<div class="line">cout &lt;&lt; <span class="stringliteral">&quot;The matrix T in the tridiagonal decomposition of A is: &quot;</span> &lt;&lt; endl;</div>
<div class="line">cout &lt;&lt; tri.matrixT() &lt;&lt; endl;</div>
<div class="line">tri.compute(2*A); <span class="comment">// re-use tri to compute eigenvalues of 2A</span></div>
<div class="line">cout &lt;&lt; <span class="stringliteral">&quot;The matrix T in the tridiagonal decomposition of 2A is: &quot;</span> &lt;&lt; endl;</div>
<div class="line">cout &lt;&lt; tri.matrixT() &lt;&lt; endl;</div>
</div><!-- fragment --><p> Output: </p>
<pre class="fragment">The matrix T in the tridiagonal decomposition of A is: 
  1.36 -0.704      0      0
-0.704 0.0147   1.71      0
     0   1.71  0.856  0.641
     0      0  0.641 -0.506
The matrix T in the tridiagonal decomposition of 2A is: 
  2.72  -1.41      0      0
 -1.41 0.0294   3.43      0
     0   3.43   1.71   1.28
     0      0   1.28  -1.01
</pre> 
<p>References <a class="el" href="classEigen_1_1PlainObjectBase.html#afbbb33d14fe7fb9683019a39ce1c659d">PlainObjectBase&lt; Derived &gt;::resize()</a>.</p>

</div>
</div>
<a class="anchor" id="ac109eefddd733d8e82841da5bb2dd8d3"></a>
<div class="memitem">
<div class="memproto">
      <table class="memname">
        <tr>
          <td class="memname"><a class="el" href="classEigen_1_1Tridiagonalization.html">Tridiagonalization</a>&lt; <a class="el" href="classEigen_1_1Tridiagonalization.html#aeb6c0eb89cc982629305f6c7e0791caf">MatrixType</a> &gt;::DiagonalReturnType diagonal </td>
          <td>(</td>
          <td class="paramname"></td><td>)</td>
          <td> const</td>
        </tr>
      </table>
</div><div class="memdoc">

<p>Returns the diagonal of the tridiagonal matrix T in the decomposition. </p>
<dl class="section return"><dt>Returns</dt><dd>expression representing the diagonal of T</dd></dl>
<dl class="section pre"><dt>Precondition</dt><dd>Either the constructor <a class="el" href="classEigen_1_1Tridiagonalization.html#aa9f9722d2cef9425e2c0da3553dfbac7" title="Constructor; computes tridiagonal decomposition of given matrix. ">Tridiagonalization(const MatrixType&amp;)</a> or the member function <a class="el" href="classEigen_1_1Tridiagonalization.html#aa69e607a4aab4fb6321ca6acbf074fc2" title="Computes tridiagonal decomposition of given matrix. ">compute(const MatrixType&amp;)</a> has been called before to compute the tridiagonal decomposition of a matrix.</dd></dl>
<p>Example: </p>
<div class="fragment"><div class="line"><a class="code" href="group__matrixtypedefs.html#gaf1d15c8c24df228ee4869535dcbfa288">MatrixXcd</a> X = <a class="code" href="classEigen_1_1DenseBase.html#a8e759dafdd9ecc446d397b7f5435f60a">MatrixXcd::Random</a>(4,4);</div>
<div class="line"><a class="code" href="group__matrixtypedefs.html#gaf1d15c8c24df228ee4869535dcbfa288">MatrixXcd</a> A = X + X.adjoint();</div>
<div class="line">cout &lt;&lt; <span class="stringliteral">&quot;Here is a random self-adjoint 4x4 matrix:&quot;</span> &lt;&lt; endl &lt;&lt; A &lt;&lt; endl &lt;&lt; endl;</div>
<div class="line"></div>
<div class="line">Tridiagonalization&lt;MatrixXcd&gt; triOfA(A);</div>
<div class="line">MatrixXd T = triOfA.matrixT();</div>
<div class="line">cout &lt;&lt; <span class="stringliteral">&quot;The tridiagonal matrix T is:&quot;</span> &lt;&lt; endl &lt;&lt; T &lt;&lt; endl &lt;&lt; endl;</div>
<div class="line"></div>
<div class="line">cout &lt;&lt; <span class="stringliteral">&quot;We can also extract the diagonals of T directly ...&quot;</span> &lt;&lt; endl;</div>
<div class="line"><a class="code" href="group__matrixtypedefs.html#ga3da45e59796fbacf67fa568297927bd1">VectorXd</a> diag = triOfA.diagonal();</div>
<div class="line">cout &lt;&lt; <span class="stringliteral">&quot;The diagonal is:&quot;</span> &lt;&lt; endl &lt;&lt; diag &lt;&lt; endl; </div>
<div class="line"><a class="code" href="group__matrixtypedefs.html#ga3da45e59796fbacf67fa568297927bd1">VectorXd</a> subdiag = triOfA.subDiagonal();</div>
<div class="line">cout &lt;&lt; <span class="stringliteral">&quot;The subdiagonal is:&quot;</span> &lt;&lt; endl &lt;&lt; subdiag &lt;&lt; endl;</div>
</div><!-- fragment --><p> Output: </p>
<pre class="fragment">Here is a random self-adjoint 4x4 matrix:
    (-0.422,0)  (0.705,-1.01) (-0.17,-0.552) (0.338,-0.357)
  (0.705,1.01)      (0.515,0) (0.241,-0.446)   (0.05,-1.64)
 (-0.17,0.552)  (0.241,0.446)      (-1.03,0)  (0.0449,1.72)
 (0.338,0.357)    (0.05,1.64) (0.0449,-1.72)       (1.36,0)

The tridiagonal matrix T is:
-0.422 -1.45     0     0
-1.45  1.01 -1.42     0
    0 -1.42   1.8  -1.2
    0     0  -1.2 -1.96

We can also extract the diagonals of T directly ...
The diagonal is:
-0.422
1.01
1.8
-1.96
The subdiagonal is:
-1.45
-1.42
-1.2
</pre><dl class="section see"><dt>See Also</dt><dd><a class="el" href="classEigen_1_1Tridiagonalization.html#aceb0f16a166f4c236a1b536b7424d292" title="Returns an expression of the tridiagonal matrix T in the decomposition. ">matrixT()</a>, <a class="el" href="classEigen_1_1Tridiagonalization.html#a8fa49216273ab7579b7bea06debb1e51" title="Returns the subdiagonal of the tridiagonal matrix T in the decomposition. ">subDiagonal()</a> </dd></dl>

</div>
</div>
<a class="anchor" id="a2ab889c75460c178d941ee24e371b206"></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_1Matrix.html">CoeffVectorType</a> householderCoefficients </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 Householder coefficients. </p>
<dl class="section return"><dt>Returns</dt><dd>a const reference to the vector of Householder coefficients</dd></dl>
<dl class="section pre"><dt>Precondition</dt><dd>Either the constructor <a class="el" href="classEigen_1_1Tridiagonalization.html#aa9f9722d2cef9425e2c0da3553dfbac7" title="Constructor; computes tridiagonal decomposition of given matrix. ">Tridiagonalization(const MatrixType&amp;)</a> or the member function <a class="el" href="classEigen_1_1Tridiagonalization.html#aa69e607a4aab4fb6321ca6acbf074fc2" title="Computes tridiagonal decomposition of given matrix. ">compute(const MatrixType&amp;)</a> has been called before to compute the tridiagonal decomposition of a matrix.</dd></dl>
<p>The Householder coefficients allow the reconstruction of the matrix <img class="formulaInl" alt="$ Q $" src="form_73.png"/> in the tridiagonal decomposition from the packed data.</p>
<p>Example: </p>
<div class="fragment"><div class="line"><a class="code" href="group__matrixtypedefs.html#gacd860ff07358f6a703c2c0d4a174e920">Matrix4d</a> X = <a class="code" href="classEigen_1_1DenseBase.html#a8e759dafdd9ecc446d397b7f5435f60a">Matrix4d::Random</a>(4,4);</div>
<div class="line"><a class="code" href="group__matrixtypedefs.html#gacd860ff07358f6a703c2c0d4a174e920">Matrix4d</a> A = X + X.transpose();</div>
<div class="line">cout &lt;&lt; <span class="stringliteral">&quot;Here is a random symmetric 4x4 matrix:&quot;</span> &lt;&lt; endl &lt;&lt; A &lt;&lt; endl;</div>
<div class="line">Tridiagonalization&lt;Matrix4d&gt; triOfA(A);</div>
<div class="line"><a class="code" href="group__matrixtypedefs.html#ga2006332f6989f501762673e21f5128f5">Vector3d</a> hc = triOfA.householderCoefficients();</div>
<div class="line">cout &lt;&lt; <span class="stringliteral">&quot;The vector of Householder coefficients is:&quot;</span> &lt;&lt; endl &lt;&lt; hc &lt;&lt; endl;</div>
</div><!-- fragment --><p> Output: </p>
<pre class="fragment">Here is a random symmetric 4x4 matrix:
   1.36   0.612   0.122   0.326
  0.612   -1.21  -0.222   0.563
  0.122  -0.222 -0.0904    1.16
  0.326   0.563    1.16    1.66
The vector of Householder coefficients is:
1.87
1.24
0
</pre><dl class="section see"><dt>See Also</dt><dd><a class="el" href="classEigen_1_1Tridiagonalization.html#a66adece364b64b26b3771662de70f2df" title="Returns the internal representation of the decomposition. ">packedMatrix()</a>, <a class="el" href="group__Householder__Module.html">Householder module</a> </dd></dl>

</div>
</div>
<a class="anchor" id="ad13845d7490115664924b3dc208ec369"></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_1Tridiagonalization.html#aa96bdbc1b19c647e3372c31301ea4999">HouseholderSequenceType</a> matrixQ </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 unitary matrix Q in the decomposition. </p>
<dl class="section return"><dt>Returns</dt><dd>object representing the matrix Q</dd></dl>
<dl class="section pre"><dt>Precondition</dt><dd>Either the constructor <a class="el" href="classEigen_1_1Tridiagonalization.html#aa9f9722d2cef9425e2c0da3553dfbac7" title="Constructor; computes tridiagonal decomposition of given matrix. ">Tridiagonalization(const MatrixType&amp;)</a> or the member function <a class="el" href="classEigen_1_1Tridiagonalization.html#aa69e607a4aab4fb6321ca6acbf074fc2" title="Computes tridiagonal decomposition of given matrix. ">compute(const MatrixType&amp;)</a> has been called before to compute the tridiagonal decomposition of a matrix.</dd></dl>
<p>This function returns a light-weight object of template class <a class="el" href="classEigen_1_1HouseholderSequence.html" title="Sequence of Householder reflections acting on subspaces with decreasing size. ">HouseholderSequence</a>. You can either apply it directly to a matrix or you can convert it to a matrix of type <a class="el" href="classEigen_1_1Tridiagonalization.html#aeb6c0eb89cc982629305f6c7e0791caf" title="Synonym for the template parameter _MatrixType. ">MatrixType</a>.</p>
<dl class="section see"><dt>See Also</dt><dd><a class="el" href="classEigen_1_1Tridiagonalization.html#aa9f9722d2cef9425e2c0da3553dfbac7" title="Constructor; computes tridiagonal decomposition of given matrix. ">Tridiagonalization(const MatrixType&amp;)</a> for an example, <a class="el" href="classEigen_1_1Tridiagonalization.html#aceb0f16a166f4c236a1b536b7424d292" title="Returns an expression of the tridiagonal matrix T in the decomposition. ">matrixT()</a>, class <a class="el" href="classEigen_1_1HouseholderSequence.html" title="Sequence of Householder reflections acting on subspaces with decreasing size. ">HouseholderSequence</a> </dd></dl>

</div>
</div>
<a class="anchor" id="aceb0f16a166f4c236a1b536b7424d292"></a>
<div class="memitem">
<div class="memproto">
<table class="mlabels">
  <tr>
  <td class="mlabels-left">
      <table class="memname">
        <tr>
          <td class="memname">MatrixTReturnType matrixT </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 an expression of the tridiagonal matrix T in the decomposition. </p>
<dl class="section return"><dt>Returns</dt><dd>expression object representing the matrix T</dd></dl>
<dl class="section pre"><dt>Precondition</dt><dd>Either the constructor <a class="el" href="classEigen_1_1Tridiagonalization.html#aa9f9722d2cef9425e2c0da3553dfbac7" title="Constructor; computes tridiagonal decomposition of given matrix. ">Tridiagonalization(const MatrixType&amp;)</a> or the member function <a class="el" href="classEigen_1_1Tridiagonalization.html#aa69e607a4aab4fb6321ca6acbf074fc2" title="Computes tridiagonal decomposition of given matrix. ">compute(const MatrixType&amp;)</a> has been called before to compute the tridiagonal decomposition of a matrix.</dd></dl>
<p>Currently, this function can be used to extract the matrix T from internal data and copy it to a dense matrix object. In most cases, it may be sufficient to directly use the packed matrix or the vector expressions returned by <a class="el" href="classEigen_1_1Tridiagonalization.html#ac109eefddd733d8e82841da5bb2dd8d3" title="Returns the diagonal of the tridiagonal matrix T in the decomposition. ">diagonal()</a> and <a class="el" href="classEigen_1_1Tridiagonalization.html#a8fa49216273ab7579b7bea06debb1e51" title="Returns the subdiagonal of the tridiagonal matrix T in the decomposition. ">subDiagonal()</a> instead of creating a new dense copy matrix with this function.</p>
<dl class="section see"><dt>See Also</dt><dd><a class="el" href="classEigen_1_1Tridiagonalization.html#aa9f9722d2cef9425e2c0da3553dfbac7" title="Constructor; computes tridiagonal decomposition of given matrix. ">Tridiagonalization(const MatrixType&amp;)</a> for an example, <a class="el" href="classEigen_1_1Tridiagonalization.html#ad13845d7490115664924b3dc208ec369" title="Returns the unitary matrix Q in the decomposition. ">matrixQ()</a>, <a class="el" href="classEigen_1_1Tridiagonalization.html#a66adece364b64b26b3771662de70f2df" title="Returns the internal representation of the decomposition. ">packedMatrix()</a>, <a class="el" href="classEigen_1_1Tridiagonalization.html#ac109eefddd733d8e82841da5bb2dd8d3" title="Returns the diagonal of the tridiagonal matrix T in the decomposition. ">diagonal()</a>, <a class="el" href="classEigen_1_1Tridiagonalization.html#a8fa49216273ab7579b7bea06debb1e51" title="Returns the subdiagonal of the tridiagonal matrix T in the decomposition. ">subDiagonal()</a> </dd></dl>

</div>
</div>
<a class="anchor" id="a66adece364b64b26b3771662de70f2df"></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_1Tridiagonalization.html#aeb6c0eb89cc982629305f6c7e0791caf">MatrixType</a>&amp; packedMatrix </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 internal representation of the decomposition. </p>
<dl class="section return"><dt>Returns</dt><dd>a const reference to a matrix with the internal representation of the decomposition.</dd></dl>
<dl class="section pre"><dt>Precondition</dt><dd>Either the constructor <a class="el" href="classEigen_1_1Tridiagonalization.html#aa9f9722d2cef9425e2c0da3553dfbac7" title="Constructor; computes tridiagonal decomposition of given matrix. ">Tridiagonalization(const MatrixType&amp;)</a> or the member function <a class="el" href="classEigen_1_1Tridiagonalization.html#aa69e607a4aab4fb6321ca6acbf074fc2" title="Computes tridiagonal decomposition of given matrix. ">compute(const MatrixType&amp;)</a> has been called before to compute the tridiagonal decomposition of a matrix.</dd></dl>
<p>The returned matrix contains the following information:</p>
<ul>
<li>the strict upper triangular part is equal to the input matrix A.</li>
<li>the diagonal and lower sub-diagonal represent the real tridiagonal symmetric matrix T.</li>
<li>the rest of the lower part contains the Householder vectors that, combined with Householder coefficients returned by <a class="el" href="classEigen_1_1Tridiagonalization.html#a2ab889c75460c178d941ee24e371b206" title="Returns the Householder coefficients. ">householderCoefficients()</a>, allows to reconstruct the matrix Q as <img class="formulaInl" alt="$ Q = H_{N-1} \ldots H_1 H_0 $" src="form_80.png"/>. Here, the matrices <img class="formulaInl" alt="$ H_i $" src="form_81.png"/> are the Householder transformations <img class="formulaInl" alt="$ H_i = (I - h_i v_i v_i^T) $" src="form_82.png"/> where <img class="formulaInl" alt="$ h_i $" src="form_83.png"/> is the <img class="formulaInl" alt="$ i $" src="form_84.png"/>th Householder coefficient and <img class="formulaInl" alt="$ v_i $" src="form_85.png"/> is the Householder vector defined by <img class="formulaInl" alt="$ v_i = [ 0, \ldots, 0, 1, M(i+2,i), \ldots, M(N-1,i) ]^T $" src="form_86.png"/> with M the matrix returned by this function.</li>
</ul>
<p>See LAPACK for further details on this packed storage.</p>
<p>Example: </p>
<div class="fragment"><div class="line"><a class="code" href="group__matrixtypedefs.html#gacd860ff07358f6a703c2c0d4a174e920">Matrix4d</a> X = <a class="code" href="classEigen_1_1DenseBase.html#a8e759dafdd9ecc446d397b7f5435f60a">Matrix4d::Random</a>(4,4);</div>
<div class="line"><a class="code" href="group__matrixtypedefs.html#gacd860ff07358f6a703c2c0d4a174e920">Matrix4d</a> A = X + X.transpose();</div>
<div class="line">cout &lt;&lt; <span class="stringliteral">&quot;Here is a random symmetric 4x4 matrix:&quot;</span> &lt;&lt; endl &lt;&lt; A &lt;&lt; endl;</div>
<div class="line">Tridiagonalization&lt;Matrix4d&gt; triOfA(A);</div>
<div class="line"><a class="code" href="group__matrixtypedefs.html#gacd860ff07358f6a703c2c0d4a174e920">Matrix4d</a> pm = triOfA.packedMatrix();</div>
<div class="line">cout &lt;&lt; <span class="stringliteral">&quot;The packed matrix M is:&quot;</span> &lt;&lt; endl &lt;&lt; pm &lt;&lt; endl;</div>
<div class="line">cout &lt;&lt; <span class="stringliteral">&quot;The diagonal and subdiagonal corresponds to the matrix T, which is:&quot;</span> </div>
<div class="line">     &lt;&lt; endl &lt;&lt; triOfA.matrixT() &lt;&lt; endl;</div>
</div><!-- fragment --><p> Output: </p>
<pre class="fragment">Here is a random symmetric 4x4 matrix:
   1.36   0.612   0.122   0.326
  0.612   -1.21  -0.222   0.563
  0.122  -0.222 -0.0904    1.16
  0.326   0.563    1.16    1.66
The packed matrix M is:
  1.36  0.612  0.122  0.326
-0.704 0.0147 -0.222  0.563
0.0925   1.71  0.856   1.16
 0.248  0.785  0.641 -0.506
The diagonal and subdiagonal corresponds to the matrix T, which is:
  1.36 -0.704      0      0
-0.704 0.0147   1.71      0
     0   1.71  0.856  0.641
     0      0  0.641 -0.506
</pre><dl class="section see"><dt>See Also</dt><dd><a class="el" href="classEigen_1_1Tridiagonalization.html#a2ab889c75460c178d941ee24e371b206" title="Returns the Householder coefficients. ">householderCoefficients()</a> </dd></dl>

</div>
</div>
<a class="anchor" id="a8fa49216273ab7579b7bea06debb1e51"></a>
<div class="memitem">
<div class="memproto">
      <table class="memname">
        <tr>
          <td class="memname"><a class="el" href="classEigen_1_1Tridiagonalization.html">Tridiagonalization</a>&lt; <a class="el" href="classEigen_1_1Tridiagonalization.html#aeb6c0eb89cc982629305f6c7e0791caf">MatrixType</a> &gt;::SubDiagonalReturnType subDiagonal </td>
          <td>(</td>
          <td class="paramname"></td><td>)</td>
          <td> const</td>
        </tr>
      </table>
</div><div class="memdoc">

<p>Returns the subdiagonal of the tridiagonal matrix T in the decomposition. </p>
<dl class="section return"><dt>Returns</dt><dd>expression representing the subdiagonal of T</dd></dl>
<dl class="section pre"><dt>Precondition</dt><dd>Either the constructor <a class="el" href="classEigen_1_1Tridiagonalization.html#aa9f9722d2cef9425e2c0da3553dfbac7" title="Constructor; computes tridiagonal decomposition of given matrix. ">Tridiagonalization(const MatrixType&amp;)</a> or the member function <a class="el" href="classEigen_1_1Tridiagonalization.html#aa69e607a4aab4fb6321ca6acbf074fc2" title="Computes tridiagonal decomposition of given matrix. ">compute(const MatrixType&amp;)</a> has been called before to compute the tridiagonal decomposition of a matrix.</dd></dl>
<dl class="section see"><dt>See Also</dt><dd><a class="el" href="classEigen_1_1Tridiagonalization.html#ac109eefddd733d8e82841da5bb2dd8d3" title="Returns the diagonal of the tridiagonal matrix T in the decomposition. ">diagonal()</a> for an example, <a class="el" href="classEigen_1_1Tridiagonalization.html#aceb0f16a166f4c236a1b536b7424d292" title="Returns an expression of the tridiagonal matrix T in the decomposition. ">matrixT()</a> </dd></dl>

</div>
</div>
<hr/>The documentation for this class was generated from the following file:<ul>
<li><a class="el" href="Tridiagonalization_8h_source.html">Tridiagonalization.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_1Tridiagonalization.html">Tridiagonalization</a></li>
    <li class="footer">Generated on Mon Oct 28 2013 11:04:31 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>