<!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: SparseLU.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('SparseLU_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">SparseLU.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 Désiré Nuentsa-Wakam <desire.nuentsa_wakam@inria.fr></span></div> <div class="line"><a name="l00005"></a><span class="lineno"> 5</span> <span class="comment">// Copyright (C) 2012 Gael Guennebaud <gael.guennebaud@inria.fr></span></div> <div class="line"><a name="l00006"></a><span class="lineno"> 6</span> <span class="comment">//</span></div> <div class="line"><a name="l00007"></a><span class="lineno"> 7</span> <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> <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> <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> </div> <div class="line"><a name="l00011"></a><span class="lineno"> 11</span> </div> <div class="line"><a name="l00012"></a><span class="lineno"> 12</span> <span class="preprocessor">#ifndef EIGEN_SPARSE_LU_H</span></div> <div class="line"><a name="l00013"></a><span class="lineno"> 13</span> <span class="preprocessor"></span><span class="preprocessor">#define EIGEN_SPARSE_LU_H</span></div> <div class="line"><a name="l00014"></a><span class="lineno"> 14</span> <span class="preprocessor"></span></div> <div class="line"><a name="l00015"></a><span class="lineno"> 15</span> <span class="keyword">namespace </span>Eigen {</div> <div class="line"><a name="l00016"></a><span class="lineno"> 16</span> </div> <div class="line"><a name="l00017"></a><span class="lineno"><a class="line" href="classEigen_1_1SparseLU.html"> 17</a></span> <span class="keyword">template</span> <<span class="keyword">typename</span> _MatrixType, <span class="keyword">typename</span> _OrderingType = COLAMDOrdering<<span class="keyword">typename</span> _MatrixType::Index> > <span class="keyword">class </span><a class="code" href="classEigen_1_1SparseLU.html">SparseLU</a>;</div> <div class="line"><a name="l00018"></a><span class="lineno"> 18</span> <span class="keyword">template</span> <<span class="keyword">typename</span> MappedSparseMatrixType> <span class="keyword">struct </span>SparseLUMatrixLReturnType;</div> <div class="line"><a name="l00019"></a><span class="lineno"> 19</span> <span class="keyword">template</span> <<span class="keyword">typename</span> MatrixLType, <span class="keyword">typename</span> MatrixUType> <span class="keyword">struct </span>SparseLUMatrixUReturnType;</div> <div class="line"><a name="l00020"></a><span class="lineno"> 20</span> </div> <div class="line"><a name="l00072"></a><span class="lineno"> 72</span> <span class="keyword">template</span> <<span class="keyword">typename</span> _MatrixType, <span class="keyword">typename</span> _OrderingType></div> <div class="line"><a name="l00073"></a><span class="lineno"> 73</span> <span class="keyword">class </span><a class="code" href="classEigen_1_1SparseLU.html">SparseLU</a> : <span class="keyword">public</span> internal::SparseLUImpl<typename _MatrixType::Scalar, typename _MatrixType::Index></div> <div class="line"><a name="l00074"></a><span class="lineno"> 74</span> {</div> <div class="line"><a name="l00075"></a><span class="lineno"> 75</span>  <span class="keyword">public</span>:</div> <div class="line"><a name="l00076"></a><span class="lineno"> 76</span>  <span class="keyword">typedef</span> _MatrixType MatrixType; </div> <div class="line"><a name="l00077"></a><span class="lineno"> 77</span>  <span class="keyword">typedef</span> _OrderingType OrderingType;</div> <div class="line"><a name="l00078"></a><span class="lineno"> 78</span>  <span class="keyword">typedef</span> <span class="keyword">typename</span> MatrixType::Scalar Scalar; </div> <div class="line"><a name="l00079"></a><span class="lineno"> 79</span>  <span class="keyword">typedef</span> <span class="keyword">typename</span> MatrixType::RealScalar RealScalar; </div> <div class="line"><a name="l00080"></a><span class="lineno"> 80</span>  <span class="keyword">typedef</span> <span class="keyword">typename</span> MatrixType::Index Index; </div> <div class="line"><a name="l00081"></a><span class="lineno"> 81</span>  <span class="keyword">typedef</span> <a class="code" href="classEigen_1_1SparseMatrix.html">SparseMatrix<Scalar,ColMajor,Index></a> <a class="code" href="classEigen_1_1SparseMatrix.html">NCMatrix</a>;</div> <div class="line"><a name="l00082"></a><span class="lineno"> 82</span>  <span class="keyword">typedef</span> internal::MappedSuperNodalMatrix<Scalar, Index> SCMatrix; </div> <div class="line"><a name="l00083"></a><span class="lineno"> 83</span>  <span class="keyword">typedef</span> <a class="code" href="classEigen_1_1Matrix.html">Matrix<Scalar,Dynamic,1></a> <a class="code" href="classEigen_1_1Matrix.html">ScalarVector</a>;</div> <div class="line"><a name="l00084"></a><span class="lineno"> 84</span>  <span class="keyword">typedef</span> <a class="code" href="classEigen_1_1Matrix.html">Matrix<Index,Dynamic,1></a> <a class="code" href="classEigen_1_1Matrix.html">IndexVector</a>;</div> <div class="line"><a name="l00085"></a><span class="lineno"> 85</span>  <span class="keyword">typedef</span> <a class="code" href="classEigen_1_1PermutationMatrix.html">PermutationMatrix<Dynamic, Dynamic, Index></a> <a class="code" href="classEigen_1_1PermutationMatrix.html">PermutationType</a>;</div> <div class="line"><a name="l00086"></a><span class="lineno"> 86</span>  <span class="keyword">typedef</span> internal::SparseLUImpl<Scalar, Index> Base;</div> <div class="line"><a name="l00087"></a><span class="lineno"> 87</span>  </div> <div class="line"><a name="l00088"></a><span class="lineno"> 88</span>  <span class="keyword">public</span>:</div> <div class="line"><a name="l00089"></a><span class="lineno"> 89</span>  <a class="code" href="classEigen_1_1SparseLU.html">SparseLU</a>():m_isInitialized(true),m_lastError(<span class="stringliteral">""</span>),m_Ustore(0,0,0,0,0,0),m_symmetricmode(false),m_diagpivotthresh(1.0),m_detPermR(1)</div> <div class="line"><a name="l00090"></a><span class="lineno"> 90</span>  {</div> <div class="line"><a name="l00091"></a><span class="lineno"> 91</span>  initperfvalues(); </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>  <a class="code" href="classEigen_1_1SparseLU.html">SparseLU</a>(<span class="keyword">const</span> MatrixType& matrix):m_isInitialized(true),m_lastError(<span class="stringliteral">""</span>),m_Ustore(0,0,0,0,0,0),m_symmetricmode(false),m_diagpivotthresh(1.0),m_detPermR(1)</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>  initperfvalues(); </div> <div class="line"><a name="l00096"></a><span class="lineno"> 96</span>  <a class="code" href="classEigen_1_1SparseLU.html#a05a0eac1c146e9e62df9ec70f5a6f69d">compute</a>(matrix);</div> <div class="line"><a name="l00097"></a><span class="lineno"> 97</span>  }</div> <div class="line"><a name="l00098"></a><span class="lineno"> 98</span>  </div> <div class="line"><a name="l00099"></a><span class="lineno"> 99</span>  ~SparseLU()</div> <div class="line"><a name="l00100"></a><span class="lineno"> 100</span>  {</div> <div class="line"><a name="l00101"></a><span class="lineno"> 101</span>  <span class="comment">// Free all explicit dynamic pointers </span></div> <div class="line"><a name="l00102"></a><span class="lineno"> 102</span>  }</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>  <span class="keywordtype">void</span> <a class="code" href="classEigen_1_1SparseLU.html#a6f2a135bd741b6b2f2558e6a404581ff">analyzePattern</a> (<span class="keyword">const</span> MatrixType& matrix);</div> <div class="line"><a name="l00105"></a><span class="lineno"> 105</span>  <span class="keywordtype">void</span> <a class="code" href="classEigen_1_1SparseLU.html#a33e2421e468033a883c041f940537a7c">factorize</a> (<span class="keyword">const</span> MatrixType& matrix);</div> <div class="line"><a name="l00106"></a><span class="lineno"> 106</span>  <span class="keywordtype">void</span> simplicialfactorize(<span class="keyword">const</span> MatrixType& matrix);</div> <div class="line"><a name="l00107"></a><span class="lineno"> 107</span>  </div> <div class="line"><a name="l00112"></a><span class="lineno"><a class="line" href="classEigen_1_1SparseLU.html#a05a0eac1c146e9e62df9ec70f5a6f69d"> 112</a></span>  <span class="keywordtype">void</span> <a class="code" href="classEigen_1_1SparseLU.html#a05a0eac1c146e9e62df9ec70f5a6f69d">compute</a> (<span class="keyword">const</span> MatrixType& matrix)</div> <div class="line"><a name="l00113"></a><span class="lineno"> 113</span>  {</div> <div class="line"><a name="l00114"></a><span class="lineno"> 114</span>  <span class="comment">// Analyze </span></div> <div class="line"><a name="l00115"></a><span class="lineno"> 115</span>  <a class="code" href="classEigen_1_1SparseLU.html#a6f2a135bd741b6b2f2558e6a404581ff">analyzePattern</a>(matrix); </div> <div class="line"><a name="l00116"></a><span class="lineno"> 116</span>  <span class="comment">//Factorize</span></div> <div class="line"><a name="l00117"></a><span class="lineno"> 117</span>  <a class="code" href="classEigen_1_1SparseLU.html#a33e2421e468033a883c041f940537a7c">factorize</a>(matrix);</div> <div class="line"><a name="l00118"></a><span class="lineno"> 118</span>  } </div> <div class="line"><a name="l00119"></a><span class="lineno"> 119</span>  </div> <div class="line"><a name="l00120"></a><span class="lineno"> 120</span>  <span class="keyword">inline</span> Index rows()<span class="keyword"> const </span>{ <span class="keywordflow">return</span> m_mat.<a class="code" href="classEigen_1_1SparseMatrix.html#a5552abd83dbd03c85cea6d61fd8875a5">rows</a>(); }</div> <div class="line"><a name="l00121"></a><span class="lineno"> 121</span>  <span class="keyword">inline</span> Index cols()<span class="keyword"> const </span>{ <span class="keywordflow">return</span> m_mat.<a class="code" href="classEigen_1_1SparseMatrix.html#aaca1908a5ec508a25ff0a8bca803e5f3">cols</a>(); }</div> <div class="line"><a name="l00123"></a><span class="lineno"><a class="line" href="classEigen_1_1SparseLU.html#acef7759add41200c1c817a6255e21dd4"> 123</a></span>  <span class="keywordtype">void</span> <a class="code" href="classEigen_1_1SparseLU.html#acef7759add41200c1c817a6255e21dd4">isSymmetric</a>(<span class="keywordtype">bool</span> sym)</div> <div class="line"><a name="l00124"></a><span class="lineno"> 124</span>  {</div> <div class="line"><a name="l00125"></a><span class="lineno"> 125</span>  m_symmetricmode = sym;</div> <div class="line"><a name="l00126"></a><span class="lineno"> 126</span>  }</div> <div class="line"><a name="l00127"></a><span class="lineno"> 127</span>  </div> <div class="line"><a name="l00134"></a><span class="lineno"><a class="line" href="classEigen_1_1SparseLU.html#a828ca6373b429d81f41f41429a7fd4ad"> 134</a></span>  SparseLUMatrixLReturnType<SCMatrix> <a class="code" href="classEigen_1_1SparseLU.html#a828ca6373b429d81f41f41429a7fd4ad">matrixL</a>()<span class="keyword"> const</span></div> <div class="line"><a name="l00135"></a><span class="lineno"> 135</span> <span class="keyword"> </span>{</div> <div class="line"><a name="l00136"></a><span class="lineno"> 136</span>  <span class="keywordflow">return</span> SparseLUMatrixLReturnType<SCMatrix>(m_Lstore);</div> <div class="line"><a name="l00137"></a><span class="lineno"> 137</span>  }</div> <div class="line"><a name="l00144"></a><span class="lineno"><a class="line" href="classEigen_1_1SparseLU.html#a6001a3f3f7ad66583d9d1ce3bb3de262"> 144</a></span>  SparseLUMatrixUReturnType<SCMatrix,MappedSparseMatrix<Scalar,ColMajor,Index> > <a class="code" href="classEigen_1_1SparseLU.html#a6001a3f3f7ad66583d9d1ce3bb3de262">matrixU</a>()<span class="keyword"> const</span></div> <div class="line"><a name="l00145"></a><span class="lineno"> 145</span> <span class="keyword"> </span>{</div> <div class="line"><a name="l00146"></a><span class="lineno"> 146</span>  <span class="keywordflow">return</span> SparseLUMatrixUReturnType<SCMatrix, MappedSparseMatrix<Scalar,ColMajor,Index> >(m_Lstore, m_Ustore);</div> <div class="line"><a name="l00147"></a><span class="lineno"> 147</span>  }</div> <div class="line"><a name="l00148"></a><span class="lineno"> 148</span> </div> <div class="line"><a name="l00153"></a><span class="lineno"><a class="line" href="classEigen_1_1SparseLU.html#ad856e0133763838bd626bfcd0cd1bf33"> 153</a></span>  <span class="keyword">inline</span> <span class="keyword">const</span> <a class="code" href="classEigen_1_1PermutationMatrix.html">PermutationType</a>& <a class="code" href="classEigen_1_1SparseLU.html#ad856e0133763838bd626bfcd0cd1bf33">rowsPermutation</a>()<span class="keyword"> const</span></div> <div class="line"><a name="l00154"></a><span class="lineno"> 154</span> <span class="keyword"> </span>{</div> <div class="line"><a name="l00155"></a><span class="lineno"> 155</span>  <span class="keywordflow">return</span> m_perm_r;</div> <div class="line"><a name="l00156"></a><span class="lineno"> 156</span>  }</div> <div class="line"><a name="l00161"></a><span class="lineno"><a class="line" href="classEigen_1_1SparseLU.html#a1229744f4e8554ca6e96fe32ac359924"> 161</a></span>  <span class="keyword">inline</span> <span class="keyword">const</span> <a class="code" href="classEigen_1_1PermutationMatrix.html">PermutationType</a>& <a class="code" href="classEigen_1_1SparseLU.html#a1229744f4e8554ca6e96fe32ac359924">colsPermutation</a>()<span class="keyword"> const</span></div> <div class="line"><a name="l00162"></a><span class="lineno"> 162</span> <span class="keyword"> </span>{</div> <div class="line"><a name="l00163"></a><span class="lineno"> 163</span>  <span class="keywordflow">return</span> m_perm_c;</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_1SparseLU.html#a20467be6f23d8b39574f3e097583d767"> 166</a></span>  <span class="keywordtype">void</span> <a class="code" href="classEigen_1_1SparseLU.html#a20467be6f23d8b39574f3e097583d767">setPivotThreshold</a>(<span class="keyword">const</span> RealScalar& thresh)</div> <div class="line"><a name="l00167"></a><span class="lineno"> 167</span>  {</div> <div class="line"><a name="l00168"></a><span class="lineno"> 168</span>  m_diagpivotthresh = thresh; </div> <div class="line"><a name="l00169"></a><span class="lineno"> 169</span>  }</div> <div class="line"><a name="l00170"></a><span class="lineno"> 170</span> </div> <div class="line"><a name="l00177"></a><span class="lineno"> 177</span>  <span class="keyword">template</span><<span class="keyword">typename</span> Rhs></div> <div class="line"><a name="l00178"></a><span class="lineno"><a class="line" href="classEigen_1_1SparseLU.html#a8166735beae3bd23e7c3a1be9f0f5ef9"> 178</a></span>  <span class="keyword">inline</span> <span class="keyword">const</span> internal::solve_retval<SparseLU, Rhs> <a class="code" href="classEigen_1_1SparseLU.html#a8166735beae3bd23e7c3a1be9f0f5ef9">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="l00179"></a><span class="lineno"> 179</span> <span class="keyword"> </span>{</div> <div class="line"><a name="l00180"></a><span class="lineno"> 180</span>  eigen_assert(m_factorizationIsOk && <span class="stringliteral">"SparseLU is not initialized."</span>); </div> <div class="line"><a name="l00181"></a><span class="lineno"> 181</span>  eigen_assert(rows()==B.rows()</div> <div class="line"><a name="l00182"></a><span class="lineno"> 182</span>  && <span class="stringliteral">"SparseLU::solve(): invalid number of rows of the right hand side matrix B"</span>);</div> <div class="line"><a name="l00183"></a><span class="lineno"> 183</span>  <span class="keywordflow">return</span> internal::solve_retval<SparseLU, Rhs>(*<span class="keyword">this</span>, B.derived());</div> <div class="line"><a name="l00184"></a><span class="lineno"> 184</span>  }</div> <div class="line"><a name="l00185"></a><span class="lineno"> 185</span> </div> <div class="line"><a name="l00190"></a><span class="lineno"> 190</span>  <span class="keyword">template</span><<span class="keyword">typename</span> Rhs></div> <div class="line"><a name="l00191"></a><span class="lineno"><a class="line" href="classEigen_1_1SparseLU.html#a39062d90ea647b4e7668452647b04b17"> 191</a></span>  <span class="keyword">inline</span> <span class="keyword">const</span> internal::sparse_solve_retval<SparseLU, Rhs> <a class="code" href="classEigen_1_1SparseLU.html#a39062d90ea647b4e7668452647b04b17">solve</a>(<span class="keyword">const</span> <a class="code" href="classEigen_1_1SparseMatrixBase.html">SparseMatrixBase<Rhs></a>& B)<span class="keyword"> const </span></div> <div class="line"><a name="l00192"></a><span class="lineno"> 192</span> <span class="keyword"> </span>{</div> <div class="line"><a name="l00193"></a><span class="lineno"> 193</span>  eigen_assert(m_factorizationIsOk && <span class="stringliteral">"SparseLU is not initialized."</span>); </div> <div class="line"><a name="l00194"></a><span class="lineno"> 194</span>  eigen_assert(rows()==B.<a class="code" href="classEigen_1_1SparseMatrixBase.html#a5552abd83dbd03c85cea6d61fd8875a5">rows</a>()</div> <div class="line"><a name="l00195"></a><span class="lineno"> 195</span>  && <span class="stringliteral">"SparseLU::solve(): invalid number of rows of the right hand side matrix B"</span>);</div> <div class="line"><a name="l00196"></a><span class="lineno"> 196</span>  <span class="keywordflow">return</span> internal::sparse_solve_retval<SparseLU, Rhs>(*<span class="keyword">this</span>, B.<a class="code" href="structEigen_1_1EigenBase.html#aa84222add803ad7c9db07dd4dd91d5d9">derived</a>());</div> <div class="line"><a name="l00197"></a><span class="lineno"> 197</span>  }</div> <div class="line"><a name="l00198"></a><span class="lineno"> 198</span>  </div> <div class="line"><a name="l00207"></a><span class="lineno"><a class="line" href="classEigen_1_1SparseLU.html#a0c06d5c2034ebb329c54235369643ad2"> 207</a></span>  <a class="code" href="group__enums.html#ga51bc1ac16f26ebe51eae1abb77bd037b">ComputationInfo</a> <a class="code" href="classEigen_1_1SparseLU.html#a0c06d5c2034ebb329c54235369643ad2">info</a>()<span class="keyword"> const</span></div> <div class="line"><a name="l00208"></a><span class="lineno"> 208</span> <span class="keyword"> </span>{</div> <div class="line"><a name="l00209"></a><span class="lineno"> 209</span>  eigen_assert(m_isInitialized && <span class="stringliteral">"Decomposition is not initialized."</span>);</div> <div class="line"><a name="l00210"></a><span class="lineno"> 210</span>  <span class="keywordflow">return</span> m_info;</div> <div class="line"><a name="l00211"></a><span class="lineno"> 211</span>  }</div> <div class="line"><a name="l00212"></a><span class="lineno"> 212</span>  </div> <div class="line"><a name="l00216"></a><span class="lineno"><a class="line" href="classEigen_1_1SparseLU.html#a040b16815cde46c7f2f6d26a2298ca4b"> 216</a></span>  std::string <a class="code" href="classEigen_1_1SparseLU.html#a040b16815cde46c7f2f6d26a2298ca4b">lastErrorMessage</a>()<span class="keyword"> const</span></div> <div class="line"><a name="l00217"></a><span class="lineno"> 217</span> <span class="keyword"> </span>{</div> <div class="line"><a name="l00218"></a><span class="lineno"> 218</span>  <span class="keywordflow">return</span> m_lastError; </div> <div class="line"><a name="l00219"></a><span class="lineno"> 219</span>  }</div> <div class="line"><a name="l00220"></a><span class="lineno"> 220</span> </div> <div class="line"><a name="l00221"></a><span class="lineno"> 221</span>  <span class="keyword">template</span><<span class="keyword">typename</span> Rhs, <span class="keyword">typename</span> Dest></div> <div class="line"><a name="l00222"></a><span class="lineno"> 222</span>  <span class="keywordtype">bool</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> &_X)<span class="keyword"> const</span></div> <div class="line"><a name="l00223"></a><span class="lineno"> 223</span> <span class="keyword"> </span>{</div> <div class="line"><a name="l00224"></a><span class="lineno"> 224</span>  Dest& X(_X.derived());</div> <div class="line"><a name="l00225"></a><span class="lineno"> 225</span>  eigen_assert(m_factorizationIsOk && <span class="stringliteral">"The matrix should be factorized first"</span>);</div> <div class="line"><a name="l00226"></a><span class="lineno"> 226</span>  EIGEN_STATIC_ASSERT((Dest::Flags&<a class="code" href="group__flags.html#ga7bd49e7b260e869e10fb9dc4fd081a85">RowMajorBit</a>)==0,</div> <div class="line"><a name="l00227"></a><span class="lineno"> 227</span>  THIS_METHOD_IS_ONLY_FOR_COLUMN_MAJOR_MATRICES);</div> <div class="line"><a name="l00228"></a><span class="lineno"> 228</span>  </div> <div class="line"><a name="l00229"></a><span class="lineno"> 229</span>  <span class="comment">// Permute the right hand side to form X = Pr*B</span></div> <div class="line"><a name="l00230"></a><span class="lineno"> 230</span>  <span class="comment">// on return, X is overwritten by the computed solution</span></div> <div class="line"><a name="l00231"></a><span class="lineno"> 231</span>  X.resize(B.rows(),B.cols());</div> <div class="line"><a name="l00232"></a><span class="lineno"> 232</span>  <span class="keywordflow">for</span>(Index j = 0; j < B.cols(); ++j)</div> <div class="line"><a name="l00233"></a><span class="lineno"> 233</span>  X.col(j) = <a class="code" href="classEigen_1_1SparseLU.html#ad856e0133763838bd626bfcd0cd1bf33">rowsPermutation</a>() * B.<a class="code" href="classEigen_1_1DenseBase.html#a58c77695de3b33405f01f2fdf3dc389d">col</a>(j);</div> <div class="line"><a name="l00234"></a><span class="lineno"> 234</span>  </div> <div class="line"><a name="l00235"></a><span class="lineno"> 235</span>  <span class="comment">//Forward substitution with L</span></div> <div class="line"><a name="l00236"></a><span class="lineno"> 236</span>  this-><a class="code" href="classEigen_1_1SparseLU.html#a828ca6373b429d81f41f41429a7fd4ad">matrixL</a>().solveInPlace(X);</div> <div class="line"><a name="l00237"></a><span class="lineno"> 237</span>  this-><a class="code" href="classEigen_1_1SparseLU.html#a6001a3f3f7ad66583d9d1ce3bb3de262">matrixU</a>().solveInPlace(X);</div> <div class="line"><a name="l00238"></a><span class="lineno"> 238</span>  </div> <div class="line"><a name="l00239"></a><span class="lineno"> 239</span>  <span class="comment">// Permute back the solution </span></div> <div class="line"><a name="l00240"></a><span class="lineno"> 240</span>  <span class="keywordflow">for</span> (Index j = 0; j < B.cols(); ++j)</div> <div class="line"><a name="l00241"></a><span class="lineno"> 241</span>  X.col(j) = <a class="code" href="classEigen_1_1SparseLU.html#a1229744f4e8554ca6e96fe32ac359924">colsPermutation</a>().<a class="code" href="classEigen_1_1PermutationBase.html#a594deabd0f368a746801f5dd14a0db2a">inverse</a>() * X.col(j);</div> <div class="line"><a name="l00242"></a><span class="lineno"> 242</span>  </div> <div class="line"><a name="l00243"></a><span class="lineno"> 243</span>  <span class="keywordflow">return</span> <span class="keyword">true</span>; </div> <div class="line"><a name="l00244"></a><span class="lineno"> 244</span>  }</div> <div class="line"><a name="l00245"></a><span class="lineno"> 245</span>  </div> <div class="line"><a name="l00256"></a><span class="lineno"><a class="line" href="classEigen_1_1SparseLU.html#a4c6b35f80ff9f6f0a3f4e74dfc121349"> 256</a></span>  Scalar <a class="code" href="classEigen_1_1SparseLU.html#a4c6b35f80ff9f6f0a3f4e74dfc121349">absDeterminant</a>()</div> <div class="line"><a name="l00257"></a><span class="lineno"> 257</span>  {</div> <div class="line"><a name="l00258"></a><span class="lineno"> 258</span>  eigen_assert(m_factorizationIsOk && <span class="stringliteral">"The matrix should be factorized first."</span>);</div> <div class="line"><a name="l00259"></a><span class="lineno"> 259</span>  <span class="comment">// Initialize with the determinant of the row matrix</span></div> <div class="line"><a name="l00260"></a><span class="lineno"> 260</span>  Scalar det = Scalar(1.);</div> <div class="line"><a name="l00261"></a><span class="lineno"> 261</span>  <span class="comment">//Note that the diagonal blocks of U are stored in supernodes,</span></div> <div class="line"><a name="l00262"></a><span class="lineno"> 262</span>  <span class="comment">// which are available in the L part :)</span></div> <div class="line"><a name="l00263"></a><span class="lineno"> 263</span>  <span class="keywordflow">for</span> (Index j = 0; j < this->cols(); ++j)</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>  <span class="keywordflow">for</span> (<span class="keyword">typename</span> SCMatrix::InnerIterator it(m_Lstore, j); it; ++it)</div> <div class="line"><a name="l00266"></a><span class="lineno"> 266</span>  {</div> <div class="line"><a name="l00267"></a><span class="lineno"> 267</span>  <span class="keywordflow">if</span>(it.row() < j) <span class="keywordflow">continue</span>;</div> <div class="line"><a name="l00268"></a><span class="lineno"> 268</span>  <span class="keywordflow">if</span>(it.row() == j)</div> <div class="line"><a name="l00269"></a><span class="lineno"> 269</span>  {</div> <div class="line"><a name="l00270"></a><span class="lineno"> 270</span>  det *= (std::abs)(it.value());</div> <div class="line"><a name="l00271"></a><span class="lineno"> 271</span>  <span class="keywordflow">break</span>;</div> <div class="line"><a name="l00272"></a><span class="lineno"> 272</span>  }</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>  }</div> <div class="line"><a name="l00275"></a><span class="lineno"> 275</span>  <span class="keywordflow">return</span> det;</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> </div> <div class="line"><a name="l00286"></a><span class="lineno"><a class="line" href="classEigen_1_1SparseLU.html#adb61fce0a1343a45b65b6a50b18408dd"> 286</a></span>  Scalar <a class="code" href="classEigen_1_1SparseLU.html#adb61fce0a1343a45b65b6a50b18408dd">logAbsDeterminant</a>()<span class="keyword"> const</span></div> <div class="line"><a name="l00287"></a><span class="lineno"> 287</span> <span class="keyword"> </span>{</div> <div class="line"><a name="l00288"></a><span class="lineno"> 288</span>  eigen_assert(m_factorizationIsOk && <span class="stringliteral">"The matrix should be factorized first."</span>);</div> <div class="line"><a name="l00289"></a><span class="lineno"> 289</span>  Scalar det = Scalar(0.);</div> <div class="line"><a name="l00290"></a><span class="lineno"> 290</span>  <span class="keywordflow">for</span> (Index j = 0; j < this->cols(); ++j)</div> <div class="line"><a name="l00291"></a><span class="lineno"> 291</span>  {</div> <div class="line"><a name="l00292"></a><span class="lineno"> 292</span>  <span class="keywordflow">for</span> (<span class="keyword">typename</span> SCMatrix::InnerIterator it(m_Lstore, j); it; ++it)</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="keywordflow">if</span>(it.row() < j) <span class="keywordflow">continue</span>;</div> <div class="line"><a name="l00295"></a><span class="lineno"> 295</span>  <span class="keywordflow">if</span>(it.row() == j)</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>  det += (std::log)((std::abs)(it.value()));</div> <div class="line"><a name="l00298"></a><span class="lineno"> 298</span>  <span class="keywordflow">break</span>;</div> <div class="line"><a name="l00299"></a><span class="lineno"> 299</span>  }</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>  <span class="keywordflow">return</span> det;</div> <div class="line"><a name="l00303"></a><span class="lineno"> 303</span>  }</div> <div class="line"><a name="l00304"></a><span class="lineno"> 304</span> </div> <div class="line"><a name="l00309"></a><span class="lineno"><a class="line" href="classEigen_1_1SparseLU.html#a7150c533f7f2f6bdb1e9c82727f426f5"> 309</a></span>  Scalar <a class="code" href="classEigen_1_1SparseLU.html#a7150c533f7f2f6bdb1e9c82727f426f5">signDeterminant</a>()</div> <div class="line"><a name="l00310"></a><span class="lineno"> 310</span>  {</div> <div class="line"><a name="l00311"></a><span class="lineno"> 311</span>  eigen_assert(m_factorizationIsOk && <span class="stringliteral">"The matrix should be factorized first."</span>);</div> <div class="line"><a name="l00312"></a><span class="lineno"> 312</span>  <span class="keywordflow">return</span> Scalar(m_detPermR);</div> <div class="line"><a name="l00313"></a><span class="lineno"> 313</span>  }</div> <div class="line"><a name="l00314"></a><span class="lineno"> 314</span> </div> <div class="line"><a name="l00315"></a><span class="lineno"> 315</span>  <span class="keyword">protected</span>:</div> <div class="line"><a name="l00316"></a><span class="lineno"> 316</span>  <span class="comment">// Functions </span></div> <div class="line"><a name="l00317"></a><span class="lineno"> 317</span>  <span class="keywordtype">void</span> initperfvalues()</div> <div class="line"><a name="l00318"></a><span class="lineno"> 318</span>  {</div> <div class="line"><a name="l00319"></a><span class="lineno"> 319</span>  m_perfv.panel_size = 1;</div> <div class="line"><a name="l00320"></a><span class="lineno"> 320</span>  m_perfv.relax = 1; </div> <div class="line"><a name="l00321"></a><span class="lineno"> 321</span>  m_perfv.maxsuper = 128; </div> <div class="line"><a name="l00322"></a><span class="lineno"> 322</span>  m_perfv.rowblk = 16; </div> <div class="line"><a name="l00323"></a><span class="lineno"> 323</span>  m_perfv.colblk = 8; </div> <div class="line"><a name="l00324"></a><span class="lineno"> 324</span>  m_perfv.fillfactor = 20; </div> <div class="line"><a name="l00325"></a><span class="lineno"> 325</span>  }</div> <div class="line"><a name="l00326"></a><span class="lineno"> 326</span>  </div> <div class="line"><a name="l00327"></a><span class="lineno"> 327</span>  <span class="comment">// Variables </span></div> <div class="line"><a name="l00328"></a><span class="lineno"> 328</span>  <span class="keyword">mutable</span> <a class="code" href="group__enums.html#ga51bc1ac16f26ebe51eae1abb77bd037b">ComputationInfo</a> m_info;</div> <div class="line"><a name="l00329"></a><span class="lineno"> 329</span>  <span class="keywordtype">bool</span> m_isInitialized;</div> <div class="line"><a name="l00330"></a><span class="lineno"> 330</span>  <span class="keywordtype">bool</span> m_factorizationIsOk;</div> <div class="line"><a name="l00331"></a><span class="lineno"> 331</span>  <span class="keywordtype">bool</span> m_analysisIsOk;</div> <div class="line"><a name="l00332"></a><span class="lineno"> 332</span>  std::string m_lastError;</div> <div class="line"><a name="l00333"></a><span class="lineno"> 333</span>  NCMatrix m_mat; <span class="comment">// The input (permuted ) matrix </span></div> <div class="line"><a name="l00334"></a><span class="lineno"> 334</span>  SCMatrix m_Lstore; <span class="comment">// The lower triangular matrix (supernodal)</span></div> <div class="line"><a name="l00335"></a><span class="lineno"> 335</span>  MappedSparseMatrix<Scalar,ColMajor,Index> m_Ustore; <span class="comment">// The upper triangular matrix</span></div> <div class="line"><a name="l00336"></a><span class="lineno"> 336</span>  PermutationType m_perm_c; <span class="comment">// Column permutation </span></div> <div class="line"><a name="l00337"></a><span class="lineno"> 337</span>  PermutationType m_perm_r ; <span class="comment">// Row permutation</span></div> <div class="line"><a name="l00338"></a><span class="lineno"> 338</span>  IndexVector m_etree; <span class="comment">// Column elimination tree </span></div> <div class="line"><a name="l00339"></a><span class="lineno"> 339</span>  </div> <div class="line"><a name="l00340"></a><span class="lineno"> 340</span>  <span class="keyword">typename</span> Base::GlobalLU_t m_glu; </div> <div class="line"><a name="l00341"></a><span class="lineno"> 341</span>  </div> <div class="line"><a name="l00342"></a><span class="lineno"> 342</span>  <span class="comment">// SparseLU options </span></div> <div class="line"><a name="l00343"></a><span class="lineno"> 343</span>  <span class="keywordtype">bool</span> m_symmetricmode;</div> <div class="line"><a name="l00344"></a><span class="lineno"> 344</span>  <span class="comment">// values for performance </span></div> <div class="line"><a name="l00345"></a><span class="lineno"> 345</span>  internal::perfvalues<Index> m_perfv; </div> <div class="line"><a name="l00346"></a><span class="lineno"> 346</span>  RealScalar m_diagpivotthresh; <span class="comment">// Specifies the threshold used for a diagonal entry to be an acceptable pivot</span></div> <div class="line"><a name="l00347"></a><span class="lineno"> 347</span>  Index m_nnzL, m_nnzU; <span class="comment">// Nonzeros in L and U factors </span></div> <div class="line"><a name="l00348"></a><span class="lineno"> 348</span>  Index m_detPermR; <span class="comment">// Determinant of the coefficient matrix</span></div> <div class="line"><a name="l00349"></a><span class="lineno"> 349</span>  <span class="keyword">private</span>:</div> <div class="line"><a name="l00350"></a><span class="lineno"> 350</span>  <span class="comment">// Disable copy constructor </span></div> <div class="line"><a name="l00351"></a><span class="lineno"> 351</span>  SparseLU (<span class="keyword">const</span> SparseLU& );</div> <div class="line"><a name="l00352"></a><span class="lineno"> 352</span>  </div> <div class="line"><a name="l00353"></a><span class="lineno"> 353</span> }; <span class="comment">// End class SparseLU</span></div> <div class="line"><a name="l00354"></a><span class="lineno"> 354</span> </div> <div class="line"><a name="l00355"></a><span class="lineno"> 355</span> </div> <div class="line"><a name="l00356"></a><span class="lineno"> 356</span> </div> <div class="line"><a name="l00357"></a><span class="lineno"> 357</span> <span class="comment">// Functions needed by the anaysis phase</span></div> <div class="line"><a name="l00368"></a><span class="lineno"> 368</span> <span class="comment"></span><span class="keyword">template</span> <<span class="keyword">typename</span> MatrixType, <span class="keyword">typename</span> OrderingType></div> <div class="line"><a name="l00369"></a><span class="lineno"><a class="line" href="classEigen_1_1SparseLU.html#a6f2a135bd741b6b2f2558e6a404581ff"> 369</a></span> <span class="keywordtype">void</span> <a class="code" href="classEigen_1_1SparseLU.html#a6f2a135bd741b6b2f2558e6a404581ff">SparseLU<MatrixType, OrderingType>::analyzePattern</a>(<span class="keyword">const</span> MatrixType& mat)</div> <div class="line"><a name="l00370"></a><span class="lineno"> 370</span> {</div> <div class="line"><a name="l00371"></a><span class="lineno"> 371</span>  </div> <div class="line"><a name="l00372"></a><span class="lineno"> 372</span>  <span class="comment">//TODO It is possible as in SuperLU to compute row and columns scaling vectors to equilibrate the matrix mat.</span></div> <div class="line"><a name="l00373"></a><span class="lineno"> 373</span>  </div> <div class="line"><a name="l00374"></a><span class="lineno"> 374</span>  OrderingType ord; </div> <div class="line"><a name="l00375"></a><span class="lineno"> 375</span>  ord(mat,m_perm_c);</div> <div class="line"><a name="l00376"></a><span class="lineno"> 376</span>  </div> <div class="line"><a name="l00377"></a><span class="lineno"> 377</span>  <span class="comment">// Apply the permutation to the column of the input matrix</span></div> <div class="line"><a name="l00378"></a><span class="lineno"> 378</span>  <span class="comment">//First copy the whole input matrix. </span></div> <div class="line"><a name="l00379"></a><span class="lineno"> 379</span>  m_mat = mat;</div> <div class="line"><a name="l00380"></a><span class="lineno"> 380</span>  <span class="keywordflow">if</span> (m_perm_c.size()) {</div> <div class="line"><a name="l00381"></a><span class="lineno"> 381</span>  m_mat.uncompress(); <span class="comment">//NOTE: The effect of this command is only to create the InnerNonzeros pointers. FIXME : This vector is filled but not subsequently used. </span></div> <div class="line"><a name="l00382"></a><span class="lineno"> 382</span>  <span class="comment">//Then, permute only the column pointers</span></div> <div class="line"><a name="l00383"></a><span class="lineno"> 383</span>  <span class="keyword">const</span> Index * outerIndexPtr;</div> <div class="line"><a name="l00384"></a><span class="lineno"> 384</span>  <span class="keywordflow">if</span> (mat.isCompressed()) outerIndexPtr = mat.outerIndexPtr();</div> <div class="line"><a name="l00385"></a><span class="lineno"> 385</span>  <span class="keywordflow">else</span></div> <div class="line"><a name="l00386"></a><span class="lineno"> 386</span>  {</div> <div class="line"><a name="l00387"></a><span class="lineno"> 387</span>  Index *outerIndexPtr_t = <span class="keyword">new</span> Index[mat.cols()+1];</div> <div class="line"><a name="l00388"></a><span class="lineno"> 388</span>  <span class="keywordflow">for</span>(Index i = 0; i <= mat.cols(); i++) outerIndexPtr_t[i] = m_mat.outerIndexPtr()[i];</div> <div class="line"><a name="l00389"></a><span class="lineno"> 389</span>  outerIndexPtr = outerIndexPtr_t;</div> <div class="line"><a name="l00390"></a><span class="lineno"> 390</span>  }</div> <div class="line"><a name="l00391"></a><span class="lineno"> 391</span>  <span class="keywordflow">for</span> (Index i = 0; i < mat.cols(); i++)</div> <div class="line"><a name="l00392"></a><span class="lineno"> 392</span>  {</div> <div class="line"><a name="l00393"></a><span class="lineno"> 393</span>  m_mat.outerIndexPtr()[m_perm_c.indices()(i)] = outerIndexPtr[i];</div> <div class="line"><a name="l00394"></a><span class="lineno"> 394</span>  m_mat.innerNonZeroPtr()[m_perm_c.indices()(i)] = outerIndexPtr[i+1] - outerIndexPtr[i];</div> <div class="line"><a name="l00395"></a><span class="lineno"> 395</span>  }</div> <div class="line"><a name="l00396"></a><span class="lineno"> 396</span>  <span class="keywordflow">if</span>(!mat.isCompressed()) <span class="keyword">delete</span>[] outerIndexPtr;</div> <div class="line"><a name="l00397"></a><span class="lineno"> 397</span>  }</div> <div class="line"><a name="l00398"></a><span class="lineno"> 398</span>  <span class="comment">// Compute the column elimination tree of the permuted matrix </span></div> <div class="line"><a name="l00399"></a><span class="lineno"> 399</span>  <a class="code" href="classEigen_1_1Matrix.html">IndexVector</a> firstRowElt;</div> <div class="line"><a name="l00400"></a><span class="lineno"> 400</span>  <a class="code" href="namespaceEigen_1_1internal.html#a22b184331e5fb270a37d1305e95cb064">internal::coletree</a>(m_mat, m_etree,firstRowElt); </div> <div class="line"><a name="l00401"></a><span class="lineno"> 401</span>  </div> <div class="line"><a name="l00402"></a><span class="lineno"> 402</span>  <span class="comment">// In symmetric mode, do not do postorder here</span></div> <div class="line"><a name="l00403"></a><span class="lineno"> 403</span>  <span class="keywordflow">if</span> (!m_symmetricmode) {</div> <div class="line"><a name="l00404"></a><span class="lineno"> 404</span>  <a class="code" href="classEigen_1_1Matrix.html">IndexVector</a> post, iwork; </div> <div class="line"><a name="l00405"></a><span class="lineno"> 405</span>  <span class="comment">// Post order etree</span></div> <div class="line"><a name="l00406"></a><span class="lineno"> 406</span>  <a class="code" href="namespaceEigen_1_1internal.html#a6afa7df34de3dce07b21194168118456">internal::treePostorder</a>(m_mat.cols(), m_etree, post); </div> <div class="line"><a name="l00407"></a><span class="lineno"> 407</span>  </div> <div class="line"><a name="l00408"></a><span class="lineno"> 408</span>  </div> <div class="line"><a name="l00409"></a><span class="lineno"> 409</span>  <span class="comment">// Renumber etree in postorder </span></div> <div class="line"><a name="l00410"></a><span class="lineno"> 410</span>  Index m = m_mat.cols(); </div> <div class="line"><a name="l00411"></a><span class="lineno"> 411</span>  iwork.<a class="code" href="classEigen_1_1PlainObjectBase.html#afbbb33d14fe7fb9683019a39ce1c659d">resize</a>(m+1);</div> <div class="line"><a name="l00412"></a><span class="lineno"> 412</span>  <span class="keywordflow">for</span> (Index i = 0; i < m; ++i) iwork(post(i)) = post(m_etree(i));</div> <div class="line"><a name="l00413"></a><span class="lineno"> 413</span>  m_etree = iwork;</div> <div class="line"><a name="l00414"></a><span class="lineno"> 414</span>  </div> <div class="line"><a name="l00415"></a><span class="lineno"> 415</span>  <span class="comment">// Postmultiply A*Pc by post, i.e reorder the matrix according to the postorder of the etree</span></div> <div class="line"><a name="l00416"></a><span class="lineno"> 416</span>  <a class="code" href="classEigen_1_1PermutationMatrix.html">PermutationType</a> post_perm(m); </div> <div class="line"><a name="l00417"></a><span class="lineno"> 417</span>  <span class="keywordflow">for</span> (Index i = 0; i < m; i++) </div> <div class="line"><a name="l00418"></a><span class="lineno"> 418</span>  post_perm.<a class="code" href="classEigen_1_1PermutationMatrix.html#a8876d615d17aad77b054a8f58b699e7d">indices</a>()(i) = post(i); </div> <div class="line"><a name="l00419"></a><span class="lineno"> 419</span>  </div> <div class="line"><a name="l00420"></a><span class="lineno"> 420</span>  <span class="comment">// Combine the two permutations : postorder the permutation for future use</span></div> <div class="line"><a name="l00421"></a><span class="lineno"> 421</span>  <span class="keywordflow">if</span>(m_perm_c.size()) {</div> <div class="line"><a name="l00422"></a><span class="lineno"> 422</span>  m_perm_c = post_perm * m_perm_c;</div> <div class="line"><a name="l00423"></a><span class="lineno"> 423</span>  }</div> <div class="line"><a name="l00424"></a><span class="lineno"> 424</span>  </div> <div class="line"><a name="l00425"></a><span class="lineno"> 425</span>  } <span class="comment">// end postordering </span></div> <div class="line"><a name="l00426"></a><span class="lineno"> 426</span>  </div> <div class="line"><a name="l00427"></a><span class="lineno"> 427</span>  m_analysisIsOk = <span class="keyword">true</span>; </div> <div class="line"><a name="l00428"></a><span class="lineno"> 428</span> }</div> <div class="line"><a name="l00429"></a><span class="lineno"> 429</span> </div> <div class="line"><a name="l00430"></a><span class="lineno"> 430</span> <span class="comment">// Functions needed by the numerical factorization phase</span></div> <div class="line"><a name="l00431"></a><span class="lineno"> 431</span> </div> <div class="line"><a name="l00432"></a><span class="lineno"> 432</span> </div> <div class="line"><a name="l00451"></a><span class="lineno"> 451</span> <span class="keyword">template</span> <<span class="keyword">typename</span> MatrixType, <span class="keyword">typename</span> OrderingType></div> <div class="line"><a name="l00452"></a><span class="lineno"><a class="line" href="classEigen_1_1SparseLU.html#a33e2421e468033a883c041f940537a7c"> 452</a></span> <span class="keywordtype">void</span> <a class="code" href="classEigen_1_1SparseLU.html#a33e2421e468033a883c041f940537a7c">SparseLU<MatrixType, OrderingType>::factorize</a>(<span class="keyword">const</span> MatrixType& matrix)</div> <div class="line"><a name="l00453"></a><span class="lineno"> 453</span> {</div> <div class="line"><a name="l00454"></a><span class="lineno"> 454</span>  <span class="keyword">using</span> internal::emptyIdxLU;</div> <div class="line"><a name="l00455"></a><span class="lineno"> 455</span>  eigen_assert(m_analysisIsOk && <span class="stringliteral">"analyzePattern() should be called first"</span>); </div> <div class="line"><a name="l00456"></a><span class="lineno"> 456</span>  eigen_assert((matrix.rows() == matrix.cols()) && <span class="stringliteral">"Only for squared matrices"</span>);</div> <div class="line"><a name="l00457"></a><span class="lineno"> 457</span>  </div> <div class="line"><a name="l00458"></a><span class="lineno"> 458</span>  <span class="keyword">typedef</span> <span class="keyword">typename</span> IndexVector::Scalar Index; </div> <div class="line"><a name="l00459"></a><span class="lineno"> 459</span>  </div> <div class="line"><a name="l00460"></a><span class="lineno"> 460</span>  </div> <div class="line"><a name="l00461"></a><span class="lineno"> 461</span>  <span class="comment">// Apply the column permutation computed in analyzepattern()</span></div> <div class="line"><a name="l00462"></a><span class="lineno"> 462</span>  <span class="comment">// m_mat = matrix * m_perm_c.inverse(); </span></div> <div class="line"><a name="l00463"></a><span class="lineno"> 463</span>  m_mat = matrix;</div> <div class="line"><a name="l00464"></a><span class="lineno"> 464</span>  <span class="keywordflow">if</span> (m_perm_c.size()) </div> <div class="line"><a name="l00465"></a><span class="lineno"> 465</span>  {</div> <div class="line"><a name="l00466"></a><span class="lineno"> 466</span>  m_mat.uncompress(); <span class="comment">//NOTE: The effect of this command is only to create the InnerNonzeros pointers.</span></div> <div class="line"><a name="l00467"></a><span class="lineno"> 467</span>  <span class="comment">//Then, permute only the column pointers</span></div> <div class="line"><a name="l00468"></a><span class="lineno"> 468</span>  <span class="keyword">const</span> Index * outerIndexPtr;</div> <div class="line"><a name="l00469"></a><span class="lineno"> 469</span>  <span class="keywordflow">if</span> (matrix.isCompressed()) outerIndexPtr = matrix.outerIndexPtr();</div> <div class="line"><a name="l00470"></a><span class="lineno"> 470</span>  <span class="keywordflow">else</span></div> <div class="line"><a name="l00471"></a><span class="lineno"> 471</span>  {</div> <div class="line"><a name="l00472"></a><span class="lineno"> 472</span>  Index* outerIndexPtr_t = <span class="keyword">new</span> Index[matrix.cols()+1];</div> <div class="line"><a name="l00473"></a><span class="lineno"> 473</span>  <span class="keywordflow">for</span>(Index i = 0; i <= matrix.cols(); i++) outerIndexPtr_t[i] = m_mat.outerIndexPtr()[i];</div> <div class="line"><a name="l00474"></a><span class="lineno"> 474</span>  outerIndexPtr = outerIndexPtr_t;</div> <div class="line"><a name="l00475"></a><span class="lineno"> 475</span>  }</div> <div class="line"><a name="l00476"></a><span class="lineno"> 476</span>  <span class="keywordflow">for</span> (Index i = 0; i < matrix.cols(); i++)</div> <div class="line"><a name="l00477"></a><span class="lineno"> 477</span>  {</div> <div class="line"><a name="l00478"></a><span class="lineno"> 478</span>  m_mat.outerIndexPtr()[m_perm_c.indices()(i)] = outerIndexPtr[i];</div> <div class="line"><a name="l00479"></a><span class="lineno"> 479</span>  m_mat.innerNonZeroPtr()[m_perm_c.indices()(i)] = outerIndexPtr[i+1] - outerIndexPtr[i];</div> <div class="line"><a name="l00480"></a><span class="lineno"> 480</span>  }</div> <div class="line"><a name="l00481"></a><span class="lineno"> 481</span>  <span class="keywordflow">if</span>(!matrix.isCompressed()) <span class="keyword">delete</span>[] outerIndexPtr;</div> <div class="line"><a name="l00482"></a><span class="lineno"> 482</span>  } </div> <div class="line"><a name="l00483"></a><span class="lineno"> 483</span>  <span class="keywordflow">else</span> </div> <div class="line"><a name="l00484"></a><span class="lineno"> 484</span>  { <span class="comment">//FIXME This should not be needed if the empty permutation is handled transparently</span></div> <div class="line"><a name="l00485"></a><span class="lineno"> 485</span>  m_perm_c.resize(matrix.cols());</div> <div class="line"><a name="l00486"></a><span class="lineno"> 486</span>  <span class="keywordflow">for</span>(Index i = 0; i < matrix.cols(); ++i) m_perm_c.indices()(i) = i;</div> <div class="line"><a name="l00487"></a><span class="lineno"> 487</span>  }</div> <div class="line"><a name="l00488"></a><span class="lineno"> 488</span>  </div> <div class="line"><a name="l00489"></a><span class="lineno"> 489</span>  Index m = m_mat.rows();</div> <div class="line"><a name="l00490"></a><span class="lineno"> 490</span>  Index n = m_mat.cols();</div> <div class="line"><a name="l00491"></a><span class="lineno"> 491</span>  Index nnz = m_mat.nonZeros();</div> <div class="line"><a name="l00492"></a><span class="lineno"> 492</span>  Index maxpanel = m_perfv.panel_size * m;</div> <div class="line"><a name="l00493"></a><span class="lineno"> 493</span>  <span class="comment">// Allocate working storage common to the factor routines</span></div> <div class="line"><a name="l00494"></a><span class="lineno"> 494</span>  Index lwork = 0;</div> <div class="line"><a name="l00495"></a><span class="lineno"> 495</span>  Index info = Base::memInit(m, n, nnz, lwork, m_perfv.fillfactor, m_perfv.panel_size, m_glu); </div> <div class="line"><a name="l00496"></a><span class="lineno"> 496</span>  <span class="keywordflow">if</span> (info) </div> <div class="line"><a name="l00497"></a><span class="lineno"> 497</span>  {</div> <div class="line"><a name="l00498"></a><span class="lineno"> 498</span>  m_lastError = <span class="stringliteral">"UNABLE TO ALLOCATE WORKING MEMORY\n\n"</span> ;</div> <div class="line"><a name="l00499"></a><span class="lineno"> 499</span>  m_factorizationIsOk = <span class="keyword">false</span>;</div> <div class="line"><a name="l00500"></a><span class="lineno"> 500</span>  return ; </div> <div class="line"><a name="l00501"></a><span class="lineno"> 501</span>  }</div> <div class="line"><a name="l00502"></a><span class="lineno"> 502</span>  </div> <div class="line"><a name="l00503"></a><span class="lineno"> 503</span>  <span class="comment">// Set up pointers for integer working arrays </span></div> <div class="line"><a name="l00504"></a><span class="lineno"> 504</span>  <a class="code" href="classEigen_1_1Matrix.html">IndexVector</a> segrep(m); segrep.<a class="code" href="classEigen_1_1PlainObjectBase.html#a04abe84a9a894de335a232681d9a0722">setZero</a>();</div> <div class="line"><a name="l00505"></a><span class="lineno"> 505</span>  <a class="code" href="classEigen_1_1Matrix.html">IndexVector</a> parent(m); parent.<a class="code" href="classEigen_1_1PlainObjectBase.html#a04abe84a9a894de335a232681d9a0722">setZero</a>();</div> <div class="line"><a name="l00506"></a><span class="lineno"> 506</span>  <a class="code" href="classEigen_1_1Matrix.html">IndexVector</a> xplore(m); xplore.<a class="code" href="classEigen_1_1PlainObjectBase.html#a04abe84a9a894de335a232681d9a0722">setZero</a>();</div> <div class="line"><a name="l00507"></a><span class="lineno"> 507</span>  <a class="code" href="classEigen_1_1Matrix.html">IndexVector</a> repfnz(maxpanel);</div> <div class="line"><a name="l00508"></a><span class="lineno"> 508</span>  <a class="code" href="classEigen_1_1Matrix.html">IndexVector</a> panel_lsub(maxpanel);</div> <div class="line"><a name="l00509"></a><span class="lineno"> 509</span>  <a class="code" href="classEigen_1_1Matrix.html">IndexVector</a> xprune(n); xprune.<a class="code" href="classEigen_1_1PlainObjectBase.html#a04abe84a9a894de335a232681d9a0722">setZero</a>();</div> <div class="line"><a name="l00510"></a><span class="lineno"> 510</span>  <a class="code" href="classEigen_1_1Matrix.html">IndexVector</a> marker(m*internal::LUNoMarker); marker.<a class="code" href="classEigen_1_1PlainObjectBase.html#a04abe84a9a894de335a232681d9a0722">setZero</a>();</div> <div class="line"><a name="l00511"></a><span class="lineno"> 511</span>  </div> <div class="line"><a name="l00512"></a><span class="lineno"> 512</span>  repfnz.<a class="code" href="classEigen_1_1PlainObjectBase.html#aabcf7b6f4543f6255253b9ccee3309b8">setConstant</a>(-1); </div> <div class="line"><a name="l00513"></a><span class="lineno"> 513</span>  panel_lsub.<a class="code" href="classEigen_1_1PlainObjectBase.html#aabcf7b6f4543f6255253b9ccee3309b8">setConstant</a>(-1);</div> <div class="line"><a name="l00514"></a><span class="lineno"> 514</span>  </div> <div class="line"><a name="l00515"></a><span class="lineno"> 515</span>  <span class="comment">// Set up pointers for scalar working arrays </span></div> <div class="line"><a name="l00516"></a><span class="lineno"> 516</span>  <a class="code" href="classEigen_1_1Matrix.html">ScalarVector</a> dense; </div> <div class="line"><a name="l00517"></a><span class="lineno"> 517</span>  dense.<a class="code" href="classEigen_1_1PlainObjectBase.html#a04abe84a9a894de335a232681d9a0722">setZero</a>(maxpanel);</div> <div class="line"><a name="l00518"></a><span class="lineno"> 518</span>  <a class="code" href="classEigen_1_1Matrix.html">ScalarVector</a> tempv; </div> <div class="line"><a name="l00519"></a><span class="lineno"> 519</span>  tempv.<a class="code" href="classEigen_1_1PlainObjectBase.html#a04abe84a9a894de335a232681d9a0722">setZero</a>(internal::LUnumTempV(m, m_perfv.panel_size, m_perfv.maxsuper, <span class="comment">/*m_perfv.rowblk*/</span>m) );</div> <div class="line"><a name="l00520"></a><span class="lineno"> 520</span>  </div> <div class="line"><a name="l00521"></a><span class="lineno"> 521</span>  <span class="comment">// Compute the inverse of perm_c</span></div> <div class="line"><a name="l00522"></a><span class="lineno"> 522</span>  <a class="code" href="classEigen_1_1PermutationMatrix.html">PermutationType</a> iperm_c(m_perm_c.inverse()); </div> <div class="line"><a name="l00523"></a><span class="lineno"> 523</span>  </div> <div class="line"><a name="l00524"></a><span class="lineno"> 524</span>  <span class="comment">// Identify initial relaxed snodes</span></div> <div class="line"><a name="l00525"></a><span class="lineno"> 525</span>  <a class="code" href="classEigen_1_1Matrix.html">IndexVector</a> relax_end(n);</div> <div class="line"><a name="l00526"></a><span class="lineno"> 526</span>  <span class="keywordflow">if</span> ( m_symmetricmode == <span class="keyword">true</span> ) </div> <div class="line"><a name="l00527"></a><span class="lineno"> 527</span>  Base::heap_relax_snode(n, m_etree, m_perfv.relax, marker, relax_end);</div> <div class="line"><a name="l00528"></a><span class="lineno"> 528</span>  <span class="keywordflow">else</span></div> <div class="line"><a name="l00529"></a><span class="lineno"> 529</span>  Base::relax_snode(n, m_etree, m_perfv.relax, marker, relax_end);</div> <div class="line"><a name="l00530"></a><span class="lineno"> 530</span>  </div> <div class="line"><a name="l00531"></a><span class="lineno"> 531</span>  </div> <div class="line"><a name="l00532"></a><span class="lineno"> 532</span>  m_perm_r.resize(m); </div> <div class="line"><a name="l00533"></a><span class="lineno"> 533</span>  m_perm_r.indices().setConstant(-1);</div> <div class="line"><a name="l00534"></a><span class="lineno"> 534</span>  marker.<a class="code" href="classEigen_1_1PlainObjectBase.html#aabcf7b6f4543f6255253b9ccee3309b8">setConstant</a>(-1);</div> <div class="line"><a name="l00535"></a><span class="lineno"> 535</span>  m_detPermR = 1; <span class="comment">// Record the determinant of the row permutation</span></div> <div class="line"><a name="l00536"></a><span class="lineno"> 536</span>  </div> <div class="line"><a name="l00537"></a><span class="lineno"> 537</span>  m_glu.supno(0) = emptyIdxLU; m_glu.xsup.setConstant(0);</div> <div class="line"><a name="l00538"></a><span class="lineno"> 538</span>  m_glu.xsup(0) = m_glu.xlsub(0) = m_glu.xusub(0) = m_glu.xlusup(0) = Index(0);</div> <div class="line"><a name="l00539"></a><span class="lineno"> 539</span>  </div> <div class="line"><a name="l00540"></a><span class="lineno"> 540</span>  <span class="comment">// Work on one 'panel' at a time. A panel is one of the following :</span></div> <div class="line"><a name="l00541"></a><span class="lineno"> 541</span>  <span class="comment">// (a) a relaxed supernode at the bottom of the etree, or</span></div> <div class="line"><a name="l00542"></a><span class="lineno"> 542</span>  <span class="comment">// (b) panel_size contiguous columns, <panel_size> defined by the user</span></div> <div class="line"><a name="l00543"></a><span class="lineno"> 543</span>  Index jcol; </div> <div class="line"><a name="l00544"></a><span class="lineno"> 544</span>  <a class="code" href="classEigen_1_1Matrix.html">IndexVector</a> panel_histo(n);</div> <div class="line"><a name="l00545"></a><span class="lineno"> 545</span>  Index pivrow; <span class="comment">// Pivotal row number in the original row matrix</span></div> <div class="line"><a name="l00546"></a><span class="lineno"> 546</span>  Index nseg1; <span class="comment">// Number of segments in U-column above panel row jcol</span></div> <div class="line"><a name="l00547"></a><span class="lineno"> 547</span>  Index nseg; <span class="comment">// Number of segments in each U-column </span></div> <div class="line"><a name="l00548"></a><span class="lineno"> 548</span>  Index irep; </div> <div class="line"><a name="l00549"></a><span class="lineno"> 549</span>  Index i, k, jj; </div> <div class="line"><a name="l00550"></a><span class="lineno"> 550</span>  <span class="keywordflow">for</span> (jcol = 0; jcol < n; )</div> <div class="line"><a name="l00551"></a><span class="lineno"> 551</span>  {</div> <div class="line"><a name="l00552"></a><span class="lineno"> 552</span>  <span class="comment">// Adjust panel size so that a panel won't overlap with the next relaxed snode. </span></div> <div class="line"><a name="l00553"></a><span class="lineno"> 553</span>  Index panel_size = m_perfv.panel_size; <span class="comment">// upper bound on panel width</span></div> <div class="line"><a name="l00554"></a><span class="lineno"> 554</span>  <span class="keywordflow">for</span> (k = jcol + 1; k < (std::min)(jcol+panel_size, n); k++)</div> <div class="line"><a name="l00555"></a><span class="lineno"> 555</span>  {</div> <div class="line"><a name="l00556"></a><span class="lineno"> 556</span>  <span class="keywordflow">if</span> (relax_end(k) != emptyIdxLU) </div> <div class="line"><a name="l00557"></a><span class="lineno"> 557</span>  {</div> <div class="line"><a name="l00558"></a><span class="lineno"> 558</span>  panel_size = k - jcol; </div> <div class="line"><a name="l00559"></a><span class="lineno"> 559</span>  <span class="keywordflow">break</span>; </div> <div class="line"><a name="l00560"></a><span class="lineno"> 560</span>  }</div> <div class="line"><a name="l00561"></a><span class="lineno"> 561</span>  }</div> <div class="line"><a name="l00562"></a><span class="lineno"> 562</span>  <span class="keywordflow">if</span> (k == n) </div> <div class="line"><a name="l00563"></a><span class="lineno"> 563</span>  panel_size = n - jcol; </div> <div class="line"><a name="l00564"></a><span class="lineno"> 564</span>  </div> <div class="line"><a name="l00565"></a><span class="lineno"> 565</span>  <span class="comment">// Symbolic outer factorization on a panel of columns </span></div> <div class="line"><a name="l00566"></a><span class="lineno"> 566</span>  Base::panel_dfs(m, panel_size, jcol, m_mat, m_perm_r.indices(), nseg1, dense, panel_lsub, segrep, repfnz, xprune, marker, parent, xplore, m_glu); </div> <div class="line"><a name="l00567"></a><span class="lineno"> 567</span>  </div> <div class="line"><a name="l00568"></a><span class="lineno"> 568</span>  <span class="comment">// Numeric sup-panel updates in topological order </span></div> <div class="line"><a name="l00569"></a><span class="lineno"> 569</span>  Base::panel_bmod(m, panel_size, jcol, nseg1, dense, tempv, segrep, repfnz, m_glu); </div> <div class="line"><a name="l00570"></a><span class="lineno"> 570</span>  </div> <div class="line"><a name="l00571"></a><span class="lineno"> 571</span>  <span class="comment">// Sparse LU within the panel, and below the panel diagonal </span></div> <div class="line"><a name="l00572"></a><span class="lineno"> 572</span>  <span class="keywordflow">for</span> ( jj = jcol; jj< jcol + panel_size; jj++) </div> <div class="line"><a name="l00573"></a><span class="lineno"> 573</span>  {</div> <div class="line"><a name="l00574"></a><span class="lineno"> 574</span>  k = (jj - jcol) * m; <span class="comment">// Column index for w-wide arrays </span></div> <div class="line"><a name="l00575"></a><span class="lineno"> 575</span>  </div> <div class="line"><a name="l00576"></a><span class="lineno"> 576</span>  nseg = nseg1; <span class="comment">// begin after all the panel segments</span></div> <div class="line"><a name="l00577"></a><span class="lineno"> 577</span>  <span class="comment">//Depth-first-search for the current column</span></div> <div class="line"><a name="l00578"></a><span class="lineno"> 578</span>  <a class="code" href="classEigen_1_1VectorBlock.html">VectorBlock<IndexVector></a> panel_lsubk(panel_lsub, k, m);</div> <div class="line"><a name="l00579"></a><span class="lineno"> 579</span>  <a class="code" href="classEigen_1_1VectorBlock.html">VectorBlock<IndexVector></a> repfnz_k(repfnz, k, m); </div> <div class="line"><a name="l00580"></a><span class="lineno"> 580</span>  info = Base::column_dfs(m, jj, m_perm_r.indices(), m_perfv.maxsuper, nseg, panel_lsubk, segrep, repfnz_k, xprune, marker, parent, xplore, m_glu); </div> <div class="line"><a name="l00581"></a><span class="lineno"> 581</span>  <span class="keywordflow">if</span> ( info ) </div> <div class="line"><a name="l00582"></a><span class="lineno"> 582</span>  {</div> <div class="line"><a name="l00583"></a><span class="lineno"> 583</span>  m_lastError = <span class="stringliteral">"UNABLE TO EXPAND MEMORY IN COLUMN_DFS() "</span>;</div> <div class="line"><a name="l00584"></a><span class="lineno"> 584</span>  m_info = <a class="code" href="group__enums.html#gga51bc1ac16f26ebe51eae1abb77bd037ba710fff14e8fc77846d4b75d8f4cc2d5c">NumericalIssue</a>; </div> <div class="line"><a name="l00585"></a><span class="lineno"> 585</span>  m_factorizationIsOk = <span class="keyword">false</span>; </div> <div class="line"><a name="l00586"></a><span class="lineno"> 586</span>  <span class="keywordflow">return</span>; </div> <div class="line"><a name="l00587"></a><span class="lineno"> 587</span>  }</div> <div class="line"><a name="l00588"></a><span class="lineno"> 588</span>  <span class="comment">// Numeric updates to this column </span></div> <div class="line"><a name="l00589"></a><span class="lineno"> 589</span>  <a class="code" href="classEigen_1_1VectorBlock.html">VectorBlock<ScalarVector></a> dense_k(dense, k, m); </div> <div class="line"><a name="l00590"></a><span class="lineno"> 590</span>  <a class="code" href="classEigen_1_1VectorBlock.html">VectorBlock<IndexVector></a> segrep_k(segrep, nseg1, m-nseg1); </div> <div class="line"><a name="l00591"></a><span class="lineno"> 591</span>  info = Base::column_bmod(jj, (nseg - nseg1), dense_k, tempv, segrep_k, repfnz_k, jcol, m_glu); </div> <div class="line"><a name="l00592"></a><span class="lineno"> 592</span>  <span class="keywordflow">if</span> ( info ) </div> <div class="line"><a name="l00593"></a><span class="lineno"> 593</span>  {</div> <div class="line"><a name="l00594"></a><span class="lineno"> 594</span>  m_lastError = <span class="stringliteral">"UNABLE TO EXPAND MEMORY IN COLUMN_BMOD() "</span>;</div> <div class="line"><a name="l00595"></a><span class="lineno"> 595</span>  m_info = <a class="code" href="group__enums.html#gga51bc1ac16f26ebe51eae1abb77bd037ba710fff14e8fc77846d4b75d8f4cc2d5c">NumericalIssue</a>; </div> <div class="line"><a name="l00596"></a><span class="lineno"> 596</span>  m_factorizationIsOk = <span class="keyword">false</span>; </div> <div class="line"><a name="l00597"></a><span class="lineno"> 597</span>  <span class="keywordflow">return</span>; </div> <div class="line"><a name="l00598"></a><span class="lineno"> 598</span>  }</div> <div class="line"><a name="l00599"></a><span class="lineno"> 599</span>  </div> <div class="line"><a name="l00600"></a><span class="lineno"> 600</span>  <span class="comment">// Copy the U-segments to ucol(*)</span></div> <div class="line"><a name="l00601"></a><span class="lineno"> 601</span>  info = Base::copy_to_ucol(jj, nseg, segrep, repfnz_k ,m_perm_r.indices(), dense_k, m_glu); </div> <div class="line"><a name="l00602"></a><span class="lineno"> 602</span>  <span class="keywordflow">if</span> ( info ) </div> <div class="line"><a name="l00603"></a><span class="lineno"> 603</span>  {</div> <div class="line"><a name="l00604"></a><span class="lineno"> 604</span>  m_lastError = <span class="stringliteral">"UNABLE TO EXPAND MEMORY IN COPY_TO_UCOL() "</span>;</div> <div class="line"><a name="l00605"></a><span class="lineno"> 605</span>  m_info = <a class="code" href="group__enums.html#gga51bc1ac16f26ebe51eae1abb77bd037ba710fff14e8fc77846d4b75d8f4cc2d5c">NumericalIssue</a>; </div> <div class="line"><a name="l00606"></a><span class="lineno"> 606</span>  m_factorizationIsOk = <span class="keyword">false</span>; </div> <div class="line"><a name="l00607"></a><span class="lineno"> 607</span>  <span class="keywordflow">return</span>; </div> <div class="line"><a name="l00608"></a><span class="lineno"> 608</span>  }</div> <div class="line"><a name="l00609"></a><span class="lineno"> 609</span>  </div> <div class="line"><a name="l00610"></a><span class="lineno"> 610</span>  <span class="comment">// Form the L-segment </span></div> <div class="line"><a name="l00611"></a><span class="lineno"> 611</span>  info = Base::pivotL(jj, m_diagpivotthresh, m_perm_r.indices(), iperm_c.indices(), pivrow, m_glu);</div> <div class="line"><a name="l00612"></a><span class="lineno"> 612</span>  <span class="keywordflow">if</span> ( info ) </div> <div class="line"><a name="l00613"></a><span class="lineno"> 613</span>  {</div> <div class="line"><a name="l00614"></a><span class="lineno"> 614</span>  m_lastError = <span class="stringliteral">"THE MATRIX IS STRUCTURALLY SINGULAR ... ZERO COLUMN AT "</span>;</div> <div class="line"><a name="l00615"></a><span class="lineno"> 615</span>  std::ostringstream returnInfo;</div> <div class="line"><a name="l00616"></a><span class="lineno"> 616</span>  returnInfo << info; </div> <div class="line"><a name="l00617"></a><span class="lineno"> 617</span>  m_lastError += returnInfo.str();</div> <div class="line"><a name="l00618"></a><span class="lineno"> 618</span>  m_info = <a class="code" href="group__enums.html#gga51bc1ac16f26ebe51eae1abb77bd037ba710fff14e8fc77846d4b75d8f4cc2d5c">NumericalIssue</a>; </div> <div class="line"><a name="l00619"></a><span class="lineno"> 619</span>  m_factorizationIsOk = <span class="keyword">false</span>; </div> <div class="line"><a name="l00620"></a><span class="lineno"> 620</span>  <span class="keywordflow">return</span>; </div> <div class="line"><a name="l00621"></a><span class="lineno"> 621</span>  }</div> <div class="line"><a name="l00622"></a><span class="lineno"> 622</span>  </div> <div class="line"><a name="l00623"></a><span class="lineno"> 623</span>  <span class="comment">// Update the determinant of the row permutation matrix</span></div> <div class="line"><a name="l00624"></a><span class="lineno"> 624</span>  <span class="keywordflow">if</span> (pivrow != jj) m_detPermR *= -1;</div> <div class="line"><a name="l00625"></a><span class="lineno"> 625</span> </div> <div class="line"><a name="l00626"></a><span class="lineno"> 626</span>  <span class="comment">// Prune columns (0:jj-1) using column jj</span></div> <div class="line"><a name="l00627"></a><span class="lineno"> 627</span>  Base::pruneL(jj, m_perm_r.indices(), pivrow, nseg, segrep, repfnz_k, xprune, m_glu); </div> <div class="line"><a name="l00628"></a><span class="lineno"> 628</span>  </div> <div class="line"><a name="l00629"></a><span class="lineno"> 629</span>  <span class="comment">// Reset repfnz for this column </span></div> <div class="line"><a name="l00630"></a><span class="lineno"> 630</span>  <span class="keywordflow">for</span> (i = 0; i < nseg; i++)</div> <div class="line"><a name="l00631"></a><span class="lineno"> 631</span>  {</div> <div class="line"><a name="l00632"></a><span class="lineno"> 632</span>  irep = segrep(i); </div> <div class="line"><a name="l00633"></a><span class="lineno"> 633</span>  repfnz_k(irep) = emptyIdxLU; </div> <div class="line"><a name="l00634"></a><span class="lineno"> 634</span>  }</div> <div class="line"><a name="l00635"></a><span class="lineno"> 635</span>  } <span class="comment">// end SparseLU within the panel </span></div> <div class="line"><a name="l00636"></a><span class="lineno"> 636</span>  jcol += panel_size; <span class="comment">// Move to the next panel</span></div> <div class="line"><a name="l00637"></a><span class="lineno"> 637</span>  } <span class="comment">// end for -- end elimination </span></div> <div class="line"><a name="l00638"></a><span class="lineno"> 638</span>  </div> <div class="line"><a name="l00639"></a><span class="lineno"> 639</span>  <span class="comment">// Count the number of nonzeros in factors </span></div> <div class="line"><a name="l00640"></a><span class="lineno"> 640</span>  Base::countnz(n, m_nnzL, m_nnzU, m_glu); </div> <div class="line"><a name="l00641"></a><span class="lineno"> 641</span>  <span class="comment">// Apply permutation to the L subscripts </span></div> <div class="line"><a name="l00642"></a><span class="lineno"> 642</span>  Base::fixupL(n, m_perm_r.indices(), m_glu); </div> <div class="line"><a name="l00643"></a><span class="lineno"> 643</span>  </div> <div class="line"><a name="l00644"></a><span class="lineno"> 644</span>  <span class="comment">// Create supernode matrix L </span></div> <div class="line"><a name="l00645"></a><span class="lineno"> 645</span>  m_Lstore.setInfos(m, n, m_glu.lusup, m_glu.xlusup, m_glu.lsub, m_glu.xlsub, m_glu.supno, m_glu.xsup); </div> <div class="line"><a name="l00646"></a><span class="lineno"> 646</span>  <span class="comment">// Create the column major upper sparse matrix U; </span></div> <div class="line"><a name="l00647"></a><span class="lineno"> 647</span>  <span class="keyword">new</span> (&m_Ustore) <a class="code" href="classEigen_1_1MappedSparseMatrix.html">MappedSparseMatrix<Scalar, ColMajor, Index></a> ( m, n, m_nnzU, m_glu.xusub.data(), m_glu.usub.data(), m_glu.ucol.data() ); </div> <div class="line"><a name="l00648"></a><span class="lineno"> 648</span>  </div> <div class="line"><a name="l00649"></a><span class="lineno"> 649</span>  m_info = <a class="code" href="group__enums.html#gga51bc1ac16f26ebe51eae1abb77bd037bafdfbdf3247bd36a1f17270d5cec74c9c">Success</a>;</div> <div class="line"><a name="l00650"></a><span class="lineno"> 650</span>  m_factorizationIsOk = <span class="keyword">true</span>;</div> <div class="line"><a name="l00651"></a><span class="lineno"> 651</span> }</div> <div class="line"><a name="l00652"></a><span class="lineno"> 652</span> </div> <div class="line"><a name="l00653"></a><span class="lineno"> 653</span> <span class="keyword">template</span><<span class="keyword">typename</span> MappedSupernodalType></div> <div class="line"><a name="l00654"></a><span class="lineno"> 654</span> <span class="keyword">struct </span>SparseLUMatrixLReturnType : internal::no_assignment_operator</div> <div class="line"><a name="l00655"></a><span class="lineno"> 655</span> {</div> <div class="line"><a name="l00656"></a><span class="lineno"> 656</span>  <span class="keyword">typedef</span> <span class="keyword">typename</span> MappedSupernodalType::Index Index;</div> <div class="line"><a name="l00657"></a><span class="lineno"> 657</span>  <span class="keyword">typedef</span> <span class="keyword">typename</span> MappedSupernodalType::Scalar Scalar;</div> <div class="line"><a name="l00658"></a><span class="lineno"> 658</span>  SparseLUMatrixLReturnType(<span class="keyword">const</span> MappedSupernodalType& mapL) : m_mapL(mapL)</div> <div class="line"><a name="l00659"></a><span class="lineno"> 659</span>  { }</div> <div class="line"><a name="l00660"></a><span class="lineno"> 660</span>  Index rows() { <span class="keywordflow">return</span> m_mapL.rows(); }</div> <div class="line"><a name="l00661"></a><span class="lineno"> 661</span>  Index cols() { <span class="keywordflow">return</span> m_mapL.cols(); }</div> <div class="line"><a name="l00662"></a><span class="lineno"> 662</span>  <span class="keyword">template</span><<span class="keyword">typename</span> Dest></div> <div class="line"><a name="l00663"></a><span class="lineno"> 663</span>  <span class="keywordtype">void</span> solveInPlace( MatrixBase<Dest> &X)<span class="keyword"> const</span></div> <div class="line"><a name="l00664"></a><span class="lineno"> 664</span> <span class="keyword"> </span>{</div> <div class="line"><a name="l00665"></a><span class="lineno"> 665</span>  m_mapL.solveInPlace(X);</div> <div class="line"><a name="l00666"></a><span class="lineno"> 666</span>  }</div> <div class="line"><a name="l00667"></a><span class="lineno"> 667</span>  <span class="keyword">const</span> MappedSupernodalType& m_mapL;</div> <div class="line"><a name="l00668"></a><span class="lineno"> 668</span> };</div> <div class="line"><a name="l00669"></a><span class="lineno"> 669</span> </div> <div class="line"><a name="l00670"></a><span class="lineno"> 670</span> <span class="keyword">template</span><<span class="keyword">typename</span> MatrixLType, <span class="keyword">typename</span> MatrixUType></div> <div class="line"><a name="l00671"></a><span class="lineno"> 671</span> <span class="keyword">struct </span>SparseLUMatrixUReturnType : internal::no_assignment_operator</div> <div class="line"><a name="l00672"></a><span class="lineno"> 672</span> {</div> <div class="line"><a name="l00673"></a><span class="lineno"> 673</span>  <span class="keyword">typedef</span> <span class="keyword">typename</span> MatrixLType::Index Index;</div> <div class="line"><a name="l00674"></a><span class="lineno"> 674</span>  <span class="keyword">typedef</span> <span class="keyword">typename</span> MatrixLType::Scalar Scalar;</div> <div class="line"><a name="l00675"></a><span class="lineno"> 675</span>  SparseLUMatrixUReturnType(<span class="keyword">const</span> MatrixLType& mapL, <span class="keyword">const</span> MatrixUType& mapU)</div> <div class="line"><a name="l00676"></a><span class="lineno"> 676</span>  : m_mapL(mapL),m_mapU(mapU)</div> <div class="line"><a name="l00677"></a><span class="lineno"> 677</span>  { }</div> <div class="line"><a name="l00678"></a><span class="lineno"> 678</span>  Index rows() { <span class="keywordflow">return</span> m_mapL.rows(); }</div> <div class="line"><a name="l00679"></a><span class="lineno"> 679</span>  Index cols() { <span class="keywordflow">return</span> m_mapL.cols(); }</div> <div class="line"><a name="l00680"></a><span class="lineno"> 680</span> </div> <div class="line"><a name="l00681"></a><span class="lineno"> 681</span>  <span class="keyword">template</span><<span class="keyword">typename</span> Dest> <span class="keywordtype">void</span> solveInPlace(MatrixBase<Dest> &X)<span class="keyword"> const</span></div> <div class="line"><a name="l00682"></a><span class="lineno"> 682</span> <span class="keyword"> </span>{</div> <div class="line"><a name="l00683"></a><span class="lineno"> 683</span>  Index nrhs = X.cols();</div> <div class="line"><a name="l00684"></a><span class="lineno"> 684</span>  Index n = X.rows();</div> <div class="line"><a name="l00685"></a><span class="lineno"> 685</span>  <span class="comment">// Backward solve with U</span></div> <div class="line"><a name="l00686"></a><span class="lineno"> 686</span>  <span class="keywordflow">for</span> (Index k = m_mapL.nsuper(); k >= 0; k--)</div> <div class="line"><a name="l00687"></a><span class="lineno"> 687</span>  {</div> <div class="line"><a name="l00688"></a><span class="lineno"> 688</span>  Index fsupc = m_mapL.supToCol()[k];</div> <div class="line"><a name="l00689"></a><span class="lineno"> 689</span>  Index lda = m_mapL.colIndexPtr()[fsupc+1] - m_mapL.colIndexPtr()[fsupc]; <span class="comment">// leading dimension</span></div> <div class="line"><a name="l00690"></a><span class="lineno"> 690</span>  Index nsupc = m_mapL.supToCol()[k+1] - fsupc;</div> <div class="line"><a name="l00691"></a><span class="lineno"> 691</span>  Index luptr = m_mapL.colIndexPtr()[fsupc];</div> <div class="line"><a name="l00692"></a><span class="lineno"> 692</span> </div> <div class="line"><a name="l00693"></a><span class="lineno"> 693</span>  <span class="keywordflow">if</span> (nsupc == 1)</div> <div class="line"><a name="l00694"></a><span class="lineno"> 694</span>  {</div> <div class="line"><a name="l00695"></a><span class="lineno"> 695</span>  <span class="keywordflow">for</span> (Index j = 0; j < nrhs; j++)</div> <div class="line"><a name="l00696"></a><span class="lineno"> 696</span>  {</div> <div class="line"><a name="l00697"></a><span class="lineno"> 697</span>  X(fsupc, j) /= m_mapL.valuePtr()[luptr];</div> <div class="line"><a name="l00698"></a><span class="lineno"> 698</span>  }</div> <div class="line"><a name="l00699"></a><span class="lineno"> 699</span>  }</div> <div class="line"><a name="l00700"></a><span class="lineno"> 700</span>  <span class="keywordflow">else</span></div> <div class="line"><a name="l00701"></a><span class="lineno"> 701</span>  {</div> <div class="line"><a name="l00702"></a><span class="lineno"> 702</span>  Map<const Matrix<Scalar,Dynamic,Dynamic>, 0, OuterStride<> > A( &(m_mapL.valuePtr()[luptr]), nsupc, nsupc, OuterStride<>(lda) );</div> <div class="line"><a name="l00703"></a><span class="lineno"> 703</span>  Map< Matrix<Scalar,Dynamic,Dynamic>, 0, OuterStride<> > U (&(X(fsupc,0)), nsupc, nrhs, OuterStride<>(n) );</div> <div class="line"><a name="l00704"></a><span class="lineno"> 704</span>  U = A.template triangularView<Upper>().solve(U);</div> <div class="line"><a name="l00705"></a><span class="lineno"> 705</span>  }</div> <div class="line"><a name="l00706"></a><span class="lineno"> 706</span> </div> <div class="line"><a name="l00707"></a><span class="lineno"> 707</span>  <span class="keywordflow">for</span> (Index j = 0; j < nrhs; ++j)</div> <div class="line"><a name="l00708"></a><span class="lineno"> 708</span>  {</div> <div class="line"><a name="l00709"></a><span class="lineno"> 709</span>  <span class="keywordflow">for</span> (Index jcol = fsupc; jcol < fsupc + nsupc; jcol++)</div> <div class="line"><a name="l00710"></a><span class="lineno"> 710</span>  {</div> <div class="line"><a name="l00711"></a><span class="lineno"> 711</span>  <span class="keyword">typename</span> MatrixUType::InnerIterator it(m_mapU, jcol);</div> <div class="line"><a name="l00712"></a><span class="lineno"> 712</span>  <span class="keywordflow">for</span> ( ; it; ++it)</div> <div class="line"><a name="l00713"></a><span class="lineno"> 713</span>  {</div> <div class="line"><a name="l00714"></a><span class="lineno"> 714</span>  Index irow = it.index();</div> <div class="line"><a name="l00715"></a><span class="lineno"> 715</span>  X(irow, j) -= X(jcol, j) * it.value();</div> <div class="line"><a name="l00716"></a><span class="lineno"> 716</span>  }</div> <div class="line"><a name="l00717"></a><span class="lineno"> 717</span>  }</div> <div class="line"><a name="l00718"></a><span class="lineno"> 718</span>  }</div> <div class="line"><a name="l00719"></a><span class="lineno"> 719</span>  } <span class="comment">// End For U-solve</span></div> <div class="line"><a name="l00720"></a><span class="lineno"> 720</span>  }</div> <div class="line"><a name="l00721"></a><span class="lineno"> 721</span>  <span class="keyword">const</span> MatrixLType& m_mapL;</div> <div class="line"><a name="l00722"></a><span class="lineno"> 722</span>  <span class="keyword">const</span> MatrixUType& m_mapU;</div> <div class="line"><a name="l00723"></a><span class="lineno"> 723</span> };</div> <div class="line"><a name="l00724"></a><span class="lineno"> 724</span> </div> <div class="line"><a name="l00725"></a><span class="lineno"> 725</span> <span class="keyword">namespace </span>internal {</div> <div class="line"><a name="l00726"></a><span class="lineno"> 726</span>  </div> <div class="line"><a name="l00727"></a><span class="lineno"> 727</span> <span class="keyword">template</span><<span class="keyword">typename</span> _MatrixType, <span class="keyword">typename</span> Derived, <span class="keyword">typename</span> Rhs></div> <div class="line"><a name="l00728"></a><span class="lineno"> 728</span> <span class="keyword">struct </span>solve_retval<SparseLU<_MatrixType,Derived>, Rhs></div> <div class="line"><a name="l00729"></a><span class="lineno"> 729</span>  : solve_retval_base<SparseLU<_MatrixType,Derived>, Rhs></div> <div class="line"><a name="l00730"></a><span class="lineno"> 730</span> {</div> <div class="line"><a name="l00731"></a><span class="lineno"> 731</span>  <span class="keyword">typedef</span> SparseLU<_MatrixType,Derived> Dec;</div> <div class="line"><a name="l00732"></a><span class="lineno"> 732</span>  EIGEN_MAKE_SOLVE_HELPERS(Dec,Rhs)</div> <div class="line"><a name="l00733"></a><span class="lineno"> 733</span> </div> <div class="line"><a name="l00734"></a><span class="lineno"> 734</span>  template<typename Dest> <span class="keywordtype">void</span> evalTo(Dest& dst)<span class="keyword"> const</span></div> <div class="line"><a name="l00735"></a><span class="lineno"> 735</span> <span class="keyword"> </span>{</div> <div class="line"><a name="l00736"></a><span class="lineno"> 736</span>  dec()._solve(rhs(),dst);</div> <div class="line"><a name="l00737"></a><span class="lineno"> 737</span>  }</div> <div class="line"><a name="l00738"></a><span class="lineno"> 738</span> };</div> <div class="line"><a name="l00739"></a><span class="lineno"> 739</span> </div> <div class="line"><a name="l00740"></a><span class="lineno"> 740</span> <span class="keyword">template</span><<span class="keyword">typename</span> _MatrixType, <span class="keyword">typename</span> Derived, <span class="keyword">typename</span> Rhs></div> <div class="line"><a name="l00741"></a><span class="lineno"> 741</span> <span class="keyword">struct </span>sparse_solve_retval<SparseLU<_MatrixType,Derived>, Rhs></div> <div class="line"><a name="l00742"></a><span class="lineno"> 742</span>  : sparse_solve_retval_base<SparseLU<_MatrixType,Derived>, Rhs></div> <div class="line"><a name="l00743"></a><span class="lineno"> 743</span> {</div> <div class="line"><a name="l00744"></a><span class="lineno"> 744</span>  <span class="keyword">typedef</span> SparseLU<_MatrixType,Derived> Dec;</div> <div class="line"><a name="l00745"></a><span class="lineno"> 745</span>  EIGEN_MAKE_SPARSE_SOLVE_HELPERS(Dec,Rhs)</div> <div class="line"><a name="l00746"></a><span class="lineno"> 746</span> </div> <div class="line"><a name="l00747"></a><span class="lineno"> 747</span>  template<typename Dest> <span class="keywordtype">void</span> evalTo(Dest& dst)<span class="keyword"> const</span></div> <div class="line"><a name="l00748"></a><span class="lineno"> 748</span> <span class="keyword"> </span>{</div> <div class="line"><a name="l00749"></a><span class="lineno"> 749</span>  this->defaultEvalTo(dst);</div> <div class="line"><a name="l00750"></a><span class="lineno"> 750</span>  }</div> <div class="line"><a name="l00751"></a><span class="lineno"> 751</span> };</div> <div class="line"><a name="l00752"></a><span class="lineno"> 752</span> } <span class="comment">// end namespace internal</span></div> <div class="line"><a name="l00753"></a><span class="lineno"> 753</span> </div> <div class="line"><a name="l00754"></a><span class="lineno"> 754</span> } <span class="comment">// End namespace Eigen </span></div> <div class="line"><a name="l00755"></a><span class="lineno"> 755</span> </div> <div class="line"><a name="l00756"></a><span class="lineno"> 756</span> <span class="preprocessor">#endif</span></div> <div class="ttc" id="classEigen_1_1SparseLU_html_a828ca6373b429d81f41f41429a7fd4ad"><div class="ttname"><a href="classEigen_1_1SparseLU.html#a828ca6373b429d81f41f41429a7fd4ad">Eigen::SparseLU::matrixL</a></div><div class="ttdeci">SparseLUMatrixLReturnType< SCMatrix > matrixL() const </div><div class="ttdef"><b>Definition:</b> SparseLU.h:134</div></div> <div class="ttc" id="classEigen_1_1SparseMatrix_html_a5552abd83dbd03c85cea6d61fd8875a5"><div class="ttname"><a href="classEigen_1_1SparseMatrix.html#a5552abd83dbd03c85cea6d61fd8875a5">Eigen::SparseMatrix::rows</a></div><div class="ttdeci">Index rows() const </div><div class="ttdef"><b>Definition:</b> SparseMatrix.h:119</div></div> <div class="ttc" id="classEigen_1_1SparseLU_html_a6f2a135bd741b6b2f2558e6a404581ff"><div class="ttname"><a href="classEigen_1_1SparseLU.html#a6f2a135bd741b6b2f2558e6a404581ff">Eigen::SparseLU::analyzePattern</a></div><div class="ttdeci">void analyzePattern(const MatrixType &matrix)</div><div class="ttdef"><b>Definition:</b> SparseLU.h:369</div></div> <div class="ttc" id="classEigen_1_1SparseMatrix_html_aaca1908a5ec508a25ff0a8bca803e5f3"><div class="ttname"><a href="classEigen_1_1SparseMatrix.html#aaca1908a5ec508a25ff0a8bca803e5f3">Eigen::SparseMatrix::cols</a></div><div class="ttdeci">Index cols() const </div><div class="ttdef"><b>Definition:</b> SparseMatrix.h:121</div></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_1PermutationBase_html_a594deabd0f368a746801f5dd14a0db2a"><div class="ttname"><a href="classEigen_1_1PermutationBase.html#a594deabd0f368a746801f5dd14a0db2a">Eigen::PermutationBase::inverse</a></div><div class="ttdeci">Transpose< PermutationBase > inverse() const </div><div class="ttdef"><b>Definition:</b> PermutationMatrix.h:201</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_1SparseLU_html_a39062d90ea647b4e7668452647b04b17"><div class="ttname"><a href="classEigen_1_1SparseLU.html#a39062d90ea647b4e7668452647b04b17">Eigen::SparseLU::solve</a></div><div class="ttdeci">const internal::sparse_solve_retval< SparseLU, Rhs > solve(const SparseMatrixBase< Rhs > &B) const </div><div class="ttdef"><b>Definition:</b> SparseLU.h:191</div></div> <div class="ttc" id="classEigen_1_1SparseLU_html_a4c6b35f80ff9f6f0a3f4e74dfc121349"><div class="ttname"><a href="classEigen_1_1SparseLU.html#a4c6b35f80ff9f6f0a3f4e74dfc121349">Eigen::SparseLU::absDeterminant</a></div><div class="ttdeci">Scalar absDeterminant()</div><div class="ttdef"><b>Definition:</b> SparseLU.h:256</div></div> <div class="ttc" id="classEigen_1_1SparseLU_html_a33e2421e468033a883c041f940537a7c"><div class="ttname"><a href="classEigen_1_1SparseLU.html#a33e2421e468033a883c041f940537a7c">Eigen::SparseLU::factorize</a></div><div class="ttdeci">void factorize(const MatrixType &matrix)</div><div class="ttdef"><b>Definition:</b> SparseLU.h:452</div></div> <div class="ttc" id="classEigen_1_1SparseLU_html_adb61fce0a1343a45b65b6a50b18408dd"><div class="ttname"><a href="classEigen_1_1SparseLU.html#adb61fce0a1343a45b65b6a50b18408dd">Eigen::SparseLU::logAbsDeterminant</a></div><div class="ttdeci">Scalar logAbsDeterminant() const </div><div class="ttdef"><b>Definition:</b> SparseLU.h:286</div></div> <div class="ttc" id="classEigen_1_1DenseBase_html_a58c77695de3b33405f01f2fdf3dc389d"><div class="ttname"><a href="classEigen_1_1DenseBase.html#a58c77695de3b33405f01f2fdf3dc389d">Eigen::DenseBase::col</a></div><div class="ttdeci">ColXpr col(Index i)</div><div class="ttdef"><b>Definition:</b> DenseBase.h:709</div></div> <div class="ttc" id="classEigen_1_1SparseLU_html_ad856e0133763838bd626bfcd0cd1bf33"><div class="ttname"><a href="classEigen_1_1SparseLU.html#ad856e0133763838bd626bfcd0cd1bf33">Eigen::SparseLU::rowsPermutation</a></div><div class="ttdeci">const PermutationType & rowsPermutation() const </div><div class="ttdef"><b>Definition:</b> SparseLU.h:153</div></div> <div class="ttc" id="classEigen_1_1SparseLU_html"><div class="ttname"><a href="classEigen_1_1SparseLU.html">Eigen::SparseLU</a></div><div class="ttdoc">Sparse supernodal LU factorization for general matrices. </div><div class="ttdef"><b>Definition:</b> SparseLU.h:17</div></div> <div class="ttc" id="classEigen_1_1SparseLU_html_a1229744f4e8554ca6e96fe32ac359924"><div class="ttname"><a href="classEigen_1_1SparseLU.html#a1229744f4e8554ca6e96fe32ac359924">Eigen::SparseLU::colsPermutation</a></div><div class="ttdeci">const PermutationType & colsPermutation() const </div><div class="ttdef"><b>Definition:</b> SparseLU.h:161</div></div> <div class="ttc" id="namespaceEigen_1_1internal_html_a22b184331e5fb270a37d1305e95cb064"><div class="ttname"><a href="namespaceEigen_1_1internal.html#a22b184331e5fb270a37d1305e95cb064">Eigen::internal::coletree</a></div><div class="ttdeci">int coletree(const MatrixType &mat, IndexVector &parent, IndexVector &firstRowElt, typename MatrixType::Index *perm=0)</div><div class="ttdef"><b>Definition:</b> SparseColEtree.h:61</div></div> <div class="ttc" id="classEigen_1_1PermutationMatrix_html"><div class="ttname"><a href="classEigen_1_1PermutationMatrix.html">Eigen::PermutationMatrix< Dynamic, Dynamic, Index ></a></div></div> <div class="ttc" id="classEigen_1_1VectorBlock_html"><div class="ttname"><a href="classEigen_1_1VectorBlock.html">Eigen::VectorBlock</a></div><div class="ttdoc">Expression of a fixed-size or dynamic-size sub-vector. </div><div class="ttdef"><b>Definition:</b> ForwardDeclarations.h:83</div></div> <div class="ttc" id="classEigen_1_1SparseMatrixBase_html"><div class="ttname"><a href="classEigen_1_1SparseMatrixBase.html">Eigen::SparseMatrixBase</a></div><div class="ttdoc">Base class of any sparse matrices or sparse expressions. </div><div class="ttdef"><b>Definition:</b> SparseMatrixBase.h:26</div></div> <div class="ttc" id="structEigen_1_1EigenBase_html_aa84222add803ad7c9db07dd4dd91d5d9"><div class="ttname"><a href="structEigen_1_1EigenBase.html#aa84222add803ad7c9db07dd4dd91d5d9">Eigen::EigenBase::derived</a></div><div class="ttdeci">Derived & derived()</div><div class="ttdef"><b>Definition:</b> EigenBase.h:34</div></div> <div class="ttc" id="classEigen_1_1SparseLU_html_a8166735beae3bd23e7c3a1be9f0f5ef9"><div class="ttname"><a href="classEigen_1_1SparseLU.html#a8166735beae3bd23e7c3a1be9f0f5ef9">Eigen::SparseLU::solve</a></div><div class="ttdeci">const internal::solve_retval< SparseLU, Rhs > solve(const MatrixBase< Rhs > &B) const </div><div class="ttdef"><b>Definition:</b> SparseLU.h:178</div></div> <div class="ttc" id="classEigen_1_1SparseLU_html_acef7759add41200c1c817a6255e21dd4"><div class="ttname"><a href="classEigen_1_1SparseLU.html#acef7759add41200c1c817a6255e21dd4">Eigen::SparseLU::isSymmetric</a></div><div class="ttdeci">void isSymmetric(bool sym)</div><div class="ttdef"><b>Definition:</b> SparseLU.h:123</div></div> <div class="ttc" id="classEigen_1_1SparseLU_html_a05a0eac1c146e9e62df9ec70f5a6f69d"><div class="ttname"><a href="classEigen_1_1SparseLU.html#a05a0eac1c146e9e62df9ec70f5a6f69d">Eigen::SparseLU::compute</a></div><div class="ttdeci">void compute(const MatrixType &matrix)</div><div class="ttdef"><b>Definition:</b> SparseLU.h:112</div></div> <div class="ttc" id="classEigen_1_1PlainObjectBase_html_aabcf7b6f4543f6255253b9ccee3309b8"><div class="ttname"><a href="classEigen_1_1PlainObjectBase.html#aabcf7b6f4543f6255253b9ccee3309b8">Eigen::PlainObjectBase::setConstant</a></div><div class="ttdeci">Derived & setConstant(Index size, const Scalar &value)</div><div class="ttdef"><b>Definition:</b> CwiseNullaryOp.h:348</div></div> <div class="ttc" id="classEigen_1_1SparseLU_html_a6001a3f3f7ad66583d9d1ce3bb3de262"><div class="ttname"><a href="classEigen_1_1SparseLU.html#a6001a3f3f7ad66583d9d1ce3bb3de262">Eigen::SparseLU::matrixU</a></div><div class="ttdeci">SparseLUMatrixUReturnType< SCMatrix, MappedSparseMatrix< Scalar, ColMajor, Index > > matrixU() const </div><div class="ttdef"><b>Definition:</b> SparseLU.h:144</div></div> <div class="ttc" id="classEigen_1_1SparseLU_html_a0c06d5c2034ebb329c54235369643ad2"><div class="ttname"><a href="classEigen_1_1SparseLU.html#a0c06d5c2034ebb329c54235369643ad2">Eigen::SparseLU::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> SparseLU.h:207</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="group__flags_html_ga7bd49e7b260e869e10fb9dc4fd081a85"><div class="ttname"><a href="group__flags.html#ga7bd49e7b260e869e10fb9dc4fd081a85">Eigen::RowMajorBit</a></div><div class="ttdeci">const unsigned int RowMajorBit</div><div class="ttdef"><b>Definition:</b> Constants.h:53</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_1SparseLU_html_a040b16815cde46c7f2f6d26a2298ca4b"><div class="ttname"><a href="classEigen_1_1SparseLU.html#a040b16815cde46c7f2f6d26a2298ca4b">Eigen::SparseLU::lastErrorMessage</a></div><div class="ttdeci">std::string lastErrorMessage() const </div><div class="ttdef"><b>Definition:</b> SparseLU.h:216</div></div> <div class="ttc" id="classEigen_1_1Matrix_html"><div class="ttname"><a href="classEigen_1_1Matrix.html">Eigen::Matrix< Scalar, Dynamic, 1 ></a></div></div> <div class="ttc" id="classEigen_1_1SparseMatrixBase_html_a5552abd83dbd03c85cea6d61fd8875a5"><div class="ttname"><a href="classEigen_1_1SparseMatrixBase.html#a5552abd83dbd03c85cea6d61fd8875a5">Eigen::SparseMatrixBase::rows</a></div><div class="ttdeci">Index rows() const </div><div class="ttdef"><b>Definition:</b> SparseMatrixBase.h:150</div></div> <div class="ttc" id="classEigen_1_1SparseLU_html_a20467be6f23d8b39574f3e097583d767"><div class="ttname"><a href="classEigen_1_1SparseLU.html#a20467be6f23d8b39574f3e097583d767">Eigen::SparseLU::setPivotThreshold</a></div><div class="ttdeci">void setPivotThreshold(const RealScalar &thresh)</div><div class="ttdef"><b>Definition:</b> SparseLU.h:166</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_1MappedSparseMatrix_html"><div class="ttname"><a href="classEigen_1_1MappedSparseMatrix.html">Eigen::MappedSparseMatrix< Scalar, ColMajor, Index ></a></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_1SparseLU_html_a7150c533f7f2f6bdb1e9c82727f426f5"><div class="ttname"><a href="classEigen_1_1SparseLU.html#a7150c533f7f2f6bdb1e9c82727f426f5">Eigen::SparseLU::signDeterminant</a></div><div class="ttdeci">Scalar signDeterminant()</div><div class="ttdef"><b>Definition:</b> SparseLU.h:309</div></div> <div class="ttc" id="namespaceEigen_1_1internal_html_a6afa7df34de3dce07b21194168118456"><div class="ttname"><a href="namespaceEigen_1_1internal.html#a6afa7df34de3dce07b21194168118456">Eigen::internal::treePostorder</a></div><div class="ttdeci">void treePostorder(Index n, IndexVector &parent, IndexVector &post)</div><div class="ttdoc">Post order a tree. </div><div class="ttdef"><b>Definition:</b> SparseColEtree.h:176</div></div> <div class="ttc" id="classEigen_1_1PlainObjectBase_html_a04abe84a9a894de335a232681d9a0722"><div class="ttname"><a href="classEigen_1_1PlainObjectBase.html#a04abe84a9a894de335a232681d9a0722">Eigen::PlainObjectBase::setZero</a></div><div class="ttdeci">Derived & setZero(Index size)</div><div class="ttdef"><b>Definition:</b> CwiseNullaryOp.h:515</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_43b7fa11304fca39592f14c296504e53.html">SparseLU</a></li><li class="navelem"><b>SparseLU.h</b></li> <li class="footer">Generated on Mon Oct 28 2013 11:04:25 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>