<!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: SuiteSparseQRSupport.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>  <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('SuiteSparseQRSupport_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"> </span>All</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(1)"><span class="SelectionMark"> </span>Classes</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(2)"><span class="SelectionMark"> </span>Namespaces</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(3)"><span class="SelectionMark"> </span>Functions</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(4)"><span class="SelectionMark"> </span>Variables</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(5)"><span class="SelectionMark"> </span>Typedefs</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(6)"><span class="SelectionMark"> </span>Enumerations</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(7)"><span class="SelectionMark"> </span>Enumerator</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(8)"><span class="SelectionMark"> </span>Friends</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(9)"><span class="SelectionMark"> </span>Groups</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(10)"><span class="SelectionMark"> </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">SuiteSparseQRSupport.h</div> </div> </div><!--header--> <div class="contents"> <div class="fragment"><div class="line"><a name="l00001"></a><span class="lineno"> 1</span> <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> <span class="comment">// for linear algebra.</span></div> <div class="line"><a name="l00003"></a><span class="lineno"> 3</span> <span class="comment">//</span></div> <div class="line"><a name="l00004"></a><span class="lineno"> 4</span> <span class="comment">// Copyright (C) 2012 Desire Nuentsa <desire.nuentsa_wakam@inria.fr></span></div> <div class="line"><a name="l00005"></a><span class="lineno"> 5</span> <span class="comment">//</span></div> <div class="line"><a name="l00006"></a><span class="lineno"> 6</span> <span class="comment">// This Source Code Form is subject to the terms of the Mozilla</span></div> <div class="line"><a name="l00007"></a><span class="lineno"> 7</span> <span class="comment">// Public License v. 2.0. If a copy of the MPL was not distributed</span></div> <div class="line"><a name="l00008"></a><span class="lineno"> 8</span> <span class="comment">// with this file, You can obtain one at http://mozilla.org/MPL/2.0/.</span></div> <div class="line"><a name="l00009"></a><span class="lineno"> 9</span> </div> <div class="line"><a name="l00010"></a><span class="lineno"> 10</span> <span class="preprocessor">#ifndef EIGEN_SUITESPARSEQRSUPPORT_H</span></div> <div class="line"><a name="l00011"></a><span class="lineno"> 11</span> <span class="preprocessor"></span><span class="preprocessor">#define EIGEN_SUITESPARSEQRSUPPORT_H</span></div> <div class="line"><a name="l00012"></a><span class="lineno"> 12</span> <span class="preprocessor"></span></div> <div class="line"><a name="l00013"></a><span class="lineno"> 13</span> <span class="keyword">namespace </span>Eigen {</div> <div class="line"><a name="l00014"></a><span class="lineno"> 14</span>  </div> <div class="line"><a name="l00015"></a><span class="lineno"><a class="line" href="classEigen_1_1SPQR.html"> 15</a></span>  <span class="keyword">template</span><<span class="keyword">typename</span> MatrixType> <span class="keyword">class </span><a class="code" href="classEigen_1_1SPQR.html">SPQR</a>; </div> <div class="line"><a name="l00016"></a><span class="lineno"> 16</span>  <span class="keyword">template</span><<span class="keyword">typename</span> SPQRType> <span class="keyword">struct </span>SPQRMatrixQReturnType; </div> <div class="line"><a name="l00017"></a><span class="lineno"> 17</span>  <span class="keyword">template</span><<span class="keyword">typename</span> SPQRType> <span class="keyword">struct </span>SPQRMatrixQTransposeReturnType; </div> <div class="line"><a name="l00018"></a><span class="lineno"> 18</span>  <span class="keyword">template</span> <<span class="keyword">typename</span> SPQRType, <span class="keyword">typename</span> Derived> <span class="keyword">struct </span>SPQR_QProduct;</div> <div class="line"><a name="l00019"></a><span class="lineno"> 19</span>  <span class="keyword">namespace </span>internal {</div> <div class="line"><a name="l00020"></a><span class="lineno"> 20</span>  <span class="keyword">template</span> <<span class="keyword">typename</span> SPQRType> <span class="keyword">struct </span>traits<SPQRMatrixQReturnType<SPQRType> ></div> <div class="line"><a name="l00021"></a><span class="lineno"> 21</span>  {</div> <div class="line"><a name="l00022"></a><span class="lineno"> 22</span>  <span class="keyword">typedef</span> <span class="keyword">typename</span> SPQRType::MatrixType ReturnType;</div> <div class="line"><a name="l00023"></a><span class="lineno"> 23</span>  };</div> <div class="line"><a name="l00024"></a><span class="lineno"> 24</span>  <span class="keyword">template</span> <<span class="keyword">typename</span> SPQRType> <span class="keyword">struct </span>traits<SPQRMatrixQTransposeReturnType<SPQRType> ></div> <div class="line"><a name="l00025"></a><span class="lineno"> 25</span>  {</div> <div class="line"><a name="l00026"></a><span class="lineno"> 26</span>  <span class="keyword">typedef</span> <span class="keyword">typename</span> SPQRType::MatrixType ReturnType;</div> <div class="line"><a name="l00027"></a><span class="lineno"> 27</span>  };</div> <div class="line"><a name="l00028"></a><span class="lineno"> 28</span>  <span class="keyword">template</span> <<span class="keyword">typename</span> SPQRType, <span class="keyword">typename</span> Derived> <span class="keyword">struct </span>traits<SPQR_QProduct<SPQRType, Derived> ></div> <div class="line"><a name="l00029"></a><span class="lineno"> 29</span>  {</div> <div class="line"><a name="l00030"></a><span class="lineno"> 30</span>  <span class="keyword">typedef</span> <span class="keyword">typename</span> Derived::PlainObject ReturnType;</div> <div class="line"><a name="l00031"></a><span class="lineno"> 31</span>  };</div> <div class="line"><a name="l00032"></a><span class="lineno"> 32</span>  } <span class="comment">// End namespace internal</span></div> <div class="line"><a name="l00033"></a><span class="lineno"> 33</span>  </div> <div class="line"><a name="l00056"></a><span class="lineno"> 56</span> <span class="keyword">template</span><<span class="keyword">typename</span> _MatrixType></div> <div class="line"><a name="l00057"></a><span class="lineno"> 57</span> <span class="keyword">class </span>SPQR</div> <div class="line"><a name="l00058"></a><span class="lineno"> 58</span> {</div> <div class="line"><a name="l00059"></a><span class="lineno"> 59</span>  <span class="keyword">public</span>:</div> <div class="line"><a name="l00060"></a><span class="lineno"> 60</span>  <span class="keyword">typedef</span> <span class="keyword">typename</span> _MatrixType::Scalar Scalar;</div> <div class="line"><a name="l00061"></a><span class="lineno"> 61</span>  <span class="keyword">typedef</span> <span class="keyword">typename</span> _MatrixType::RealScalar RealScalar;</div> <div class="line"><a name="l00062"></a><span class="lineno"> 62</span>  <span class="keyword">typedef</span> UF_long Index ; </div> <div class="line"><a name="l00063"></a><span class="lineno"> 63</span>  <span class="keyword">typedef</span> SparseMatrix<Scalar, ColMajor, Index> MatrixType;</div> <div class="line"><a name="l00064"></a><span class="lineno"> 64</span>  <span class="keyword">typedef</span> PermutationMatrix<Dynamic, Dynamic> PermutationType;</div> <div class="line"><a name="l00065"></a><span class="lineno"> 65</span>  <span class="keyword">public</span>:</div> <div class="line"><a name="l00066"></a><span class="lineno"> 66</span>  SPQR() </div> <div class="line"><a name="l00067"></a><span class="lineno"> 67</span>  : m_ordering(SPQR_ORDERING_DEFAULT),</div> <div class="line"><a name="l00068"></a><span class="lineno"> 68</span>  m_allow_tol(SPQR_DEFAULT_TOL),</div> <div class="line"><a name="l00069"></a><span class="lineno"> 69</span>  m_tolerance (NumTraits<Scalar>::epsilon())</div> <div class="line"><a name="l00070"></a><span class="lineno"> 70</span>  { </div> <div class="line"><a name="l00071"></a><span class="lineno"> 71</span>  cholmod_l_start(&m_cc);</div> <div class="line"><a name="l00072"></a><span class="lineno"> 72</span>  }</div> <div class="line"><a name="l00073"></a><span class="lineno"> 73</span>  </div> <div class="line"><a name="l00074"></a><span class="lineno"> 74</span>  SPQR(<span class="keyword">const</span> _MatrixType& matrix) </div> <div class="line"><a name="l00075"></a><span class="lineno"> 75</span>  : m_ordering(SPQR_ORDERING_DEFAULT),</div> <div class="line"><a name="l00076"></a><span class="lineno"> 76</span>  m_allow_tol(SPQR_DEFAULT_TOL),</div> <div class="line"><a name="l00077"></a><span class="lineno"> 77</span>  m_tolerance (NumTraits<Scalar>::epsilon())</div> <div class="line"><a name="l00078"></a><span class="lineno"> 78</span>  {</div> <div class="line"><a name="l00079"></a><span class="lineno"> 79</span>  cholmod_l_start(&m_cc);</div> <div class="line"><a name="l00080"></a><span class="lineno"> 80</span>  compute(matrix);</div> <div class="line"><a name="l00081"></a><span class="lineno"> 81</span>  }</div> <div class="line"><a name="l00082"></a><span class="lineno"> 82</span>  </div> <div class="line"><a name="l00083"></a><span class="lineno"> 83</span>  ~SPQR()</div> <div class="line"><a name="l00084"></a><span class="lineno"> 84</span>  {</div> <div class="line"><a name="l00085"></a><span class="lineno"> 85</span>  <span class="comment">// Calls SuiteSparseQR_free()</span></div> <div class="line"><a name="l00086"></a><span class="lineno"> 86</span>  cholmod_l_free_sparse(&m_H, &m_cc);</div> <div class="line"><a name="l00087"></a><span class="lineno"> 87</span>  cholmod_l_free_sparse(&m_cR, &m_cc);</div> <div class="line"><a name="l00088"></a><span class="lineno"> 88</span>  cholmod_l_free_dense(&m_HTau, &m_cc);</div> <div class="line"><a name="l00089"></a><span class="lineno"> 89</span>  std::free(m_E);</div> <div class="line"><a name="l00090"></a><span class="lineno"> 90</span>  std::free(m_HPinv);</div> <div class="line"><a name="l00091"></a><span class="lineno"> 91</span>  cholmod_l_finish(&m_cc);</div> <div class="line"><a name="l00092"></a><span class="lineno"> 92</span>  }</div> <div class="line"><a name="l00093"></a><span class="lineno"> 93</span>  <span class="keywordtype">void</span> compute(<span class="keyword">const</span> _MatrixType& matrix)</div> <div class="line"><a name="l00094"></a><span class="lineno"> 94</span>  {</div> <div class="line"><a name="l00095"></a><span class="lineno"> 95</span>  MatrixType mat(matrix);</div> <div class="line"><a name="l00096"></a><span class="lineno"> 96</span>  cholmod_sparse A; </div> <div class="line"><a name="l00097"></a><span class="lineno"> 97</span>  A = <a class="code" href="namespaceEigen.html#ab3b27d71e21e833c5287fcbe49509478">viewAsCholmod</a>(mat);</div> <div class="line"><a name="l00098"></a><span class="lineno"> 98</span>  Index col = matrix.cols();</div> <div class="line"><a name="l00099"></a><span class="lineno"> 99</span>  m_rank = SuiteSparseQR<Scalar>(m_ordering, m_tolerance, col, &A, </div> <div class="line"><a name="l00100"></a><span class="lineno"> 100</span>  &m_cR, &m_E, &m_H, &m_HPinv, &m_HTau, &m_cc);</div> <div class="line"><a name="l00101"></a><span class="lineno"> 101</span> </div> <div class="line"><a name="l00102"></a><span class="lineno"> 102</span>  <span class="keywordflow">if</span> (!m_cR)</div> <div class="line"><a name="l00103"></a><span class="lineno"> 103</span>  {</div> <div class="line"><a name="l00104"></a><span class="lineno"> 104</span>  m_info = <a class="code" href="group__enums.html#gga51bc1ac16f26ebe51eae1abb77bd037ba710fff14e8fc77846d4b75d8f4cc2d5c">NumericalIssue</a>; </div> <div class="line"><a name="l00105"></a><span class="lineno"> 105</span>  m_isInitialized = <span class="keyword">false</span>;</div> <div class="line"><a name="l00106"></a><span class="lineno"> 106</span>  <span class="keywordflow">return</span>;</div> <div class="line"><a name="l00107"></a><span class="lineno"> 107</span>  }</div> <div class="line"><a name="l00108"></a><span class="lineno"> 108</span>  m_info = <a class="code" href="group__enums.html#gga51bc1ac16f26ebe51eae1abb77bd037bafdfbdf3247bd36a1f17270d5cec74c9c">Success</a>;</div> <div class="line"><a name="l00109"></a><span class="lineno"> 109</span>  m_isInitialized = <span class="keyword">true</span>;</div> <div class="line"><a name="l00110"></a><span class="lineno"> 110</span>  m_isRUpToDate = <span class="keyword">false</span>;</div> <div class="line"><a name="l00111"></a><span class="lineno"> 111</span>  }</div> <div class="line"><a name="l00115"></a><span class="lineno"><a class="line" href="classEigen_1_1SPQR.html#a5552abd83dbd03c85cea6d61fd8875a5"> 115</a></span>  <span class="keyword">inline</span> Index <a class="code" href="classEigen_1_1SPQR.html#a5552abd83dbd03c85cea6d61fd8875a5">rows</a>()<span class="keyword"> const </span>{<span class="keywordflow">return</span> m_H->nrow; }</div> <div class="line"><a name="l00116"></a><span class="lineno"> 116</span>  </div> <div class="line"><a name="l00120"></a><span class="lineno"><a class="line" href="classEigen_1_1SPQR.html#aaca1908a5ec508a25ff0a8bca803e5f3"> 120</a></span>  <span class="keyword">inline</span> Index <a class="code" href="classEigen_1_1SPQR.html#aaca1908a5ec508a25ff0a8bca803e5f3">cols</a>()<span class="keyword"> const </span>{ <span class="keywordflow">return</span> m_cR->ncol; }</div> <div class="line"><a name="l00121"></a><span class="lineno"> 121</span>  </div> <div class="line"><a name="l00126"></a><span class="lineno"> 126</span>  <span class="keyword">template</span><<span class="keyword">typename</span> Rhs></div> <div class="line"><a name="l00127"></a><span class="lineno"><a class="line" href="classEigen_1_1SPQR.html#a7fcaf5446687d32b11ec4fb9fc5587e2"> 127</a></span>  <span class="keyword">inline</span> <span class="keyword">const</span> internal::solve_retval<SPQR, Rhs> <a class="code" href="classEigen_1_1SPQR.html#a7fcaf5446687d32b11ec4fb9fc5587e2">solve</a>(<span class="keyword">const</span> <a class="code" href="classEigen_1_1MatrixBase.html">MatrixBase<Rhs></a>& B)<span class="keyword"> const </span></div> <div class="line"><a name="l00128"></a><span class="lineno"> 128</span> <span class="keyword"> </span>{</div> <div class="line"><a name="l00129"></a><span class="lineno"> 129</span>  eigen_assert(m_isInitialized && <span class="stringliteral">" The QR factorization should be computed first, call compute()"</span>);</div> <div class="line"><a name="l00130"></a><span class="lineno"> 130</span>  eigen_assert(this-><a class="code" href="classEigen_1_1SPQR.html#a5552abd83dbd03c85cea6d61fd8875a5">rows</a>()==B.rows()</div> <div class="line"><a name="l00131"></a><span class="lineno"> 131</span>  && <span class="stringliteral">"SPQR::solve(): invalid number of rows of the right hand side matrix B"</span>);</div> <div class="line"><a name="l00132"></a><span class="lineno"> 132</span>  <span class="keywordflow">return</span> internal::solve_retval<SPQR, Rhs>(*<span class="keyword">this</span>, B.derived());</div> <div class="line"><a name="l00133"></a><span class="lineno"> 133</span>  }</div> <div class="line"><a name="l00134"></a><span class="lineno"> 134</span>  </div> <div class="line"><a name="l00135"></a><span class="lineno"> 135</span>  <span class="keyword">template</span><<span class="keyword">typename</span> Rhs, <span class="keyword">typename</span> Dest></div> <div class="line"><a name="l00136"></a><span class="lineno"> 136</span>  <span class="keywordtype">void</span> _solve(<span class="keyword">const</span> <a class="code" href="classEigen_1_1MatrixBase.html">MatrixBase<Rhs></a> &b, <a class="code" href="classEigen_1_1MatrixBase.html">MatrixBase<Dest></a> &dest)<span class="keyword"> const</span></div> <div class="line"><a name="l00137"></a><span class="lineno"> 137</span> <span class="keyword"> </span>{</div> <div class="line"><a name="l00138"></a><span class="lineno"> 138</span>  eigen_assert(m_isInitialized && <span class="stringliteral">" The QR factorization should be computed first, call compute()"</span>);</div> <div class="line"><a name="l00139"></a><span class="lineno"> 139</span>  eigen_assert(b.cols()==1 && <span class="stringliteral">"This method is for vectors only"</span>);</div> <div class="line"><a name="l00140"></a><span class="lineno"> 140</span>  </div> <div class="line"><a name="l00141"></a><span class="lineno"> 141</span>  <span class="comment">//Compute Q^T * b</span></div> <div class="line"><a name="l00142"></a><span class="lineno"> 142</span>  Dest y; </div> <div class="line"><a name="l00143"></a><span class="lineno"> 143</span>  y = <a class="code" href="classEigen_1_1SPQR.html#a2cca4c4270406a8d173bcc4339ef9394">matrixQ</a>().transpose() * b;</div> <div class="line"><a name="l00144"></a><span class="lineno"> 144</span>  <span class="comment">// Solves with the triangular matrix R</span></div> <div class="line"><a name="l00145"></a><span class="lineno"> 145</span>  Index rk = this-><a class="code" href="classEigen_1_1SPQR.html#a363d1c09d77f09d6ea2d2789776e7be3">rank</a>();</div> <div class="line"><a name="l00146"></a><span class="lineno"> 146</span>  y.topRows(rk) = this-><a class="code" href="classEigen_1_1SPQR.html#a9e36438e57c52afd02c489ab9833a440">matrixR</a>().<a class="code" href="classEigen_1_1SparseMatrixBase.html#a6f5fc5fe9d3fb70e62d4a9b1795704a8">topLeftCorner</a>(rk, rk).template triangularView<Upper>().<a class="code" href="classEigen_1_1SPQR.html#a7fcaf5446687d32b11ec4fb9fc5587e2">solve</a>(y.topRows(rk));</div> <div class="line"><a name="l00147"></a><span class="lineno"> 147</span>  y.bottomRows(<a class="code" href="classEigen_1_1SPQR.html#aaca1908a5ec508a25ff0a8bca803e5f3">cols</a>()-rk).setZero();</div> <div class="line"><a name="l00148"></a><span class="lineno"> 148</span>  <span class="comment">// Apply the column permutation </span></div> <div class="line"><a name="l00149"></a><span class="lineno"> 149</span>  dest.<a class="code" href="classEigen_1_1DenseBase.html#afdfced3e157f74b3bc68892aad178424">topRows</a>(<a class="code" href="classEigen_1_1SPQR.html#aaca1908a5ec508a25ff0a8bca803e5f3">cols</a>()) = <a class="code" href="classEigen_1_1SPQR.html#a60c6e65462cb1cab3a89d3afb37b35a0">colsPermutation</a>() * y.topRows(<a class="code" href="classEigen_1_1SPQR.html#aaca1908a5ec508a25ff0a8bca803e5f3">cols</a>());</div> <div class="line"><a name="l00150"></a><span class="lineno"> 150</span>  </div> <div class="line"><a name="l00151"></a><span class="lineno"> 151</span>  m_info = <a class="code" href="group__enums.html#gga51bc1ac16f26ebe51eae1abb77bd037bafdfbdf3247bd36a1f17270d5cec74c9c">Success</a>;</div> <div class="line"><a name="l00152"></a><span class="lineno"> 152</span>  }</div> <div class="line"><a name="l00153"></a><span class="lineno"> 153</span>  </div> <div class="line"><a name="l00156"></a><span class="lineno"><a class="line" href="classEigen_1_1SPQR.html#a9e36438e57c52afd02c489ab9833a440"> 156</a></span>  <span class="keyword">const</span> <a class="code" href="classEigen_1_1SparseMatrix.html">MatrixType</a> <a class="code" href="classEigen_1_1SPQR.html#a9e36438e57c52afd02c489ab9833a440">matrixR</a>()<span class="keyword"> const</span></div> <div class="line"><a name="l00157"></a><span class="lineno"> 157</span> <span class="keyword"> </span>{</div> <div class="line"><a name="l00158"></a><span class="lineno"> 158</span>  eigen_assert(m_isInitialized && <span class="stringliteral">" The QR factorization should be computed first, call compute()"</span>);</div> <div class="line"><a name="l00159"></a><span class="lineno"> 159</span>  <span class="keywordflow">if</span>(!m_isRUpToDate) {</div> <div class="line"><a name="l00160"></a><span class="lineno"> 160</span>  m_R = viewAsEigen<Scalar,ColMajor, typename MatrixType::Index>(*m_cR);</div> <div class="line"><a name="l00161"></a><span class="lineno"> 161</span>  m_isRUpToDate = <span class="keyword">true</span>;</div> <div class="line"><a name="l00162"></a><span class="lineno"> 162</span>  }</div> <div class="line"><a name="l00163"></a><span class="lineno"> 163</span>  <span class="keywordflow">return</span> m_R;</div> <div class="line"><a name="l00164"></a><span class="lineno"> 164</span>  }</div> <div class="line"><a name="l00166"></a><span class="lineno"><a class="line" href="classEigen_1_1SPQR.html#a2cca4c4270406a8d173bcc4339ef9394"> 166</a></span>  SPQRMatrixQReturnType<SPQR> <a class="code" href="classEigen_1_1SPQR.html#a2cca4c4270406a8d173bcc4339ef9394">matrixQ</a>()<span class="keyword"> const</span></div> <div class="line"><a name="l00167"></a><span class="lineno"> 167</span> <span class="keyword"> </span>{</div> <div class="line"><a name="l00168"></a><span class="lineno"> 168</span>  <span class="keywordflow">return</span> SPQRMatrixQReturnType<SPQR>(*this);</div> <div class="line"><a name="l00169"></a><span class="lineno"> 169</span>  }</div> <div class="line"><a name="l00171"></a><span class="lineno"><a class="line" href="classEigen_1_1SPQR.html#a60c6e65462cb1cab3a89d3afb37b35a0"> 171</a></span>  <a class="code" href="classEigen_1_1PermutationMatrix.html">PermutationType</a> <a class="code" href="classEigen_1_1SPQR.html#a60c6e65462cb1cab3a89d3afb37b35a0">colsPermutation</a>()<span class="keyword"> const</span></div> <div class="line"><a name="l00172"></a><span class="lineno"> 172</span> <span class="keyword"> </span>{ </div> <div class="line"><a name="l00173"></a><span class="lineno"> 173</span>  eigen_assert(m_isInitialized && <span class="stringliteral">"Decomposition is not initialized."</span>);</div> <div class="line"><a name="l00174"></a><span class="lineno"> 174</span>  Index n = m_cR->ncol;</div> <div class="line"><a name="l00175"></a><span class="lineno"> 175</span>  <a class="code" href="classEigen_1_1PermutationMatrix.html">PermutationType</a> colsPerm(n);</div> <div class="line"><a name="l00176"></a><span class="lineno"> 176</span>  <span class="keywordflow">for</span>(Index j = 0; j <n; j++) colsPerm.<a class="code" href="classEigen_1_1PermutationMatrix.html#a8876d615d17aad77b054a8f58b699e7d">indices</a>()(j) = m_E[j];</div> <div class="line"><a name="l00177"></a><span class="lineno"> 177</span>  <span class="keywordflow">return</span> colsPerm; </div> <div class="line"><a name="l00178"></a><span class="lineno"> 178</span>  </div> <div class="line"><a name="l00179"></a><span class="lineno"> 179</span>  }</div> <div class="line"><a name="l00184"></a><span class="lineno"><a class="line" href="classEigen_1_1SPQR.html#a363d1c09d77f09d6ea2d2789776e7be3"> 184</a></span>  Index <a class="code" href="classEigen_1_1SPQR.html#a363d1c09d77f09d6ea2d2789776e7be3">rank</a>()<span class="keyword"> const</span></div> <div class="line"><a name="l00185"></a><span class="lineno"> 185</span> <span class="keyword"> </span>{</div> <div class="line"><a name="l00186"></a><span class="lineno"> 186</span>  eigen_assert(m_isInitialized && <span class="stringliteral">"Decomposition is not initialized."</span>);</div> <div class="line"><a name="l00187"></a><span class="lineno"> 187</span>  <span class="keywordflow">return</span> m_cc.SPQR_istat[4];</div> <div class="line"><a name="l00188"></a><span class="lineno"> 188</span>  }</div> <div class="line"><a name="l00190"></a><span class="lineno"><a class="line" href="classEigen_1_1SPQR.html#a2a2c8ca72c8bbd609c162950128ac969"> 190</a></span>  <span class="keywordtype">void</span> <a class="code" href="classEigen_1_1SPQR.html#a2a2c8ca72c8bbd609c162950128ac969">setSPQROrdering</a>(<span class="keywordtype">int</span> ord) { m_ordering = ord;}</div> <div class="line"><a name="l00192"></a><span class="lineno"><a class="line" href="classEigen_1_1SPQR.html#a5d0d1c42d19d3e538a123813f7b62c5b"> 192</a></span>  <span class="keywordtype">void</span> <a class="code" href="classEigen_1_1SPQR.html#a5d0d1c42d19d3e538a123813f7b62c5b">setPivotThreshold</a>(<span class="keyword">const</span> RealScalar& tol) { m_tolerance = tol; }</div> <div class="line"><a name="l00193"></a><span class="lineno"> 193</span>  </div> <div class="line"><a name="l00195"></a><span class="lineno"><a class="line" href="classEigen_1_1SPQR.html#a6900999956a38bc60855289918a62d82"> 195</a></span>  cholmod_common *<a class="code" href="classEigen_1_1SPQR.html#a6900999956a38bc60855289918a62d82">cholmodCommon</a>()<span class="keyword"> const </span>{ <span class="keywordflow">return</span> &m_cc; }</div> <div class="line"><a name="l00196"></a><span class="lineno"> 196</span>  </div> <div class="line"><a name="l00197"></a><span class="lineno"> 197</span>  </div> <div class="line"><a name="l00203"></a><span class="lineno"><a class="line" href="classEigen_1_1SPQR.html#a0c06d5c2034ebb329c54235369643ad2"> 203</a></span>  <a class="code" href="group__enums.html#ga51bc1ac16f26ebe51eae1abb77bd037b">ComputationInfo</a> <a class="code" href="classEigen_1_1SPQR.html#a0c06d5c2034ebb329c54235369643ad2">info</a>()<span class="keyword"> const</span></div> <div class="line"><a name="l00204"></a><span class="lineno"> 204</span> <span class="keyword"> </span>{</div> <div class="line"><a name="l00205"></a><span class="lineno"> 205</span>  eigen_assert(m_isInitialized && <span class="stringliteral">"Decomposition is not initialized."</span>);</div> <div class="line"><a name="l00206"></a><span class="lineno"> 206</span>  <span class="keywordflow">return</span> m_info;</div> <div class="line"><a name="l00207"></a><span class="lineno"> 207</span>  }</div> <div class="line"><a name="l00208"></a><span class="lineno"> 208</span>  <span class="keyword">protected</span>:</div> <div class="line"><a name="l00209"></a><span class="lineno"> 209</span>  <span class="keywordtype">bool</span> m_isInitialized;</div> <div class="line"><a name="l00210"></a><span class="lineno"> 210</span>  <span class="keywordtype">bool</span> m_analysisIsOk;</div> <div class="line"><a name="l00211"></a><span class="lineno"> 211</span>  <span class="keywordtype">bool</span> m_factorizationIsOk;</div> <div class="line"><a name="l00212"></a><span class="lineno"> 212</span>  <span class="keyword">mutable</span> <span class="keywordtype">bool</span> m_isRUpToDate;</div> <div class="line"><a name="l00213"></a><span class="lineno"> 213</span>  <span class="keyword">mutable</span> <a class="code" href="group__enums.html#ga51bc1ac16f26ebe51eae1abb77bd037b">ComputationInfo</a> m_info;</div> <div class="line"><a name="l00214"></a><span class="lineno"> 214</span>  <span class="keywordtype">int</span> m_ordering; <span class="comment">// Ordering method to use, see SPQR's manual</span></div> <div class="line"><a name="l00215"></a><span class="lineno"> 215</span>  <span class="keywordtype">int</span> m_allow_tol; <span class="comment">// Allow to use some tolerance during numerical factorization.</span></div> <div class="line"><a name="l00216"></a><span class="lineno"> 216</span>  RealScalar m_tolerance; <span class="comment">// treat columns with 2-norm below this tolerance as zero</span></div> <div class="line"><a name="l00217"></a><span class="lineno"> 217</span>  <span class="keyword">mutable</span> cholmod_sparse *m_cR; <span class="comment">// The sparse R factor in cholmod format</span></div> <div class="line"><a name="l00218"></a><span class="lineno"> 218</span>  <span class="keyword">mutable</span> MatrixType m_R; <span class="comment">// The sparse matrix R in Eigen format</span></div> <div class="line"><a name="l00219"></a><span class="lineno"> 219</span>  <span class="keyword">mutable</span> Index *m_E; <span class="comment">// The permutation applied to columns</span></div> <div class="line"><a name="l00220"></a><span class="lineno"> 220</span>  <span class="keyword">mutable</span> cholmod_sparse *m_H; <span class="comment">//The householder vectors</span></div> <div class="line"><a name="l00221"></a><span class="lineno"> 221</span>  <span class="keyword">mutable</span> Index *m_HPinv; <span class="comment">// The row permutation of H</span></div> <div class="line"><a name="l00222"></a><span class="lineno"> 222</span>  <span class="keyword">mutable</span> cholmod_dense *m_HTau; <span class="comment">// The Householder coefficients</span></div> <div class="line"><a name="l00223"></a><span class="lineno"> 223</span>  <span class="keyword">mutable</span> Index m_rank; <span class="comment">// The rank of the matrix</span></div> <div class="line"><a name="l00224"></a><span class="lineno"> 224</span>  <span class="keyword">mutable</span> cholmod_common m_cc; <span class="comment">// Workspace and parameters</span></div> <div class="line"><a name="l00225"></a><span class="lineno"> 225</span>  <span class="keyword">template</span><<span class="keyword">typename</span> ,<span class="keyword">typename</span> > <span class="keyword">friend</span> <span class="keyword">struct </span>SPQR_QProduct;</div> <div class="line"><a name="l00226"></a><span class="lineno"> 226</span> };</div> <div class="line"><a name="l00227"></a><span class="lineno"> 227</span> </div> <div class="line"><a name="l00228"></a><span class="lineno"> 228</span> <span class="keyword">template</span> <<span class="keyword">typename</span> SPQRType, <span class="keyword">typename</span> Derived></div> <div class="line"><a name="l00229"></a><span class="lineno"> 229</span> <span class="keyword">struct </span>SPQR_QProduct : ReturnByValue<SPQR_QProduct<SPQRType,Derived> ></div> <div class="line"><a name="l00230"></a><span class="lineno"> 230</span> {</div> <div class="line"><a name="l00231"></a><span class="lineno"> 231</span>  <span class="keyword">typedef</span> <span class="keyword">typename</span> SPQRType::Scalar Scalar;</div> <div class="line"><a name="l00232"></a><span class="lineno"> 232</span>  <span class="keyword">typedef</span> <span class="keyword">typename</span> SPQRType::Index Index;</div> <div class="line"><a name="l00233"></a><span class="lineno"> 233</span>  <span class="comment">//Define the constructor to get reference to argument types</span></div> <div class="line"><a name="l00234"></a><span class="lineno"> 234</span>  SPQR_QProduct(<span class="keyword">const</span> SPQRType& spqr, <span class="keyword">const</span> Derived& other, <span class="keywordtype">bool</span> transpose) : m_spqr(spqr),m_other(other),m_transpose(transpose) {}</div> <div class="line"><a name="l00235"></a><span class="lineno"> 235</span>  </div> <div class="line"><a name="l00236"></a><span class="lineno"> 236</span>  <span class="keyword">inline</span> Index rows()<span class="keyword"> const </span>{ <span class="keywordflow">return</span> m_transpose ? m_spqr.rows() : m_spqr.cols(); }</div> <div class="line"><a name="l00237"></a><span class="lineno"> 237</span>  <span class="keyword">inline</span> Index cols()<span class="keyword"> const </span>{ <span class="keywordflow">return</span> m_other.cols(); }</div> <div class="line"><a name="l00238"></a><span class="lineno"> 238</span>  <span class="comment">// Assign to a vector</span></div> <div class="line"><a name="l00239"></a><span class="lineno"> 239</span>  <span class="keyword">template</span><<span class="keyword">typename</span> ResType></div> <div class="line"><a name="l00240"></a><span class="lineno"> 240</span>  <span class="keywordtype">void</span> evalTo(ResType& res)<span class="keyword"> const</span></div> <div class="line"><a name="l00241"></a><span class="lineno"> 241</span> <span class="keyword"> </span>{</div> <div class="line"><a name="l00242"></a><span class="lineno"> 242</span>  cholmod_dense y_cd;</div> <div class="line"><a name="l00243"></a><span class="lineno"> 243</span>  cholmod_dense *x_cd; </div> <div class="line"><a name="l00244"></a><span class="lineno"> 244</span>  <span class="keywordtype">int</span> method = m_transpose ? SPQR_QTX : SPQR_QX; </div> <div class="line"><a name="l00245"></a><span class="lineno"> 245</span>  cholmod_common *cc = m_spqr.cholmodCommon();</div> <div class="line"><a name="l00246"></a><span class="lineno"> 246</span>  y_cd = <a class="code" href="namespaceEigen.html#ab3b27d71e21e833c5287fcbe49509478">viewAsCholmod</a>(m_other.const_cast_derived());</div> <div class="line"><a name="l00247"></a><span class="lineno"> 247</span>  x_cd = SuiteSparseQR_qmult<Scalar>(method, m_spqr.m_H, m_spqr.m_HTau, m_spqr.m_HPinv, &y_cd, cc);</div> <div class="line"><a name="l00248"></a><span class="lineno"> 248</span>  res = Matrix<Scalar,ResType::RowsAtCompileTime,ResType::ColsAtCompileTime>::Map(reinterpret_cast<Scalar*>(x_cd->x), x_cd->nrow, x_cd->ncol);</div> <div class="line"><a name="l00249"></a><span class="lineno"> 249</span>  cholmod_l_free_dense(&x_cd, cc);</div> <div class="line"><a name="l00250"></a><span class="lineno"> 250</span>  }</div> <div class="line"><a name="l00251"></a><span class="lineno"> 251</span>  <span class="keyword">const</span> SPQRType& m_spqr; </div> <div class="line"><a name="l00252"></a><span class="lineno"> 252</span>  <span class="keyword">const</span> Derived& m_other; </div> <div class="line"><a name="l00253"></a><span class="lineno"> 253</span>  <span class="keywordtype">bool</span> m_transpose; </div> <div class="line"><a name="l00254"></a><span class="lineno"> 254</span>  </div> <div class="line"><a name="l00255"></a><span class="lineno"> 255</span> };</div> <div class="line"><a name="l00256"></a><span class="lineno"> 256</span> <span class="keyword">template</span><<span class="keyword">typename</span> SPQRType></div> <div class="line"><a name="l00257"></a><span class="lineno"> 257</span> <span class="keyword">struct </span>SPQRMatrixQReturnType{</div> <div class="line"><a name="l00258"></a><span class="lineno"> 258</span>  </div> <div class="line"><a name="l00259"></a><span class="lineno"> 259</span>  SPQRMatrixQReturnType(<span class="keyword">const</span> SPQRType& spqr) : m_spqr(spqr) {}</div> <div class="line"><a name="l00260"></a><span class="lineno"> 260</span>  <span class="keyword">template</span><<span class="keyword">typename</span> Derived></div> <div class="line"><a name="l00261"></a><span class="lineno"> 261</span>  SPQR_QProduct<SPQRType, Derived> <a class="code" href="namespaceEigen.html#a81fb70d0dc1c6deb42e9816647607247">operator*</a>(<span class="keyword">const</span> MatrixBase<Derived>& other)</div> <div class="line"><a name="l00262"></a><span class="lineno"> 262</span>  {</div> <div class="line"><a name="l00263"></a><span class="lineno"> 263</span>  <span class="keywordflow">return</span> SPQR_QProduct<SPQRType,Derived>(m_spqr,other.derived(),<span class="keyword">false</span>);</div> <div class="line"><a name="l00264"></a><span class="lineno"> 264</span>  }</div> <div class="line"><a name="l00265"></a><span class="lineno"> 265</span>  SPQRMatrixQTransposeReturnType<SPQRType> adjoint()<span class="keyword"> const</span></div> <div class="line"><a name="l00266"></a><span class="lineno"> 266</span> <span class="keyword"> </span>{</div> <div class="line"><a name="l00267"></a><span class="lineno"> 267</span>  <span class="keywordflow">return</span> SPQRMatrixQTransposeReturnType<SPQRType>(m_spqr);</div> <div class="line"><a name="l00268"></a><span class="lineno"> 268</span>  }</div> <div class="line"><a name="l00269"></a><span class="lineno"> 269</span>  <span class="comment">// To use for operations with the transpose of Q</span></div> <div class="line"><a name="l00270"></a><span class="lineno"> 270</span>  SPQRMatrixQTransposeReturnType<SPQRType> transpose()<span class="keyword"> const</span></div> <div class="line"><a name="l00271"></a><span class="lineno"> 271</span> <span class="keyword"> </span>{</div> <div class="line"><a name="l00272"></a><span class="lineno"> 272</span>  <span class="keywordflow">return</span> SPQRMatrixQTransposeReturnType<SPQRType>(m_spqr);</div> <div class="line"><a name="l00273"></a><span class="lineno"> 273</span>  }</div> <div class="line"><a name="l00274"></a><span class="lineno"> 274</span>  <span class="keyword">const</span> SPQRType& m_spqr;</div> <div class="line"><a name="l00275"></a><span class="lineno"> 275</span> };</div> <div class="line"><a name="l00276"></a><span class="lineno"> 276</span> </div> <div class="line"><a name="l00277"></a><span class="lineno"> 277</span> <span class="keyword">template</span><<span class="keyword">typename</span> SPQRType></div> <div class="line"><a name="l00278"></a><span class="lineno"> 278</span> <span class="keyword">struct </span>SPQRMatrixQTransposeReturnType{</div> <div class="line"><a name="l00279"></a><span class="lineno"> 279</span>  SPQRMatrixQTransposeReturnType(<span class="keyword">const</span> SPQRType& spqr) : m_spqr(spqr) {}</div> <div class="line"><a name="l00280"></a><span class="lineno"> 280</span>  <span class="keyword">template</span><<span class="keyword">typename</span> Derived></div> <div class="line"><a name="l00281"></a><span class="lineno"> 281</span>  SPQR_QProduct<SPQRType,Derived> <a class="code" href="namespaceEigen.html#a81fb70d0dc1c6deb42e9816647607247">operator*</a>(<span class="keyword">const</span> MatrixBase<Derived>& other)</div> <div class="line"><a name="l00282"></a><span class="lineno"> 282</span>  {</div> <div class="line"><a name="l00283"></a><span class="lineno"> 283</span>  <span class="keywordflow">return</span> SPQR_QProduct<SPQRType,Derived>(m_spqr,other.derived(), <span class="keyword">true</span>);</div> <div class="line"><a name="l00284"></a><span class="lineno"> 284</span>  }</div> <div class="line"><a name="l00285"></a><span class="lineno"> 285</span>  <span class="keyword">const</span> SPQRType& m_spqr;</div> <div class="line"><a name="l00286"></a><span class="lineno"> 286</span> };</div> <div class="line"><a name="l00287"></a><span class="lineno"> 287</span> </div> <div class="line"><a name="l00288"></a><span class="lineno"> 288</span> <span class="keyword">namespace </span>internal {</div> <div class="line"><a name="l00289"></a><span class="lineno"> 289</span>  </div> <div class="line"><a name="l00290"></a><span class="lineno"> 290</span> <span class="keyword">template</span><<span class="keyword">typename</span> _MatrixType, <span class="keyword">typename</span> Rhs></div> <div class="line"><a name="l00291"></a><span class="lineno"> 291</span> <span class="keyword">struct </span>solve_retval<SPQR<_MatrixType>, Rhs></div> <div class="line"><a name="l00292"></a><span class="lineno"> 292</span>  : solve_retval_base<SPQR<_MatrixType>, Rhs></div> <div class="line"><a name="l00293"></a><span class="lineno"> 293</span> {</div> <div class="line"><a name="l00294"></a><span class="lineno"> 294</span>  <span class="keyword">typedef</span> SPQR<_MatrixType> Dec;</div> <div class="line"><a name="l00295"></a><span class="lineno"> 295</span>  EIGEN_MAKE_SOLVE_HELPERS(Dec,Rhs)</div> <div class="line"><a name="l00296"></a><span class="lineno"> 296</span> </div> <div class="line"><a name="l00297"></a><span class="lineno"> 297</span>  template<typename Dest> <span class="keywordtype">void</span> evalTo(Dest& dst)<span class="keyword"> const</span></div> <div class="line"><a name="l00298"></a><span class="lineno"> 298</span> <span class="keyword"> </span>{</div> <div class="line"><a name="l00299"></a><span class="lineno"> 299</span>  dec()._solve(rhs(),dst);</div> <div class="line"><a name="l00300"></a><span class="lineno"> 300</span>  }</div> <div class="line"><a name="l00301"></a><span class="lineno"> 301</span> };</div> <div class="line"><a name="l00302"></a><span class="lineno"> 302</span> </div> <div class="line"><a name="l00303"></a><span class="lineno"> 303</span> } <span class="comment">// end namespace internal</span></div> <div class="line"><a name="l00304"></a><span class="lineno"> 304</span> </div> <div class="line"><a name="l00305"></a><span class="lineno"> 305</span> }<span class="comment">// End namespace Eigen</span></div> <div class="line"><a name="l00306"></a><span class="lineno"> 306</span> <span class="preprocessor">#endif</span></div> <div class="ttc" id="classEigen_1_1PermutationMatrix_html_a8876d615d17aad77b054a8f58b699e7d"><div class="ttname"><a href="classEigen_1_1PermutationMatrix.html#a8876d615d17aad77b054a8f58b699e7d">Eigen::PermutationMatrix::indices</a></div><div class="ttdeci">const IndicesType & indices() const </div><div class="ttdef"><b>Definition:</b> PermutationMatrix.h:358</div></div> <div class="ttc" id="classEigen_1_1SparseMatrix_html"><div class="ttname"><a href="classEigen_1_1SparseMatrix.html">Eigen::SparseMatrix< Scalar, ColMajor, Index ></a></div></div> <div class="ttc" id="classEigen_1_1SPQR_html_a7fcaf5446687d32b11ec4fb9fc5587e2"><div class="ttname"><a href="classEigen_1_1SPQR.html#a7fcaf5446687d32b11ec4fb9fc5587e2">Eigen::SPQR::solve</a></div><div class="ttdeci">const internal::solve_retval< SPQR, Rhs > solve(const MatrixBase< Rhs > &B) const </div><div class="ttdef"><b>Definition:</b> SuiteSparseQRSupport.h:127</div></div> <div class="ttc" id="group__enums_html_gga51bc1ac16f26ebe51eae1abb77bd037ba710fff14e8fc77846d4b75d8f4cc2d5c"><div class="ttname"><a href="group__enums.html#gga51bc1ac16f26ebe51eae1abb77bd037ba710fff14e8fc77846d4b75d8f4cc2d5c">Eigen::NumericalIssue</a></div><div class="ttdef"><b>Definition:</b> Constants.h:378</div></div> <div class="ttc" id="classEigen_1_1SPQR_html_a5552abd83dbd03c85cea6d61fd8875a5"><div class="ttname"><a href="classEigen_1_1SPQR.html#a5552abd83dbd03c85cea6d61fd8875a5">Eigen::SPQR::rows</a></div><div class="ttdeci">Index rows() const </div><div class="ttdef"><b>Definition:</b> SuiteSparseQRSupport.h:115</div></div> <div class="ttc" id="classEigen_1_1SPQR_html_a0c06d5c2034ebb329c54235369643ad2"><div class="ttname"><a href="classEigen_1_1SPQR.html#a0c06d5c2034ebb329c54235369643ad2">Eigen::SPQR::info</a></div><div class="ttdeci">ComputationInfo info() const </div><div class="ttdoc">Reports whether previous computation was successful. </div><div class="ttdef"><b>Definition:</b> SuiteSparseQRSupport.h:203</div></div> <div class="ttc" id="namespaceEigen_html_a81fb70d0dc1c6deb42e9816647607247"><div class="ttname"><a href="namespaceEigen.html#a81fb70d0dc1c6deb42e9816647607247">Eigen::operator*</a></div><div class="ttdeci">const internal::permut_matrix_product_retval< PermutationDerived, Derived, OnTheRight > operator*(const MatrixBase< Derived > &matrix, const PermutationBase< PermutationDerived > &permutation)</div><div class="ttdef"><b>Definition:</b> PermutationMatrix.h:510</div></div> <div class="ttc" id="classEigen_1_1SPQR_html_a2cca4c4270406a8d173bcc4339ef9394"><div class="ttname"><a href="classEigen_1_1SPQR.html#a2cca4c4270406a8d173bcc4339ef9394">Eigen::SPQR::matrixQ</a></div><div class="ttdeci">SPQRMatrixQReturnType< SPQR > matrixQ() const </div><div class="ttdoc">Get an expression of the matrix Q. </div><div class="ttdef"><b>Definition:</b> SuiteSparseQRSupport.h:166</div></div> <div class="ttc" id="classEigen_1_1SPQR_html_a6900999956a38bc60855289918a62d82"><div class="ttname"><a href="classEigen_1_1SPQR.html#a6900999956a38bc60855289918a62d82">Eigen::SPQR::cholmodCommon</a></div><div class="ttdeci">cholmod_common * cholmodCommon() const </div><div class="ttdef"><b>Definition:</b> SuiteSparseQRSupport.h:195</div></div> <div class="ttc" id="classEigen_1_1PermutationMatrix_html"><div class="ttname"><a href="classEigen_1_1PermutationMatrix.html">Eigen::PermutationMatrix</a></div><div class="ttdoc">Permutation matrix. </div><div class="ttdef"><b>Definition:</b> PermutationMatrix.h:283</div></div> <div class="ttc" id="namespaceEigen_html_ab3b27d71e21e833c5287fcbe49509478"><div class="ttname"><a href="namespaceEigen.html#ab3b27d71e21e833c5287fcbe49509478">Eigen::viewAsCholmod</a></div><div class="ttdeci">cholmod_sparse viewAsCholmod(SparseMatrix< _Scalar, _Options, _Index > &mat)</div><div class="ttdef"><b>Definition:</b> CholmodSupport.h:52</div></div> <div class="ttc" id="classEigen_1_1SPQR_html_a2a2c8ca72c8bbd609c162950128ac969"><div class="ttname"><a href="classEigen_1_1SPQR.html#a2a2c8ca72c8bbd609c162950128ac969">Eigen::SPQR::setSPQROrdering</a></div><div class="ttdeci">void setSPQROrdering(int ord)</div><div class="ttdoc">Set the fill-reducing ordering method to be used. </div><div class="ttdef"><b>Definition:</b> SuiteSparseQRSupport.h:190</div></div> <div class="ttc" id="classEigen_1_1SPQR_html_a5d0d1c42d19d3e538a123813f7b62c5b"><div class="ttname"><a href="classEigen_1_1SPQR.html#a5d0d1c42d19d3e538a123813f7b62c5b">Eigen::SPQR::setPivotThreshold</a></div><div class="ttdeci">void setPivotThreshold(const RealScalar &tol)</div><div class="ttdoc">Set the tolerance tol to treat columns with 2-norm &lt; =tol as zero. </div><div class="ttdef"><b>Definition:</b> SuiteSparseQRSupport.h:192</div></div> <div class="ttc" id="classEigen_1_1DenseBase_html_afdfced3e157f74b3bc68892aad178424"><div class="ttname"><a href="classEigen_1_1DenseBase.html#afdfced3e157f74b3bc68892aad178424">Eigen::DenseBase::topRows</a></div><div class="ttdeci">RowsBlockXpr topRows(Index n)</div><div class="ttdef"><b>Definition:</b> DenseBase.h:381</div></div> <div class="ttc" id="classEigen_1_1SPQR_html_aaca1908a5ec508a25ff0a8bca803e5f3"><div class="ttname"><a href="classEigen_1_1SPQR.html#aaca1908a5ec508a25ff0a8bca803e5f3">Eigen::SPQR::cols</a></div><div class="ttdeci">Index cols() const </div><div class="ttdef"><b>Definition:</b> SuiteSparseQRSupport.h:120</div></div> <div class="ttc" id="classEigen_1_1SPQR_html_a9e36438e57c52afd02c489ab9833a440"><div class="ttname"><a href="classEigen_1_1SPQR.html#a9e36438e57c52afd02c489ab9833a440">Eigen::SPQR::matrixR</a></div><div class="ttdeci">const MatrixType matrixR() const </div><div class="ttdef"><b>Definition:</b> SuiteSparseQRSupport.h:156</div></div> <div class="ttc" id="classEigen_1_1SPQR_html"><div class="ttname"><a href="classEigen_1_1SPQR.html">Eigen::SPQR</a></div><div class="ttdoc">Sparse QR factorization based on SuiteSparseQR library. </div><div class="ttdef"><b>Definition:</b> SuiteSparseQRSupport.h:15</div></div> <div class="ttc" id="classEigen_1_1SPQR_html_a363d1c09d77f09d6ea2d2789776e7be3"><div class="ttname"><a href="classEigen_1_1SPQR.html#a363d1c09d77f09d6ea2d2789776e7be3">Eigen::SPQR::rank</a></div><div class="ttdeci">Index rank() const </div><div class="ttdef"><b>Definition:</b> SuiteSparseQRSupport.h:184</div></div> <div class="ttc" id="group__enums_html_gga51bc1ac16f26ebe51eae1abb77bd037bafdfbdf3247bd36a1f17270d5cec74c9c"><div class="ttname"><a href="group__enums.html#gga51bc1ac16f26ebe51eae1abb77bd037bafdfbdf3247bd36a1f17270d5cec74c9c">Eigen::Success</a></div><div class="ttdef"><b>Definition:</b> Constants.h:376</div></div> <div class="ttc" id="classEigen_1_1SparseMatrixBase_html_a6f5fc5fe9d3fb70e62d4a9b1795704a8"><div class="ttname"><a href="classEigen_1_1SparseMatrixBase.html#a6f5fc5fe9d3fb70e62d4a9b1795704a8">Eigen::SparseMatrixBase::topLeftCorner</a></div><div class="ttdeci">Block< Derived > topLeftCorner(Index cRows, Index cCols)</div><div class="ttdef"><b>Definition:</b> SparseMatrixBase.h:157</div></div> <div class="ttc" id="group__enums_html_ga51bc1ac16f26ebe51eae1abb77bd037b"><div class="ttname"><a href="group__enums.html#ga51bc1ac16f26ebe51eae1abb77bd037b">Eigen::ComputationInfo</a></div><div class="ttdeci">ComputationInfo</div><div class="ttdef"><b>Definition:</b> Constants.h:374</div></div> <div class="ttc" id="classEigen_1_1MatrixBase_html"><div class="ttname"><a href="classEigen_1_1MatrixBase.html">Eigen::MatrixBase</a></div><div class="ttdoc">Base class for all dense matrices, vectors, and expressions. </div><div class="ttdef"><b>Definition:</b> MatrixBase.h:48</div></div> <div class="ttc" id="classEigen_1_1SPQR_html_a60c6e65462cb1cab3a89d3afb37b35a0"><div class="ttname"><a href="classEigen_1_1SPQR.html#a60c6e65462cb1cab3a89d3afb37b35a0">Eigen::SPQR::colsPermutation</a></div><div class="ttdeci">PermutationType colsPermutation() const </div><div class="ttdoc">Get the permutation that was applied to columns of A. </div><div class="ttdef"><b>Definition:</b> SuiteSparseQRSupport.h:171</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_257263196cadf2f5040b468f234e045f.html">SPQRSupport</a></li><li class="navelem"><b>SuiteSparseQRSupport.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>