<!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: LDLT< _MatrixType, _UpLo > 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_1LDLT.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_1LDLT-members.html">List of all members</a> | <a href="#pub-methods">Public Member Functions</a> </div> <div class="headertitle"> <div class="title">LDLT< _MatrixType, _UpLo > Class Template Reference<div class="ingroups"><a class="el" href="group__Cholesky__Module.html">Cholesky 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 _MatrixType, int _UpLo><br/> class Eigen::LDLT< _MatrixType, _UpLo ></h3> <p>Robust Cholesky decomposition of a matrix with pivoting. </p> <dl class="params"><dt>Parameters</dt><dd> <table class="params"> <tr><td class="paramname">MatrixType</td><td>the type of the matrix of which to compute the LDL^T Cholesky decomposition </td></tr> <tr><td class="paramname">UpLo</td><td>the triangular part that will be used for the decompositon: Lower (default) or Upper. The other triangular part won't be read.</td></tr> </table> </dd> </dl> <p>Perform a robust Cholesky decomposition of a positive semidefinite or negative semidefinite matrix <img class="formulaInl" alt="$ A $" src="form_1.png"/> such that <img class="formulaInl" alt="$ A = P^TLDL^*P $" src="form_2.png"/>, where P is a permutation matrix, L is lower triangular with a unit diagonal and D is a diagonal matrix.</p> <p>The decomposition uses pivoting to ensure stability, so that L will have zeros in the bottom right rank(A) - n submatrix. Avoiding the square root on D also stabilizes the computation.</p> <p>Remember that Cholesky decompositions are not rank-revealing. Also, do not use a Cholesky decomposition to determine whether a system of equations has a solution.</p> <dl class="section see"><dt>See Also</dt><dd><a class="el" href="classEigen_1_1MatrixBase.html#a197a04cda6b4606ec2416fd3f950371f">MatrixBase::ldlt()</a>, class <a class="el" href="classEigen_1_1LLT.html" title="Standard Cholesky decomposition (LL^T) of a matrix and associated features. ">LLT</a> </dd></dl> </div><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:a46fcf31fdacf5205d1b6e6d64161a4b9"><td class="memItemLeft" align="right" valign="top"><a class="el" href="classEigen_1_1LDLT.html">LDLT</a> & </td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1LDLT.html#a46fcf31fdacf5205d1b6e6d64161a4b9">compute</a> (const MatrixType &matrix)</td></tr> <tr class="separator:a46fcf31fdacf5205d1b6e6d64161a4b9"><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_1LDLT.html#a0c06d5c2034ebb329c54235369643ad2">info</a> () const </td></tr> <tr class="memdesc:a0c06d5c2034ebb329c54235369643ad2"><td class="mdescLeft"> </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"> </td></tr> <tr class="memitem:a13ed609c8609698a470bb8dd0e43c09d"><td class="memItemLeft" align="right" valign="top">bool </td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1LDLT.html#a13ed609c8609698a470bb8dd0e43c09d">isNegative</a> (void) const </td></tr> <tr class="separator:a13ed609c8609698a470bb8dd0e43c09d"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:a02588e810c938a215dfd59d933488ba3"><td class="memItemLeft" align="right" valign="top">bool </td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1LDLT.html#a02588e810c938a215dfd59d933488ba3">isPositive</a> () const </td></tr> <tr class="separator:a02588e810c938a215dfd59d933488ba3"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:aec4c764dd032c14d861798976367e74d"><td class="memItemLeft" align="right" valign="top"> </td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1LDLT.html#aec4c764dd032c14d861798976367e74d">LDLT</a> ()</td></tr> <tr class="memdesc:aec4c764dd032c14d861798976367e74d"><td class="mdescLeft"> </td><td class="mdescRight">Default Constructor. <a href="#aec4c764dd032c14d861798976367e74d">More...</a><br/></td></tr> <tr class="separator:aec4c764dd032c14d861798976367e74d"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:ad329b592888d8f5db088dfa01504f6ad"><td class="memItemLeft" align="right" valign="top"> </td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1LDLT.html#ad329b592888d8f5db088dfa01504f6ad">LDLT</a> (Index size)</td></tr> <tr class="memdesc:ad329b592888d8f5db088dfa01504f6ad"><td class="mdescLeft"> </td><td class="mdescRight">Default Constructor with memory preallocation. <a href="#ad329b592888d8f5db088dfa01504f6ad">More...</a><br/></td></tr> <tr class="separator:ad329b592888d8f5db088dfa01504f6ad"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:a4f9c1dcdc8070fb63772e3d940e41fb3"><td class="memItemLeft" align="right" valign="top"> </td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1LDLT.html#a4f9c1dcdc8070fb63772e3d940e41fb3">LDLT</a> (const MatrixType &matrix)</td></tr> <tr class="memdesc:a4f9c1dcdc8070fb63772e3d940e41fb3"><td class="mdescLeft"> </td><td class="mdescRight">Constructor with decomposition. <a href="#a4f9c1dcdc8070fb63772e3d940e41fb3">More...</a><br/></td></tr> <tr class="separator:a4f9c1dcdc8070fb63772e3d940e41fb3"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:a625f69b684b4434db1cf3cc434e86fe6"><td class="memItemLeft" align="right" valign="top">Traits::MatrixL </td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1LDLT.html#a625f69b684b4434db1cf3cc434e86fe6">matrixL</a> () const </td></tr> <tr class="separator:a625f69b684b4434db1cf3cc434e86fe6"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:abfb1a0b0f7ea94b420697055019a5b6a"><td class="memItemLeft" align="right" valign="top">const MatrixType & </td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1LDLT.html#abfb1a0b0f7ea94b420697055019a5b6a">matrixLDLT</a> () const </td></tr> <tr class="separator:abfb1a0b0f7ea94b420697055019a5b6a"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:adc32fdba9f5b478afb2a96d53c6eacbb"><td class="memItemLeft" align="right" valign="top">Traits::MatrixU </td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1LDLT.html#adc32fdba9f5b478afb2a96d53c6eacbb">matrixU</a> () const </td></tr> <tr class="separator:adc32fdba9f5b478afb2a96d53c6eacbb"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:a910d592df729f6d1cee0ae2abf5f9612"><td class="memTemplParams" colspan="2">template<typename Derived > </td></tr> <tr class="memitem:a910d592df729f6d1cee0ae2abf5f9612"><td class="memTemplItemLeft" align="right" valign="top"><a class="el" href="classEigen_1_1LDLT.html">LDLT</a>< MatrixType, _UpLo > & </td><td class="memTemplItemRight" valign="bottom"><a class="el" href="classEigen_1_1LDLT.html#a910d592df729f6d1cee0ae2abf5f9612">rankUpdate</a> (const <a class="el" href="classEigen_1_1MatrixBase.html">MatrixBase</a>< Derived > &w, const typename <a class="el" href="structEigen_1_1NumTraits.html">NumTraits</a>< typename MatrixType::Scalar >::Real &sigma)</td></tr> <tr class="separator:a910d592df729f6d1cee0ae2abf5f9612"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:ac36925ac693435a090efee1cb5d6d16a"><td class="memItemLeft" align="right" valign="top">MatrixType </td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1LDLT.html#ac36925ac693435a090efee1cb5d6d16a">reconstructedMatrix</a> () const </td></tr> <tr class="separator:ac36925ac693435a090efee1cb5d6d16a"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:a47affd1a10b589811fc4828c1a2e0c6d"><td class="memItemLeft" align="right" valign="top">void </td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1LDLT.html#a47affd1a10b589811fc4828c1a2e0c6d">setZero</a> ()</td></tr> <tr class="separator:a47affd1a10b589811fc4828c1a2e0c6d"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:a71f78ede18adb53bdfe988161653e812"><td class="memTemplParams" colspan="2">template<typename Rhs > </td></tr> <tr class="memitem:a71f78ede18adb53bdfe988161653e812"><td class="memTemplItemLeft" align="right" valign="top">const internal::solve_retval<br class="typebreak"/> < <a class="el" href="classEigen_1_1LDLT.html">LDLT</a>, Rhs > </td><td class="memTemplItemRight" valign="bottom"><a class="el" href="classEigen_1_1LDLT.html#a71f78ede18adb53bdfe988161653e812">solve</a> (const <a class="el" href="classEigen_1_1MatrixBase.html">MatrixBase</a>< Rhs > &b) const </td></tr> <tr class="separator:a71f78ede18adb53bdfe988161653e812"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:a9a07717c12ca93c79d792ae77f9767ae"><td class="memItemLeft" align="right" valign="top">const <a class="el" href="classEigen_1_1Transpositions.html">TranspositionType</a> & </td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1LDLT.html#a9a07717c12ca93c79d792ae77f9767ae">transpositionsP</a> () const </td></tr> <tr class="separator:a9a07717c12ca93c79d792ae77f9767ae"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:adc5b1ae5b9cbc8a64912b8818cef9b9d"><td class="memItemLeft" align="right" valign="top"><a class="el" href="classEigen_1_1Diagonal.html">Diagonal</a>< const MatrixType > </td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1LDLT.html#adc5b1ae5b9cbc8a64912b8818cef9b9d">vectorD</a> () const </td></tr> <tr class="separator:adc5b1ae5b9cbc8a64912b8818cef9b9d"><td class="memSeparator" colspan="2"> </td></tr> </table> <h2 class="groupheader">Constructor & Destructor Documentation</h2> <a class="anchor" id="aec4c764dd032c14d861798976367e74d"></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_1LDLT.html">LDLT</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> <p>The default constructor is useful in cases in which the user intends to perform decompositions via <a class="el" href="classEigen_1_1LDLT.html#a46fcf31fdacf5205d1b6e6d64161a4b9">LDLT::compute(const MatrixType&)</a>. </p> </div> </div> <a class="anchor" id="ad329b592888d8f5db088dfa01504f6ad"></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_1LDLT.html">LDLT</a> </td> <td>(</td> <td class="paramtype">Index </td> <td class="paramname"><em>size</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>Default Constructor with memory preallocation. </p> <p>Like the default constructor but with preallocation of the internal data according to the specified problem <em>size</em>. </p> <dl class="section see"><dt>See Also</dt><dd><a class="el" href="classEigen_1_1LDLT.html#aec4c764dd032c14d861798976367e74d" title="Default Constructor. ">LDLT()</a> </dd></dl> </div> </div> <a class="anchor" id="a4f9c1dcdc8070fb63772e3d940e41fb3"></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_1LDLT.html">LDLT</a> </td> <td>(</td> <td class="paramtype">const MatrixType & </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>Constructor with decomposition. </p> <p>This calculates the decomposition for the input <em>matrix</em>. </p> <dl class="section see"><dt>See Also</dt><dd><a class="el" href="classEigen_1_1LDLT.html#ad329b592888d8f5db088dfa01504f6ad" title="Default Constructor with memory preallocation. ">LDLT(Index size)</a> </dd></dl> <p>References <a class="el" href="classEigen_1_1LDLT.html#a46fcf31fdacf5205d1b6e6d64161a4b9">LDLT< _MatrixType, _UpLo >::compute()</a>.</p> </div> </div> <h2 class="groupheader">Member Function Documentation</h2> <a class="anchor" id="a46fcf31fdacf5205d1b6e6d64161a4b9"></a> <div class="memitem"> <div class="memproto"> <table class="memname"> <tr> <td class="memname"><a class="el" href="classEigen_1_1LDLT.html">LDLT</a>< MatrixType, _UpLo > & compute </td> <td>(</td> <td class="paramtype">const MatrixType & </td> <td class="paramname"><em>a</em></td><td>)</td> <td></td> </tr> </table> </div><div class="memdoc"> <p>Compute / recompute the <a class="el" href="classEigen_1_1LDLT.html" title="Robust Cholesky decomposition of a matrix with pivoting. ">LDLT</a> decomposition A = L D L^* = U^* D U of <em>matrix</em> </p> <p>Referenced by <a class="el" href="classEigen_1_1LDLT.html#a4f9c1dcdc8070fb63772e3d940e41fb3">LDLT< _MatrixType, _UpLo >::LDLT()</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 matrix.appears to be negative. </dd></dl> <p>References <a class="el" href="group__enums.html#gga51bc1ac16f26ebe51eae1abb77bd037bafdfbdf3247bd36a1f17270d5cec74c9c">Eigen::Success</a>.</p> </div> </div> <a class="anchor" id="a13ed609c8609698a470bb8dd0e43c09d"></a> <div class="memitem"> <div class="memproto"> <table class="mlabels"> <tr> <td class="mlabels-left"> <table class="memname"> <tr> <td class="memname">bool isNegative </td> <td>(</td> <td class="paramtype">void </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>true if the matrix is negative (semidefinite) </dd></dl> </div> </div> <a class="anchor" id="a02588e810c938a215dfd59d933488ba3"></a> <div class="memitem"> <div class="memproto"> <table class="mlabels"> <tr> <td class="mlabels-left"> <table class="memname"> <tr> <td class="memname">bool isPositive </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>true if the matrix is positive (semidefinite) </dd></dl> </div> </div> <a class="anchor" id="a625f69b684b4434db1cf3cc434e86fe6"></a> <div class="memitem"> <div class="memproto"> <table class="mlabels"> <tr> <td class="mlabels-left"> <table class="memname"> <tr> <td class="memname">Traits::MatrixL 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>a view of the lower triangular matrix L </dd></dl> </div> </div> <a class="anchor" id="abfb1a0b0f7ea94b420697055019a5b6a"></a> <div class="memitem"> <div class="memproto"> <table class="mlabels"> <tr> <td class="mlabels-left"> <table class="memname"> <tr> <td class="memname">const MatrixType& matrixLDLT </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 internal <a class="el" href="classEigen_1_1LDLT.html" title="Robust Cholesky decomposition of a matrix with pivoting. ">LDLT</a> decomposition matrix</dd></dl> <p>TODO: document the storage layout </p> </div> </div> <a class="anchor" id="adc32fdba9f5b478afb2a96d53c6eacbb"></a> <div class="memitem"> <div class="memproto"> <table class="mlabels"> <tr> <td class="mlabels-left"> <table class="memname"> <tr> <td class="memname">Traits::MatrixU 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>a view of the upper triangular matrix U </dd></dl> </div> </div> <a class="anchor" id="a910d592df729f6d1cee0ae2abf5f9612"></a> <div class="memitem"> <div class="memproto"> <table class="memname"> <tr> <td class="memname"><a class="el" href="classEigen_1_1LDLT.html">LDLT</a><MatrixType,_UpLo>& rankUpdate </td> <td>(</td> <td class="paramtype">const <a class="el" href="classEigen_1_1MatrixBase.html">MatrixBase</a>< Derived > & </td> <td class="paramname"><em>w</em>, </td> </tr> <tr> <td class="paramkey"></td> <td></td> <td class="paramtype">const typename <a class="el" href="structEigen_1_1NumTraits.html">NumTraits</a>< typename MatrixType::Scalar >::Real & </td> <td class="paramname"><em>sigma</em> </td> </tr> <tr> <td></td> <td>)</td> <td></td><td></td> </tr> </table> </div><div class="memdoc"> <p>Update the <a class="el" href="classEigen_1_1LDLT.html" title="Robust Cholesky decomposition of a matrix with pivoting. ">LDLT</a> decomposition: given A = L D L^T, efficiently compute the decomposition of A + sigma w w^T. </p> <dl class="params"><dt>Parameters</dt><dd> <table class="params"> <tr><td class="paramname">w</td><td>a vector to be incorporated into the decomposition. </td></tr> <tr><td class="paramname">sigma</td><td>a scalar, +1 for updates and -1 for "downdates," which correspond to removing previously-added column vectors. Optional; default value is +1. </td></tr> </table> </dd> </dl> <dl class="section see"><dt>See Also</dt><dd><a class="el" href="classEigen_1_1LDLT.html#a47affd1a10b589811fc4828c1a2e0c6d">setZero()</a> </dd></dl> </div> </div> <a class="anchor" id="ac36925ac693435a090efee1cb5d6d16a"></a> <div class="memitem"> <div class="memproto"> <table class="memname"> <tr> <td class="memname">MatrixType reconstructedMatrix </td> <td>(</td> <td class="paramname"></td><td>)</td> <td> const</td> </tr> </table> </div><div class="memdoc"> <dl class="section return"><dt>Returns</dt><dd>the matrix represented by the decomposition, i.e., it returns the product: P^T L D L^* P. This function is provided for debug purpose. </dd></dl> </div> </div> <a class="anchor" id="a47affd1a10b589811fc4828c1a2e0c6d"></a> <div class="memitem"> <div class="memproto"> <table class="mlabels"> <tr> <td class="mlabels-left"> <table class="memname"> <tr> <td class="memname">void setZero </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>Clear any existing decomposition </p> <dl class="section see"><dt>See Also</dt><dd>rankUpdate(w,sigma) </dd></dl> </div> </div> <a class="anchor" id="a71f78ede18adb53bdfe988161653e812"></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<<a class="el" href="classEigen_1_1LDLT.html">LDLT</a>, 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>a solution x of <img class="formulaInl" alt="$ A x = b $" src="form_3.png"/> using the current decomposition of A.</dd></dl> <p>This function also supports in-place solves using the syntax <code>x = decompositionObject.solve(x)</code> .</p> <p>This method just tries to find as good a solution as possible. If you want to check whether a solution exists or if it is accurate, just call this function to get a result and then compute the error of this result, or use <a class="el" href="classEigen_1_1DenseBase.html#a158c2184951e6e415c2e9b98db8e8966">MatrixBase::isApprox()</a> directly, for instance like this:</p> <div class="fragment"><div class="line"><span class="keywordtype">bool</span> a_solution_exists = (A*result).isApprox(b, precision); </div> </div><!-- fragment --><p> This method avoids dividing by zero, so that the non-existence of a solution doesn't by itself mean that you'll get <code>inf</code> or <code>nan</code> values.</p> <p>More precisely, this method solves <img class="formulaInl" alt="$ A x = b $" src="form_3.png"/> using the decomposition <img class="formulaInl" alt="$ A = P^T L D L^* P $" src="form_4.png"/> by solving the systems <img class="formulaInl" alt="$ P^T y_1 = b $" src="form_5.png"/>, <img class="formulaInl" alt="$ L y_2 = y_1 $" src="form_6.png"/>, <img class="formulaInl" alt="$ D y_3 = y_2 $" src="form_7.png"/>, <img class="formulaInl" alt="$ L^* y_4 = y_3 $" src="form_8.png"/> and <img class="formulaInl" alt="$ P x = y_4 $" src="form_9.png"/> in succession. If the matrix <img class="formulaInl" alt="$ A $" src="form_1.png"/> is singular, then <img class="formulaInl" alt="$ D $" src="form_10.png"/> will also be singular (all the other matrices are invertible). In that case, the least-square solution of <img class="formulaInl" alt="$ D y_3 = y_2 $" src="form_7.png"/> is computed. This does not mean that this function computes the least-square solution of <img class="formulaInl" alt="$ A x = b $" src="form_3.png"/> is <img class="formulaInl" alt="$ A $" src="form_1.png"/> is singular.</p> <dl class="section see"><dt>See Also</dt><dd><a class="el" href="classEigen_1_1MatrixBase.html#a197a04cda6b4606ec2416fd3f950371f">MatrixBase::ldlt()</a> </dd></dl> </div> </div> <a class="anchor" id="a9a07717c12ca93c79d792ae77f9767ae"></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_1Transpositions.html">TranspositionType</a>& transpositionsP </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 permutation matrix P as a transposition sequence. </dd></dl> </div> </div> <a class="anchor" id="adc5b1ae5b9cbc8a64912b8818cef9b9d"></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_1Diagonal.html">Diagonal</a><const MatrixType> vectorD </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 coefficients of the diagonal matrix D </dd></dl> </div> </div> <hr/>The documentation for this class was generated from the following file:<ul> <li><a class="el" href="LDLT_8h_source.html">LDLT.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_1LDLT.html">LDLT</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>