Sophie

Sophie

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

eigen3-doc-3.2-3.mga4.noarch.rpm

<!DOCTYPE html PUBLIC "-//W3C//DTD XHTML 1.0 Transitional//EN" "http://www.w3.org/TR/xhtml1/DTD/xhtml1-transitional.dtd">
<html xmlns="http://www.w3.org/1999/xhtml">
<head>
<meta http-equiv="Content-Type" content="text/xhtml;charset=UTF-8"/>
<meta http-equiv="X-UA-Compatible" content="IE=9"/>
<meta name="generator" content="Doxygen 1.8.5"/>
<title>Eigen: SparseLU&lt; _MatrixType, _OrderingType &gt; Class Template Reference</title>
<link href="tabs.css" rel="stylesheet" type="text/css"/>
<script type="text/javascript" src="jquery.js"></script>
<script type="text/javascript" src="dynsections.js"></script>
<link href="navtree.css" rel="stylesheet" type="text/css"/>
<script type="text/javascript" src="resize.js"></script>
<script type="text/javascript" src="navtree.js"></script>
<script type="text/javascript">
  $(document).ready(initResizable);
  $(window).load(resizeHeight);
</script>
<link href="search/search.css" rel="stylesheet" type="text/css"/>
<script type="text/javascript" src="search/search.js"></script>
<script type="text/javascript">
  $(document).ready(function() { searchBox.OnSelectItem(0); });
</script>
<link href="doxygen.css"   rel="stylesheet" type="text/css" />
<link href="eigendoxy.css" rel="stylesheet" type="text/css">
<!--  -->
<script type="text/javascript" src="eigen_navtree_hacks.js"></script>
<!-- <script type="text/javascript"> -->
<!-- </script> -->
</head>
<body>
<div id="top"><!-- do not remove this div, it is closed by doxygen! -->
<!-- <a name="top"></a> -->
<div id="titlearea">
<table cellspacing="0" cellpadding="0">
 <tbody>
 <tr style="height: 56px;">
  <td id="projectlogo"><img alt="Logo" src="Eigen_Silly_Professor_64x64.png"/></td>
  <td style="padding-left: 0.5em;">
   <div id="projectname"><a href="http://eigen.tuxfamily.org">Eigen</a>
   &#160;<span id="projectnumber">3.2.0</span>
   </div>
  </td>
   <td>        <div id="MSearchBox" class="MSearchBoxInactive">
        <span class="left">
          <img id="MSearchSelect" src="search/mag_sel.png"
               onmouseover="return searchBox.OnSearchSelectShow()"
               onmouseout="return searchBox.OnSearchSelectHide()"
               alt=""/>
          <input type="text" id="MSearchField" value="Search" accesskey="S"
               onfocus="searchBox.OnSearchFieldFocus(true)" 
               onblur="searchBox.OnSearchFieldFocus(false)" 
               onkeyup="searchBox.OnSearchFieldChange(event)"/>
          </span><span class="right">
            <a id="MSearchClose" href="javascript:searchBox.CloseResultsWindow()"><img id="MSearchCloseImg" border="0" src="search/close.png" alt=""/></a>
          </span>
        </div>
</td>
 </tr>
 </tbody>
</table>
</div>
<!-- end header part -->
<!-- Generated by Doxygen 1.8.5 -->
<script type="text/javascript">
var searchBox = new SearchBox("searchBox", "search",false,'Search');
</script>
</div><!-- top -->
<div id="side-nav" class="ui-resizable side-nav-resizable">
  <div id="nav-tree">
    <div id="nav-tree-contents">
      <div id="nav-sync" class="sync"></div>
    </div>
  </div>
  <div id="splitbar" style="-moz-user-select:none;" 
       class="ui-resizable-handle">
  </div>
</div>
<script type="text/javascript">
$(document).ready(function(){initNavTree('classEigen_1_1SparseLU.html','');});
</script>
<div id="doc-content">
<!-- window showing the filter options -->
<div id="MSearchSelectWindow"
     onmouseover="return searchBox.OnSearchSelectShow()"
     onmouseout="return searchBox.OnSearchSelectHide()"
     onkeydown="return searchBox.OnSearchSelectKey(event)">
<a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(0)"><span class="SelectionMark">&#160;</span>All</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(1)"><span class="SelectionMark">&#160;</span>Classes</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(2)"><span class="SelectionMark">&#160;</span>Namespaces</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(3)"><span class="SelectionMark">&#160;</span>Functions</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(4)"><span class="SelectionMark">&#160;</span>Variables</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(5)"><span class="SelectionMark">&#160;</span>Typedefs</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(6)"><span class="SelectionMark">&#160;</span>Enumerations</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(7)"><span class="SelectionMark">&#160;</span>Enumerator</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(8)"><span class="SelectionMark">&#160;</span>Friends</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(9)"><span class="SelectionMark">&#160;</span>Groups</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(10)"><span class="SelectionMark">&#160;</span>Pages</a></div>

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

<div class="header">
  <div class="summary">
