<!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: SelfAdjointView< 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_1SelfAdjointView.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_1SelfAdjointView-members.html">List of all members</a> | <a href="#pub-types">Public Types</a> | <a href="#pub-methods">Public Member Functions</a> | <a href="#friends">Friends</a> </div> <div class="headertitle"> <div class="title">SelfAdjointView< MatrixType, UpLo > Class Template Reference<div class="ingroups"><a class="el" href="group__Core__Module.html">Core 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, unsigned int UpLo><br/> class Eigen::SelfAdjointView< MatrixType, UpLo ></h3> <p>Expression of a selfadjoint matrix from a triangular part of a dense matrix. </p> <dl class="params"><dt>Parameters</dt><dd> <table class="params"> <tr><td class="paramname">MatrixType</td><td>the type of the dense matrix storing the coefficients </td></tr> <tr><td class="paramname">TriangularPart</td><td>can be either <code><a class="el" href="group__enums.html#ggab59c1bec446b10af208f977a871d910baf886b397626076218462d53d50eb96bc">Lower</a></code> or <code><a class="el" href="group__enums.html#ggab59c1bec446b10af208f977a871d910bae70afef0d3ff7aca74e17e85ff6c9f2e">Upper</a></code> </td></tr> </table> </dd> </dl> <p>This class is an expression of a sefladjoint matrix from a triangular part of a matrix with given dense storage of the coefficients. It is the return type of MatrixBase::selfadjointView() and most of the time this is the only way that it is used.</p> <dl class="section see"><dt>See Also</dt><dd>class TriangularBase, MatrixBase::selfadjointView() </dd></dl> </div> <p>Inherits TriangularBase< Derived >.</p> <table class="memberdecls"> <tr class="heading"><td colspan="2"><h2 class="groupheader"><a name="pub-types"></a> Public Types</h2></td></tr> <tr class="memitem:af7cc63140d17b562998188356f30da6f"><td class="memItemLeft" align="right" valign="top">typedef <a class="el" href="classEigen_1_1Matrix.html">Matrix</a>< <a class="el" href="classEigen_1_1SelfAdjointView.html#acb5c3dc237f99cf17167e8a629f01b43">RealScalar</a>, <br class="typebreak"/> internal::traits< MatrixType ><br class="typebreak"/> ::ColsAtCompileTime, 1 > </td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1SelfAdjointView.html#af7cc63140d17b562998188356f30da6f">EigenvaluesReturnType</a></td></tr> <tr class="separator:af7cc63140d17b562998188356f30da6f"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:acb5c3dc237f99cf17167e8a629f01b43"><td class="memItemLeft" align="right" valign="top">typedef <a class="el" href="structEigen_1_1NumTraits.html">NumTraits</a>< <a class="el" href="classEigen_1_1SelfAdjointView.html#a9ef555098df604af49b0a00d389d823b">Scalar</a> >::Real </td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1SelfAdjointView.html#acb5c3dc237f99cf17167e8a629f01b43">RealScalar</a></td></tr> <tr class="separator:acb5c3dc237f99cf17167e8a629f01b43"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:a9ef555098df604af49b0a00d389d823b"><td class="memItemLeft" align="right" valign="top"><a class="anchor" id="a9ef555098df604af49b0a00d389d823b"></a> typedef internal::traits<br class="typebreak"/> < <a class="el" href="classEigen_1_1SelfAdjointView.html">SelfAdjointView</a> >::<a class="el" href="classEigen_1_1SelfAdjointView.html#a9ef555098df604af49b0a00d389d823b">Scalar</a> </td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1SelfAdjointView.html#a9ef555098df604af49b0a00d389d823b">Scalar</a></td></tr> <tr class="memdesc:a9ef555098df604af49b0a00d389d823b"><td class="mdescLeft"> </td><td class="mdescRight">The type of coefficients in this matrix. <br/></td></tr> <tr class="separator:a9ef555098df604af49b0a00d389d823b"><td class="memSeparator" colspan="2"> </td></tr> </table><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:afacbe8d905dd8b00eb52e72c2e938aed"><td class="memItemLeft" align="right" valign="top"><a class="el" href="classEigen_1_1SelfAdjointView.html#a9ef555098df604af49b0a00d389d823b">Scalar</a> </td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1SelfAdjointView.html#afacbe8d905dd8b00eb52e72c2e938aed">coeff</a> (Index row, Index col) const </td></tr> <tr class="separator:afacbe8d905dd8b00eb52e72c2e938aed"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:a4386da578a119f7b500d0b22cab2de0e"><td class="memItemLeft" align="right" valign="top"><a class="el" href="classEigen_1_1SelfAdjointView.html#a9ef555098df604af49b0a00d389d823b">Scalar</a> & </td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1SelfAdjointView.html#a4386da578a119f7b500d0b22cab2de0e">coeffRef</a> (Index row, Index col)</td></tr> <tr class="separator:a4386da578a119f7b500d0b22cab2de0e"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:af7177136c93ee65b7c41c1798e23292d"><td class="memTemplParams" colspan="2">template<typename Other > </td></tr> <tr class="memitem:af7177136c93ee65b7c41c1798e23292d"><td class="memTemplItemLeft" align="right" valign="top">void </td><td class="memTemplItemRight" valign="bottom"><b>copyCoeff</b> (Index row, Index col, Other &other)</td></tr> <tr class="separator:af7177136c93ee65b7c41c1798e23292d"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:aa84222add803ad7c9db07dd4dd91d5d9"><td class="memItemLeft" align="right" valign="top">Derived & </td><td class="memItemRight" valign="bottom"><a class="el" href="structEigen_1_1EigenBase.html#aa84222add803ad7c9db07dd4dd91d5d9">derived</a> ()</td></tr> <tr class="separator:aa84222add803ad7c9db07dd4dd91d5d9"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:a2850de52161dfe2bb6e380f1892276b8"><td class="memItemLeft" align="right" valign="top">const Derived & </td><td class="memItemRight" valign="bottom"><a class="el" href="structEigen_1_1EigenBase.html#a2850de52161dfe2bb6e380f1892276b8">derived</a> () const </td></tr> <tr class="separator:a2850de52161dfe2bb6e380f1892276b8"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:a847d0db612cd20e1b796c6a70d0f0d54"><td class="memItemLeft" align="right" valign="top"><a class="el" href="classEigen_1_1SelfAdjointView.html#af7cc63140d17b562998188356f30da6f">EigenvaluesReturnType</a> </td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1SelfAdjointView.html#a847d0db612cd20e1b796c6a70d0f0d54">eigenvalues</a> () const </td></tr> <tr class="memdesc:a847d0db612cd20e1b796c6a70d0f0d54"><td class="mdescLeft"> </td><td class="mdescRight">Computes the eigenvalues of a matrix. <a href="#a847d0db612cd20e1b796c6a70d0f0d54">More...</a><br/></td></tr> <tr class="separator:a847d0db612cd20e1b796c6a70d0f0d54"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:a66ac5c9b6c04321ab12364bba06fee66"><td class="memTemplParams" colspan="2">template<typename DenseDerived > </td></tr> <tr class="memitem:a66ac5c9b6c04321ab12364bba06fee66"><td class="memTemplItemLeft" align="right" valign="top">void </td><td class="memTemplItemRight" valign="bottom"><b>evalTo</b> (<a class="el" href="classEigen_1_1MatrixBase.html">MatrixBase</a>< DenseDerived > &other) const </td></tr> <tr class="separator:a66ac5c9b6c04321ab12364bba06fee66"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:affcc5fe7ffb333e0b99a5532d086cd2f"><td class="memTemplParams" colspan="2">template<typename DenseDerived > </td></tr> <tr class="memitem:affcc5fe7ffb333e0b99a5532d086cd2f"><td class="memTemplItemLeft" align="right" valign="top">void </td><td class="memTemplItemRight" valign="bottom"><b>evalToLazy</b> (<a class="el" href="classEigen_1_1MatrixBase.html">MatrixBase</a>< DenseDerived > &other) const </td></tr> <tr class="separator:affcc5fe7ffb333e0b99a5532d086cd2f"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:a245df21a519cd34f5eb6f724cf5bf8e7"><td class="memItemLeft" align="right" valign="top">const <a class="el" href="classEigen_1_1LDLT.html">LDLT</a>< PlainObject, UpLo > </td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1SelfAdjointView.html#a245df21a519cd34f5eb6f724cf5bf8e7">ldlt</a> () const </td></tr> <tr class="separator:a245df21a519cd34f5eb6f724cf5bf8e7"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:ab6d4081f604fb2bc15ad285cf6bf54a2"><td class="memItemLeft" align="right" valign="top">const <a class="el" href="classEigen_1_1LLT.html">LLT</a>< PlainObject, UpLo > </td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1SelfAdjointView.html#ab6d4081f604fb2bc15ad285cf6bf54a2">llt</a> () const </td></tr> <tr class="separator:ab6d4081f604fb2bc15ad285cf6bf54a2"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:a8b615fb312e03e69a42aa86f50c448f1"><td class="memTemplParams" colspan="2">template<typename OtherDerived > </td></tr> <tr class="memitem:a8b615fb312e03e69a42aa86f50c448f1"><td class="memTemplItemLeft" align="right" valign="top">SelfadjointProductMatrix<br class="typebreak"/> < MatrixType, Mode, false, <br class="typebreak"/> OtherDerived, <br class="typebreak"/> 0, OtherDerived::IsVectorAtCompileTime > </td><td class="memTemplItemRight" valign="bottom"><a class="el" href="classEigen_1_1SelfAdjointView.html#a8b615fb312e03e69a42aa86f50c448f1">operator*</a> (const <a class="el" href="classEigen_1_1MatrixBase.html">MatrixBase</a>< OtherDerived > &rhs) const </td></tr> <tr class="separator:a8b615fb312e03e69a42aa86f50c448f1"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:a67279bc5b39fba3763acabf1d52ecad8"><td class="memItemLeft" align="right" valign="top"><a class="el" href="classEigen_1_1SelfAdjointView.html#acb5c3dc237f99cf17167e8a629f01b43">RealScalar</a> </td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1SelfAdjointView.html#a67279bc5b39fba3763acabf1d52ecad8">operatorNorm</a> () const </td></tr> <tr class="memdesc:a67279bc5b39fba3763acabf1d52ecad8"><td class="mdescLeft"> </td><td class="mdescRight">Computes the L2 operator norm. <a href="#a67279bc5b39fba3763acabf1d52ecad8">More...</a><br/></td></tr> <tr class="separator:a67279bc5b39fba3763acabf1d52ecad8"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:a024a26e49b7865e23e39479de6dc56fc"><td class="memTemplParams" colspan="2">template<typename DerivedU , typename DerivedV > </td></tr> <tr class="memitem:a024a26e49b7865e23e39479de6dc56fc"><td class="memTemplItemLeft" align="right" valign="top"><a class="el" href="classEigen_1_1SelfAdjointView.html">SelfAdjointView</a> & </td><td class="memTemplItemRight" valign="bottom"><a class="el" href="classEigen_1_1SelfAdjointView.html#a024a26e49b7865e23e39479de6dc56fc">rankUpdate</a> (const <a class="el" href="classEigen_1_1MatrixBase.html">MatrixBase</a>< DerivedU > &u, const <a class="el" href="classEigen_1_1MatrixBase.html">MatrixBase</a>< DerivedV > &v, const <a class="el" href="classEigen_1_1SelfAdjointView.html#a9ef555098df604af49b0a00d389d823b">Scalar</a> &alpha=<a class="el" href="classEigen_1_1SelfAdjointView.html#a9ef555098df604af49b0a00d389d823b">Scalar</a>(1))</td></tr> <tr class="separator:a024a26e49b7865e23e39479de6dc56fc"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:aad06bacc2e9b8345ae92279acc362584"><td class="memTemplParams" colspan="2">template<typename DerivedU > </td></tr> <tr class="memitem:aad06bacc2e9b8345ae92279acc362584"><td class="memTemplItemLeft" align="right" valign="top"><a class="el" href="classEigen_1_1SelfAdjointView.html">SelfAdjointView</a> & </td><td class="memTemplItemRight" valign="bottom"><a class="el" href="classEigen_1_1SelfAdjointView.html#aad06bacc2e9b8345ae92279acc362584">rankUpdate</a> (const <a class="el" href="classEigen_1_1MatrixBase.html">MatrixBase</a>< DerivedU > &u, const <a class="el" href="classEigen_1_1SelfAdjointView.html#a9ef555098df604af49b0a00d389d823b">Scalar</a> &alpha=<a class="el" href="classEigen_1_1SelfAdjointView.html#a9ef555098df604af49b0a00d389d823b">Scalar</a>(1))</td></tr> <tr class="separator:aad06bacc2e9b8345ae92279acc362584"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:a76f5bc8a03ec105ab4be1f2b91b7b5d5"><td class="memItemLeft" align="right" valign="top">Index </td><td class="memItemRight" valign="bottom"><a class="el" href="structEigen_1_1EigenBase.html#a76f5bc8a03ec105ab4be1f2b91b7b5d5">size</a> () const </td></tr> <tr class="separator:a76f5bc8a03ec105ab4be1f2b91b7b5d5"><td class="memSeparator" colspan="2"> </td></tr> </table><table class="memberdecls"> <tr class="heading"><td colspan="2"><h2 class="groupheader"><a name="friends"></a> Friends</h2></td></tr> <tr class="memitem:a1db4d59f32248eaf6a85abc51ba5ca6d"><td class="memTemplParams" colspan="2">template<typename OtherDerived > </td></tr> <tr class="memitem:a1db4d59f32248eaf6a85abc51ba5ca6d"><td class="memTemplItemLeft" align="right" valign="top">SelfadjointProductMatrix<br class="typebreak"/> < OtherDerived, <br class="typebreak"/> 0, OtherDerived::IsVectorAtCompileTime, <br class="typebreak"/> MatrixType, Mode, false > </td><td class="memTemplItemRight" valign="bottom"><a class="el" href="classEigen_1_1SelfAdjointView.html#a1db4d59f32248eaf6a85abc51ba5ca6d">operator*</a> (const <a class="el" href="classEigen_1_1MatrixBase.html">MatrixBase</a>< OtherDerived > &lhs, const <a class="el" href="classEigen_1_1SelfAdjointView.html">SelfAdjointView</a> &rhs)</td></tr> <tr class="separator:a1db4d59f32248eaf6a85abc51ba5ca6d"><td class="memSeparator" colspan="2"> </td></tr> </table> <h2 class="groupheader">Member Typedef Documentation</h2> <a class="anchor" id="af7cc63140d17b562998188356f30da6f"></a> <div class="memitem"> <div class="memproto"> <table class="memname"> <tr> <td class="memname">typedef <a class="el" href="classEigen_1_1Matrix.html">Matrix</a><<a class="el" href="classEigen_1_1SelfAdjointView.html#acb5c3dc237f99cf17167e8a629f01b43">RealScalar</a>, internal::traits<MatrixType>::ColsAtCompileTime, 1> <a class="el" href="classEigen_1_1SelfAdjointView.html#af7cc63140d17b562998188356f30da6f">EigenvaluesReturnType</a></td> </tr> </table> </div><div class="memdoc"> <p>Return type of <a class="el" href="classEigen_1_1SelfAdjointView.html#a847d0db612cd20e1b796c6a70d0f0d54" title="Computes the eigenvalues of a matrix. ">eigenvalues()</a> </p> </div> </div> <a class="anchor" id="acb5c3dc237f99cf17167e8a629f01b43"></a> <div class="memitem"> <div class="memproto"> <table class="memname"> <tr> <td class="memname">typedef <a class="el" href="structEigen_1_1NumTraits.html">NumTraits</a><<a class="el" href="classEigen_1_1SelfAdjointView.html#a9ef555098df604af49b0a00d389d823b">Scalar</a>>::Real <a class="el" href="classEigen_1_1SelfAdjointView.html#acb5c3dc237f99cf17167e8a629f01b43">RealScalar</a></td> </tr> </table> </div><div class="memdoc"> <p>Real part of <a class="el" href="classEigen_1_1SelfAdjointView.html#a9ef555098df604af49b0a00d389d823b" title="The type of coefficients in this matrix. ">Scalar</a> </p> </div> </div> <h2 class="groupheader">Member Function Documentation</h2> <a class="anchor" id="afacbe8d905dd8b00eb52e72c2e938aed"></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_1SelfAdjointView.html#a9ef555098df604af49b0a00d389d823b">Scalar</a> coeff </td> <td>(</td> <td class="paramtype">Index </td> <td class="paramname"><em>row</em>, </td> </tr> <tr> <td class="paramkey"></td> <td></td> <td class="paramtype">Index </td> <td class="paramname"><em>col</em> </td> </tr> <tr> <td></td> <td>)</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 see"><dt>See Also</dt><dd>MatrixBase::coeff() </dd></dl> <dl class="section warning"><dt>Warning</dt><dd>the coordinates must fit into the referenced triangular part </dd></dl> </div> </div> <a class="anchor" id="a4386da578a119f7b500d0b22cab2de0e"></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_1SelfAdjointView.html#a9ef555098df604af49b0a00d389d823b">Scalar</a>& coeffRef </td> <td>(</td> <td class="paramtype">Index </td> <td class="paramname"><em>row</em>, </td> </tr> <tr> <td class="paramkey"></td> <td></td> <td class="paramtype">Index </td> <td class="paramname"><em>col</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"> <dl class="section see"><dt>See Also</dt><dd>MatrixBase::coeffRef() </dd></dl> <dl class="section warning"><dt>Warning</dt><dd>the coordinates must fit into the referenced triangular part </dd></dl> </div> </div> <a class="anchor" id="af7177136c93ee65b7c41c1798e23292d"></a> <div class="memitem"> <div class="memproto"> <table class="mlabels"> <tr> <td class="mlabels-left"> <table class="memname"> <tr> <td class="memname">void copyCoeff </td> <td>(</td> <td class="paramtype">Index </td> <td class="paramname"><em>row</em>, </td> </tr> <tr> <td class="paramkey"></td> <td></td> <td class="paramtype">Index </td> <td class="paramname"><em>col</em>, </td> </tr> <tr> <td class="paramkey"></td> <td></td> <td class="paramtype">Other & </td> <td class="paramname"><em>other</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 class="mlabel">inherited</span></span> </td> </tr> </table> </div><div class="memdoc"> <dl class="section see"><dt>See Also</dt><dd>MatrixBase::copyCoeff(row,col) </dd></dl> </div> </div> <a class="anchor" id="aa84222add803ad7c9db07dd4dd91d5d9"></a> <div class="memitem"> <div class="memproto"> <table class="mlabels"> <tr> <td class="mlabels-left"> <table class="memname"> <tr> <td class="memname">Derived& derived </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 class="mlabel">inherited</span></span> </td> </tr> </table> </div><div class="memdoc"> <dl class="section return"><dt>Returns</dt><dd>a reference to the derived object </dd></dl> <p>Referenced by <a class="el" href="classEigen_1_1MatrixBase.html#a9aaee4fba54baa7d0ca508e094f2d6f9">MatrixBase< Derived >::applyOnTheLeft()</a>, <a class="el" href="classEigen_1_1MatrixBase.html#a6e5b1e71a1916ddd01e7c324831c40dc">MatrixBase< Derived >::applyOnTheRight()</a>, <a class="el" href="classEigen_1_1PermutationBase.html#a975548d0fd92d4b2a2cb09c2adf10ce6">PermutationBase< PermutationMatrix< SizeAtCompileTime, MaxSizeAtCompileTime, IndexType > >::applyTranspositionOnTheLeft()</a>, <a class="el" href="classEigen_1_1PermutationBase.html#a3f93571373e20799e196513109f5b26e">PermutationBase< PermutationMatrix< SizeAtCompileTime, MaxSizeAtCompileTime, IndexType > >::applyTranspositionOnTheRight()</a>, <a class="el" href="classEigen_1_1SparseMatrixBase.html#ae762210a1ab89bed1673edc09a6ec8b7">SparseMatrixBase< CwiseBinaryOp< BinaryOp, Lhs, Rhs > >::binaryExpr()</a>, <a class="el" href="classEigen_1_1SparseMatrixBase.html#a878151caf2ff25ebbc76186485f81207">SparseMatrixBase< CwiseBinaryOp< BinaryOp, Lhs, Rhs > >::bottomLeftCorner()</a>, <a class="el" href="structEigen_1_1EigenBase.html#aaca1908a5ec508a25ff0a8bca803e5f3">EigenBase< SparseSymmetricPermutationProduct< MatrixType, UpLo > >::cols()</a>, <a class="el" href="classEigen_1_1SparseMatrixBase.html#aaca1908a5ec508a25ff0a8bca803e5f3">SparseMatrixBase< CwiseBinaryOp< BinaryOp, Lhs, Rhs > >::cols()</a>, <a class="el" href="classEigen_1_1SparseMatrixBase.html#a571c36d7211bd03dd1a06ac46c002fe5">SparseMatrixBase< CwiseBinaryOp< BinaryOp, Lhs, Rhs > >::eval()</a>, <a class="el" href="classEigen_1_1PermutationBase.html#a8876d615d17aad77b054a8f58b699e7d">PermutationBase< PermutationMatrix< SizeAtCompileTime, MaxSizeAtCompileTime, IndexType > >::indices()</a>, <a class="el" href="classEigen_1_1PermutationBase.html#a594deabd0f368a746801f5dd14a0db2a">PermutationBase< PermutationMatrix< SizeAtCompileTime, MaxSizeAtCompileTime, IndexType > >::inverse()</a>, <a class="el" href="classEigen_1_1RotationBase.html#acd4f0475d084e54c0dde2ef3f78822c0">RotationBase< Derived, 3 >::operator*()</a>, <a class="el" href="classEigen_1_1SparseSelfAdjointView.html#a6becf66c0294e516cbe73141c4276452">SparseSelfAdjointView< MatrixType, UpLo >::operator*()</a>, <a class="el" href="classEigen_1_1PermutationBase.html#aa06f5ef5169ce3e671abbc8776d2c339">PermutationBase< PermutationMatrix< SizeAtCompileTime, MaxSizeAtCompileTime, IndexType > >::operator*()</a>, <a class="el" href="classEigen_1_1SparseMatrixBase.html#ae7d7083917427445e50f3c4b5bf8ba6a">SparseMatrixBase< Derived >::operator*()</a>, <a class="el" href="classEigen_1_1Transform.html#a3846f7a4b892ed530881f52a5460ddf6">Transform< Scalar, Dim, Mode, _Options >::operator*()</a>, <a class="el" href="namespaceEigen.html#a81fb70d0dc1c6deb42e9816647607247">Eigen::operator*()</a>, <a class="el" href="classEigen_1_1MatrixBase.html#ac9beaa3954c989bbbc59b361d9423e78">MatrixBase< Derived >::operator*=()</a>, <a class="el" href="classEigen_1_1PermutationBase.html#ac96d811c408d8e282ad0036b1b9ed42f">PermutationBase< PermutationMatrix< SizeAtCompileTime, MaxSizeAtCompileTime, IndexType > >::operator=()</a>, <a class="el" href="classEigen_1_1DenseBase.html#acde54311cdc2b5731a6d27f7a5f7ad1b">DenseBase< Derived >::operator=()</a>, <a class="el" href="classEigen_1_1Transform.html#a746fdf41f0377f203f231a1705cf6012">Transform< Scalar, Dim, Mode, _Options >::operator=()</a>, <a class="el" href="classEigen_1_1PlainObjectBase.html#a65a58af12990c2e395a7f0b4012d3b94">PlainObjectBase< Matrix< _Scalar, _Rows, _Cols, _Options, _MaxRows, _MaxCols > >::operator=()</a>, <a class="el" href="classEigen_1_1PlainObjectBase.html#ab5207d4c7f9d20e5fc74bfc9b54620b4">PlainObjectBase< Matrix< _Scalar, _Rows, _Cols, _Options, _MaxRows, _MaxCols > >::PlainObjectBase()</a>, <a class="el" href="classEigen_1_1PlainObjectBase.html#abe7b7417203825a6e434449cf0ac6529">PlainObjectBase< Matrix< _Scalar, _Rows, _Cols, _Options, _MaxRows, _MaxCols > >::resizeLike()</a>, <a class="el" href="structEigen_1_1EigenBase.html#a5552abd83dbd03c85cea6d61fd8875a5">EigenBase< SparseSymmetricPermutationProduct< MatrixType, UpLo > >::rows()</a>, <a class="el" href="classEigen_1_1SparseMatrixBase.html#a5552abd83dbd03c85cea6d61fd8875a5">SparseMatrixBase< CwiseBinaryOp< BinaryOp, Lhs, Rhs > >::rows()</a>, <a class="el" href="classEigen_1_1SimplicialCholeskyBase.html#abe2318376d9cb41fa612b793035da397">SimplicialCholeskyBase< SimplicialLDLT< _MatrixType, _UpLo > >::solve()</a>, <a class="el" href="classEigen_1_1IterativeSolverBase.html#a7892f6fbc52bfa658e1324a63d376d1c">IterativeSolverBase< ConjugateGradient< _MatrixType, _UpLo, _Preconditioner > >::solve()</a>, <a class="el" href="classEigen_1_1SparseLU.html#a39062d90ea647b4e7668452647b04b17">SparseLU< _MatrixType, _OrderingType >::solve()</a>, <a class="el" href="classEigen_1_1UmfPackLU.html#a6c606e11055ed3d815183705c57b22b4">UmfPackLU< _MatrixType >::solve()</a>, <a class="el" href="classEigen_1_1CholmodBase.html#a81db76105fedf365cdbc200774911847">CholmodBase< _MatrixType, _UpLo, CholmodSimplicialLLT< _MatrixType, _UpLo > >::solve()</a>, <a class="el" href="classEigen_1_1SuperLUBase.html#a0bf5dc2154e2f03bf871eca997c6e6e6">SuperLUBase< _MatrixType, SuperILU< _MatrixType > >::solve()</a>, <a class="el" href="classEigen_1_1SparseMatrix.html#ab42390876ea7ef3fe4845b936c281c99">SparseMatrix< Scalar, RowMajor >::SparseMatrix()</a>, <a class="el" href="classEigen_1_1PermutationBase.html#ac4bd143bb38762fd91e546dbda01069b">PermutationBase< PermutationMatrix< SizeAtCompileTime, MaxSizeAtCompileTime, IndexType > >::toDenseMatrix()</a>, <a class="el" href="classEigen_1_1SparseMatrixBase.html#ae06d0a5c4008014fb717866aec8d30c7">SparseMatrixBase< CwiseBinaryOp< BinaryOp, Lhs, Rhs > >::topLeftCorner()</a>, <a class="el" href="classEigen_1_1Transform.html#ae80ecf3abb4b22d409c5adbd19850aa8">Transform< Scalar, Dim, Mode, _Options >::Transform()</a>, <a class="el" href="classEigen_1_1PermutationBase.html#a71ed95486aeea9fb396e85461ad5b73a">PermutationBase< PermutationMatrix< SizeAtCompileTime, MaxSizeAtCompileTime, IndexType > >::transpose()</a>, and <a class="el" href="classEigen_1_1SparseMatrixBase.html#aba9ae3923f6caa962ef3418d6872c369">SparseMatrixBase< CwiseBinaryOp< BinaryOp, Lhs, Rhs > >::unaryViewExpr()</a>.</p> </div> </div> <a class="anchor" id="a2850de52161dfe2bb6e380f1892276b8"></a> <div class="memitem"> <div class="memproto"> <table class="mlabels"> <tr> <td class="mlabels-left"> <table class="memname"> <tr> <td class="memname">const Derived& derived </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 class="mlabel">inherited</span></span> </td> </tr> </table> </div><div class="memdoc"> <dl class="section return"><dt>Returns</dt><dd>a const reference to the derived object </dd></dl> </div> </div> <a class="anchor" id="a847d0db612cd20e1b796c6a70d0f0d54"></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_1SelfAdjointView.html">SelfAdjointView</a>< MatrixType, UpLo >::<a class="el" href="classEigen_1_1SelfAdjointView.html#af7cc63140d17b562998188356f30da6f">EigenvaluesReturnType</a> eigenvalues </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>Computes the eigenvalues of a matrix. </p> <dl class="section return"><dt>Returns</dt><dd>Column vector containing the eigenvalues.</dd></dl> <p>This is defined in the Eigenvalues module.</p> <div class="fragment"><div class="line"><span class="preprocessor">#include <Eigen/Eigenvalues></span> </div> </div><!-- fragment --><p> This function computes the eigenvalues with the help of the <a class="el" href="classEigen_1_1SelfAdjointEigenSolver.html" title="Computes eigenvalues and eigenvectors of selfadjoint matrices. ">SelfAdjointEigenSolver</a> class. The eigenvalues are repeated according to their algebraic multiplicity, so there are as many eigenvalues as rows in the matrix.</p> <p>Example: </p> <div class="fragment"><div class="line">MatrixXd ones = <a class="code" href="classEigen_1_1DenseBase.html#a2278addf9a3c977d40322571a0df8ac9">MatrixXd::Ones</a>(3,3);</div> <div class="line"><a class="code" href="group__matrixtypedefs.html#ga3da45e59796fbacf67fa568297927bd1">VectorXd</a> eivals = ones.selfadjointView<<a class="code" href="group__enums.html#ggab59c1bec446b10af208f977a871d910baf886b397626076218462d53d50eb96bc">Lower</a>>().<a class="code" href="classEigen_1_1SelfAdjointView.html#a847d0db612cd20e1b796c6a70d0f0d54">eigenvalues</a>();</div> <div class="line">cout << <span class="stringliteral">"The eigenvalues of the 3x3 matrix of ones are:"</span> << endl << eivals << endl;</div> </div><!-- fragment --><p> Output: </p> <pre class="fragment">The eigenvalues of the 3x3 matrix of ones are: -3.09e-16 0 3 </pre><dl class="section see"><dt>See Also</dt><dd><a class="el" href="classEigen_1_1SelfAdjointEigenSolver.html#af54b25fe7d2a3f578269381e9e5592a2" title="Returns the eigenvalues of given matrix. ">SelfAdjointEigenSolver::eigenvalues()</a>, <a class="el" href="classEigen_1_1MatrixBase.html#a0ffa061371b1bd1b9f14ecef94b4502e" title="Computes the eigenvalues of a matrix. ">MatrixBase::eigenvalues()</a> </dd></dl> </div> </div> <a class="anchor" id="a66ac5c9b6c04321ab12364bba06fee66"></a> <div class="memitem"> <div class="memproto"> <table class="mlabels"> <tr> <td class="mlabels-left"> <table class="memname"> <tr> <td class="memname">void evalTo </td> <td>(</td> <td class="paramtype"><a class="el" href="classEigen_1_1MatrixBase.html">MatrixBase</a>< DenseDerived > & </td> <td class="paramname"><em>other</em></td><td>)</td> <td> const</td> </tr> </table> </td> <td class="mlabels-right"> <span class="mlabels"><span class="mlabel">inherited</span></span> </td> </tr> </table> </div><div class="memdoc"> <p>Assigns a triangular or selfadjoint matrix to a dense matrix. If the matrix is triangular, the opposite part is set to zero. </p> </div> </div> <a class="anchor" id="affcc5fe7ffb333e0b99a5532d086cd2f"></a> <div class="memitem"> <div class="memproto"> <table class="mlabels"> <tr> <td class="mlabels-left"> <table class="memname"> <tr> <td class="memname">void evalToLazy </td> <td>(</td> <td class="paramtype"><a class="el" href="classEigen_1_1MatrixBase.html">MatrixBase</a>< DenseDerived > & </td> <td class="paramname"><em>other</em></td><td>)</td> <td> const</td> </tr> </table> </td> <td class="mlabels-right"> <span class="mlabels"><span class="mlabel">inherited</span></span> </td> </tr> </table> </div><div class="memdoc"> <p>Assigns a triangular or selfadjoint matrix to a dense matrix. If the matrix is triangular, the opposite part is set to zero. </p> </div> </div> <a class="anchor" id="a245df21a519cd34f5eb6f724cf5bf8e7"></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_1LDLT.html">LDLT</a>< typename <a class="el" href="classEigen_1_1SelfAdjointView.html">SelfAdjointView</a>< MatrixType, UpLo >::PlainObject, UpLo > ldlt </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>This is defined in the Cholesky module.</p> <div class="fragment"><div class="line"><span class="preprocessor">#include <Eigen/Cholesky></span> </div> </div><!-- fragment --> <dl class="section return"><dt>Returns</dt><dd>the Cholesky decomposition with full pivoting without square root of <code>*this</code> </dd></dl> </div> </div> <a class="anchor" id="ab6d4081f604fb2bc15ad285cf6bf54a2"></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_1LLT.html">LLT</a>< typename <a class="el" href="classEigen_1_1SelfAdjointView.html">SelfAdjointView</a>< MatrixType, UpLo >::PlainObject, UpLo > llt </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>This is defined in the Cholesky module.</p> <div class="fragment"><div class="line"><span class="preprocessor">#include <Eigen/Cholesky></span> </div> </div><!-- fragment --> <dl class="section return"><dt>Returns</dt><dd>the <a class="el" href="classEigen_1_1LLT.html" title="Standard Cholesky decomposition (LL^T) of a matrix and associated features. ">LLT</a> decomposition of <code>*this</code> </dd></dl> </div> </div> <a class="anchor" id="a8b615fb312e03e69a42aa86f50c448f1"></a> <div class="memitem"> <div class="memproto"> <table class="mlabels"> <tr> <td class="mlabels-left"> <table class="memname"> <tr> <td class="memname">SelfadjointProductMatrix<MatrixType,Mode,false,OtherDerived,0,OtherDerived::IsVectorAtCompileTime> operator* </td> <td>(</td> <td class="paramtype">const <a class="el" href="classEigen_1_1MatrixBase.html">MatrixBase</a>< OtherDerived > & </td> <td class="paramname"><em>rhs</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>Efficient self-adjoint matrix times vector/matrix product </p> </div> </div> <a class="anchor" id="a67279bc5b39fba3763acabf1d52ecad8"></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_1SelfAdjointView.html">SelfAdjointView</a>< MatrixType, UpLo >::<a class="el" href="classEigen_1_1SelfAdjointView.html#acb5c3dc237f99cf17167e8a629f01b43">RealScalar</a> operatorNorm </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>Computes the L2 operator norm. </p> <dl class="section return"><dt>Returns</dt><dd>Operator norm of the matrix.</dd></dl> <p>This is defined in the Eigenvalues module.</p> <div class="fragment"><div class="line"><span class="preprocessor">#include <Eigen/Eigenvalues></span> </div> </div><!-- fragment --><p> This function computes the L2 operator norm of a self-adjoint matrix. For a self-adjoint matrix, the operator norm is the largest eigenvalue.</p> <p>The current implementation uses the eigenvalues of the matrix, as computed by <a class="el" href="classEigen_1_1SelfAdjointView.html#a847d0db612cd20e1b796c6a70d0f0d54" title="Computes the eigenvalues of a matrix. ">eigenvalues()</a>, to compute the operator norm of the matrix.</p> <p>Example: </p> <div class="fragment"><div class="line">MatrixXd ones = <a class="code" href="classEigen_1_1DenseBase.html#a2278addf9a3c977d40322571a0df8ac9">MatrixXd::Ones</a>(3,3);</div> <div class="line">cout << <span class="stringliteral">"The operator norm of the 3x3 matrix of ones is "</span></div> <div class="line"> << ones.selfadjointView<<a class="code" href="group__enums.html#ggab59c1bec446b10af208f977a871d910baf886b397626076218462d53d50eb96bc">Lower</a>>().<a class="code" href="classEigen_1_1SelfAdjointView.html#a67279bc5b39fba3763acabf1d52ecad8">operatorNorm</a>() << endl;</div> </div><!-- fragment --><p> Output: </p> <pre class="fragment">The operator norm of the 3x3 matrix of ones is 3 </pre><dl class="section see"><dt>See Also</dt><dd><a class="el" href="classEigen_1_1SelfAdjointView.html#a847d0db612cd20e1b796c6a70d0f0d54" title="Computes the eigenvalues of a matrix. ">eigenvalues()</a>, <a class="el" href="classEigen_1_1MatrixBase.html#a0650a523c77a498a88b4998809d0bd14" title="Computes the L2 operator norm. ">MatrixBase::operatorNorm()</a> </dd></dl> </div> </div> <a class="anchor" id="a024a26e49b7865e23e39479de6dc56fc"></a> <div class="memitem"> <div class="memproto"> <table class="memname"> <tr> <td class="memname"><a class="el" href="classEigen_1_1SelfAdjointView.html">SelfAdjointView</a>& rankUpdate </td> <td>(</td> <td class="paramtype">const <a class="el" href="classEigen_1_1MatrixBase.html">MatrixBase</a>< DerivedU > & </td> <td class="paramname"><em>u</em>, </td> </tr> <tr> <td class="paramkey"></td> <td></td> <td class="paramtype">const <a class="el" href="classEigen_1_1MatrixBase.html">MatrixBase</a>< DerivedV > & </td> <td class="paramname"><em>v</em>, </td> </tr> <tr> <td class="paramkey"></td> <td></td> <td class="paramtype">const <a class="el" href="classEigen_1_1SelfAdjointView.html#a9ef555098df604af49b0a00d389d823b">Scalar</a> & </td> <td class="paramname"><em>alpha</em> = <code><a class="el" href="classEigen_1_1SelfAdjointView.html#a9ef555098df604af49b0a00d389d823b">Scalar</a>(1)</code> </td> </tr> <tr> <td></td> <td>)</td> <td></td><td></td> </tr> </table> </div><div class="memdoc"> <p>Perform a symmetric rank 2 update of the selfadjoint matrix <code>*this</code>: <img class="formulaInl" alt="$ this = this + \alpha u v^* + conj(\alpha) v u^* $" src="form_27.png"/> </p> <dl class="section return"><dt>Returns</dt><dd>a reference to <code>*this</code> </dd></dl> <p>The vectors <em>u</em> and <code>v</code> <b>must</b> be column vectors, however they can be a adjoint expression without any overhead. Only the meaningful triangular part of the matrix is updated, the rest is left unchanged.</p> <dl class="section see"><dt>See Also</dt><dd>rankUpdate(const MatrixBase<DerivedU>&, Scalar) </dd></dl> </div> </div> <a class="anchor" id="aad06bacc2e9b8345ae92279acc362584"></a> <div class="memitem"> <div class="memproto"> <table class="memname"> <tr> <td class="memname"><a class="el" href="classEigen_1_1SelfAdjointView.html">SelfAdjointView</a>& rankUpdate </td> <td>(</td> <td class="paramtype">const <a class="el" href="classEigen_1_1MatrixBase.html">MatrixBase</a>< DerivedU > & </td> <td class="paramname"><em>u</em>, </td> </tr> <tr> <td class="paramkey"></td> <td></td> <td class="paramtype">const <a class="el" href="classEigen_1_1SelfAdjointView.html#a9ef555098df604af49b0a00d389d823b">Scalar</a> & </td> <td class="paramname"><em>alpha</em> = <code><a class="el" href="classEigen_1_1SelfAdjointView.html#a9ef555098df604af49b0a00d389d823b">Scalar</a>(1)</code> </td> </tr> <tr> <td></td> <td>)</td> <td></td><td></td> </tr> </table> </div><div class="memdoc"> <p>Perform a symmetric rank K update of the selfadjoint matrix <code>*this</code>: <img class="formulaInl" alt="$ this = this + \alpha ( u u^* ) $" src="form_28.png"/> where <em>u</em> is a vector or matrix.</p> <dl class="section return"><dt>Returns</dt><dd>a reference to <code>*this</code> </dd></dl> <p>Note that to perform <img class="formulaInl" alt="$ this = this + \alpha ( u^* u ) $" src="form_29.png"/> you can simply call this function with u.adjoint().</p> <dl class="section see"><dt>See Also</dt><dd>rankUpdate(const MatrixBase<DerivedU>&, const MatrixBase<DerivedV>&, Scalar) </dd></dl> </div> </div> <a class="anchor" id="a76f5bc8a03ec105ab4be1f2b91b7b5d5"></a> <div class="memitem"> <div class="memproto"> <table class="mlabels"> <tr> <td class="mlabels-left"> <table class="memname"> <tr> <td class="memname">Index size </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 class="mlabel">inherited</span></span> </td> </tr> </table> </div><div class="memdoc"> <dl class="section return"><dt>Returns</dt><dd>the number of coefficients, which is <a class="el" href="structEigen_1_1EigenBase.html#a5552abd83dbd03c85cea6d61fd8875a5">rows()</a>*cols(). </dd></dl> <dl class="section see"><dt>See Also</dt><dd><a class="el" href="structEigen_1_1EigenBase.html#a5552abd83dbd03c85cea6d61fd8875a5">rows()</a>, <a class="el" href="structEigen_1_1EigenBase.html#aaca1908a5ec508a25ff0a8bca803e5f3">cols()</a>, SizeAtCompileTime. </dd></dl> </div> </div> <h2 class="groupheader">Friends And Related Function Documentation</h2> <a class="anchor" id="a1db4d59f32248eaf6a85abc51ba5ca6d"></a> <div class="memitem"> <div class="memproto"> <table class="mlabels"> <tr> <td class="mlabels-left"> <table class="memname"> <tr> <td class="memname">SelfadjointProductMatrix<OtherDerived,0,OtherDerived::IsVectorAtCompileTime,MatrixType,Mode,false> operator* </td> <td>(</td> <td class="paramtype">const <a class="el" href="classEigen_1_1MatrixBase.html">MatrixBase</a>< OtherDerived > & </td> <td class="paramname"><em>lhs</em>, </td> </tr> <tr> <td class="paramkey"></td> <td></td> <td class="paramtype">const <a class="el" href="classEigen_1_1SelfAdjointView.html">SelfAdjointView</a>< MatrixType, UpLo > & </td> <td class="paramname"><em>rhs</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">friend</span></span> </td> </tr> </table> </div><div class="memdoc"> <p>Efficient vector/matrix times self-adjoint matrix product </p> </div> </div> <hr/>The documentation for this class was generated from the following files:<ul> <li><a class="el" href="SelfAdjointView_8h_source.html">SelfAdjointView.h</a></li> <li><a class="el" href="LDLT_8h_source.html">LDLT.h</a></li> <li><a class="el" href="LLT_8h_source.html">LLT.h</a></li> <li><a class="el" href="SelfadjointProduct_8h_source.html">SelfadjointProduct.h</a></li> <li><a class="el" href="SelfadjointRank2Update_8h_source.html">SelfadjointRank2Update.h</a></li> <li><a class="el" href="MatrixBaseEigenvalues_8h_source.html">MatrixBaseEigenvalues.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_1SelfAdjointView.html">SelfAdjointView</a></li> <li class="footer">Generated on Mon Oct 28 2013 11:04:30 for Eigen by <a href="http://www.doxygen.org/index.html"> <img class="footer" src="doxygen.png" alt="doxygen"/></a> 1.8.5 </li> </ul> </div> <!-- Piwik --> <!-- <script type="text/javascript"> var pkBaseURL = (("https:" == document.location.protocol) ? "https://stats.sylphide-consulting.com/piwik/" : "http://stats.sylphide-consulting.com/piwik/"); document.write(unescape("%3Cscript src='" + pkBaseURL + "piwik.js' type='text/javascript'%3E%3C/script%3E")); </script><script type="text/javascript"> try { var piwikTracker = Piwik.getTracker(pkBaseURL + "piwik.php", 20); piwikTracker.trackPageView(); piwikTracker.enableLinkTracking(); } catch( err ) {} </script><noscript><p><img src="http://stats.sylphide-consulting.com/piwik/piwik.php?idsite=20" style="border:0" alt="" /></p></noscript> --> <!-- End Piwik Tracking Code --> </body> </html>