<!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: IterativeSolverBase< Derived > 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>  <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_1IterativeSolverBase.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="summary"> <a href="classEigen_1_1IterativeSolverBase-members.html">List of all members</a> | <a href="#pub-methods">Public Member Functions</a> </div> <div class="headertitle"> <div class="title">IterativeSolverBase< Derived > Class Template Reference<div class="ingroups"><a class="el" href="group__IterativeLinearSolvers__Module.html">IterativeLinearSolvers 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<typename Derived><br/> class Eigen::IterativeSolverBase< Derived ></h3> <p>Base class for linear iterative solvers. </p> <dl class="section see"><dt>See Also</dt><dd>class <a class="el" href="classEigen_1_1SimplicialCholesky.html">SimplicialCholesky</a>, <a class="el" href="classEigen_1_1DiagonalPreconditioner.html" title="A preconditioner based on the digonal entries. ">DiagonalPreconditioner</a>, <a class="el" href="classEigen_1_1IdentityPreconditioner.html" title="A naive preconditioner which approximates any matrix as the identity matrix. ">IdentityPreconditioner</a> </dd></dl> </div> <p>Inherits noncopyable.</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:ae044142449c071af13afc4fbbe633682"><td class="memItemLeft" align="right" valign="top">Derived & </td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1IterativeSolverBase.html#ae044142449c071af13afc4fbbe633682">analyzePattern</a> (const MatrixType &A)</td></tr> <tr class="separator:ae044142449c071af13afc4fbbe633682"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:a693a8142a905e7130bacd77586a6f6a0"><td class="memItemLeft" align="right" valign="top">Derived & </td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1IterativeSolverBase.html#a693a8142a905e7130bacd77586a6f6a0">compute</a> (const MatrixType &A)</td></tr> <tr class="separator:a693a8142a905e7130bacd77586a6f6a0"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:a736349b3462f6ebf37c34b450db2ed18"><td class="memItemLeft" align="right" valign="top">RealScalar </td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1IterativeSolverBase.html#a736349b3462f6ebf37c34b450db2ed18">error</a> () const </td></tr> <tr class="separator:a736349b3462f6ebf37c34b450db2ed18"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:ab837a9d99eea1f79dd4721727ef98d37"><td class="memItemLeft" align="right" valign="top">Derived & </td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1IterativeSolverBase.html#ab837a9d99eea1f79dd4721727ef98d37">factorize</a> (const MatrixType &A)</td></tr> <tr class="separator:ab837a9d99eea1f79dd4721727ef98d37"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:a0c06d5c2034ebb329c54235369643ad2"><td class="memItemLeft" align="right" valign="top"><a class="el" href="group__enums.html#ga51bc1ac16f26ebe51eae1abb77bd037b">ComputationInfo</a> </td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1IterativeSolverBase.html#a0c06d5c2034ebb329c54235369643ad2">info</a> () const </td></tr> <tr class="separator:a0c06d5c2034ebb329c54235369643ad2"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:a07c5f1fe9fa1556edc79452b52839829"><td class="memItemLeft" align="right" valign="top">int </td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1IterativeSolverBase.html#a07c5f1fe9fa1556edc79452b52839829">iterations</a> () const </td></tr> <tr class="separator:a07c5f1fe9fa1556edc79452b52839829"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:a5140d3e2be1561de7f003504ef1be45d"><td class="memItemLeft" align="right" valign="top"> </td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1IterativeSolverBase.html#a5140d3e2be1561de7f003504ef1be45d">IterativeSolverBase</a> ()</td></tr> <tr class="separator:a5140d3e2be1561de7f003504ef1be45d"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:a81b1c0c3c9d19709b15e03a47ea90ef0"><td class="memItemLeft" align="right" valign="top"> </td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1IterativeSolverBase.html#a81b1c0c3c9d19709b15e03a47ea90ef0">IterativeSolverBase</a> (const MatrixType &A)</td></tr> <tr class="separator:a81b1c0c3c9d19709b15e03a47ea90ef0"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:a9c020a2e26d3081eeff611b802962530"><td class="memItemLeft" align="right" valign="top">int </td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1IterativeSolverBase.html#a9c020a2e26d3081eeff611b802962530">maxIterations</a> () const </td></tr> <tr class="separator:a9c020a2e26d3081eeff611b802962530"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:ae4bec7fa55f17f7bb64b6b2fc4039b58"><td class="memItemLeft" align="right" valign="top">Preconditioner & </td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1IterativeSolverBase.html#ae4bec7fa55f17f7bb64b6b2fc4039b58">preconditioner</a> ()</td></tr> <tr class="separator:ae4bec7fa55f17f7bb64b6b2fc4039b58"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:a11aed1f82ecbec6d7436915c6d804d7a"><td class="memItemLeft" align="right" valign="top">const Preconditioner & </td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1IterativeSolverBase.html#a11aed1f82ecbec6d7436915c6d804d7a">preconditioner</a> () const </td></tr> <tr class="separator:a11aed1f82ecbec6d7436915c6d804d7a"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:ad9aa2f8ecfda11aba72cfd3c5eefae20"><td class="memItemLeft" align="right" valign="top">Derived & </td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1IterativeSolverBase.html#ad9aa2f8ecfda11aba72cfd3c5eefae20">setMaxIterations</a> (int maxIters)</td></tr> <tr class="separator:ad9aa2f8ecfda11aba72cfd3c5eefae20"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:a39b43ecdc61211a15d1e986e28a03513"><td class="memItemLeft" align="right" valign="top">Derived & </td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1IterativeSolverBase.html#a39b43ecdc61211a15d1e986e28a03513">setTolerance</a> (const RealScalar &<a class="el" href="classEigen_1_1IterativeSolverBase.html#aab135f12fadf495a511f20c71557e1aa">tolerance</a>)</td></tr> <tr class="separator:a39b43ecdc61211a15d1e986e28a03513"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:acae1088225b7c746b16848d02d5e1c2a"><td class="memTemplParams" colspan="2">template<typename Rhs > </td></tr> <tr class="memitem:acae1088225b7c746b16848d02d5e1c2a"><td class="memTemplItemLeft" align="right" valign="top">const internal::solve_retval<br class="typebreak"/> < Derived, Rhs > </td><td class="memTemplItemRight" valign="bottom"><a class="el" href="classEigen_1_1IterativeSolverBase.html#acae1088225b7c746b16848d02d5e1c2a">solve</a> (const <a class="el" href="classEigen_1_1MatrixBase.html">MatrixBase</a>< Rhs > &b) const </td></tr> <tr class="separator:acae1088225b7c746b16848d02d5e1c2a"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:a7892f6fbc52bfa658e1324a63d376d1c"><td class="memTemplParams" colspan="2">template<typename Rhs > </td></tr> <tr class="memitem:a7892f6fbc52bfa658e1324a63d376d1c"><td class="memTemplItemLeft" align="right" valign="top">const <br class="typebreak"/> internal::sparse_solve_retval<br class="typebreak"/> < <a class="el" href="classEigen_1_1IterativeSolverBase.html">IterativeSolverBase</a>, Rhs > </td><td class="memTemplItemRight" valign="bottom"><a class="el" href="classEigen_1_1IterativeSolverBase.html#a7892f6fbc52bfa658e1324a63d376d1c">solve</a> (const <a class="el" href="classEigen_1_1SparseMatrixBase.html">SparseMatrixBase</a>< Rhs > &b) const </td></tr> <tr class="separator:a7892f6fbc52bfa658e1324a63d376d1c"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:aab135f12fadf495a511f20c71557e1aa"><td class="memItemLeft" align="right" valign="top">RealScalar </td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1IterativeSolverBase.html#aab135f12fadf495a511f20c71557e1aa">tolerance</a> () const </td></tr> <tr class="separator:aab135f12fadf495a511f20c71557e1aa"><td class="memSeparator" colspan="2"> </td></tr> </table> <h2 class="groupheader">Constructor & Destructor Documentation</h2> <a class="anchor" id="a5140d3e2be1561de7f003504ef1be45d"></a> <div class="memitem"> <div class="memproto"> <table class="mlabels"> <tr> <td class="mlabels-left"> <table class="memname"> <tr> <td class="memname"><a class="el" href="classEigen_1_1IterativeSolverBase.html">IterativeSolverBase</a> </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"> <p>Default constructor. </p> </div> </div> <a class="anchor" id="a81b1c0c3c9d19709b15e03a47ea90ef0"></a> <div class="memitem"> <div class="memproto"> <table class="mlabels"> <tr> <td class="mlabels-left"> <table class="memname"> <tr> <td class="memname"><a class="el" href="classEigen_1_1IterativeSolverBase.html">IterativeSolverBase</a> </td> <td>(</td> <td class="paramtype">const MatrixType & </td> <td class="paramname"><em>A</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>Initialize the solver with matrix <em>A</em> for further <code>Ax=b</code> solving.</p> <p>This constructor is a shortcut for the default constructor followed by a call to <a class="el" href="classEigen_1_1IterativeSolverBase.html#a693a8142a905e7130bacd77586a6f6a0">compute()</a>.</p> <dl class="section warning"><dt>Warning</dt><dd>this class stores a reference to the matrix A as well as some precomputed values that depend on it. Therefore, if <em>A</em> is changed this class becomes invalid. Call <a class="el" href="classEigen_1_1IterativeSolverBase.html#a693a8142a905e7130bacd77586a6f6a0">compute()</a> to update it with the new matrix A, or modify a copy of A. </dd></dl> </div> </div> <h2 class="groupheader">Member Function Documentation</h2> <a class="anchor" id="ae044142449c071af13afc4fbbe633682"></a> <div class="memitem"> <div class="memproto"> <table class="mlabels"> <tr> <td class="mlabels-left"> <table class="memname"> <tr> <td class="memname">Derived& analyzePattern </td> <td>(</td> <td class="paramtype">const MatrixType & </td> <td class="paramname"><em>A</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>Initializes the iterative solver for the sparcity pattern of the matrix <em>A</em> for further solving <code>Ax=b</code> problems.</p> <p>Currently, this function mostly call analyzePattern on the preconditioner. In the future we might, for instance, implement column reodering for faster matrix vector products. </p> </div> </div> <a class="anchor" id="a693a8142a905e7130bacd77586a6f6a0"></a> <div class="memitem"> <div class="memproto"> <table class="mlabels"> <tr> <td class="mlabels-left"> <table class="memname"> <tr> <td class="memname">Derived& compute </td> <td>(</td> <td class="paramtype">const MatrixType & </td> <td class="paramname"><em>A</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>Initializes the iterative solver with the matrix <em>A</em> for further solving <code>Ax=b</code> problems.</p> <p>Currently, this function mostly initialized/compute the preconditioner. In the future we might, for instance, implement column reodering for faster matrix vector products.</p> <dl class="section warning"><dt>Warning</dt><dd>this class stores a reference to the matrix A as well as some precomputed values that depend on it. Therefore, if <em>A</em> is changed this class becomes invalid. Call <a class="el" href="classEigen_1_1IterativeSolverBase.html#a693a8142a905e7130bacd77586a6f6a0">compute()</a> to update it with the new matrix A, or modify a copy of A. </dd></dl> <p>Referenced by <a class="el" href="classEigen_1_1IterativeSolverBase.html#a81b1c0c3c9d19709b15e03a47ea90ef0">IterativeSolverBase< ConjugateGradient< _MatrixType, _UpLo, _Preconditioner > >::IterativeSolverBase()</a>.</p> </div> </div> <a class="anchor" id="a736349b3462f6ebf37c34b450db2ed18"></a> <div class="memitem"> <div class="memproto"> <table class="mlabels"> <tr> <td class="mlabels-left"> <table class="memname"> <tr> <td class="memname">RealScalar error </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 tolerance error reached during the last solve </dd></dl> </div> </div> <a class="anchor" id="ab837a9d99eea1f79dd4721727ef98d37"></a> <div class="memitem"> <div class="memproto"> <table class="mlabels"> <tr> <td class="mlabels-left"> <table class="memname"> <tr> <td class="memname">Derived& factorize </td> <td>(</td> <td class="paramtype">const MatrixType & </td> <td class="paramname"><em>A</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>Initializes the iterative solver with the numerical values of the matrix <em>A</em> for further solving <code>Ax=b</code> problems.</p> <p>Currently, this function mostly call factorize on the preconditioner.</p> <dl class="section warning"><dt>Warning</dt><dd>this class stores a reference to the matrix A as well as some precomputed values that depend on it. Therefore, if <em>A</em> is changed this class becomes invalid. Call <a class="el" href="classEigen_1_1IterativeSolverBase.html#a693a8142a905e7130bacd77586a6f6a0">compute()</a> to update it with the new matrix A, or modify a copy of A. </dd></dl> </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"> <dl class="section return"><dt>Returns</dt><dd>Success if the iterations converged, and NoConvergence otherwise. </dd></dl> </div> </div> <a class="anchor" id="a07c5f1fe9fa1556edc79452b52839829"></a> <div class="memitem"> <div class="memproto"> <table class="mlabels"> <tr> <td class="mlabels-left"> <table class="memname"> <tr> <td class="memname">int iterations </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 number of iterations performed during the last solve </dd></dl> </div> </div> <a class="anchor" id="a9c020a2e26d3081eeff611b802962530"></a> <div class="memitem"> <div class="memproto"> <table class="mlabels"> <tr> <td class="mlabels-left"> <table class="memname"> <tr> <td class="memname">int maxIterations </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 max number of iterations </dd></dl> </div> </div> <a class="anchor" id="ae4bec7fa55f17f7bb64b6b2fc4039b58"></a> <div class="memitem"> <div class="memproto"> <table class="mlabels"> <tr> <td class="mlabels-left"> <table class="memname"> <tr> <td class="memname">Preconditioner& preconditioner </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 read-write reference to the preconditioner for custom configuration. </dd></dl> </div> </div> <a class="anchor" id="a11aed1f82ecbec6d7436915c6d804d7a"></a> <div class="memitem"> <div class="memproto"> <table class="mlabels"> <tr> <td class="mlabels-left"> <table class="memname"> <tr> <td class="memname">const Preconditioner& preconditioner </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 read-only reference to the preconditioner. </dd></dl> </div> </div> <a class="anchor" id="ad9aa2f8ecfda11aba72cfd3c5eefae20"></a> <div class="memitem"> <div class="memproto"> <table class="mlabels"> <tr> <td class="mlabels-left"> <table class="memname"> <tr> <td class="memname">Derived& setMaxIterations </td> <td>(</td> <td class="paramtype">int </td> <td class="paramname"><em>maxIters</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>Sets the max number of iterations </p> </div> </div> <a class="anchor" id="a39b43ecdc61211a15d1e986e28a03513"></a> <div class="memitem"> <div class="memproto"> <table class="mlabels"> <tr> <td class="mlabels-left"> <table class="memname"> <tr> <td class="memname">Derived& setTolerance </td> <td>(</td> <td class="paramtype">const RealScalar & </td> <td class="paramname"><em>tolerance</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>Sets the tolerance threshold used by the stopping criteria </p> </div> </div> <a class="anchor" id="acae1088225b7c746b16848d02d5e1c2a"></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<Derived, Rhs> solve </td> <td>(</td> <td class="paramtype">const <a class="el" href="classEigen_1_1MatrixBase.html">MatrixBase</a>< Rhs > & </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_3.png"/> using the current decomposition of A.</dd></dl> <dl class="section see"><dt>See Also</dt><dd><a class="el" href="classEigen_1_1IterativeSolverBase.html#a693a8142a905e7130bacd77586a6f6a0">compute()</a> </dd></dl> </div> </div> <a class="anchor" id="a7892f6fbc52bfa658e1324a63d376d1c"></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<<a class="el" href="classEigen_1_1IterativeSolverBase.html">IterativeSolverBase</a>, Rhs> solve </td> <td>(</td> <td class="paramtype">const <a class="el" href="classEigen_1_1SparseMatrixBase.html">SparseMatrixBase</a>< Rhs > & </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_3.png"/> using the current decomposition of A.</dd></dl> <dl class="section see"><dt>See Also</dt><dd><a class="el" href="classEigen_1_1IterativeSolverBase.html#a693a8142a905e7130bacd77586a6f6a0">compute()</a> </dd></dl> </div> </div> <a class="anchor" id="aab135f12fadf495a511f20c71557e1aa"></a> <div class="memitem"> <div class="memproto"> <table class="mlabels"> <tr> <td class="mlabels-left"> <table class="memname"> <tr> <td class="memname">RealScalar tolerance </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 tolerance threshold used by the stopping criteria </dd></dl> <p>Referenced by <a class="el" href="classEigen_1_1IterativeSolverBase.html#a39b43ecdc61211a15d1e986e28a03513">IterativeSolverBase< ConjugateGradient< _MatrixType, _UpLo, _Preconditioner > >::setTolerance()</a>.</p> </div> </div> <hr/>The documentation for this class was generated from the following file:<ul> <li><a class="el" href="IterativeSolverBase_8h_source.html">IterativeSolverBase.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_1IterativeSolverBase.html">IterativeSolverBase</a></li> <li class="footer">Generated on Mon Oct 28 2013 11:04:29 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>