<!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: ColPivHouseholderQR< MatrixType > 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_1ColPivHouseholderQR.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_1ColPivHouseholderQR-members.html">List of all members</a> | <a href="#pub-methods">Public Member Functions</a> </div> <div class="headertitle"> <div class="title">ColPivHouseholderQR< MatrixType > Class Template Reference<div class="ingroups"><a class="el" href="group__QR__Module.html">QR 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><br/> class Eigen::ColPivHouseholderQR< MatrixType ></h3> <p>Householder rank-revealing QR decomposition of a matrix with column-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 we are computing the QR decomposition</td></tr> </table> </dd> </dl> <p>This class performs a rank-revealing QR decomposition of a matrix <b>A</b> into matrices <b>P</b>, <b>Q</b> and <b>R</b> such that </p> <p class="formulaDsp"> <img class="formulaDsp" alt="\[ \mathbf{A} \, \mathbf{P} = \mathbf{Q} \, \mathbf{R} \]" src="form_157.png"/> </p> <p> by using Householder transformations. Here, <b>P</b> is a permutation matrix, <b>Q</b> a unitary matrix and <b>R</b> an upper triangular matrix.</p> <p>This decomposition performs column pivoting in order to be rank-revealing and improve numerical stability. It is slower than <a class="el" href="classEigen_1_1HouseholderQR.html" title="Householder QR decomposition of a matrix. ">HouseholderQR</a>, and faster than <a class="el" href="classEigen_1_1FullPivHouseholderQR.html" title="Householder rank-revealing QR decomposition of a matrix with full pivoting. ">FullPivHouseholderQR</a>.</p> <dl class="section see"><dt>See Also</dt><dd><a class="el" href="classEigen_1_1MatrixBase.html#a05afed751d3a7277951d1918468e0872">MatrixBase::colPivHouseholderQr()</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:a2bb8fa5593d2f447c16ba2e3cfb9c432"><td class="memItemLeft" align="right" valign="top">MatrixType::RealScalar </td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1ColPivHouseholderQR.html#a2bb8fa5593d2f447c16ba2e3cfb9c432">absDeterminant</a> () const </td></tr> <tr class="separator:a2bb8fa5593d2f447c16ba2e3cfb9c432"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:a72f7fb8c1359fa8d9743068e5699228e"><td class="memItemLeft" align="right" valign="top"> </td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1ColPivHouseholderQR.html#a72f7fb8c1359fa8d9743068e5699228e">ColPivHouseholderQR</a> ()</td></tr> <tr class="memdesc:a72f7fb8c1359fa8d9743068e5699228e"><td class="mdescLeft"> </td><td class="mdescRight">Default Constructor. <a href="#a72f7fb8c1359fa8d9743068e5699228e">More...</a><br/></td></tr> <tr class="separator:a72f7fb8c1359fa8d9743068e5699228e"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:a4a15f3e5e33478b9969d8c71aabd09ea"><td class="memItemLeft" align="right" valign="top"> </td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1ColPivHouseholderQR.html#a4a15f3e5e33478b9969d8c71aabd09ea">ColPivHouseholderQR</a> (Index rows, Index cols)</td></tr> <tr class="memdesc:a4a15f3e5e33478b9969d8c71aabd09ea"><td class="mdescLeft"> </td><td class="mdescRight">Default Constructor with memory preallocation. <a href="#a4a15f3e5e33478b9969d8c71aabd09ea">More...</a><br/></td></tr> <tr class="separator:a4a15f3e5e33478b9969d8c71aabd09ea"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:ac93cc6a4fd4247486c2c8eab79b4c0ad"><td class="memItemLeft" align="right" valign="top"> </td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1ColPivHouseholderQR.html#ac93cc6a4fd4247486c2c8eab79b4c0ad">ColPivHouseholderQR</a> (const MatrixType &matrix)</td></tr> <tr class="memdesc:ac93cc6a4fd4247486c2c8eab79b4c0ad"><td class="mdescLeft"> </td><td class="mdescRight">Constructs a QR factorization from a given matrix. <a href="#ac93cc6a4fd4247486c2c8eab79b4c0ad">More...</a><br/></td></tr> <tr class="separator:ac93cc6a4fd4247486c2c8eab79b4c0ad"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:a1229744f4e8554ca6e96fe32ac359924"><td class="memItemLeft" align="right" valign="top">const <a class="el" href="classEigen_1_1PermutationMatrix.html">PermutationType</a> & </td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1ColPivHouseholderQR.html#a1229744f4e8554ca6e96fe32ac359924">colsPermutation</a> () const </td></tr> <tr class="separator:a1229744f4e8554ca6e96fe32ac359924"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:a812205d31cae2005b7779c7a3a442f1b"><td class="memItemLeft" align="right" valign="top"><a class="el" href="classEigen_1_1ColPivHouseholderQR.html">ColPivHouseholderQR</a> & </td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1ColPivHouseholderQR.html#a812205d31cae2005b7779c7a3a442f1b">compute</a> (const MatrixType &matrix)</td></tr> <tr class="separator:a812205d31cae2005b7779c7a3a442f1b"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:a7c6323871c4f080fc6e2d3ad7fc607fc"><td class="memItemLeft" align="right" valign="top">Index </td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1ColPivHouseholderQR.html#a7c6323871c4f080fc6e2d3ad7fc607fc">dimensionOfKernel</a> () const </td></tr> <tr class="separator:a7c6323871c4f080fc6e2d3ad7fc607fc"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:a321e96844544eacde2272c3a487048e9"><td class="memItemLeft" align="right" valign="top">const HCoeffsType & </td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1ColPivHouseholderQR.html#a321e96844544eacde2272c3a487048e9">hCoeffs</a> () const </td></tr> <tr class="separator:a321e96844544eacde2272c3a487048e9"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:a97150efba7e2da32f48a7f823f047a89"><td class="memItemLeft" align="right" valign="top"><a class="el" href="classEigen_1_1HouseholderSequence.html">HouseholderSequenceType</a> </td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1ColPivHouseholderQR.html#a97150efba7e2da32f48a7f823f047a89">householderQ</a> (void) const </td></tr> <tr class="separator:a97150efba7e2da32f48a7f823f047a89"><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_1ColPivHouseholderQR.html#a0c06d5c2034ebb329c54235369643ad2">info</a> () const </td></tr> <tr class="memdesc:a0c06d5c2034ebb329c54235369643ad2"><td class="mdescLeft"> </td><td class="mdescRight">Reports whether the QR factorization was succesful. <a href="#a0c06d5c2034ebb329c54235369643ad2">More...</a><br/></td></tr> <tr class="separator:a0c06d5c2034ebb329c54235369643ad2"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:a3021b1285c868ef94b26f387766bdf92"><td class="memItemLeft" align="right" valign="top">const internal::solve_retval<br class="typebreak"/> < <a class="el" href="classEigen_1_1ColPivHouseholderQR.html">ColPivHouseholderQR</a>, <br class="typebreak"/> typename <br class="typebreak"/> MatrixType::IdentityReturnType > </td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1ColPivHouseholderQR.html#a3021b1285c868ef94b26f387766bdf92">inverse</a> () const </td></tr> <tr class="separator:a3021b1285c868ef94b26f387766bdf92"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:a1e119085e53eca65e9ba15451c102d40"><td class="memItemLeft" align="right" valign="top">bool </td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1ColPivHouseholderQR.html#a1e119085e53eca65e9ba15451c102d40">isInjective</a> () const </td></tr> <tr class="separator:a1e119085e53eca65e9ba15451c102d40"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:ab60c7d993c9eba31668fb8886d621094"><td class="memItemLeft" align="right" valign="top">bool </td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1ColPivHouseholderQR.html#ab60c7d993c9eba31668fb8886d621094">isInvertible</a> () const </td></tr> <tr class="separator:ab60c7d993c9eba31668fb8886d621094"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:a5968b9ca46303b3cc7250e7b120ab7e6"><td class="memItemLeft" align="right" valign="top">bool </td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1ColPivHouseholderQR.html#a5968b9ca46303b3cc7250e7b120ab7e6">isSurjective</a> () const </td></tr> <tr class="separator:a5968b9ca46303b3cc7250e7b120ab7e6"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:a82a3f15f0cda1f4f874d50951344b5a9"><td class="memItemLeft" align="right" valign="top">MatrixType::RealScalar </td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1ColPivHouseholderQR.html#a82a3f15f0cda1f4f874d50951344b5a9">logAbsDeterminant</a> () const </td></tr> <tr class="separator:a82a3f15f0cda1f4f874d50951344b5a9"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:a7b142db8924dd7fad99acbdd0ff4bdd1"><td class="memItemLeft" align="right" valign="top">const MatrixType & </td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1ColPivHouseholderQR.html#a7b142db8924dd7fad99acbdd0ff4bdd1">matrixQR</a> () const </td></tr> <tr class="separator:a7b142db8924dd7fad99acbdd0ff4bdd1"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:a84cd3aba44220008bdbd81506703de37"><td class="memItemLeft" align="right" valign="top">const MatrixType & </td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1ColPivHouseholderQR.html#a84cd3aba44220008bdbd81506703de37">matrixR</a> () const </td></tr> <tr class="separator:a84cd3aba44220008bdbd81506703de37"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:a067e9d4143ce0558fc684b736128a4ed"><td class="memItemLeft" align="right" valign="top">RealScalar </td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1ColPivHouseholderQR.html#a067e9d4143ce0558fc684b736128a4ed">maxPivot</a> () const </td></tr> <tr class="separator:a067e9d4143ce0558fc684b736128a4ed"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:a6dda9285f13dec9f49e9c17229a89988"><td class="memItemLeft" align="right" valign="top">Index </td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1ColPivHouseholderQR.html#a6dda9285f13dec9f49e9c17229a89988">nonzeroPivots</a> () const </td></tr> <tr class="separator:a6dda9285f13dec9f49e9c17229a89988"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:a363d1c09d77f09d6ea2d2789776e7be3"><td class="memItemLeft" align="right" valign="top">Index </td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1ColPivHouseholderQR.html#a363d1c09d77f09d6ea2d2789776e7be3">rank</a> () const </td></tr> <tr class="separator:a363d1c09d77f09d6ea2d2789776e7be3"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:a5631914aef48f4f719789c823edacb2c"><td class="memItemLeft" align="right" valign="top"><a class="el" href="classEigen_1_1ColPivHouseholderQR.html">ColPivHouseholderQR</a> & </td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1ColPivHouseholderQR.html#a5631914aef48f4f719789c823edacb2c">setThreshold</a> (const RealScalar &<a class="el" href="classEigen_1_1ColPivHouseholderQR.html#aa5a87faaa92a3081045d1f934e292ef0">threshold</a>)</td></tr> <tr class="separator:a5631914aef48f4f719789c823edacb2c"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:a446f775405186b238f20bd092fd107a6"><td class="memItemLeft" align="right" valign="top"><a class="el" href="classEigen_1_1ColPivHouseholderQR.html">ColPivHouseholderQR</a> & </td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1ColPivHouseholderQR.html#a446f775405186b238f20bd092fd107a6">setThreshold</a> (Default_t)</td></tr> <tr class="separator:a446f775405186b238f20bd092fd107a6"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:afbe1cd1202964011ae7e7411577749a0"><td class="memTemplParams" colspan="2">template<typename Rhs > </td></tr> <tr class="memitem:afbe1cd1202964011ae7e7411577749a0"><td class="memTemplItemLeft" align="right" valign="top">const internal::solve_retval<br class="typebreak"/> < <a class="el" href="classEigen_1_1ColPivHouseholderQR.html">ColPivHouseholderQR</a>, Rhs > </td><td class="memTemplItemRight" valign="bottom"><a class="el" href="classEigen_1_1ColPivHouseholderQR.html#afbe1cd1202964011ae7e7411577749a0">solve</a> (const <a class="el" href="classEigen_1_1MatrixBase.html">MatrixBase</a>< Rhs > &b) const </td></tr> <tr class="separator:afbe1cd1202964011ae7e7411577749a0"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:aa5a87faaa92a3081045d1f934e292ef0"><td class="memItemLeft" align="right" valign="top">RealScalar </td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1ColPivHouseholderQR.html#aa5a87faaa92a3081045d1f934e292ef0">threshold</a> () const </td></tr> <tr class="separator:aa5a87faaa92a3081045d1f934e292ef0"><td class="memSeparator" colspan="2"> </td></tr> </table> <h2 class="groupheader">Constructor & Destructor Documentation</h2> <a class="anchor" id="a72f7fb8c1359fa8d9743068e5699228e"></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_1ColPivHouseholderQR.html">ColPivHouseholderQR</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_1ColPivHouseholderQR.html#a812205d31cae2005b7779c7a3a442f1b">ColPivHouseholderQR::compute(const MatrixType&)</a>. </p> </div> </div> <a class="anchor" id="a4a15f3e5e33478b9969d8c71aabd09ea"></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_1ColPivHouseholderQR.html">ColPivHouseholderQR</a> </td> <td>(</td> <td class="paramtype">Index </td> <td class="paramname"><em>rows</em>, </td> </tr> <tr> <td class="paramkey"></td> <td></td> <td class="paramtype">Index </td> <td class="paramname"><em>cols</em> </td> </tr> <tr> <td></td> <td>)</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_1ColPivHouseholderQR.html#a72f7fb8c1359fa8d9743068e5699228e" title="Default Constructor. ">ColPivHouseholderQR()</a> </dd></dl> </div> </div> <a class="anchor" id="ac93cc6a4fd4247486c2c8eab79b4c0ad"></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_1ColPivHouseholderQR.html">ColPivHouseholderQR</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>Constructs a QR factorization from a given matrix. </p> <p>This constructor computes the QR factorization of the matrix <em>matrix</em> by calling the method <a class="el" href="classEigen_1_1ColPivHouseholderQR.html#a812205d31cae2005b7779c7a3a442f1b">compute()</a>. It is a short cut for:</p> <div class="fragment"><div class="line">* ColPivHouseholderQR<MatrixType> qr(matrix.rows(), matrix.cols());</div> <div class="line">* qr.compute(matrix);</div> <div class="line">* </div> </div><!-- fragment --><dl class="section see"><dt>See Also</dt><dd><a class="el" href="classEigen_1_1ColPivHouseholderQR.html#a812205d31cae2005b7779c7a3a442f1b">compute()</a> </dd></dl> <p>References <a class="el" href="classEigen_1_1ColPivHouseholderQR.html#a812205d31cae2005b7779c7a3a442f1b">ColPivHouseholderQR< MatrixType >::compute()</a>.</p> </div> </div> <h2 class="groupheader">Member Function Documentation</h2> <a class="anchor" id="a2bb8fa5593d2f447c16ba2e3cfb9c432"></a> <div class="memitem"> <div class="memproto"> <table class="memname"> <tr> <td class="memname">MatrixType::RealScalar absDeterminant </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 absolute value of the determinant of the matrix of which *this is the QR decomposition. It has only linear complexity (that is, O(n) where n is the dimension of the square matrix) as the QR decomposition has already been computed.</dd></dl> <dl class="section note"><dt>Note</dt><dd>This is only for square matrices.</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_1ColPivHouseholderQR.html#a82a3f15f0cda1f4f874d50951344b5a9">logAbsDeterminant()</a> instead.</dd></dl> <dl class="section see"><dt>See Also</dt><dd><a class="el" href="classEigen_1_1ColPivHouseholderQR.html#a82a3f15f0cda1f4f874d50951344b5a9">logAbsDeterminant()</a>, <a class="el" href="classEigen_1_1MatrixBase.html#ad63cea11a4bf220298dce4489a1704c7">MatrixBase::determinant()</a> </dd></dl> </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>& 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 const reference to the column permutation matrix </dd></dl> </div> </div> <a class="anchor" id="a812205d31cae2005b7779c7a3a442f1b"></a> <div class="memitem"> <div class="memproto"> <table class="memname"> <tr> <td class="memname"><a class="el" href="classEigen_1_1ColPivHouseholderQR.html">ColPivHouseholderQR</a>< MatrixType > & compute </td> <td>(</td> <td class="paramtype">const MatrixType & </td> <td class="paramname"><em>matrix</em></td><td>)</td> <td></td> </tr> </table> </div><div class="memdoc"> <p>Performs the QR factorization of the given matrix <em>matrix</em>. The result of the factorization is stored into <code>*this</code>, and a reference to <code>*this</code> is returned.</p> <dl class="section see"><dt>See Also</dt><dd>class <a class="el" href="classEigen_1_1ColPivHouseholderQR.html" title="Householder rank-revealing QR decomposition of a matrix with column-pivoting. ">ColPivHouseholderQR</a>, <a class="el" href="classEigen_1_1ColPivHouseholderQR.html#ac93cc6a4fd4247486c2c8eab79b4c0ad" title="Constructs a QR factorization from a given matrix. ">ColPivHouseholderQR(const MatrixType&)</a> </dd></dl> <p>Referenced by <a class="el" href="classEigen_1_1ColPivHouseholderQR.html#ac93cc6a4fd4247486c2c8eab79b4c0ad">ColPivHouseholderQR< MatrixType >::ColPivHouseholderQR()</a>.</p> </div> </div> <a class="anchor" id="a7c6323871c4f080fc6e2d3ad7fc607fc"></a> <div class="memitem"> <div class="memproto"> <table class="mlabels"> <tr> <td class="mlabels-left"> <table class="memname"> <tr> <td class="memname">Index dimensionOfKernel </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 dimension of the kernel of the matrix of which *this is the QR decomposition.</dd></dl> <dl class="section note"><dt>Note</dt><dd>This method has to determine which pivots should be considered nonzero. For that, it uses the threshold value that you can control by calling <a class="el" href="classEigen_1_1ColPivHouseholderQR.html#a5631914aef48f4f719789c823edacb2c">setThreshold(const RealScalar&)</a>. </dd></dl> <p>References <a class="el" href="classEigen_1_1ColPivHouseholderQR.html#a363d1c09d77f09d6ea2d2789776e7be3">ColPivHouseholderQR< MatrixType >::rank()</a>.</p> </div> </div> <a class="anchor" id="a321e96844544eacde2272c3a487048e9"></a> <div class="memitem"> <div class="memproto"> <table class="mlabels"> <tr> <td class="mlabels-left"> <table class="memname"> <tr> <td class="memname">const HCoeffsType& hCoeffs </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 const reference to the vector of Householder coefficients used to represent the factor <code>Q</code>.</dd></dl> <p>For advanced uses only. </p> </div> </div> <a class="anchor" id="a97150efba7e2da32f48a7f823f047a89"></a> <div class="memitem"> <div class="memproto"> <table class="memname"> <tr> <td class="memname"><a class="el" href="classEigen_1_1ColPivHouseholderQR.html">ColPivHouseholderQR</a>< MatrixType >::<a class="el" href="classEigen_1_1HouseholderSequence.html">HouseholderSequenceType</a> householderQ </td> <td>(</td> <td class="paramtype">void </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 Q as a sequence of householder transformations </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"> <p>Reports whether the QR factorization was succesful. </p> <dl class="section note"><dt>Note</dt><dd>This function always returns <code>Success</code>. It is provided for compatibility with other factorization routines. </dd></dl> <dl class="section return"><dt>Returns</dt><dd><code>Success</code> </dd></dl> <p>References <a class="el" href="group__enums.html#gga51bc1ac16f26ebe51eae1abb77bd037bafdfbdf3247bd36a1f17270d5cec74c9c">Eigen::Success</a>.</p> </div> </div> <a class="anchor" id="a3021b1285c868ef94b26f387766bdf92"></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_1ColPivHouseholderQR.html">ColPivHouseholderQR</a>, typename MatrixType::IdentityReturnType> inverse </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 inverse of the matrix of which *this is the QR decomposition.</dd></dl> <dl class="section note"><dt>Note</dt><dd>If this matrix is not invertible, the returned matrix has undefined coefficients. Use <a class="el" href="classEigen_1_1ColPivHouseholderQR.html#ab60c7d993c9eba31668fb8886d621094">isInvertible()</a> to first determine whether this matrix is invertible. </dd></dl> </div> </div> <a class="anchor" id="a1e119085e53eca65e9ba15451c102d40"></a> <div class="memitem"> <div class="memproto"> <table class="mlabels"> <tr> <td class="mlabels-left"> <table class="memname"> <tr> <td class="memname">bool isInjective </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 of which *this is the QR decomposition represents an injective linear map, i.e. has trivial kernel; false otherwise.</dd></dl> <dl class="section note"><dt>Note</dt><dd>This method has to determine which pivots should be considered nonzero. For that, it uses the threshold value that you can control by calling <a class="el" href="classEigen_1_1ColPivHouseholderQR.html#a5631914aef48f4f719789c823edacb2c">setThreshold(const RealScalar&)</a>. </dd></dl> <p>References <a class="el" href="classEigen_1_1ColPivHouseholderQR.html#a363d1c09d77f09d6ea2d2789776e7be3">ColPivHouseholderQR< MatrixType >::rank()</a>.</p> <p>Referenced by <a class="el" href="classEigen_1_1ColPivHouseholderQR.html#ab60c7d993c9eba31668fb8886d621094">ColPivHouseholderQR< MatrixType >::isInvertible()</a>.</p> </div> </div> <a class="anchor" id="ab60c7d993c9eba31668fb8886d621094"></a> <div class="memitem"> <div class="memproto"> <table class="mlabels"> <tr> <td class="mlabels-left"> <table class="memname"> <tr> <td class="memname">bool isInvertible </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 of which *this is the QR decomposition is invertible.</dd></dl> <dl class="section note"><dt>Note</dt><dd>This method has to determine which pivots should be considered nonzero. For that, it uses the threshold value that you can control by calling <a class="el" href="classEigen_1_1ColPivHouseholderQR.html#a5631914aef48f4f719789c823edacb2c">setThreshold(const RealScalar&)</a>. </dd></dl> <p>References <a class="el" href="classEigen_1_1ColPivHouseholderQR.html#a1e119085e53eca65e9ba15451c102d40">ColPivHouseholderQR< MatrixType >::isInjective()</a>, and <a class="el" href="classEigen_1_1ColPivHouseholderQR.html#a5968b9ca46303b3cc7250e7b120ab7e6">ColPivHouseholderQR< MatrixType >::isSurjective()</a>.</p> </div> </div> <a class="anchor" id="a5968b9ca46303b3cc7250e7b120ab7e6"></a> <div class="memitem"> <div class="memproto"> <table class="mlabels"> <tr> <td class="mlabels-left"> <table class="memname"> <tr> <td class="memname">bool isSurjective </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 of which *this is the QR decomposition represents a surjective linear map; false otherwise.</dd></dl> <dl class="section note"><dt>Note</dt><dd>This method has to determine which pivots should be considered nonzero. For that, it uses the threshold value that you can control by calling <a class="el" href="classEigen_1_1ColPivHouseholderQR.html#a5631914aef48f4f719789c823edacb2c">setThreshold(const RealScalar&)</a>. </dd></dl> <p>References <a class="el" href="classEigen_1_1ColPivHouseholderQR.html#a363d1c09d77f09d6ea2d2789776e7be3">ColPivHouseholderQR< MatrixType >::rank()</a>.</p> <p>Referenced by <a class="el" href="classEigen_1_1ColPivHouseholderQR.html#ab60c7d993c9eba31668fb8886d621094">ColPivHouseholderQR< MatrixType >::isInvertible()</a>.</p> </div> </div> <a class="anchor" id="a82a3f15f0cda1f4f874d50951344b5a9"></a> <div class="memitem"> <div class="memproto"> <table class="memname"> <tr> <td class="memname">MatrixType::RealScalar logAbsDeterminant </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 natural log of the absolute value of the determinant of the matrix of which *this is the QR decomposition. It has only linear complexity (that is, O(n) where n is the dimension of the square matrix) as the QR decomposition has already been computed.</dd></dl> <dl class="section note"><dt>Note</dt><dd>This is only for square matrices.</dd> <dd> This method is useful to work around the risk of overflow/underflow that's inherent to determinant computation.</dd></dl> <dl class="section see"><dt>See Also</dt><dd><a class="el" href="classEigen_1_1ColPivHouseholderQR.html#a2bb8fa5593d2f447c16ba2e3cfb9c432">absDeterminant()</a>, <a class="el" href="classEigen_1_1MatrixBase.html#ad63cea11a4bf220298dce4489a1704c7">MatrixBase::determinant()</a> </dd></dl> </div> </div> <a class="anchor" id="a7b142db8924dd7fad99acbdd0ff4bdd1"></a> <div class="memitem"> <div class="memproto"> <table class="mlabels"> <tr> <td class="mlabels-left"> <table class="memname"> <tr> <td class="memname">const MatrixType& matrixQR </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 matrix where the Householder QR decomposition is stored </dd></dl> </div> </div> <a class="anchor" id="a84cd3aba44220008bdbd81506703de37"></a> <div class="memitem"> <div class="memproto"> <table class="mlabels"> <tr> <td class="mlabels-left"> <table class="memname"> <tr> <td class="memname">const MatrixType& matrixR </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 matrix where the result Householder QR is stored </dd></dl> <dl class="section warning"><dt>Warning</dt><dd>The strict lower part of this matrix contains internal values. Only the upper triangular part should be referenced. To get it, use <div class="fragment"><div class="line"><a class="code" href="classEigen_1_1ColPivHouseholderQR.html#a84cd3aba44220008bdbd81506703de37">matrixR</a>().template triangularView<Upper>() </div> </div><!-- fragment --> For rank-deficient matrices, use <div class="fragment"><div class="line">* <a class="code" href="classEigen_1_1ColPivHouseholderQR.html#a84cd3aba44220008bdbd81506703de37">matrixR</a>().topLeftCorner(<a class="code" href="classEigen_1_1ColPivHouseholderQR.html#a363d1c09d77f09d6ea2d2789776e7be3">rank</a>(), <a class="code" href="classEigen_1_1ColPivHouseholderQR.html#a363d1c09d77f09d6ea2d2789776e7be3">rank</a>()).template triangularView<Upper>() </div> <div class="line">* </div> </div><!-- fragment --> </dd></dl> </div> </div> <a class="anchor" id="a067e9d4143ce0558fc684b736128a4ed"></a> <div class="memitem"> <div class="memproto"> <table class="mlabels"> <tr> <td class="mlabels-left"> <table class="memname"> <tr> <td class="memname">RealScalar maxPivot </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 absolute value of the biggest pivot, i.e. the biggest diagonal coefficient of R. </dd></dl> </div> </div> <a class="anchor" id="a6dda9285f13dec9f49e9c17229a89988"></a> <div class="memitem"> <div class="memproto"> <table class="mlabels"> <tr> <td class="mlabels-left"> <table class="memname"> <tr> <td class="memname">Index nonzeroPivots </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 nonzero pivots in the QR decomposition. Here nonzero is meant in the exact sense, not in a fuzzy sense. So that notion isn't really intrinsically interesting, but it is still useful when implementing algorithms.</dd></dl> <dl class="section see"><dt>See Also</dt><dd><a class="el" href="classEigen_1_1ColPivHouseholderQR.html#a363d1c09d77f09d6ea2d2789776e7be3">rank()</a> </dd></dl> </div> </div> <a class="anchor" id="a363d1c09d77f09d6ea2d2789776e7be3"></a> <div class="memitem"> <div class="memproto"> <table class="mlabels"> <tr> <td class="mlabels-left"> <table class="memname"> <tr> <td class="memname">Index rank </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 rank of the matrix of which *this is the QR decomposition.</dd></dl> <dl class="section note"><dt>Note</dt><dd>This method has to determine which pivots should be considered nonzero. For that, it uses the threshold value that you can control by calling <a class="el" href="classEigen_1_1ColPivHouseholderQR.html#a5631914aef48f4f719789c823edacb2c">setThreshold(const RealScalar&)</a>. </dd></dl> <p>References <a class="el" href="classEigen_1_1ColPivHouseholderQR.html#aa5a87faaa92a3081045d1f934e292ef0">ColPivHouseholderQR< MatrixType >::threshold()</a>.</p> <p>Referenced by <a class="el" href="classEigen_1_1ColPivHouseholderQR.html#a7c6323871c4f080fc6e2d3ad7fc607fc">ColPivHouseholderQR< MatrixType >::dimensionOfKernel()</a>, <a class="el" href="classEigen_1_1ColPivHouseholderQR.html#a1e119085e53eca65e9ba15451c102d40">ColPivHouseholderQR< MatrixType >::isInjective()</a>, and <a class="el" href="classEigen_1_1ColPivHouseholderQR.html#a5968b9ca46303b3cc7250e7b120ab7e6">ColPivHouseholderQR< MatrixType >::isSurjective()</a>.</p> </div> </div> <a class="anchor" id="a5631914aef48f4f719789c823edacb2c"></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_1ColPivHouseholderQR.html">ColPivHouseholderQR</a>& setThreshold </td> <td>(</td> <td class="paramtype">const RealScalar & </td> <td class="paramname"><em>threshold</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>Allows to prescribe a threshold to be used by certain methods, such as <a class="el" href="classEigen_1_1ColPivHouseholderQR.html#a363d1c09d77f09d6ea2d2789776e7be3">rank()</a>, who need to determine when pivots are to be considered nonzero. This is not used for the QR decomposition itself.</p> <p>When it needs to get the threshold value, <a class="el" href="namespaceEigen.html" title="Namespace containing all symbols from the Eigen library. ">Eigen</a> calls <a class="el" href="classEigen_1_1ColPivHouseholderQR.html#aa5a87faaa92a3081045d1f934e292ef0">threshold()</a>. By default, this uses a formula to automatically determine a reasonable threshold. Once you have called the present method <a class="el" href="classEigen_1_1ColPivHouseholderQR.html#a5631914aef48f4f719789c823edacb2c">setThreshold(const RealScalar&)</a>, your value is used instead.</p> <dl class="params"><dt>Parameters</dt><dd> <table class="params"> <tr><td class="paramname">threshold</td><td>The new value to use as the threshold.</td></tr> </table> </dd> </dl> <p>A pivot will be considered nonzero if its absolute value is strictly greater than <img class="formulaInl" alt="$ \vert pivot \vert \leqslant threshold \times \vert maxpivot \vert $" src="form_156.png"/> where maxpivot is the biggest pivot.</p> <p>If you want to come back to the default behavior, call <a class="el" href="classEigen_1_1ColPivHouseholderQR.html#a446f775405186b238f20bd092fd107a6">setThreshold(Default_t)</a> </p> <p>References <a class="el" href="classEigen_1_1ColPivHouseholderQR.html#aa5a87faaa92a3081045d1f934e292ef0">ColPivHouseholderQR< MatrixType >::threshold()</a>.</p> </div> </div> <a class="anchor" id="a446f775405186b238f20bd092fd107a6"></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_1ColPivHouseholderQR.html">ColPivHouseholderQR</a>& setThreshold </td> <td>(</td> <td class="paramtype">Default_t </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>Allows to come back to the default behavior, letting <a class="el" href="namespaceEigen.html" title="Namespace containing all symbols from the Eigen library. ">Eigen</a> use its default formula for determining the threshold.</p> <p>You should pass the special object Eigen::Default as parameter here. </p> <div class="fragment"><div class="line">qr.setThreshold(Eigen::Default); </div> </div><!-- fragment --><p>See the documentation of <a class="el" href="classEigen_1_1ColPivHouseholderQR.html#a5631914aef48f4f719789c823edacb2c">setThreshold(const RealScalar&)</a>. </p> </div> </div> <a class="anchor" id="afbe1cd1202964011ae7e7411577749a0"></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_1ColPivHouseholderQR.html">ColPivHouseholderQR</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"> <p>This method finds a solution x to the equation Ax=b, where A is the matrix of which *this is the QR decomposition, if any exists.</p> <dl class="params"><dt>Parameters</dt><dd> <table class="params"> <tr><td class="paramname">b</td><td>the right-hand-side of the equation to solve.</td></tr> </table> </dd> </dl> <dl class="section return"><dt>Returns</dt><dd>a solution.</dd></dl> <dl class="section note"><dt>Note</dt><dd>The case where b is a matrix is not yet implemented. Also, this code is space inefficient.</dd></dl> <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>If there exists more than one solution, this method will arbitrarily choose one.</p> <p>Example: </p> <div class="fragment"><div class="line"><a class="code" href="group__matrixtypedefs.html#ga2b787393a7fc870d99aa634f60b2510c">Matrix3f</a> m = <a class="code" href="classEigen_1_1DenseBase.html#a8e759dafdd9ecc446d397b7f5435f60a">Matrix3f::Random</a>();</div> <div class="line"><a class="code" href="group__matrixtypedefs.html#ga2b787393a7fc870d99aa634f60b2510c">Matrix3f</a> y = <a class="code" href="classEigen_1_1DenseBase.html#a8e759dafdd9ecc446d397b7f5435f60a">Matrix3f::Random</a>();</div> <div class="line">cout << <span class="stringliteral">"Here is the matrix m:"</span> << endl << m << endl;</div> <div class="line">cout << <span class="stringliteral">"Here is the matrix y:"</span> << endl << y << endl;</div> <div class="line"><a class="code" href="group__matrixtypedefs.html#ga2b787393a7fc870d99aa634f60b2510c">Matrix3f</a> x;</div> <div class="line">x = m.colPivHouseholderQr().solve(y);</div> <div class="line">assert(y.isApprox(m*x));</div> <div class="line">cout << <span class="stringliteral">"Here is a solution x to the equation mx=y:"</span> << endl << x << endl;</div> </div><!-- fragment --><p> Output: </p> <pre class="fragment">Here is the matrix m: 0.68 0.597 -0.33 -0.211 0.823 0.536 0.566 -0.605 -0.444 Here is the matrix y: 0.108 -0.27 0.832 -0.0452 0.0268 0.271 0.258 0.904 0.435 Here is a solution x to the equation mx=y: 0.609 2.68 1.67 -0.231 -1.57 0.0713 0.51 3.51 1.05 </pre> </div> </div> <a class="anchor" id="aa5a87faaa92a3081045d1f934e292ef0"></a> <div class="memitem"> <div class="memproto"> <table class="mlabels"> <tr> <td class="mlabels-left"> <table class="memname"> <tr> <td class="memname">RealScalar threshold </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>Returns the threshold that will be used by certain methods such as <a class="el" href="classEigen_1_1ColPivHouseholderQR.html#a363d1c09d77f09d6ea2d2789776e7be3">rank()</a>.</p> <p>See the documentation of <a class="el" href="classEigen_1_1ColPivHouseholderQR.html#a5631914aef48f4f719789c823edacb2c">setThreshold(const RealScalar&)</a>. </p> <p>Referenced by <a class="el" href="classEigen_1_1ColPivHouseholderQR.html#a363d1c09d77f09d6ea2d2789776e7be3">ColPivHouseholderQR< MatrixType >::rank()</a>, and <a class="el" href="classEigen_1_1ColPivHouseholderQR.html#a5631914aef48f4f719789c823edacb2c">ColPivHouseholderQR< MatrixType >::setThreshold()</a>.</p> </div> </div> <hr/>The documentation for this class was generated from the following files:<ul> <li><a class="el" href="ForwardDeclarations_8h_source.html">ForwardDeclarations.h</a></li> <li><a class="el" href="ColPivHouseholderQR_8h_source.html">ColPivHouseholderQR.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_1ColPivHouseholderQR.html">ColPivHouseholderQR</a></li> <li class="footer">Generated on Mon Oct 28 2013 11:04:28 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>