<a href="classEigen_1_1SparseLU-members.html">List of all members</a> &#124;
<a href="#pub-methods">Public Member Functions</a>  </div>
  <div class="headertitle">
<div class="title">SparseLU&lt; _MatrixType, _OrderingType &gt; Class Template Reference<div class="ingroups"><a class="el" href="group__SparseLU__Module.html">SparseLU module</a></div></div>  </div>
</div><!--header-->
<div class="contents">
<a name="details" id="details"></a><h2 class="groupheader">Detailed Description</h2>
<div class="textblock"><h3>template&lt;typename _MatrixType, typename _OrderingType = COLAMDOrdering&lt;typename _MatrixType::Index&gt;&gt;<br/>
class Eigen::SparseLU&lt; _MatrixType, _OrderingType &gt;</h3>

<p><a class="el" href="structEigen_1_1Sparse.html">Sparse</a> supernodal LU factorization for general matrices. </p>
<p>This class implements the supernodal LU factorization for general matrices. It uses the main techniques from the sequential <a class="el" href="classEigen_1_1SuperLU.html" title="A sparse direct LU factorization and solver based on the SuperLU library. ">SuperLU</a> package (<a href="http://crd-legacy.lbl.gov/~xiaoye/SuperLU/">http://crd-legacy.lbl.gov/~xiaoye/SuperLU/</a>). It handles transparently real and complex arithmetics with single and double precision, depending on the scalar type of your input matrix. The code has been optimized to provide BLAS-3 operations during supernode-panel updates. It benefits directly from the built-in high-performant <a class="el" href="namespaceEigen.html" title="Namespace containing all symbols from the Eigen library. ">Eigen</a> BLAS routines. Moreover, when the size of a supernode is very small, the BLAS calls are avoided to enable a better optimization from the compiler. For best performance, you should compile it with NDEBUG flag to avoid the numerous bounds checking on vectors.</p>
<p>An important parameter of this class is the ordering method. It is used to reorder the columns (and eventually the rows) of the matrix to reduce the number of new elements that are created during numerical factorization. The cheapest method available is COLAMD. See <a class="el" href="group__OrderingMethods__Module.html">the OrderingMethods module </a> for the list of built-in and external ordering methods.</p>
<p>Simple example with key steps </p>
<div class="fragment"><div class="line">* <a class="code" href="group__matrixtypedefs.html#ga3da45e59796fbacf67fa568297927bd1">VectorXd</a> x(n), b(n);</div>
<div class="line">* SparseMatrix&lt;double, ColMajor&gt; A;</div>
<div class="line">* SparseLU&lt;SparseMatrix&lt;scalar, ColMajor&gt;, COLAMDOrdering&lt;Index&gt; &gt;   solver;</div>
<div class="line">* <span class="comment">// fill A and b;</span></div>
<div class="line">* <span class="comment">// Compute the ordering permutation vector from the structural pattern of A</span></div>
<div class="line">* solver.analyzePattern(A); </div>
<div class="line">* <span class="comment">// Compute the numerical factorization </span></div>
<div class="line">* solver.factorize(A); </div>
<div class="line">* <span class="comment">//Use the factors to solve the linear system </span></div>
<div class="line">* x = solver.solve(b); </div>
<div class="line">* </div>
</div><!-- fragment --><dl class="section warning"><dt>Warning</dt><dd>The input matrix A should be in a <b>compressed</b> and <b>column-major</b> form. Otherwise an expensive copy will be made. You can call the inexpensive makeCompressed() to get a compressed matrix.</dd></dl>
<dl class="section note"><dt>Note</dt><dd>Unlike the initial <a class="el" href="classEigen_1_1SuperLU.html" title="A sparse direct LU factorization and solver based on the SuperLU library. ">SuperLU</a> implementation, there is no step to equilibrate the matrix. For badly scaled matrices, this step can be useful to reduce the pivoting during factorization. If this is the case for your matrices, you can try the basic scaling method at "unsupported/Eigen/src/IterativeSolvers/Scaling.h"</dd></dl>
<dl class="tparams"><dt>Template Parameters</dt><dd>
  <table class="tparams">
    <tr><td class="paramname">_MatrixType</td><td>The type of the sparse matrix. It must be a column-major SparseMatrix&lt;&gt; </td></tr>
    <tr><td class="paramname">_OrderingType</td><td>The ordering method to use, either AMD, COLAMD or METIS. Default is COLMAD</td></tr>
  </table>
  </dd>
</dl>
<dl class="section see"><dt>See Also</dt><dd><a class="el" href="group__TopicSparseSystems.html#TutorialSparseDirectSolvers">Sparse solvers</a> </dd>
<dd>
<a class="el" href="group__OrderingMethods__Module.html">OrderingMethods module</a> </dd></dl>
</div>
<p>Inherits SparseLUImpl&lt; _MatrixType::Scalar, _MatrixType::Index &gt;.</p>
<table class="memberdecls">
<tr class="heading"><td colspan="2"><h2 class="groupheader"><a name="pub-methods"></a>
Public Member Functions</h2></td></tr>
<tr class="memitem:a4c6b35f80ff9f6f0a3f4e74dfc121349"><td class="memItemLeft" align="right" valign="top">Scalar&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1SparseLU.html#a4c6b35f80ff9f6f0a3f4e74dfc121349">absDeterminant</a> ()</td></tr>
<tr class="separator:a4c6b35f80ff9f6f0a3f4e74dfc121349"><td class="memSeparator" colspan="2">&#160;</td></tr>
<tr class="memitem:a6f2a135bd741b6b2f2558e6a404581ff"><td class="memItemLeft" align="right" valign="top">void&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1SparseLU.html#a6f2a135bd741b6b2f2558e6a404581ff">analyzePattern</a> (const MatrixType &amp;matrix)</td></tr>
<tr class="separator:a6f2a135bd741b6b2f2558e6a404581ff"><td class="memSeparator" colspan="2">&#160;</td></tr>
<tr class="memitem:a1229744f4e8554ca6e96fe32ac359924"><td class="memItemLeft" align="right" valign="top">const <a class="el" href="classEigen_1_1PermutationMatrix.html">PermutationType</a> &amp;&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1SparseLU.html#a1229744f4e8554ca6e96fe32ac359924">colsPermutation</a> () const </td></tr>
<tr class="separator:a1229744f4e8554ca6e96fe32ac359924"><td class="memSeparator" colspan="2">&#160;</td></tr>
<tr class="memitem:a05a0eac1c146e9e62df9ec70f5a6f69d"><td class="memItemLeft" align="right" valign="top">void&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1SparseLU.html#a05a0eac1c146e9e62df9ec70f5a6f69d">compute</a> (const MatrixType &amp;matrix)</td></tr>
<tr class="separator:a05a0eac1c146e9e62df9ec70f5a6f69d"><td class="memSeparator" colspan="2">&#160;</td></tr>
<tr class="memitem:a33e2421e468033a883c041f940537a7c"><td class="memItemLeft" align="right" valign="top">void&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1SparseLU.html#a33e2421e468033a883c041f940537a7c">factorize</a> (const MatrixType &amp;matrix)</td></tr>
<tr class="separator:a33e2421e468033a883c041f940537a7c"><td class="memSeparator" colspan="2">&#160;</td></tr>
<tr class="memitem:a0c06d5c2034ebb329c54235369643ad2"><td class="memItemLeft" align="right" valign="top"><a class="el" href="group__enums.html#ga51bc1ac16f26ebe51eae1abb77bd037b">ComputationInfo</a>&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1SparseLU.html#a0c06d5c2034ebb329c54235369643ad2">info</a> () const </td></tr>
<tr class="memdesc:a0c06d5c2034ebb329c54235369643ad2"><td class="mdescLeft">&#160;</td><td class="mdescRight">Reports whether previous computation was successful.  <a href="#a0c06d5c2034ebb329c54235369643ad2">More...</a><br/></td></tr>
<tr class="separator:a0c06d5c2034ebb329c54235369643ad2"><td class="memSeparator" colspan="2">&#160;</td></tr>
<tr class="memitem:acef7759add41200c1c817a6255e21dd4"><td class="memItemLeft" align="right" valign="top">void&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1SparseLU.html#acef7759add41200c1c817a6255e21dd4">isSymmetric</a> (bool sym)</td></tr>
<tr class="separator:acef7759add41200c1c817a6255e21dd4"><td class="memSeparator" colspan="2">&#160;</td></tr>
<tr class="memitem:a040b16815cde46c7f2f6d26a2298ca4b"><td class="memItemLeft" align="right" valign="top">std::string&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1SparseLU.html#a040b16815cde46c7f2f6d26a2298ca4b">lastErrorMessage</a> () const </td></tr>
<tr class="separator:a040b16815cde46c7f2f6d26a2298ca4b"><td class="memSeparator" colspan="2">&#160;</td></tr>
<tr class="memitem:adb61fce0a1343a45b65b6a50b18408dd"><td class="memItemLeft" align="right" valign="top">Scalar&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1SparseLU.html#adb61fce0a1343a45b65b6a50b18408dd">logAbsDeterminant</a> () const </td></tr>
<tr class="separator:adb61fce0a1343a45b65b6a50b18408dd"><td class="memSeparator" colspan="2">&#160;</td></tr>
<tr class="memitem:a828ca6373b429d81f41f41429a7fd4ad"><td class="memItemLeft" align="right" valign="top">SparseLUMatrixLReturnType<br class="typebreak"/>
&lt; SCMatrix &gt;&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1SparseLU.html#a828ca6373b429d81f41f41429a7fd4ad">matrixL</a> () const </td></tr>
<tr class="separator:a828ca6373b429d81f41f41429a7fd4ad"><td class="memSeparator" colspan="2">&#160;</td></tr>
<tr class="memitem:a6001a3f3f7ad66583d9d1ce3bb3de262"><td class="memItemLeft" align="right" valign="top">SparseLUMatrixUReturnType<br class="typebreak"/>
&lt; SCMatrix, <a class="el" href="classEigen_1_1MappedSparseMatrix.html">MappedSparseMatrix</a><br class="typebreak"/>
&lt; Scalar, <a class="el" href="group__enums.html#gga0c5bde183ecefe103f70b49ad9740bcdac86184b0e3be936fbfd20249a057a0bf">ColMajor</a>, Index &gt; &gt;&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1SparseLU.html#a6001a3f3f7ad66583d9d1ce3bb3de262">matrixU</a> () const </td></tr>
<tr class="separator:a6001a3f3f7ad66583d9d1ce3bb3de262"><td class="memSeparator" colspan="2">&#160;</td></tr>
<tr class="memitem:ad856e0133763838bd626bfcd0cd1bf33"><td class="memItemLeft" align="right" valign="top">const <a class="el" href="classEigen_1_1PermutationMatrix.html">PermutationType</a> &amp;&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1SparseLU.html#ad856e0133763838bd626bfcd0cd1bf33">rowsPermutation</a> () const </td></tr>
<tr class="separator:ad856e0133763838bd626bfcd0cd1bf33"><td class="memSeparator" colspan="2">&#160;</td></tr>
<tr class="memitem:a20467be6f23d8b39574f3e097583d767"><td class="memItemLeft" align="right" valign="top">void&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1SparseLU.html#a20467be6f23d8b39574f3e097583d767">setPivotThreshold</a> (const RealScalar &amp;thresh)</td></tr>
<tr class="separator:a20467be6f23d8b39574f3e097583d767"><td class="memSeparator" colspan="2">&#160;</td></tr>
<tr class="memitem:a7150c533f7f2f6bdb1e9c82727f426f5"><td class="memItemLeft" align="right" valign="top">Scalar&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1SparseLU.html#a7150c533f7f2f6bdb1e9c82727f426f5">signDeterminant</a> ()</td></tr>
<tr class="separator:a7150c533f7f2f6bdb1e9c82727f426f5"><td class="memSeparator" colspan="2">&#160;</td></tr>
<tr class="memitem:a8166735beae3bd23e7c3a1be9f0f5ef9"><td class="memTemplParams" colspan="2">template&lt;typename Rhs &gt; </td></tr>
<tr class="memitem:a8166735beae3bd23e7c3a1be9f0f5ef9"><td class="memTemplItemLeft" align="right" valign="top">const internal::solve_retval<br class="typebreak"/>
&lt; <a class="el" href="classEigen_1_1SparseLU.html">SparseLU</a>, Rhs &gt;&#160;</td><td class="memTemplItemRight" valign="bottom"><a class="el" href="classEigen_1_1SparseLU.html#a8166735beae3bd23e7c3a1be9f0f5ef9">solve</a> (const <a class="el" href="classEigen_1_1MatrixBase.html">MatrixBase</a>&lt; Rhs &gt; &amp;B) const </td></tr>
<tr class="separator:a8166735beae3bd23e7c3a1be9f0f5ef9"><td class="memSeparator" colspan="2">&#160;</td></tr>
<tr class="memitem:a39062d90ea647b4e7668452647b04b17"><td class="memTemplParams" colspan="2">template&lt;typename Rhs &gt; </td></tr>
<tr class="memitem:a39062d90ea647b4e7668452647b04b17"><td class="memTemplItemLeft" align="right" valign="top">const <br class="typebreak"/>
internal::sparse_solve_retval<br class="typebreak"/>
&lt; <a class="el" href="classEigen_1_1SparseLU.html">SparseLU</a>, Rhs &gt;&#160;</td><td class="memTemplItemRight" valign="bottom"><a class="el" href="classEigen_1_1SparseLU.html#a39062d90ea647b4e7668452647b04b17">solve</a> (const <a class="el" href="classEigen_1_1SparseMatrixBase.html">SparseMatrixBase</a>&lt; Rhs &gt; &amp;B) const </td></tr>
<tr class="separator:a39062d90ea647b4e7668452647b04b17"><td class="memSeparator" colspan="2">&#160;</td></tr>
</table>
<h2 class="groupheader">Member Function Documentation</h2>
<a class="anchor" id="a4c6b35f80ff9f6f0a3f4e74dfc121349"></a>
<div class="memitem">
<div class="memproto">
<table class="mlabels">
  <tr>
  <td class="mlabels-left">
      <table class="memname">
        <tr>
          <td class="memname">Scalar absDeterminant </td>
          <td>(</td>
          <td class="paramname"></td><td>)</td>
          <td></td>
        </tr>
      </table>
  </td>
  <td class="mlabels-right">
<span class="mlabels"><span class="mlabel">inline</span></span>  </td>
  </tr>
</table>
</div><div class="memdoc">
<dl class="section return"><dt>Returns</dt><dd>the absolute value of the determinant of the matrix of which *this is the QR decomposition.</dd></dl>
<dl class="section warning"><dt>Warning</dt><dd>a determinant can be very big or small, so for matrices of large enough dimension, there is a risk of overflow/underflow. One way to work around that is to use <a class="el" href="classEigen_1_1SparseLU.html#adb61fce0a1343a45b65b6a50b18408dd">logAbsDeterminant()</a> instead.</dd></dl>
<dl class="section see"><dt>See Also</dt><dd><a class="el" href="classEigen_1_1SparseLU.html#adb61fce0a1343a45b65b6a50b18408dd">logAbsDeterminant()</a>, <a class="el" href="classEigen_1_1SparseLU.html#a7150c533f7f2f6bdb1e9c82727f426f5">signDeterminant()</a> </dd></dl>

</div>
</div>
<a class="anchor" id="a6f2a135bd741b6b2f2558e6a404581ff"></a>
<div class="memitem">
<div class="memproto">
      <table class="memname">
        <tr>
          <td class="memname">void analyzePattern </td>
          <td>(</td>
          <td class="paramtype">const MatrixType &amp;&#160;</td>
          <td class="paramname"><em>mat</em></td><td>)</td>
          <td></td>
        </tr>
      </table>
</div><div class="memdoc">
<p>Compute the column permutation to minimize the fill-in</p>
<ul>
<li>Apply this permutation to the input matrix -</li>
<li>Compute the column elimination tree on the permuted matrix</li>
<li>Postorder the elimination tree and the column permutation </li>
</ul>

<p>References <a class="el" href="namespaceEigen_1_1internal.html#a22b184331e5fb270a37d1305e95cb064">Eigen::internal::coletree()</a>, <a class="el" href="classEigen_1_1PermutationMatrix.html#a8876d615d17aad77b054a8f58b699e7d">PermutationMatrix&lt; SizeAtCompileTime, MaxSizeAtCompileTime, IndexType &gt;::indices()</a>, <a class="el" href="classEigen_1_1PlainObjectBase.html#afbbb33d14fe7fb9683019a39ce1c659d">PlainObjectBase&lt; Derived &gt;::resize()</a>, and <a class="el" href="namespaceEigen_1_1internal.html#a6afa7df34de3dce07b21194168118456">Eigen::internal::treePostorder()</a>.</p>

<p>Referenced by <a class="el" href="classEigen_1_1SparseLU.html#a05a0eac1c146e9e62df9ec70f5a6f69d">SparseLU&lt; _MatrixType, _OrderingType &gt;::compute()</a>.</p>

</div>
</div>
<a class="anchor" id="a1229744f4e8554ca6e96fe32ac359924"></a>
<div class="memitem">
<div class="memproto">
<table class="mlabels">
  <tr>
  <td class="mlabels-left">
      <table class="memname">
        <tr>
          <td class="memname">const <a class="el" href="classEigen_1_1PermutationMatrix.html">PermutationType</a>&amp; colsPermutation </td>
          <td>(</td>
          <td class="paramname"></td><td>)</td>
          <td> const</td>
        </tr>
      </table>
  </td>
  <td class="mlabels-right">
<span class="mlabels"><span class="mlabel">inline</span></span>  </td>
  </tr>
</table>
</div><div class="memdoc">
<dl class="section return"><dt>Returns</dt><dd>a reference to the column matrix permutation <img class="formulaInl" alt="$ P_c^T $" src="form_161.png"/> such that <img class="formulaInl" alt="$P_r A P_c^T = L U$" src="form_160.png"/> </dd></dl>
<dl class="section see"><dt>See Also</dt><dd><a class="el" href="classEigen_1_1SparseLU.html#ad856e0133763838bd626bfcd0cd1bf33">rowsPermutation()</a> </dd></dl>

</div>
</div>
<a class="anchor" id="a05a0eac1c146e9e62df9ec70f5a6f69d"></a>
<div class="memitem">
<div class="memproto">
<table class="mlabels">
  <tr>
  <td class="mlabels-left">
      <table class="memname">
        <tr>
          <td class="memname">void compute </td>
          <td>(</td>
          <td class="paramtype">const MatrixType &amp;&#160;</td>
          <td class="paramname"><em>matrix</em></td><td>)</td>
          <td></td>
        </tr>
      </table>
  </td>
  <td class="mlabels-right">
<span class="mlabels"><span class="mlabel">inline</span></span>  </td>
  </tr>
</table>
</div><div class="memdoc">
<p>Compute the symbolic and numeric factorization of the input sparse matrix. The input matrix should be in column-major storage. </p>

<p>References <a class="el" href="classEigen_1_1SparseLU.html#a6f2a135bd741b6b2f2558e6a404581ff">SparseLU&lt; _MatrixType, _OrderingType &gt;::analyzePattern()</a>, and <a class="el" href="classEigen_1_1SparseLU.html#a33e2421e468033a883c041f940537a7c">SparseLU&lt; _MatrixType, _OrderingType &gt;::factorize()</a>.</p>

</div>
</div>
<a class="anchor" id="a33e2421e468033a883c041f940537a7c"></a>
<div class="memitem">
<div class="memproto">
      <table class="memname">
        <tr>
          <td class="memname">void factorize </td>
          <td>(</td>
          <td class="paramtype">const MatrixType &amp;&#160;</td>
          <td class="paramname"><em>matrix</em></td><td>)</td>
          <td></td>
        </tr>
      </table>
</div><div class="memdoc">
<ul>
<li>Numerical factorization</li>
<li><p class="startli">Interleaved with the symbolic factorization On exit, info is</p>
<p class="startli">= 0: successful factorization</p>
</li>
</ul>
<blockquote class="doxtable">
<p>0: if info = i, and i is</p>
<p></p>
</blockquote>
<pre class="fragment">  &lt;= A-&gt;ncol: U(i,i) is exactly zero. The factorization has
     been completed, but the factor U is exactly singular,
     and division by zero will occur if it is used to solve a
     system of equations.

  &gt; A-&gt;ncol: number of bytes allocated when memory allocation
    failure occurred, plus A-&gt;ncol. If lwork = -1, it is
    the estimated amount of space needed, plus A-&gt;ncol.  </pre> 
<p>References <a class="el" href="group__enums.html#gga51bc1ac16f26ebe51eae1abb77bd037ba710fff14e8fc77846d4b75d8f4cc2d5c">Eigen::NumericalIssue</a>, <a class="el" href="classEigen_1_1PlainObjectBase.html#aabcf7b6f4543f6255253b9ccee3309b8">PlainObjectBase&lt; Derived &gt;::setConstant()</a>, <a class="el" href="classEigen_1_1PlainObjectBase.html#a04abe84a9a894de335a232681d9a0722">PlainObjectBase&lt; Derived &gt;::setZero()</a>, and <a class="el" href="group__enums.html#gga51bc1ac16f26ebe51eae1abb77bd037bafdfbdf3247bd36a1f17270d5cec74c9c">Eigen::Success</a>.</p>

<p>Referenced by <a class="el" href="classEigen_1_1SparseLU.html#a05a0eac1c146e9e62df9ec70f5a6f69d">SparseLU&lt; _MatrixType, _OrderingType &gt;::compute()</a>.</p>

</div>
</div>
<a class="anchor" id="a0c06d5c2034ebb329c54235369643ad2"></a>
<div class="memitem">
<div class="memproto">
<table class="mlabels">
  <tr>
  <td class="mlabels-left">
      <table class="memname">
        <tr>
          <td class="memname"><a class="el" href="group__enums.html#ga51bc1ac16f26ebe51eae1abb77bd037b">ComputationInfo</a> info </td>
          <td>(</td>
          <td class="paramname"></td><td>)</td>
          <td> const</td>
        </tr>
      </table>
  </td>
  <td class="mlabels-right">
<span class="mlabels"><span class="mlabel">inline</span></span>  </td>
  </tr>
</table>
</div><div class="memdoc">

<p>Reports whether previous computation was successful. </p>
<dl class="section return"><dt>Returns</dt><dd><code>Success</code> if computation was succesful, <code>NumericalIssue</code> if the LU factorization reports a problem, zero diagonal for instance <code>InvalidInput</code> if the input matrix is invalid</dd></dl>
<dl class="section see"><dt>See Also</dt><dd>iparm() </dd></dl>

</div>
</div>
<a class="anchor" id="acef7759add41200c1c817a6255e21dd4"></a>
<div class="memitem">
<div class="memproto">
<table class="mlabels">
  <tr>
  <td class="mlabels-left">
      <table class="memname">
        <tr>
          <td class="memname">void isSymmetric </td>
          <td>(</td>
          <td class="paramtype">bool&#160;</td>
          <td class="paramname"><em>sym</em></td><td>)</td>
          <td></td>
        </tr>
      </table>
  </td>
  <td class="mlabels-right">
<span class="mlabels"><span class="mlabel">inline</span></span>  </td>
  </tr>
</table>
</div><div class="memdoc">
<p>Indicate that the pattern of the input matrix is symmetric </p>

</div>
</div>
<a class="anchor" id="a040b16815cde46c7f2f6d26a2298ca4b"></a>
<div class="memitem">
<div class="memproto">
<table class="mlabels">
  <tr>
  <td class="mlabels-left">
      <table class="memname">
        <tr>
          <td class="memname">std::string lastErrorMessage </td>
          <td>(</td>
          <td class="paramname"></td><td>)</td>
          <td> const</td>
        </tr>
      </table>
  </td>
  <td class="mlabels-right">
<span class="mlabels"><span class="mlabel">inline</span></span>  </td>
  </tr>
</table>
</div><div class="memdoc">
<dl class="section return"><dt>Returns</dt><dd>A string describing the type of error </dd></dl>

</div>
</div>
<a class="anchor" id="adb61fce0a1343a45b65b6a50b18408dd"></a>
<div class="memitem">
<div class="memproto">
<table class="mlabels">
  <tr>
  <td class="mlabels-left">
      <table class="memname">
        <tr>
          <td class="memname">Scalar logAbsDeterminant </td>
          <td>(</td>
          <td class="paramname"></td><td>)</td>
          <td> const</td>
        </tr>
      </table>
  </td>
  <td class="mlabels-right">
<span class="mlabels"><span class="mlabel">inline</span></span>  </td>
  </tr>
</table>
</div><div class="memdoc">
<dl class="section return"><dt>Returns</dt><dd>the natural log of the absolute value of the determinant of the matrix of which **this is the QR decomposition</dd></dl>
<dl class="section note"><dt>Note</dt><dd>This method is useful to work around the risk of overflow/underflow that's inherent to the determinant computation.</dd></dl>
<dl class="section see"><dt>See Also</dt><dd><a class="el" href="classEigen_1_1SparseLU.html#a4c6b35f80ff9f6f0a3f4e74dfc121349">absDeterminant()</a>, <a class="el" href="classEigen_1_1SparseLU.html#a7150c533f7f2f6bdb1e9c82727f426f5">signDeterminant()</a> </dd></dl>

</div>
</div>
<a class="anchor" id="a828ca6373b429d81f41f41429a7fd4ad"></a>
<div class="memitem">
<div class="memproto">
<table class="mlabels">
  <tr>
  <td class="mlabels-left">
      <table class="memname">
        <tr>
          <td class="memname">SparseLUMatrixLReturnType&lt;SCMatrix&gt; matrixL </td>
          <td>(</td>
          <td class="paramname"></td><td>)</td>
          <td> const</td>
        </tr>
      </table>
  </td>
  <td class="mlabels-right">
<span class="mlabels"><span class="mlabel">inline</span></span>  </td>
  </tr>
</table>
</div><div class="memdoc">
<dl class="section return"><dt>Returns</dt><dd>an expression of the matrix L, internally stored as supernodes The only operation available with this expression is the triangular solve <div class="fragment"><div class="line">* y = b; <a class="code" href="classEigen_1_1SparseLU.html#a828ca6373b429d81f41f41429a7fd4ad">matrixL</a>().solveInPlace(y);</div>
<div class="line">* </div>
</div><!-- fragment --> </dd></dl>

</div>
</div>
<a class="anchor" id="a6001a3f3f7ad66583d9d1ce3bb3de262"></a>
<div class="memitem">
<div class="memproto">
<table class="mlabels">
  <tr>
  <td class="mlabels-left">
      <table class="memname">
        <tr>
          <td class="memname">SparseLUMatrixUReturnType&lt;SCMatrix,<a class="el" href="classEigen_1_1MappedSparseMatrix.html">MappedSparseMatrix</a>&lt;Scalar,<a class="el" href="group__enums.html#gga0c5bde183ecefe103f70b49ad9740bcdac86184b0e3be936fbfd20249a057a0bf">ColMajor</a>,Index&gt; &gt; matrixU </td>
          <td>(</td>
          <td class="paramname"></td><td>)</td>
          <td> const</td>
        </tr>
      </table>
  </td>
  <td class="mlabels-right">
<span class="mlabels"><span class="mlabel">inline</span></span>  </td>
  </tr>
</table>
</div><div class="memdoc">
<dl class="section return"><dt>Returns</dt><dd>an expression of the matrix U, The only operation available with this expression is the triangular solve <div class="fragment"><div class="line">* y = b; <a class="code" href="classEigen_1_1SparseLU.html#a6001a3f3f7ad66583d9d1ce3bb3de262">matrixU</a>().solveInPlace(y);</div>
<div class="line">* </div>
</div><!-- fragment --> </dd></dl>

</div>
</div>
<a class="anchor" id="ad856e0133763838bd626bfcd0cd1bf33"></a>
<div class="memitem">
<div class="memproto">
<table class="mlabels">
  <tr>
  <td class="mlabels-left">
      <table class="memname">
        <tr>
          <td class="memname">const <a class="el" href="classEigen_1_1PermutationMatrix.html">PermutationType</a>&amp; rowsPermutation </td>
          <td>(</td>
          <td class="paramname"></td><td>)</td>
          <td> const</td>
        </tr>
      </table>
  </td>
  <td class="mlabels-right">
<span class="mlabels"><span class="mlabel">inline</span></span>  </td>
  </tr>
</table>
</div><div class="memdoc">
<dl class="section return"><dt>Returns</dt><dd>a reference to the row matrix permutation <img class="formulaInl" alt="$ P_r $" src="form_159.png"/> such that <img class="formulaInl" alt="$P_r A P_c^T = L U$" src="form_160.png"/> </dd></dl>
<dl class="section see"><dt>See Also</dt><dd><a class="el" href="classEigen_1_1SparseLU.html#a1229744f4e8554ca6e96fe32ac359924">colsPermutation()</a> </dd></dl>

</div>
</div>
<a class="anchor" id="a20467be6f23d8b39574f3e097583d767"></a>
<div class="memitem">
<div class="memproto">
<table class="mlabels">
  <tr>
  <td class="mlabels-left">
      <table class="memname">
        <tr>
          <td class="memname">void setPivotThreshold </td>
          <td>(</td>
          <td class="paramtype">const RealScalar &amp;&#160;</td>
          <td class="paramname"><em>thresh</em></td><td>)</td>
          <td></td>
        </tr>
      </table>
  </td>
  <td class="mlabels-right">
<span class="mlabels"><span class="mlabel">inline</span></span>  </td>
  </tr>
</table>
</div><div class="memdoc">
<p>Set the threshold used for a diagonal entry to be an acceptable pivot. </p>

</div>
</div>
<a class="anchor" id="a7150c533f7f2f6bdb1e9c82727f426f5"></a>
<div class="memitem">
<div class="memproto">
<table class="mlabels">
  <tr>
  <td class="mlabels-left">
      <table class="memname">
        <tr>
          <td class="memname">Scalar signDeterminant </td>
          <td>(</td>
          <td class="paramname"></td><td>)</td>
          <td></td>
        </tr>
      </table>
  </td>
  <td class="mlabels-right">
<span class="mlabels"><span class="mlabel">inline</span></span>  </td>
  </tr>
</table>
</div><div class="memdoc">
<dl class="section return"><dt>Returns</dt><dd>A number representing the sign of the determinant</dd></dl>
<dl class="section see"><dt>See Also</dt><dd><a class="el" href="classEigen_1_1SparseLU.html#a4c6b35f80ff9f6f0a3f4e74dfc121349">absDeterminant()</a>, <a class="el" href="classEigen_1_1SparseLU.html#adb61fce0a1343a45b65b6a50b18408dd">logAbsDeterminant()</a> </dd></dl>

</div>
</div>
<a class="anchor" id="a8166735beae3bd23e7c3a1be9f0f5ef9"></a>
<div class="memitem">
<div class="memproto">
<table class="mlabels">
  <tr>
  <td class="mlabels-left">
      <table class="memname">
        <tr>
          <td class="memname">const internal::solve_retval&lt;<a class="el" href="classEigen_1_1SparseLU.html">SparseLU</a>, Rhs&gt; solve </td>
          <td>(</td>
          <td class="paramtype">const <a class="el" href="classEigen_1_1MatrixBase.html">MatrixBase</a>&lt; Rhs &gt; &amp;&#160;</td>
          <td class="paramname"><em>B</em></td><td>)</td>
          <td> const</td>
        </tr>
      </table>
  </td>
  <td class="mlabels-right">
<span class="mlabels"><span class="mlabel">inline</span></span>  </td>
  </tr>
</table>
</div><div class="memdoc">
<dl class="section return"><dt>Returns</dt><dd>the solution X of <img class="formulaInl" alt="$ A X = B $" src="form_162.png"/> using the current decomposition of A.</dd></dl>
<dl class="section warning"><dt>Warning</dt><dd>the destination matrix X in X = this-&gt;solve(B) must be colmun-major.</dd></dl>
<dl class="section see"><dt>See Also</dt><dd><a class="el" href="classEigen_1_1SparseLU.html#a05a0eac1c146e9e62df9ec70f5a6f69d">compute()</a> </dd></dl>

</div>
</div>
<a class="anchor" id="a39062d90ea647b4e7668452647b04b17"></a>
<div class="memitem">
<div class="memproto">
<table class="mlabels">
  <tr>
  <td class="mlabels-left">
      <table class="memname">
        <tr>
          <td class="memname">const internal::sparse_solve_retval&lt;<a class="el" href="classEigen_1_1SparseLU.html">SparseLU</a>, Rhs&gt; solve </td>
          <td>(</td>
          <td class="paramtype">const <a class="el" href="classEigen_1_1SparseMatrixBase.html">SparseMatrixBase</a>&lt; Rhs &gt; &amp;&#160;</td>
          <td class="paramname"><em>B</em></td><td>)</td>
          <td> const</td>
        </tr>
      </table>
  </td>
  <td class="mlabels-right">
<span class="mlabels"><span class="mlabel">inline</span></span>  </td>
  </tr>
</table>
</div><div class="memdoc">
<dl class="section return"><dt>Returns</dt><dd>the solution X of <img class="formulaInl" alt="$ A X = B $" src="form_162.png"/> using the current decomposition of A.</dd></dl>
<dl class="section see"><dt>See Also</dt><dd><a class="el" href="classEigen_1_1SparseLU.html#a05a0eac1c146e9e62df9ec70f5a6f69d">compute()</a> </dd></dl>

<p>References <a class="el" href="structEigen_1_1EigenBase.html#aa84222add803ad7c9db07dd4dd91d5d9">EigenBase&lt; Derived &gt;::derived()</a>, and <a class="el" href="classEigen_1_1SparseMatrixBase.html#a5552abd83dbd03c85cea6d61fd8875a5">SparseMatrixBase&lt; Derived &gt;::rows()</a>.</p>

</div>
</div>
<hr/>The documentation for this class was generated from the following file:<ul>
<li><a class="el" href="SparseLU_8h_source.html">SparseLU.h</a></li>
</ul>
</div><!-- contents -->
</div><!-- doc-content -->
<!-- start footer part -->
<div id="nav-path" class="navpath"><!-- id is needed for treeview function! -->
  <ul>
    <li class="navelem"><a class="el" href="namespaceEigen.html">Eigen</a></li><li class="navelem"><a class="el" href="classEigen_1_1SparseLU.html">SparseLU</a></li>
    <li class="footer">Generated on Mon Oct 28 2013 11:04:30 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>