Sophie

Sophie

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

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.h Source File</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('Tridiagonalization_8h_source.html','');});
</script>
<div id="doc-content">
<!-- window showing the filter options -->
<div id="MSearchSelectWindow"
     onmouseover="return searchBox.OnSearchSelectShow()"
     onmouseout="return searchBox.OnSearchSelectHide()"
     onkeydown="return searchBox.OnSearchSelectKey(event)">
<a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(0)"><span class="SelectionMark">&#160;</span>All</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(1)"><span class="SelectionMark">&#160;</span>Classes</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(2)"><span class="SelectionMark">&#160;</span>Namespaces</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(3)"><span class="SelectionMark">&#160;</span>Functions</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(4)"><span class="SelectionMark">&#160;</span>Variables</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(5)"><span class="SelectionMark">&#160;</span>Typedefs</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(6)"><span class="SelectionMark">&#160;</span>Enumerations</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(7)"><span class="SelectionMark">&#160;</span>Enumerator</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(8)"><span class="SelectionMark">&#160;</span>Friends</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(9)"><span class="SelectionMark">&#160;</span>Groups</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(10)"><span class="SelectionMark">&#160;</span>Pages</a></div>

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

<div class="header">
  <div class="headertitle">
<div class="title">Tridiagonalization.h</div>  </div>
</div><!--header-->
<div class="contents">
<div class="fragment"><div class="line"><a name="l00001"></a><span class="lineno">    1</span>&#160;<span class="comment">// This file is part of Eigen, a lightweight C++ template library</span></div>
<div class="line"><a name="l00002"></a><span class="lineno">    2</span>&#160;<span class="comment">// for linear algebra.</span></div>
<div class="line"><a name="l00003"></a><span class="lineno">    3</span>&#160;<span class="comment">//</span></div>
<div class="line"><a name="l00004"></a><span class="lineno">    4</span>&#160;<span class="comment">// Copyright (C) 2008 Gael Guennebaud &lt;gael.guennebaud@inria.fr&gt;</span></div>
<div class="line"><a name="l00005"></a><span class="lineno">    5</span>&#160;<span class="comment">// Copyright (C) 2010 Jitse Niesen &lt;jitse@maths.leeds.ac.uk&gt;</span></div>
<div class="line"><a name="l00006"></a><span class="lineno">    6</span>&#160;<span class="comment">//</span></div>
<div class="line"><a name="l00007"></a><span class="lineno">    7</span>&#160;<span class="comment">// This Source Code Form is subject to the terms of the Mozilla</span></div>
<div class="line"><a name="l00008"></a><span class="lineno">    8</span>&#160;<span class="comment">// Public License v. 2.0. If a copy of the MPL was not distributed</span></div>
<div class="line"><a name="l00009"></a><span class="lineno">    9</span>&#160;<span class="comment">// with this file, You can obtain one at http://mozilla.org/MPL/2.0/.</span></div>
<div class="line"><a name="l00010"></a><span class="lineno">   10</span>&#160;</div>
<div class="line"><a name="l00011"></a><span class="lineno">   11</span>&#160;<span class="preprocessor">#ifndef EIGEN_TRIDIAGONALIZATION_H</span></div>
<div class="line"><a name="l00012"></a><span class="lineno">   12</span>&#160;<span class="preprocessor"></span><span class="preprocessor">#define EIGEN_TRIDIAGONALIZATION_H</span></div>
<div class="line"><a name="l00013"></a><span class="lineno">   13</span>&#160;<span class="preprocessor"></span></div>
<div class="line"><a name="l00014"></a><span class="lineno">   14</span>&#160;<span class="keyword">namespace </span>Eigen { </div>
<div class="line"><a name="l00015"></a><span class="lineno">   15</span>&#160;</div>
<div class="line"><a name="l00016"></a><span class="lineno">   16</span>&#160;<span class="keyword">namespace </span>internal {</div>
<div class="line"><a name="l00017"></a><span class="lineno">   17</span>&#160;  </div>
<div class="line"><a name="l00018"></a><span class="lineno">   18</span>&#160;<span class="keyword">template</span>&lt;<span class="keyword">typename</span> MatrixType&gt; <span class="keyword">struct </span>TridiagonalizationMatrixTReturnType;</div>
<div class="line"><a name="l00019"></a><span class="lineno">   19</span>&#160;<span class="keyword">template</span>&lt;<span class="keyword">typename</span> MatrixType&gt;</div>
<div class="line"><a name="l00020"></a><span class="lineno">   20</span>&#160;<span class="keyword">struct </span>traits&lt;TridiagonalizationMatrixTReturnType&lt;MatrixType&gt; &gt;</div>
<div class="line"><a name="l00021"></a><span class="lineno">   21</span>&#160;{</div>
<div class="line"><a name="l00022"></a><span class="lineno">   22</span>&#160;  <span class="keyword">typedef</span> <span class="keyword">typename</span> MatrixType::PlainObject ReturnType;</div>
<div class="line"><a name="l00023"></a><span class="lineno">   23</span>&#160;};</div>
<div class="line"><a name="l00024"></a><span class="lineno">   24</span>&#160;</div>
<div class="line"><a name="l00025"></a><span class="lineno">   25</span>&#160;<span class="keyword">template</span>&lt;<span class="keyword">typename</span> MatrixType, <span class="keyword">typename</span> CoeffVectorType&gt;</div>
<div class="line"><a name="l00026"></a><span class="lineno">   26</span>&#160;<span class="keywordtype">void</span> tridiagonalization_inplace(MatrixType&amp; matA, CoeffVectorType&amp; hCoeffs);</div>
<div class="line"><a name="l00027"></a><span class="lineno">   27</span>&#160;}</div>
<div class="line"><a name="l00028"></a><span class="lineno">   28</span>&#160;</div>
<div class="line"><a name="l00061"></a><span class="lineno"><a class="line" href="classEigen_1_1Tridiagonalization.html">   61</a></span>&#160;<span class="keyword">template</span>&lt;<span class="keyword">typename</span> _MatrixType&gt; <span class="keyword">class </span><a class="code" href="classEigen_1_1Tridiagonalization.html">Tridiagonalization</a></div>
<div class="line"><a name="l00062"></a><span class="lineno">   62</span>&#160;{</div>
<div class="line"><a name="l00063"></a><span class="lineno">   63</span>&#160;  <span class="keyword">public</span>:</div>
<div class="line"><a name="l00064"></a><span class="lineno">   64</span>&#160;</div>
<div class="line"><a name="l00066"></a><span class="lineno"><a class="line" href="classEigen_1_1Tridiagonalization.html#aeb6c0eb89cc982629305f6c7e0791caf">   66</a></span>&#160;    <span class="keyword">typedef</span> _MatrixType <a class="code" href="classEigen_1_1Tridiagonalization.html#aeb6c0eb89cc982629305f6c7e0791caf">MatrixType</a>;</div>
<div class="line"><a name="l00067"></a><span class="lineno">   67</span>&#160;</div>
<div class="line"><a name="l00068"></a><span class="lineno">   68</span>&#160;    <span class="keyword">typedef</span> <span class="keyword">typename</span> MatrixType::Scalar Scalar;</div>
<div class="line"><a name="l00069"></a><span class="lineno">   69</span>&#160;    <span class="keyword">typedef</span> <span class="keyword">typename</span> <a class="code" href="structEigen_1_1NumTraits.html">NumTraits&lt;Scalar&gt;::Real</a> RealScalar;</div>
<div class="line"><a name="l00070"></a><span class="lineno">   70</span>&#160;    <span class="keyword">typedef</span> <span class="keyword">typename</span> MatrixType::Index Index;</div>
<div class="line"><a name="l00071"></a><span class="lineno">   71</span>&#160;</div>
<div class="line"><a name="l00072"></a><span class="lineno">   72</span>&#160;    <span class="keyword">enum</span> {</div>
<div class="line"><a name="l00073"></a><span class="lineno">   73</span>&#160;      Size = MatrixType::RowsAtCompileTime,</div>
<div class="line"><a name="l00074"></a><span class="lineno">   74</span>&#160;      SizeMinusOne = Size == <a class="code" href="namespaceEigen.html#adc9da5be31bdce40c25a92c27999c0e3">Dynamic</a> ? <a class="code" href="namespaceEigen.html#adc9da5be31bdce40c25a92c27999c0e3">Dynamic</a> : (Size &gt; 1 ? Size - 1 : 1),</div>
<div class="line"><a name="l00075"></a><span class="lineno">   75</span>&#160;      Options = MatrixType::Options,</div>
<div class="line"><a name="l00076"></a><span class="lineno">   76</span>&#160;      MaxSize = MatrixType::MaxRowsAtCompileTime,</div>
<div class="line"><a name="l00077"></a><span class="lineno">   77</span>&#160;      MaxSizeMinusOne = MaxSize == <a class="code" href="namespaceEigen.html#adc9da5be31bdce40c25a92c27999c0e3">Dynamic</a> ? <a class="code" href="namespaceEigen.html#adc9da5be31bdce40c25a92c27999c0e3">Dynamic</a> : (MaxSize &gt; 1 ? MaxSize - 1 : 1)</div>
<div class="line"><a name="l00078"></a><span class="lineno">   78</span>&#160;    };</div>
<div class="line"><a name="l00079"></a><span class="lineno">   79</span>&#160;</div>
<div class="line"><a name="l00080"></a><span class="lineno">   80</span>&#160;    <span class="keyword">typedef</span> <a class="code" href="classEigen_1_1Matrix.html">Matrix&lt;Scalar, SizeMinusOne, 1, Options &amp; ~RowMajor, MaxSizeMinusOne, 1&gt;</a> CoeffVectorType;</div>
<div class="line"><a name="l00081"></a><span class="lineno">   81</span>&#160;    <span class="keyword">typedef</span> <span class="keyword">typename</span> internal::plain_col_type&lt;MatrixType, RealScalar&gt;::type DiagonalType;</div>
<div class="line"><a name="l00082"></a><span class="lineno">   82</span>&#160;    <span class="keyword">typedef</span> <a class="code" href="classEigen_1_1Matrix.html">Matrix&lt;RealScalar, SizeMinusOne, 1, Options &amp; ~RowMajor, MaxSizeMinusOne, 1&gt;</a> SubDiagonalType;</div>
<div class="line"><a name="l00083"></a><span class="lineno">   83</span>&#160;    <span class="keyword">typedef</span> <span class="keyword">typename</span> internal::remove_all&lt;typename MatrixType::RealReturnType&gt;::type MatrixTypeRealView;</div>
<div class="line"><a name="l00084"></a><span class="lineno">   84</span>&#160;    <span class="keyword">typedef</span> internal::TridiagonalizationMatrixTReturnType&lt;MatrixTypeRealView&gt; MatrixTReturnType;</div>
<div class="line"><a name="l00085"></a><span class="lineno">   85</span>&#160;</div>
<div class="line"><a name="l00086"></a><span class="lineno">   86</span>&#160;    <span class="keyword">typedef</span> <span class="keyword">typename</span> internal::conditional&lt;NumTraits&lt;Scalar&gt;::IsComplex,</div>
<div class="line"><a name="l00087"></a><span class="lineno">   87</span>&#160;              <span class="keyword">typename</span> internal::add_const_on_value_type&lt;typename Diagonal&lt;const MatrixType&gt;::RealReturnType&gt;::type,</div>
<div class="line"><a name="l00088"></a><span class="lineno">   88</span>&#160;              <span class="keyword">const</span> <a class="code" href="classEigen_1_1Diagonal.html">Diagonal&lt;const MatrixType&gt;</a></div>
<div class="line"><a name="l00089"></a><span class="lineno">   89</span>&#160;            &gt;::type DiagonalReturnType;</div>
<div class="line"><a name="l00090"></a><span class="lineno">   90</span>&#160;</div>
<div class="line"><a name="l00091"></a><span class="lineno">   91</span>&#160;    <span class="keyword">typedef</span> <span class="keyword">typename</span> internal::conditional&lt;NumTraits&lt;Scalar&gt;::IsComplex,</div>
<div class="line"><a name="l00092"></a><span class="lineno">   92</span>&#160;              <span class="keyword">typename</span> internal::add_const_on_value_type&lt;<span class="keyword">typename</span> <a class="code" href="classEigen_1_1Diagonal.html">Diagonal</a>&lt;</div>
<div class="line"><a name="l00093"></a><span class="lineno">   93</span>&#160;                <a class="code" href="classEigen_1_1Block.html">Block&lt;const MatrixType,SizeMinusOne,SizeMinusOne&gt;</a> &gt;::RealReturnType&gt;::type,</div>
<div class="line"><a name="l00094"></a><span class="lineno">   94</span>&#160;              <span class="keyword">const</span> <a class="code" href="classEigen_1_1Diagonal.html">Diagonal</a>&lt;</div>
<div class="line"><a name="l00095"></a><span class="lineno">   95</span>&#160;                <a class="code" href="classEigen_1_1Block.html">Block&lt;const MatrixType,SizeMinusOne,SizeMinusOne&gt;</a> &gt;</div>
<div class="line"><a name="l00096"></a><span class="lineno">   96</span>&#160;            &gt;::type SubDiagonalReturnType;</div>
<div class="line"><a name="l00097"></a><span class="lineno">   97</span>&#160;</div>
<div class="line"><a name="l00099"></a><span class="lineno"><a class="line" href="classEigen_1_1Tridiagonalization.html#aa96bdbc1b19c647e3372c31301ea4999">   99</a></span>&#160;    <span class="keyword">typedef</span> <a class="code" href="classEigen_1_1HouseholderSequence.html">HouseholderSequence&lt;MatrixType,typename internal::remove_all&lt;typename CoeffVectorType::ConjugateReturnType&gt;::type</a>&gt; <a class="code" href="classEigen_1_1Tridiagonalization.html#aa96bdbc1b19c647e3372c31301ea4999">HouseholderSequenceType</a>;</div>
<div class="line"><a name="l00100"></a><span class="lineno">  100</span>&#160;</div>
<div class="line"><a name="l00113"></a><span class="lineno"><a class="line" href="classEigen_1_1Tridiagonalization.html#a0698ae78b0ab6f239c475b73b9c6bbee">  113</a></span>&#160;    <a class="code" href="classEigen_1_1Tridiagonalization.html#a0698ae78b0ab6f239c475b73b9c6bbee">Tridiagonalization</a>(Index size = Size==<a class="code" href="namespaceEigen.html#adc9da5be31bdce40c25a92c27999c0e3">Dynamic</a> ? 2 : Size)</div>
<div class="line"><a name="l00114"></a><span class="lineno">  114</span>&#160;      : m_matrix(size,size),</div>
<div class="line"><a name="l00115"></a><span class="lineno">  115</span>&#160;        m_hCoeffs(size &gt; 1 ? size-1 : 1),</div>
<div class="line"><a name="l00116"></a><span class="lineno">  116</span>&#160;        m_isInitialized(false)</div>
<div class="line"><a name="l00117"></a><span class="lineno">  117</span>&#160;    {}</div>
<div class="line"><a name="l00118"></a><span class="lineno">  118</span>&#160;</div>
<div class="line"><a name="l00129"></a><span class="lineno"><a class="line" href="classEigen_1_1Tridiagonalization.html#aa9f9722d2cef9425e2c0da3553dfbac7">  129</a></span>&#160;    <a class="code" href="classEigen_1_1Tridiagonalization.html#aa9f9722d2cef9425e2c0da3553dfbac7">Tridiagonalization</a>(<span class="keyword">const</span> <a class="code" href="classEigen_1_1Tridiagonalization.html#aeb6c0eb89cc982629305f6c7e0791caf">MatrixType</a>&amp; matrix)</div>
<div class="line"><a name="l00130"></a><span class="lineno">  130</span>&#160;      : m_matrix(matrix),</div>
<div class="line"><a name="l00131"></a><span class="lineno">  131</span>&#160;        m_hCoeffs(matrix.cols() &gt; 1 ? matrix.cols()-1 : 1),</div>
<div class="line"><a name="l00132"></a><span class="lineno">  132</span>&#160;        m_isInitialized(false)</div>
<div class="line"><a name="l00133"></a><span class="lineno">  133</span>&#160;    {</div>
<div class="line"><a name="l00134"></a><span class="lineno">  134</span>&#160;      internal::tridiagonalization_inplace(m_matrix, m_hCoeffs);</div>
<div class="line"><a name="l00135"></a><span class="lineno">  135</span>&#160;      m_isInitialized = <span class="keyword">true</span>;</div>
<div class="line"><a name="l00136"></a><span class="lineno">  136</span>&#160;    }</div>
<div class="line"><a name="l00137"></a><span class="lineno">  137</span>&#160;</div>
<div class="line"><a name="l00155"></a><span class="lineno"><a class="line" href="classEigen_1_1Tridiagonalization.html#aa69e607a4aab4fb6321ca6acbf074fc2">  155</a></span>&#160;    <a class="code" href="classEigen_1_1Tridiagonalization.html">Tridiagonalization</a>&amp; <a class="code" href="classEigen_1_1Tridiagonalization.html#aa69e607a4aab4fb6321ca6acbf074fc2">compute</a>(<span class="keyword">const</span> <a class="code" href="classEigen_1_1Tridiagonalization.html#aeb6c0eb89cc982629305f6c7e0791caf">MatrixType</a>&amp; matrix)</div>
<div class="line"><a name="l00156"></a><span class="lineno">  156</span>&#160;    {</div>
<div class="line"><a name="l00157"></a><span class="lineno">  157</span>&#160;      m_matrix = matrix;</div>
<div class="line"><a name="l00158"></a><span class="lineno">  158</span>&#160;      m_hCoeffs.<a class="code" href="classEigen_1_1PlainObjectBase.html#afbbb33d14fe7fb9683019a39ce1c659d">resize</a>(matrix.rows()-1, 1);</div>
<div class="line"><a name="l00159"></a><span class="lineno">  159</span>&#160;      internal::tridiagonalization_inplace(m_matrix, m_hCoeffs);</div>
<div class="line"><a name="l00160"></a><span class="lineno">  160</span>&#160;      m_isInitialized = <span class="keyword">true</span>;</div>
<div class="line"><a name="l00161"></a><span class="lineno">  161</span>&#160;      <span class="keywordflow">return</span> *<span class="keyword">this</span>;</div>
<div class="line"><a name="l00162"></a><span class="lineno">  162</span>&#160;    }</div>
<div class="line"><a name="l00163"></a><span class="lineno">  163</span>&#160;</div>
<div class="line"><a name="l00180"></a><span class="lineno"><a class="line" href="classEigen_1_1Tridiagonalization.html#a2ab889c75460c178d941ee24e371b206">  180</a></span>&#160;    <span class="keyword">inline</span> <a class="code" href="classEigen_1_1Matrix.html">CoeffVectorType</a> <a class="code" href="classEigen_1_1Tridiagonalization.html#a2ab889c75460c178d941ee24e371b206">householderCoefficients</a>()<span class="keyword"> const</span></div>
<div class="line"><a name="l00181"></a><span class="lineno">  181</span>&#160;<span class="keyword">    </span>{</div>
<div class="line"><a name="l00182"></a><span class="lineno">  182</span>&#160;      eigen_assert(m_isInitialized &amp;&amp; <span class="stringliteral">&quot;Tridiagonalization is not initialized.&quot;</span>);</div>
<div class="line"><a name="l00183"></a><span class="lineno">  183</span>&#160;      <span class="keywordflow">return</span> m_hCoeffs;</div>
<div class="line"><a name="l00184"></a><span class="lineno">  184</span>&#160;    }</div>
<div class="line"><a name="l00185"></a><span class="lineno">  185</span>&#160;</div>
<div class="line"><a name="l00217"></a><span class="lineno"><a class="line" href="classEigen_1_1Tridiagonalization.html#a66adece364b64b26b3771662de70f2df">  217</a></span>&#160;    <span class="keyword">inline</span> <span class="keyword">const</span> <a class="code" href="classEigen_1_1Tridiagonalization.html#aeb6c0eb89cc982629305f6c7e0791caf">MatrixType</a>&amp; <a class="code" href="classEigen_1_1Tridiagonalization.html#a66adece364b64b26b3771662de70f2df">packedMatrix</a>()<span class="keyword"> const</span></div>
<div class="line"><a name="l00218"></a><span class="lineno">  218</span>&#160;<span class="keyword">    </span>{</div>
<div class="line"><a name="l00219"></a><span class="lineno">  219</span>&#160;      eigen_assert(m_isInitialized &amp;&amp; <span class="stringliteral">&quot;Tridiagonalization is not initialized.&quot;</span>);</div>
<div class="line"><a name="l00220"></a><span class="lineno">  220</span>&#160;      <span class="keywordflow">return</span> m_matrix;</div>
<div class="line"><a name="l00221"></a><span class="lineno">  221</span>&#160;    }</div>
<div class="line"><a name="l00222"></a><span class="lineno">  222</span>&#160;</div>
<div class="line"><a name="l00238"></a><span class="lineno"><a class="line" href="classEigen_1_1Tridiagonalization.html#ad13845d7490115664924b3dc208ec369">  238</a></span>&#160;    <a class="code" href="classEigen_1_1HouseholderSequence.html">HouseholderSequenceType</a> <a class="code" href="classEigen_1_1Tridiagonalization.html#ad13845d7490115664924b3dc208ec369">matrixQ</a>()<span class="keyword"> const</span></div>
<div class="line"><a name="l00239"></a><span class="lineno">  239</span>&#160;<span class="keyword">    </span>{</div>
<div class="line"><a name="l00240"></a><span class="lineno">  240</span>&#160;      eigen_assert(m_isInitialized &amp;&amp; <span class="stringliteral">&quot;Tridiagonalization is not initialized.&quot;</span>);</div>
<div class="line"><a name="l00241"></a><span class="lineno">  241</span>&#160;      <span class="keywordflow">return</span> <a class="code" href="classEigen_1_1Tridiagonalization.html#aa96bdbc1b19c647e3372c31301ea4999">HouseholderSequenceType</a>(m_matrix, m_hCoeffs.conjugate())</div>
<div class="line"><a name="l00242"></a><span class="lineno">  242</span>&#160;             .setLength(m_matrix.rows() - 1)</div>
<div class="line"><a name="l00243"></a><span class="lineno">  243</span>&#160;             .setShift(1);</div>
<div class="line"><a name="l00244"></a><span class="lineno">  244</span>&#160;    }</div>
<div class="line"><a name="l00245"></a><span class="lineno">  245</span>&#160;</div>
<div class="line"><a name="l00263"></a><span class="lineno"><a class="line" href="classEigen_1_1Tridiagonalization.html#aceb0f16a166f4c236a1b536b7424d292">  263</a></span>&#160;    MatrixTReturnType <a class="code" href="classEigen_1_1Tridiagonalization.html#aceb0f16a166f4c236a1b536b7424d292">matrixT</a>()<span class="keyword"> const</span></div>
<div class="line"><a name="l00264"></a><span class="lineno">  264</span>&#160;<span class="keyword">    </span>{</div>
<div class="line"><a name="l00265"></a><span class="lineno">  265</span>&#160;      eigen_assert(m_isInitialized &amp;&amp; <span class="stringliteral">&quot;Tridiagonalization is not initialized.&quot;</span>);</div>
<div class="line"><a name="l00266"></a><span class="lineno">  266</span>&#160;      <span class="keywordflow">return</span> MatrixTReturnType(m_matrix.real());</div>
<div class="line"><a name="l00267"></a><span class="lineno">  267</span>&#160;    }</div>
<div class="line"><a name="l00268"></a><span class="lineno">  268</span>&#160;</div>
<div class="line"><a name="l00282"></a><span class="lineno">  282</span>&#160;    DiagonalReturnType <a class="code" href="classEigen_1_1Tridiagonalization.html#ac109eefddd733d8e82841da5bb2dd8d3">diagonal</a>() <span class="keyword">const</span>;</div>
<div class="line"><a name="l00283"></a><span class="lineno">  283</span>&#160;</div>
<div class="line"><a name="l00294"></a><span class="lineno">  294</span>&#160;    SubDiagonalReturnType <a class="code" href="classEigen_1_1Tridiagonalization.html#a8fa49216273ab7579b7bea06debb1e51">subDiagonal</a>() <span class="keyword">const</span>;</div>
<div class="line"><a name="l00295"></a><span class="lineno">  295</span>&#160;</div>
<div class="line"><a name="l00296"></a><span class="lineno">  296</span>&#160;  <span class="keyword">protected</span>:</div>
<div class="line"><a name="l00297"></a><span class="lineno">  297</span>&#160;</div>
<div class="line"><a name="l00298"></a><span class="lineno">  298</span>&#160;    <a class="code" href="classEigen_1_1Tridiagonalization.html#aeb6c0eb89cc982629305f6c7e0791caf">MatrixType</a> m_matrix;</div>
<div class="line"><a name="l00299"></a><span class="lineno">  299</span>&#160;    CoeffVectorType m_hCoeffs;</div>
<div class="line"><a name="l00300"></a><span class="lineno">  300</span>&#160;    <span class="keywordtype">bool</span> m_isInitialized;</div>
<div class="line"><a name="l00301"></a><span class="lineno">  301</span>&#160;};</div>
<div class="line"><a name="l00302"></a><span class="lineno">  302</span>&#160;</div>
<div class="line"><a name="l00303"></a><span class="lineno">  303</span>&#160;<span class="keyword">template</span>&lt;<span class="keyword">typename</span> MatrixType&gt;</div>
<div class="line"><a name="l00304"></a><span class="lineno">  304</span>&#160;<span class="keyword">typename</span> Tridiagonalization&lt;MatrixType&gt;::DiagonalReturnType</div>
<div class="line"><a name="l00305"></a><span class="lineno"><a class="line" href="classEigen_1_1Tridiagonalization.html#ac109eefddd733d8e82841da5bb2dd8d3">  305</a></span>&#160;<a class="code" href="classEigen_1_1Tridiagonalization.html#ac109eefddd733d8e82841da5bb2dd8d3">Tridiagonalization&lt;MatrixType&gt;::diagonal</a>()<span class="keyword"> const</span></div>
<div class="line"><a name="l00306"></a><span class="lineno">  306</span>&#160;<span class="keyword"></span>{</div>
<div class="line"><a name="l00307"></a><span class="lineno">  307</span>&#160;  eigen_assert(m_isInitialized &amp;&amp; <span class="stringliteral">&quot;Tridiagonalization is not initialized.&quot;</span>);</div>
<div class="line"><a name="l00308"></a><span class="lineno">  308</span>&#160;  <span class="keywordflow">return</span> m_matrix.diagonal();</div>
<div class="line"><a name="l00309"></a><span class="lineno">  309</span>&#160;}</div>
<div class="line"><a name="l00310"></a><span class="lineno">  310</span>&#160;</div>
<div class="line"><a name="l00311"></a><span class="lineno">  311</span>&#160;<span class="keyword">template</span>&lt;<span class="keyword">typename</span> MatrixType&gt;</div>
<div class="line"><a name="l00312"></a><span class="lineno">  312</span>&#160;<span class="keyword">typename</span> Tridiagonalization&lt;MatrixType&gt;::SubDiagonalReturnType</div>
<div class="line"><a name="l00313"></a><span class="lineno"><a class="line" href="classEigen_1_1Tridiagonalization.html#a8fa49216273ab7579b7bea06debb1e51">  313</a></span>&#160;<a class="code" href="classEigen_1_1Tridiagonalization.html#a8fa49216273ab7579b7bea06debb1e51">Tridiagonalization&lt;MatrixType&gt;::subDiagonal</a>()<span class="keyword"> const</span></div>
<div class="line"><a name="l00314"></a><span class="lineno">  314</span>&#160;<span class="keyword"></span>{</div>
<div class="line"><a name="l00315"></a><span class="lineno">  315</span>&#160;  eigen_assert(m_isInitialized &amp;&amp; <span class="stringliteral">&quot;Tridiagonalization is not initialized.&quot;</span>);</div>
<div class="line"><a name="l00316"></a><span class="lineno">  316</span>&#160;  Index n = m_matrix.rows();</div>
<div class="line"><a name="l00317"></a><span class="lineno">  317</span>&#160;  <span class="keywordflow">return</span> <a class="code" href="classEigen_1_1Block.html">Block&lt;const MatrixType,SizeMinusOne,SizeMinusOne&gt;</a>(m_matrix, 1, 0, n-1,n-1).diagonal();</div>
<div class="line"><a name="l00318"></a><span class="lineno">  318</span>&#160;}</div>
<div class="line"><a name="l00319"></a><span class="lineno">  319</span>&#160;</div>
<div class="line"><a name="l00320"></a><span class="lineno">  320</span>&#160;<span class="keyword">namespace </span>internal {</div>
<div class="line"><a name="l00321"></a><span class="lineno">  321</span>&#160;</div>
<div class="line"><a name="l00345"></a><span class="lineno">  345</span>&#160;<span class="keyword">template</span>&lt;<span class="keyword">typename</span> MatrixType, <span class="keyword">typename</span> CoeffVectorType&gt;</div>
<div class="line"><a name="l00346"></a><span class="lineno">  346</span>&#160;<span class="keywordtype">void</span> tridiagonalization_inplace(MatrixType&amp; matA, CoeffVectorType&amp; hCoeffs)</div>
<div class="line"><a name="l00347"></a><span class="lineno">  347</span>&#160;{</div>
<div class="line"><a name="l00348"></a><span class="lineno">  348</span>&#160;  <span class="keyword">using</span> numext::conj;</div>
<div class="line"><a name="l00349"></a><span class="lineno">  349</span>&#160;  <span class="keyword">typedef</span> <span class="keyword">typename</span> MatrixType::Index Index;</div>
<div class="line"><a name="l00350"></a><span class="lineno">  350</span>&#160;  <span class="keyword">typedef</span> <span class="keyword">typename</span> MatrixType::Scalar Scalar;</div>
<div class="line"><a name="l00351"></a><span class="lineno">  351</span>&#160;  <span class="keyword">typedef</span> <span class="keyword">typename</span> MatrixType::RealScalar RealScalar;</div>
<div class="line"><a name="l00352"></a><span class="lineno">  352</span>&#160;  Index n = matA.rows();</div>
<div class="line"><a name="l00353"></a><span class="lineno">  353</span>&#160;  eigen_assert(n==matA.cols());</div>
<div class="line"><a name="l00354"></a><span class="lineno">  354</span>&#160;  eigen_assert(n==hCoeffs.size()+1 || n==1);</div>
<div class="line"><a name="l00355"></a><span class="lineno">  355</span>&#160;  </div>
<div class="line"><a name="l00356"></a><span class="lineno">  356</span>&#160;  <span class="keywordflow">for</span> (Index i = 0; i&lt;n-1; ++i)</div>
<div class="line"><a name="l00357"></a><span class="lineno">  357</span>&#160;  {</div>
<div class="line"><a name="l00358"></a><span class="lineno">  358</span>&#160;    Index remainingSize = n-i-1;</div>
<div class="line"><a name="l00359"></a><span class="lineno">  359</span>&#160;    RealScalar beta;</div>
<div class="line"><a name="l00360"></a><span class="lineno">  360</span>&#160;    Scalar h;</div>
<div class="line"><a name="l00361"></a><span class="lineno">  361</span>&#160;    matA.col(i).tail(remainingSize).makeHouseholderInPlace(h, beta);</div>
<div class="line"><a name="l00362"></a><span class="lineno">  362</span>&#160;</div>
<div class="line"><a name="l00363"></a><span class="lineno">  363</span>&#160;    <span class="comment">// Apply similarity transformation to remaining columns,</span></div>
<div class="line"><a name="l00364"></a><span class="lineno">  364</span>&#160;    <span class="comment">// i.e., A = H A H&#39; where H = I - h v v&#39; and v = matA.col(i).tail(n-i-1)</span></div>
<div class="line"><a name="l00365"></a><span class="lineno">  365</span>&#160;    matA.col(i).coeffRef(i+1) = 1;</div>
<div class="line"><a name="l00366"></a><span class="lineno">  366</span>&#160;</div>
<div class="line"><a name="l00367"></a><span class="lineno">  367</span>&#160;    hCoeffs.tail(n-i-1).noalias() = (matA.bottomRightCorner(remainingSize,remainingSize).template selfadjointView&lt;Lower&gt;()</div>
<div class="line"><a name="l00368"></a><span class="lineno">  368</span>&#160;                                  * (conj(h) * matA.col(i).tail(remainingSize)));</div>
<div class="line"><a name="l00369"></a><span class="lineno">  369</span>&#160;</div>
<div class="line"><a name="l00370"></a><span class="lineno">  370</span>&#160;    hCoeffs.tail(n-i-1) += (conj(h)*Scalar(-0.5)*(hCoeffs.tail(remainingSize).dot(matA.col(i).tail(remainingSize)))) * matA.col(i).tail(n-i-1);</div>
<div class="line"><a name="l00371"></a><span class="lineno">  371</span>&#160;</div>
<div class="line"><a name="l00372"></a><span class="lineno">  372</span>&#160;    matA.bottomRightCorner(remainingSize, remainingSize).template selfadjointView&lt;Lower&gt;()</div>
<div class="line"><a name="l00373"></a><span class="lineno">  373</span>&#160;      .rankUpdate(matA.col(i).tail(remainingSize), hCoeffs.tail(remainingSize), -1);</div>
<div class="line"><a name="l00374"></a><span class="lineno">  374</span>&#160;</div>
<div class="line"><a name="l00375"></a><span class="lineno">  375</span>&#160;    matA.col(i).coeffRef(i+1) = beta;</div>
<div class="line"><a name="l00376"></a><span class="lineno">  376</span>&#160;    hCoeffs.coeffRef(i) = h;</div>
<div class="line"><a name="l00377"></a><span class="lineno">  377</span>&#160;  }</div>
<div class="line"><a name="l00378"></a><span class="lineno">  378</span>&#160;}</div>
<div class="line"><a name="l00379"></a><span class="lineno">  379</span>&#160;</div>
<div class="line"><a name="l00380"></a><span class="lineno">  380</span>&#160;<span class="comment">// forward declaration, implementation at the end of this file</span></div>
<div class="line"><a name="l00381"></a><span class="lineno">  381</span>&#160;<span class="keyword">template</span>&lt;<span class="keyword">typename</span> MatrixType,</div>
<div class="line"><a name="l00382"></a><span class="lineno">  382</span>&#160;         <span class="keywordtype">int</span> Size=MatrixType::ColsAtCompileTime,</div>
<div class="line"><a name="l00383"></a><span class="lineno">  383</span>&#160;         <span class="keywordtype">bool</span> IsComplex=NumTraits&lt;typename MatrixType::Scalar&gt;::IsComplex&gt;</div>
<div class="line"><a name="l00384"></a><span class="lineno">  384</span>&#160;<span class="keyword">struct </span>tridiagonalization_inplace_selector;</div>
<div class="line"><a name="l00385"></a><span class="lineno">  385</span>&#160;</div>
<div class="line"><a name="l00426"></a><span class="lineno">  426</span>&#160;<span class="keyword">template</span>&lt;<span class="keyword">typename</span> MatrixType, <span class="keyword">typename</span> DiagonalType, <span class="keyword">typename</span> SubDiagonalType&gt;</div>
<div class="line"><a name="l00427"></a><span class="lineno"><a class="line" href="namespaceEigen_1_1internal.html#aa53570cf2e676b41631f08397658ca0f">  427</a></span>&#160;<span class="keywordtype">void</span> tridiagonalization_inplace(MatrixType&amp; mat, DiagonalType&amp; diag, SubDiagonalType&amp; subdiag, <span class="keywordtype">bool</span> extractQ)</div>
<div class="line"><a name="l00428"></a><span class="lineno">  428</span>&#160;{</div>
<div class="line"><a name="l00429"></a><span class="lineno">  429</span>&#160;  eigen_assert(mat.cols()==mat.rows() &amp;&amp; diag.size()==mat.rows() &amp;&amp; subdiag.size()==mat.rows()-1);</div>
<div class="line"><a name="l00430"></a><span class="lineno">  430</span>&#160;  tridiagonalization_inplace_selector&lt;MatrixType&gt;::run(mat, diag, subdiag, extractQ);</div>
<div class="line"><a name="l00431"></a><span class="lineno">  431</span>&#160;}</div>
<div class="line"><a name="l00432"></a><span class="lineno">  432</span>&#160;</div>
<div class="line"><a name="l00436"></a><span class="lineno">  436</span>&#160;<span class="keyword">template</span>&lt;<span class="keyword">typename</span> MatrixType, <span class="keywordtype">int</span> Size, <span class="keywordtype">bool</span> IsComplex&gt;</div>
<div class="line"><a name="l00437"></a><span class="lineno">  437</span>&#160;<span class="keyword">struct </span>tridiagonalization_inplace_selector</div>
<div class="line"><a name="l00438"></a><span class="lineno">  438</span>&#160;{</div>
<div class="line"><a name="l00439"></a><span class="lineno">  439</span>&#160;  <span class="keyword">typedef</span> <span class="keyword">typename</span> <a class="code" href="classEigen_1_1Matrix.html">Tridiagonalization&lt;MatrixType&gt;::CoeffVectorType</a> CoeffVectorType;</div>
<div class="line"><a name="l00440"></a><span class="lineno">  440</span>&#160;  <span class="keyword">typedef</span> <span class="keyword">typename</span> <a class="code" href="classEigen_1_1HouseholderSequence.html">Tridiagonalization&lt;MatrixType&gt;::HouseholderSequenceType</a> HouseholderSequenceType;</div>
<div class="line"><a name="l00441"></a><span class="lineno">  441</span>&#160;  <span class="keyword">typedef</span> <span class="keyword">typename</span> MatrixType::Index Index;</div>
<div class="line"><a name="l00442"></a><span class="lineno">  442</span>&#160;  <span class="keyword">template</span>&lt;<span class="keyword">typename</span> DiagonalType, <span class="keyword">typename</span> SubDiagonalType&gt;</div>
<div class="line"><a name="l00443"></a><span class="lineno">  443</span>&#160;  <span class="keyword">static</span> <span class="keywordtype">void</span> run(MatrixType&amp; mat, DiagonalType&amp; diag, SubDiagonalType&amp; subdiag, <span class="keywordtype">bool</span> extractQ)</div>
<div class="line"><a name="l00444"></a><span class="lineno">  444</span>&#160;  {</div>
<div class="line"><a name="l00445"></a><span class="lineno">  445</span>&#160;    CoeffVectorType hCoeffs(mat.cols()-1);</div>
<div class="line"><a name="l00446"></a><span class="lineno">  446</span>&#160;    tridiagonalization_inplace(mat,hCoeffs);</div>
<div class="line"><a name="l00447"></a><span class="lineno">  447</span>&#160;    diag = mat.diagonal().real();</div>
<div class="line"><a name="l00448"></a><span class="lineno">  448</span>&#160;    subdiag = mat.template diagonal&lt;-1&gt;().real();</div>
<div class="line"><a name="l00449"></a><span class="lineno">  449</span>&#160;    <span class="keywordflow">if</span>(extractQ)</div>
<div class="line"><a name="l00450"></a><span class="lineno">  450</span>&#160;      mat = HouseholderSequenceType(mat, hCoeffs.conjugate())</div>
<div class="line"><a name="l00451"></a><span class="lineno">  451</span>&#160;            .setLength(mat.rows() - 1)</div>
<div class="line"><a name="l00452"></a><span class="lineno">  452</span>&#160;            .setShift(1);</div>
<div class="line"><a name="l00453"></a><span class="lineno">  453</span>&#160;  }</div>
<div class="line"><a name="l00454"></a><span class="lineno">  454</span>&#160;};</div>
<div class="line"><a name="l00455"></a><span class="lineno">  455</span>&#160;</div>
<div class="line"><a name="l00460"></a><span class="lineno">  460</span>&#160;<span class="keyword">template</span>&lt;<span class="keyword">typename</span> MatrixType&gt;</div>
<div class="line"><a name="l00461"></a><span class="lineno">  461</span>&#160;<span class="keyword">struct </span>tridiagonalization_inplace_selector&lt;MatrixType,3,false&gt;</div>
<div class="line"><a name="l00462"></a><span class="lineno">  462</span>&#160;{</div>
<div class="line"><a name="l00463"></a><span class="lineno">  463</span>&#160;  <span class="keyword">typedef</span> <span class="keyword">typename</span> MatrixType::Scalar Scalar;</div>
<div class="line"><a name="l00464"></a><span class="lineno">  464</span>&#160;  <span class="keyword">typedef</span> <span class="keyword">typename</span> MatrixType::RealScalar RealScalar;</div>
<div class="line"><a name="l00465"></a><span class="lineno">  465</span>&#160;</div>
<div class="line"><a name="l00466"></a><span class="lineno">  466</span>&#160;  <span class="keyword">template</span>&lt;<span class="keyword">typename</span> DiagonalType, <span class="keyword">typename</span> SubDiagonalType&gt;</div>
<div class="line"><a name="l00467"></a><span class="lineno">  467</span>&#160;  <span class="keyword">static</span> <span class="keywordtype">void</span> run(MatrixType&amp; mat, DiagonalType&amp; diag, SubDiagonalType&amp; subdiag, <span class="keywordtype">bool</span> extractQ)</div>
<div class="line"><a name="l00468"></a><span class="lineno">  468</span>&#160;  {</div>
<div class="line"><a name="l00469"></a><span class="lineno">  469</span>&#160;    <span class="keyword">using</span> std::sqrt;</div>
<div class="line"><a name="l00470"></a><span class="lineno">  470</span>&#160;    diag[0] = mat(0,0);</div>
<div class="line"><a name="l00471"></a><span class="lineno">  471</span>&#160;    RealScalar v1norm2 = numext::abs2(mat(2,0));</div>
<div class="line"><a name="l00472"></a><span class="lineno">  472</span>&#160;    <span class="keywordflow">if</span>(v1norm2 == RealScalar(0))</div>
<div class="line"><a name="l00473"></a><span class="lineno">  473</span>&#160;    {</div>
<div class="line"><a name="l00474"></a><span class="lineno">  474</span>&#160;      diag[1] = mat(1,1);</div>
<div class="line"><a name="l00475"></a><span class="lineno">  475</span>&#160;      diag[2] = mat(2,2);</div>
<div class="line"><a name="l00476"></a><span class="lineno">  476</span>&#160;      subdiag[0] = mat(1,0);</div>
<div class="line"><a name="l00477"></a><span class="lineno">  477</span>&#160;      subdiag[1] = mat(2,1);</div>
<div class="line"><a name="l00478"></a><span class="lineno">  478</span>&#160;      <span class="keywordflow">if</span> (extractQ)</div>
<div class="line"><a name="l00479"></a><span class="lineno">  479</span>&#160;        mat.setIdentity();</div>
<div class="line"><a name="l00480"></a><span class="lineno">  480</span>&#160;    }</div>
<div class="line"><a name="l00481"></a><span class="lineno">  481</span>&#160;    <span class="keywordflow">else</span></div>
<div class="line"><a name="l00482"></a><span class="lineno">  482</span>&#160;    {</div>
<div class="line"><a name="l00483"></a><span class="lineno">  483</span>&#160;      RealScalar beta = sqrt(numext::abs2(mat(1,0)) + v1norm2);</div>
<div class="line"><a name="l00484"></a><span class="lineno">  484</span>&#160;      RealScalar invBeta = RealScalar(1)/beta;</div>
<div class="line"><a name="l00485"></a><span class="lineno">  485</span>&#160;      Scalar m01 = mat(1,0) * invBeta;</div>
<div class="line"><a name="l00486"></a><span class="lineno">  486</span>&#160;      Scalar m02 = mat(2,0) * invBeta;</div>
<div class="line"><a name="l00487"></a><span class="lineno">  487</span>&#160;      Scalar q = RealScalar(2)*m01*mat(2,1) + m02*(mat(2,2) - mat(1,1));</div>
<div class="line"><a name="l00488"></a><span class="lineno">  488</span>&#160;      diag[1] = mat(1,1) + m02*q;</div>
<div class="line"><a name="l00489"></a><span class="lineno">  489</span>&#160;      diag[2] = mat(2,2) - m02*q;</div>
<div class="line"><a name="l00490"></a><span class="lineno">  490</span>&#160;      subdiag[0] = beta;</div>
<div class="line"><a name="l00491"></a><span class="lineno">  491</span>&#160;      subdiag[1] = mat(2,1) - m01 * q;</div>
<div class="line"><a name="l00492"></a><span class="lineno">  492</span>&#160;      <span class="keywordflow">if</span> (extractQ)</div>
<div class="line"><a name="l00493"></a><span class="lineno">  493</span>&#160;      {</div>
<div class="line"><a name="l00494"></a><span class="lineno">  494</span>&#160;        mat &lt;&lt; 1,   0,    0,</div>
<div class="line"><a name="l00495"></a><span class="lineno">  495</span>&#160;               0, m01,  m02,</div>
<div class="line"><a name="l00496"></a><span class="lineno">  496</span>&#160;               0, m02, -m01;</div>
<div class="line"><a name="l00497"></a><span class="lineno">  497</span>&#160;      }</div>
<div class="line"><a name="l00498"></a><span class="lineno">  498</span>&#160;    }</div>
<div class="line"><a name="l00499"></a><span class="lineno">  499</span>&#160;  }</div>
<div class="line"><a name="l00500"></a><span class="lineno">  500</span>&#160;};</div>
<div class="line"><a name="l00501"></a><span class="lineno">  501</span>&#160;</div>
<div class="line"><a name="l00505"></a><span class="lineno">  505</span>&#160;<span class="keyword">template</span>&lt;<span class="keyword">typename</span> MatrixType, <span class="keywordtype">bool</span> IsComplex&gt;</div>
<div class="line"><a name="l00506"></a><span class="lineno">  506</span>&#160;<span class="keyword">struct </span>tridiagonalization_inplace_selector&lt;MatrixType,1,IsComplex&gt;</div>
<div class="line"><a name="l00507"></a><span class="lineno">  507</span>&#160;{</div>
<div class="line"><a name="l00508"></a><span class="lineno">  508</span>&#160;  <span class="keyword">typedef</span> <span class="keyword">typename</span> MatrixType::Scalar Scalar;</div>
<div class="line"><a name="l00509"></a><span class="lineno">  509</span>&#160;</div>
<div class="line"><a name="l00510"></a><span class="lineno">  510</span>&#160;  <span class="keyword">template</span>&lt;<span class="keyword">typename</span> DiagonalType, <span class="keyword">typename</span> SubDiagonalType&gt;</div>
<div class="line"><a name="l00511"></a><span class="lineno">  511</span>&#160;  <span class="keyword">static</span> <span class="keywordtype">void</span> run(MatrixType&amp; mat, DiagonalType&amp; diag, SubDiagonalType&amp;, <span class="keywordtype">bool</span> extractQ)</div>
<div class="line"><a name="l00512"></a><span class="lineno">  512</span>&#160;  {</div>
<div class="line"><a name="l00513"></a><span class="lineno">  513</span>&#160;    diag(0,0) = numext::real(mat(0,0));</div>
<div class="line"><a name="l00514"></a><span class="lineno">  514</span>&#160;    <span class="keywordflow">if</span>(extractQ)</div>
<div class="line"><a name="l00515"></a><span class="lineno">  515</span>&#160;      mat(0,0) = Scalar(1);</div>
<div class="line"><a name="l00516"></a><span class="lineno">  516</span>&#160;  }</div>
<div class="line"><a name="l00517"></a><span class="lineno">  517</span>&#160;};</div>
<div class="line"><a name="l00518"></a><span class="lineno">  518</span>&#160;</div>
<div class="line"><a name="l00526"></a><span class="lineno">  526</span>&#160;<span class="keyword">template</span>&lt;<span class="keyword">typename</span> MatrixType&gt; <span class="keyword">struct </span>TridiagonalizationMatrixTReturnType</div>
<div class="line"><a name="l00527"></a><span class="lineno">  527</span>&#160;: <span class="keyword">public</span> ReturnByValue&lt;TridiagonalizationMatrixTReturnType&lt;MatrixType&gt; &gt;</div>
<div class="line"><a name="l00528"></a><span class="lineno">  528</span>&#160;{</div>
<div class="line"><a name="l00529"></a><span class="lineno">  529</span>&#160;    <span class="keyword">typedef</span> <span class="keyword">typename</span> MatrixType::Index Index;</div>
<div class="line"><a name="l00530"></a><span class="lineno">  530</span>&#160;  <span class="keyword">public</span>:</div>
<div class="line"><a name="l00535"></a><span class="lineno">  535</span>&#160;    TridiagonalizationMatrixTReturnType(<span class="keyword">const</span> MatrixType&amp; mat) : m_matrix(mat) { }</div>
<div class="line"><a name="l00536"></a><span class="lineno">  536</span>&#160;</div>
<div class="line"><a name="l00537"></a><span class="lineno">  537</span>&#160;    <span class="keyword">template</span> &lt;<span class="keyword">typename</span> ResultType&gt;</div>
<div class="line"><a name="l00538"></a><span class="lineno">  538</span>&#160;    <span class="keyword">inline</span> <span class="keywordtype">void</span> evalTo(ResultType&amp; result)<span class="keyword"> const</span></div>
<div class="line"><a name="l00539"></a><span class="lineno">  539</span>&#160;<span class="keyword">    </span>{</div>
<div class="line"><a name="l00540"></a><span class="lineno">  540</span>&#160;      result.setZero();</div>
<div class="line"><a name="l00541"></a><span class="lineno">  541</span>&#160;      result.template diagonal&lt;1&gt;() = m_matrix.template diagonal&lt;-1&gt;().conjugate();</div>
<div class="line"><a name="l00542"></a><span class="lineno">  542</span>&#160;      result.diagonal() = m_matrix.diagonal();</div>
<div class="line"><a name="l00543"></a><span class="lineno">  543</span>&#160;      result.template diagonal&lt;-1&gt;() = m_matrix.template diagonal&lt;-1&gt;();</div>
<div class="line"><a name="l00544"></a><span class="lineno">  544</span>&#160;    }</div>
<div class="line"><a name="l00545"></a><span class="lineno">  545</span>&#160;</div>
<div class="line"><a name="l00546"></a><span class="lineno">  546</span>&#160;    Index rows()<span class="keyword"> const </span>{ <span class="keywordflow">return</span> m_matrix.rows(); }</div>
<div class="line"><a name="l00547"></a><span class="lineno">  547</span>&#160;    Index cols()<span class="keyword"> const </span>{ <span class="keywordflow">return</span> m_matrix.cols(); }</div>
<div class="line"><a name="l00548"></a><span class="lineno">  548</span>&#160;</div>
<div class="line"><a name="l00549"></a><span class="lineno">  549</span>&#160;  <span class="keyword">protected</span>:</div>
<div class="line"><a name="l00550"></a><span class="lineno">  550</span>&#160;    <span class="keyword">typename</span> MatrixType::Nested m_matrix;</div>
<div class="line"><a name="l00551"></a><span class="lineno">  551</span>&#160;};</div>
<div class="line"><a name="l00552"></a><span class="lineno">  552</span>&#160;</div>
<div class="line"><a name="l00553"></a><span class="lineno">  553</span>&#160;} <span class="comment">// end namespace internal</span></div>
<div class="line"><a name="l00554"></a><span class="lineno">  554</span>&#160;</div>
<div class="line"><a name="l00555"></a><span class="lineno">  555</span>&#160;} <span class="comment">// end namespace Eigen</span></div>
<div class="line"><a name="l00556"></a><span class="lineno">  556</span>&#160;</div>
<div class="line"><a name="l00557"></a><span class="lineno">  557</span>&#160;<span class="preprocessor">#endif // EIGEN_TRIDIAGONALIZATION_H</span></div>
<div class="ttc" id="classEigen_1_1Tridiagonalization_html_a8fa49216273ab7579b7bea06debb1e51"><div class="ttname"><a href="classEigen_1_1Tridiagonalization.html#a8fa49216273ab7579b7bea06debb1e51">Eigen::Tridiagonalization::subDiagonal</a></div><div class="ttdeci">SubDiagonalReturnType subDiagonal() const </div><div class="ttdoc">Returns the subdiagonal of the tridiagonal matrix T in the decomposition. </div><div class="ttdef"><b>Definition:</b> Tridiagonalization.h:313</div></div>
<div class="ttc" id="structEigen_1_1NumTraits_html"><div class="ttname"><a href="structEigen_1_1NumTraits.html">Eigen::NumTraits</a></div><div class="ttdoc">Holds information about the various numeric (i.e. scalar) types allowed by Eigen. ...</div><div class="ttdef"><b>Definition:</b> NumTraits.h:88</div></div>
<div class="ttc" id="classEigen_1_1Tridiagonalization_html"><div class="ttname"><a href="classEigen_1_1Tridiagonalization.html">Eigen::Tridiagonalization</a></div><div class="ttdoc">Tridiagonal decomposition of a selfadjoint matrix. </div><div class="ttdef"><b>Definition:</b> Tridiagonalization.h:61</div></div>
<div class="ttc" id="namespaceEigen_html_adc9da5be31bdce40c25a92c27999c0e3"><div class="ttname"><a href="namespaceEigen.html#adc9da5be31bdce40c25a92c27999c0e3">Eigen::Dynamic</a></div><div class="ttdeci">const int Dynamic</div><div class="ttdef"><b>Definition:</b> Constants.h:21</div></div>
<div class="ttc" id="classEigen_1_1Tridiagonalization_html_aa96bdbc1b19c647e3372c31301ea4999"><div class="ttname"><a href="classEigen_1_1Tridiagonalization.html#aa96bdbc1b19c647e3372c31301ea4999">Eigen::Tridiagonalization::HouseholderSequenceType</a></div><div class="ttdeci">HouseholderSequence&lt; MatrixType, typename internal::remove_all&lt; typename CoeffVectorType::ConjugateReturnType &gt;::type &gt; HouseholderSequenceType</div><div class="ttdoc">Return type of matrixQ() </div><div class="ttdef"><b>Definition:</b> Tridiagonalization.h:99</div></div>
<div class="ttc" id="classEigen_1_1HouseholderSequence_html"><div class="ttname"><a href="classEigen_1_1HouseholderSequence.html">Eigen::HouseholderSequence</a></div><div class="ttdoc">Sequence of Householder reflections acting on subspaces with decreasing size. </div><div class="ttdef"><b>Definition:</b> ForwardDeclarations.h:227</div></div>
<div class="ttc" id="classEigen_1_1Tridiagonalization_html_aa9f9722d2cef9425e2c0da3553dfbac7"><div class="ttname"><a href="classEigen_1_1Tridiagonalization.html#aa9f9722d2cef9425e2c0da3553dfbac7">Eigen::Tridiagonalization::Tridiagonalization</a></div><div class="ttdeci">Tridiagonalization(const MatrixType &amp;matrix)</div><div class="ttdoc">Constructor; computes tridiagonal decomposition of given matrix. </div><div class="ttdef"><b>Definition:</b> Tridiagonalization.h:129</div></div>
<div class="ttc" id="classEigen_1_1Tridiagonalization_html_aceb0f16a166f4c236a1b536b7424d292"><div class="ttname"><a href="classEigen_1_1Tridiagonalization.html#aceb0f16a166f4c236a1b536b7424d292">Eigen::Tridiagonalization::matrixT</a></div><div class="ttdeci">MatrixTReturnType matrixT() const </div><div class="ttdoc">Returns an expression of the tridiagonal matrix T in the decomposition. </div><div class="ttdef"><b>Definition:</b> Tridiagonalization.h:263</div></div>
<div class="ttc" id="classEigen_1_1Tridiagonalization_html_aa69e607a4aab4fb6321ca6acbf074fc2"><div class="ttname"><a href="classEigen_1_1Tridiagonalization.html#aa69e607a4aab4fb6321ca6acbf074fc2">Eigen::Tridiagonalization::compute</a></div><div class="ttdeci">Tridiagonalization &amp; compute(const MatrixType &amp;matrix)</div><div class="ttdoc">Computes tridiagonal decomposition of given matrix. </div><div class="ttdef"><b>Definition:</b> Tridiagonalization.h:155</div></div>
<div class="ttc" id="classEigen_1_1Tridiagonalization_html_a0698ae78b0ab6f239c475b73b9c6bbee"><div class="ttname"><a href="classEigen_1_1Tridiagonalization.html#a0698ae78b0ab6f239c475b73b9c6bbee">Eigen::Tridiagonalization::Tridiagonalization</a></div><div class="ttdeci">Tridiagonalization(Index size=Size==Dynamic?2:Size)</div><div class="ttdoc">Default constructor. </div><div class="ttdef"><b>Definition:</b> Tridiagonalization.h:113</div></div>
<div class="ttc" id="classEigen_1_1Tridiagonalization_html_a66adece364b64b26b3771662de70f2df"><div class="ttname"><a href="classEigen_1_1Tridiagonalization.html#a66adece364b64b26b3771662de70f2df">Eigen::Tridiagonalization::packedMatrix</a></div><div class="ttdeci">const MatrixType &amp; packedMatrix() const </div><div class="ttdoc">Returns the internal representation of the decomposition. </div><div class="ttdef"><b>Definition:</b> Tridiagonalization.h:217</div></div>
<div class="ttc" id="classEigen_1_1Block_html"><div class="ttname"><a href="classEigen_1_1Block.html">Eigen::Block</a></div><div class="ttdoc">Expression of a fixed-size or dynamic-size block. </div><div class="ttdef"><b>Definition:</b> Block.h:102</div></div>
<div class="ttc" id="classEigen_1_1Tridiagonalization_html_a2ab889c75460c178d941ee24e371b206"><div class="ttname"><a href="classEigen_1_1Tridiagonalization.html#a2ab889c75460c178d941ee24e371b206">Eigen::Tridiagonalization::householderCoefficients</a></div><div class="ttdeci">CoeffVectorType householderCoefficients() const </div><div class="ttdoc">Returns the Householder coefficients. </div><div class="ttdef"><b>Definition:</b> Tridiagonalization.h:180</div></div>
<div class="ttc" id="classEigen_1_1Tridiagonalization_html_ad13845d7490115664924b3dc208ec369"><div class="ttname"><a href="classEigen_1_1Tridiagonalization.html#ad13845d7490115664924b3dc208ec369">Eigen::Tridiagonalization::matrixQ</a></div><div class="ttdeci">HouseholderSequenceType matrixQ() const </div><div class="ttdoc">Returns the unitary matrix Q in the decomposition. </div><div class="ttdef"><b>Definition:</b> Tridiagonalization.h:238</div></div>
<div class="ttc" id="classEigen_1_1Diagonal_html"><div class="ttname"><a href="classEigen_1_1Diagonal.html">Eigen::Diagonal</a></div><div class="ttdoc">Expression of a diagonal/subdiagonal/superdiagonal in a matrix. </div><div class="ttdef"><b>Definition:</b> Diagonal.h:64</div></div>
<div class="ttc" id="classEigen_1_1PlainObjectBase_html_afbbb33d14fe7fb9683019a39ce1c659d"><div class="ttname"><a href="classEigen_1_1PlainObjectBase.html#afbbb33d14fe7fb9683019a39ce1c659d">Eigen::PlainObjectBase::resize</a></div><div class="ttdeci">void resize(Index nbRows, Index nbCols)</div><div class="ttdef"><b>Definition:</b> PlainObjectBase.h:232</div></div>
<div class="ttc" id="classEigen_1_1Tridiagonalization_html_ac109eefddd733d8e82841da5bb2dd8d3"><div class="ttname"><a href="classEigen_1_1Tridiagonalization.html#ac109eefddd733d8e82841da5bb2dd8d3">Eigen::Tridiagonalization::diagonal</a></div><div class="ttdeci">DiagonalReturnType diagonal() const </div><div class="ttdoc">Returns the diagonal of the tridiagonal matrix T in the decomposition. </div><div class="ttdef"><b>Definition:</b> Tridiagonalization.h:305</div></div>
<div class="ttc" id="classEigen_1_1Matrix_html"><div class="ttname"><a href="classEigen_1_1Matrix.html">Eigen::Matrix&lt; Scalar, SizeMinusOne, 1, Options &amp;~RowMajor, MaxSizeMinusOne, 1 &gt;</a></div></div>
<div class="ttc" id="classEigen_1_1Tridiagonalization_html_aeb6c0eb89cc982629305f6c7e0791caf"><div class="ttname"><a href="classEigen_1_1Tridiagonalization.html#aeb6c0eb89cc982629305f6c7e0791caf">Eigen::Tridiagonalization::MatrixType</a></div><div class="ttdeci">_MatrixType MatrixType</div><div class="ttdoc">Synonym for the template parameter _MatrixType. </div><div class="ttdef"><b>Definition:</b> Tridiagonalization.h:66</div></div>
</div><!-- fragment --></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="dir_e49d68e3078f12dfcf157021597ad168.html">Eigen</a></li><li class="navelem"><a class="el" href="dir_64b228556dc7f9fe757d43bb57fbfc24.html">src</a></li><li class="navelem"><a class="el" href="dir_b3e8aad20632d7b7c332c11ff568cf95.html">Eigenvalues</a></li><li class="navelem"><b>Tridiagonalization.h</b></li>
    <li class="footer">Generated on Mon Oct 28 2013 11:04:26 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>