<!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: Geometry module</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('group__Geometry__Module.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="#groups">Modules</a> | <a href="#nested-classes">Classes</a> | <a href="#typedef-members">Typedefs</a> | <a href="#func-members">Functions</a> </div> <div class="headertitle"> <div class="title">Geometry module<div class="ingroups"><a class="el" href="group__Geometry__Reference.html">Reference</a></div></div> </div> </div><!--header--> <div class="contents"> <a name="details" id="details"></a><h2 class="groupheader">Detailed Description</h2> <p>This module provides support for:</p> <ul> <li>fixed-size homogeneous transformations</li> <li>translation, scaling, 2D and 3D rotations</li> <li>quaternions</li> <li>cross product</li> <li><a class="el" href="classEigen_1_1MatrixBase.html#a42fdca0c5854b85d9482853c8a085dc1">orthognal vector generation</a></li> <li>some linear components: parametrized-lines and hyperplanes</li> </ul> <div class="fragment"><div class="line">* #include <Eigen/Geometry></div> <div class="line">* </div> </div><!-- fragment --> <table class="memberdecls"> <tr class="heading"><td colspan="2"><h2 class="groupheader"><a name="groups"></a> Modules</h2></td></tr> <tr class="memitem:group__alignedboxtypedefs"><td class="memItemLeft" align="right" valign="top"> </td><td class="memItemRight" valign="bottom"><a class="el" href="group__alignedboxtypedefs.html">Global aligned box typedefs</a></td></tr> <tr class="separator:"><td class="memSeparator" colspan="2"> </td></tr> </table><table class="memberdecls"> <tr class="heading"><td colspan="2"><h2 class="groupheader"><a name="nested-classes"></a> Classes</h2></td></tr> <tr class="memitem:"><td class="memItemLeft" align="right" valign="top">class  </td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1AlignedBox.html">AlignedBox< _Scalar, _AmbientDim ></a></td></tr> <tr class="memdesc:"><td class="mdescLeft"> </td><td class="mdescRight">An axis aligned box. <a href="classEigen_1_1AlignedBox.html#details">More...</a><br/></td></tr> <tr class="separator:"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:"><td class="memItemLeft" align="right" valign="top">class  </td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1AngleAxis.html">AngleAxis< Scalar ></a></td></tr> <tr class="memdesc:"><td class="mdescLeft"> </td><td class="mdescRight">Represents a 3D rotation as a rotation angle around an arbitrary 3D axis. <a href="classEigen_1_1AngleAxis.html#details">More...</a><br/></td></tr> <tr class="separator:"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:"><td class="memItemLeft" align="right" valign="top">class  </td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1Homogeneous.html">Homogeneous< MatrixType, Direction ></a></td></tr> <tr class="memdesc:"><td class="mdescLeft"> </td><td class="mdescRight">Expression of one (or a set of) homogeneous vector(s) <a href="classEigen_1_1Homogeneous.html#details">More...</a><br/></td></tr> <tr class="separator:"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:"><td class="memItemLeft" align="right" valign="top">class  </td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1Hyperplane.html">Hyperplane< _Scalar, _AmbientDim, Options ></a></td></tr> <tr class="memdesc:"><td class="mdescLeft"> </td><td class="mdescRight">A hyperplane. <a href="classEigen_1_1Hyperplane.html#details">More...</a><br/></td></tr> <tr class="separator:"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:"><td class="memItemLeft" align="right" valign="top">class  </td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1Map_3_01const_01Quaternion_3_01__Scalar_01_4_00_01__Options_01_4.html">Map< const Quaternion< _Scalar >, _Options ></a></td></tr> <tr class="memdesc:"><td class="mdescLeft"> </td><td class="mdescRight"><a class="el" href="classEigen_1_1Quaternion.html" title="The quaternion class used to represent 3D orientations and rotations. ">Quaternion</a> expression mapping a constant memory buffer. <a href="classEigen_1_1Map_3_01const_01Quaternion_3_01__Scalar_01_4_00_01__Options_01_4.html#details">More...</a><br/></td></tr> <tr class="separator:"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:"><td class="memItemLeft" align="right" valign="top">class  </td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1Map_3_01Quaternion_3_01__Scalar_01_4_00_01__Options_01_4.html">Map< Quaternion< _Scalar >, _Options ></a></td></tr> <tr class="memdesc:"><td class="mdescLeft"> </td><td class="mdescRight">Expression of a quaternion from a memory buffer. <a href="classEigen_1_1Map_3_01Quaternion_3_01__Scalar_01_4_00_01__Options_01_4.html#details">More...</a><br/></td></tr> <tr class="separator:"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:"><td class="memItemLeft" align="right" valign="top">class  </td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1ParametrizedLine.html">ParametrizedLine< _Scalar, _AmbientDim, Options ></a></td></tr> <tr class="memdesc:"><td class="mdescLeft"> </td><td class="mdescRight">A parametrized line. <a href="classEigen_1_1ParametrizedLine.html#details">More...</a><br/></td></tr> <tr class="separator:"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:"><td class="memItemLeft" align="right" valign="top">class  </td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1Quaternion.html">Quaternion< Scalar, Options ></a></td></tr> <tr class="memdesc:"><td class="mdescLeft"> </td><td class="mdescRight">The quaternion class used to represent 3D orientations and rotations. <a href="classEigen_1_1Quaternion.html#details">More...</a><br/></td></tr> <tr class="separator:"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:"><td class="memItemLeft" align="right" valign="top">class  </td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1QuaternionBase.html">QuaternionBase< Derived ></a></td></tr> <tr class="memdesc:"><td class="mdescLeft"> </td><td class="mdescRight">Base class for quaternion expressions. <a href="classEigen_1_1QuaternionBase.html#details">More...</a><br/></td></tr> <tr class="separator:"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:"><td class="memItemLeft" align="right" valign="top">class  </td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1Rotation2D.html">Rotation2D< Scalar ></a></td></tr> <tr class="memdesc:"><td class="mdescLeft"> </td><td class="mdescRight">Represents a rotation/orientation in a 2 dimensional space. <a href="classEigen_1_1Rotation2D.html#details">More...</a><br/></td></tr> <tr class="separator:"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:"><td class="memItemLeft" align="right" valign="top">class  </td><td class="memItemRight" valign="bottom"><a class="el" href="classScaling.html">Scaling</a></td></tr> <tr class="memdesc:"><td class="mdescLeft"> </td><td class="mdescRight">Represents a generic uniform scaling transformation. <a href="classScaling.html#details">More...</a><br/></td></tr> <tr class="separator:"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:"><td class="memItemLeft" align="right" valign="top">class  </td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1Transform.html">Transform< Scalar, Dim, Mode, _Options ></a></td></tr> <tr class="memdesc:"><td class="mdescLeft"> </td><td class="mdescRight">Represents an homogeneous transformation in a N dimensional space. <a href="classEigen_1_1Transform.html#details">More...</a><br/></td></tr> <tr class="separator:"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:"><td class="memItemLeft" align="right" valign="top">class  </td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1Translation.html">Translation< Scalar, Dim ></a></td></tr> <tr class="memdesc:"><td class="mdescLeft"> </td><td class="mdescRight">Represents a translation transformation. <a href="classEigen_1_1Translation.html#details">More...</a><br/></td></tr> <tr class="separator:"><td class="memSeparator" colspan="2"> </td></tr> </table><table class="memberdecls"> <tr class="heading"><td colspan="2"><h2 class="groupheader"><a name="typedef-members"></a> Typedefs</h2></td></tr> <tr class="memitem:gaac62cb423a4af85793c656a9848282a0"><td class="memItemLeft" align="right" valign="top">typedef Transform< double, <br class="typebreak"/> 2, Affine > </td><td class="memItemRight" valign="bottom"><a class="el" href="group__Geometry__Module.html#gaac62cb423a4af85793c656a9848282a0">Affine2d</a></td></tr> <tr class="separator:gaac62cb423a4af85793c656a9848282a0"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:gaaacee014168e2286115b9acd93de9c77"><td class="memItemLeft" align="right" valign="top">typedef Transform< float, <br class="typebreak"/> 2, Affine > </td><td class="memItemRight" valign="bottom"><a class="el" href="group__Geometry__Module.html#gaaacee014168e2286115b9acd93de9c77">Affine2f</a></td></tr> <tr class="separator:gaaacee014168e2286115b9acd93de9c77"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:ga1cee3e45cad60050d3c94cf603e2ff45"><td class="memItemLeft" align="right" valign="top">typedef Transform< double, <br class="typebreak"/> 3, Affine > </td><td class="memItemRight" valign="bottom"><a class="el" href="group__Geometry__Module.html#ga1cee3e45cad60050d3c94cf603e2ff45">Affine3d</a></td></tr> <tr class="separator:ga1cee3e45cad60050d3c94cf603e2ff45"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:gaca1e57a15de08fe20c934c5600f19314"><td class="memItemLeft" align="right" valign="top">typedef Transform< float, <br class="typebreak"/> 3, Affine > </td><td class="memItemRight" valign="bottom"><a class="el" href="group__Geometry__Module.html#gaca1e57a15de08fe20c934c5600f19314">Affine3f</a></td></tr> <tr class="separator:gaca1e57a15de08fe20c934c5600f19314"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:gaa44476c8fb0e300b47e21111db7f8224"><td class="memItemLeft" align="right" valign="top">typedef Transform< double, <br class="typebreak"/> 2, AffineCompact > </td><td class="memItemRight" valign="bottom"><a class="el" href="group__Geometry__Module.html#gaa44476c8fb0e300b47e21111db7f8224">AffineCompact2d</a></td></tr> <tr class="separator:gaa44476c8fb0e300b47e21111db7f8224"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:ga2ffd4a2121dd59a3955735299e5e6364"><td class="memItemLeft" align="right" valign="top">typedef Transform< float, <br class="typebreak"/> 2, AffineCompact > </td><td class="memItemRight" valign="bottom"><a class="el" href="group__Geometry__Module.html#ga2ffd4a2121dd59a3955735299e5e6364">AffineCompact2f</a></td></tr> <tr class="separator:ga2ffd4a2121dd59a3955735299e5e6364"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:ga40d0ecae258b06c9633323da78162515"><td class="memItemLeft" align="right" valign="top">typedef Transform< double, <br class="typebreak"/> 3, AffineCompact > </td><td class="memItemRight" valign="bottom"><a class="el" href="group__Geometry__Module.html#ga40d0ecae258b06c9633323da78162515">AffineCompact3d</a></td></tr> <tr class="separator:ga40d0ecae258b06c9633323da78162515"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:ga1c999eef555fe940f147f3e8d18be6de"><td class="memItemLeft" align="right" valign="top">typedef Transform< float, <br class="typebreak"/> 3, AffineCompact > </td><td class="memItemRight" valign="bottom"><a class="el" href="group__Geometry__Module.html#ga1c999eef555fe940f147f3e8d18be6de">AffineCompact3f</a></td></tr> <tr class="separator:ga1c999eef555fe940f147f3e8d18be6de"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:gaae4796ca310c4b5b0d857333a3c99afc"><td class="memItemLeft" align="right" valign="top">typedef DiagonalMatrix< double, 2 > </td><td class="memItemRight" valign="bottom"><a class="el" href="group__Geometry__Module.html#gaae4796ca310c4b5b0d857333a3c99afc">AlignedScaling2d</a></td></tr> <tr class="separator:gaae4796ca310c4b5b0d857333a3c99afc"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:ga2d228eea1fcc6870e75d608ef2aa9565"><td class="memItemLeft" align="right" valign="top">typedef DiagonalMatrix< float, 2 > </td><td class="memItemRight" valign="bottom"><a class="el" href="group__Geometry__Module.html#ga2d228eea1fcc6870e75d608ef2aa9565">AlignedScaling2f</a></td></tr> <tr class="separator:ga2d228eea1fcc6870e75d608ef2aa9565"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:gabe609cfffef979a1b2ab553486b9f48b"><td class="memItemLeft" align="right" valign="top">typedef DiagonalMatrix< double, 3 > </td><td class="memItemRight" valign="bottom"><a class="el" href="group__Geometry__Module.html#gabe609cfffef979a1b2ab553486b9f48b">AlignedScaling3d</a></td></tr> <tr class="separator:gabe609cfffef979a1b2ab553486b9f48b"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:ga462a06cebdc958f79cf7b50df07461ff"><td class="memItemLeft" align="right" valign="top">typedef DiagonalMatrix< float, 3 > </td><td class="memItemRight" valign="bottom"><a class="el" href="group__Geometry__Module.html#ga462a06cebdc958f79cf7b50df07461ff">AlignedScaling3f</a></td></tr> <tr class="separator:ga462a06cebdc958f79cf7b50df07461ff"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:ga0db1cc067c51aaa6dedf5805ee0c53d7"><td class="memItemLeft" align="right" valign="top">typedef AngleAxis< double > </td><td class="memItemRight" valign="bottom"><a class="el" href="group__Geometry__Module.html#ga0db1cc067c51aaa6dedf5805ee0c53d7">AngleAxisd</a></td></tr> <tr class="separator:ga0db1cc067c51aaa6dedf5805ee0c53d7"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:ga811d6fdab2002723bc7a72f055ce8c1d"><td class="memItemLeft" align="right" valign="top">typedef AngleAxis< float > </td><td class="memItemRight" valign="bottom"><a class="el" href="group__Geometry__Module.html#ga811d6fdab2002723bc7a72f055ce8c1d">AngleAxisf</a></td></tr> <tr class="separator:ga811d6fdab2002723bc7a72f055ce8c1d"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:ga3802603a0c0017d3bc38e1c2483bc0c8"><td class="memItemLeft" align="right" valign="top">typedef Transform< double, <br class="typebreak"/> 2, Isometry > </td><td class="memItemRight" valign="bottom"><a class="el" href="group__Geometry__Module.html#ga3802603a0c0017d3bc38e1c2483bc0c8">Isometry2d</a></td></tr> <tr class="separator:ga3802603a0c0017d3bc38e1c2483bc0c8"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:ga4ee356e527fdc5aa9b048542ff6d4278"><td class="memItemLeft" align="right" valign="top">typedef Transform< float, <br class="typebreak"/> 2, Isometry > </td><td class="memItemRight" valign="bottom"><a class="el" href="group__Geometry__Module.html#ga4ee356e527fdc5aa9b048542ff6d4278">Isometry2f</a></td></tr> <tr class="separator:ga4ee356e527fdc5aa9b048542ff6d4278"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:gaf6a2d462b1cb12f96cb6baedc132a3f2"><td class="memItemLeft" align="right" valign="top">typedef Transform< double, <br class="typebreak"/> 3, Isometry > </td><td class="memItemRight" valign="bottom"><a class="el" href="group__Geometry__Module.html#gaf6a2d462b1cb12f96cb6baedc132a3f2">Isometry3d</a></td></tr> <tr class="separator:gaf6a2d462b1cb12f96cb6baedc132a3f2"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:gaaab0779d16c81454c038cee3324585f2"><td class="memItemLeft" align="right" valign="top">typedef Transform< float, <br class="typebreak"/> 3, Isometry > </td><td class="memItemRight" valign="bottom"><a class="el" href="group__Geometry__Module.html#gaaab0779d16c81454c038cee3324585f2">Isometry3f</a></td></tr> <tr class="separator:gaaab0779d16c81454c038cee3324585f2"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:gab2b63f5ecf343e58bee2eb482592bc49"><td class="memItemLeft" align="right" valign="top">typedef Transform< double, <br class="typebreak"/> 2, Projective > </td><td class="memItemRight" valign="bottom"><a class="el" href="group__Geometry__Module.html#gab2b63f5ecf343e58bee2eb482592bc49">Projective2d</a></td></tr> <tr class="separator:gab2b63f5ecf343e58bee2eb482592bc49"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:ga0af228e7443e5447de55bc7ba10b4d8e"><td class="memItemLeft" align="right" valign="top">typedef Transform< float, <br class="typebreak"/> 2, Projective > </td><td class="memItemRight" valign="bottom"><a class="el" href="group__Geometry__Module.html#ga0af228e7443e5447de55bc7ba10b4d8e">Projective2f</a></td></tr> <tr class="separator:ga0af228e7443e5447de55bc7ba10b4d8e"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:gad835f3f4d1783ee8b7ff85af90b9d13e"><td class="memItemLeft" align="right" valign="top">typedef Transform< double, <br class="typebreak"/> 3, Projective > </td><td class="memItemRight" valign="bottom"><a class="el" href="group__Geometry__Module.html#gad835f3f4d1783ee8b7ff85af90b9d13e">Projective3d</a></td></tr> <tr class="separator:gad835f3f4d1783ee8b7ff85af90b9d13e"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:gaaf2985f9f12b96ba19c07a5ecd47056b"><td class="memItemLeft" align="right" valign="top">typedef Transform< float, <br class="typebreak"/> 3, Projective > </td><td class="memItemRight" valign="bottom"><a class="el" href="group__Geometry__Module.html#gaaf2985f9f12b96ba19c07a5ecd47056b">Projective3f</a></td></tr> <tr class="separator:gaaf2985f9f12b96ba19c07a5ecd47056b"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:ga0d2bd45f1215359f8e7c0d7ab53c4acb"><td class="memItemLeft" align="right" valign="top">typedef Quaternion< double > </td><td class="memItemRight" valign="bottom"><a class="el" href="group__Geometry__Module.html#ga0d2bd45f1215359f8e7c0d7ab53c4acb">Quaterniond</a></td></tr> <tr class="separator:ga0d2bd45f1215359f8e7c0d7ab53c4acb"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:gaf65cf6f803890e57488d7de750bef682"><td class="memItemLeft" align="right" valign="top">typedef Quaternion< float > </td><td class="memItemRight" valign="bottom"><a class="el" href="group__Geometry__Module.html#gaf65cf6f803890e57488d7de750bef682">Quaternionf</a></td></tr> <tr class="separator:gaf65cf6f803890e57488d7de750bef682"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:gabfe832b0fc8c4c8b4dcf8b59ebe25603"><td class="memItemLeft" align="right" valign="top">typedef Map< Quaternion<br class="typebreak"/> < double >, Aligned > </td><td class="memItemRight" valign="bottom"><a class="el" href="group__Geometry__Module.html#gabfe832b0fc8c4c8b4dcf8b59ebe25603">QuaternionMapAlignedd</a></td></tr> <tr class="separator:gabfe832b0fc8c4c8b4dcf8b59ebe25603"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:ga987bcae344d68a98892d0e9404138404"><td class="memItemLeft" align="right" valign="top">typedef Map< Quaternion< float ><br class="typebreak"/> , Aligned > </td><td class="memItemRight" valign="bottom"><a class="el" href="group__Geometry__Module.html#ga987bcae344d68a98892d0e9404138404">QuaternionMapAlignedf</a></td></tr> <tr class="separator:ga987bcae344d68a98892d0e9404138404"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:gaca6aed0c662d8272c53663c0093aacaa"><td class="memItemLeft" align="right" valign="top">typedef Map< Quaternion<br class="typebreak"/> < double >, 0 > </td><td class="memItemRight" valign="bottom"><a class="el" href="group__Geometry__Module.html#gaca6aed0c662d8272c53663c0093aacaa">QuaternionMapd</a></td></tr> <tr class="separator:gaca6aed0c662d8272c53663c0093aacaa"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:ga8d125241c02e63c656b75b38b2382816"><td class="memItemLeft" align="right" valign="top">typedef Map< Quaternion< float >, 0 > </td><td class="memItemRight" valign="bottom"><a class="el" href="group__Geometry__Module.html#ga8d125241c02e63c656b75b38b2382816">QuaternionMapf</a></td></tr> <tr class="separator:ga8d125241c02e63c656b75b38b2382816"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:ga197b898c965d08135ebfb21bf41b23e2"><td class="memItemLeft" align="right" valign="top">typedef Rotation2D< double > </td><td class="memItemRight" valign="bottom"><a class="el" href="group__Geometry__Module.html#ga197b898c965d08135ebfb21bf41b23e2">Rotation2Dd</a></td></tr> <tr class="separator:ga197b898c965d08135ebfb21bf41b23e2"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:ga9a907f59280ce5650106a74904bab16d"><td class="memItemLeft" align="right" valign="top">typedef Rotation2D< float > </td><td class="memItemRight" valign="bottom"><a class="el" href="group__Geometry__Module.html#ga9a907f59280ce5650106a74904bab16d">Rotation2Df</a></td></tr> <tr class="separator:ga9a907f59280ce5650106a74904bab16d"><td class="memSeparator" colspan="2"> </td></tr> </table><table class="memberdecls"> <tr class="heading"><td colspan="2"><h2 class="groupheader"><a name="func-members"></a> Functions</h2></td></tr> <tr class="memitem:gad118fececd448d7485ffea4858775e5a"><td class="memItemLeft" align="right" valign="top">Matrix< Scalar, 3, 1 > </td><td class="memItemRight" valign="bottom"><a class="el" href="group__Geometry__Module.html#gad118fececd448d7485ffea4858775e5a">eulerAngles</a> (Index a0, Index a1, Index a2) const </td></tr> <tr class="separator:gad118fececd448d7485ffea4858775e5a"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:gab3f5a82a24490b936f8694cf8fef8e60"><td class="memTemplParams" colspan="2">template<typename Derived , typename OtherDerived > </td></tr> <tr class="memitem:gab3f5a82a24490b936f8694cf8fef8e60"><td class="memTemplItemLeft" align="right" valign="top">internal::umeyama_transform_matrix_type<br class="typebreak"/> < Derived, OtherDerived ><br class="typebreak"/> ::type </td><td class="memTemplItemRight" valign="bottom"><a class="el" href="group__Geometry__Module.html#gab3f5a82a24490b936f8694cf8fef8e60">umeyama</a> (const MatrixBase< Derived > &src, const MatrixBase< OtherDerived > &dst, bool with_scaling=true)</td></tr> <tr class="memdesc:gab3f5a82a24490b936f8694cf8fef8e60"><td class="mdescLeft"> </td><td class="mdescRight">Returns the transformation between two point sets. <a href="#gab3f5a82a24490b936f8694cf8fef8e60">More...</a><br/></td></tr> <tr class="separator:gab3f5a82a24490b936f8694cf8fef8e60"><td class="memSeparator" colspan="2"> </td></tr> </table> <h2 class="groupheader">Typedef Documentation</h2> <a class="anchor" id="gaac62cb423a4af85793c656a9848282a0"></a> <div class="memitem"> <div class="memproto"> <table class="memname"> <tr> <td class="memname">typedef Transform<double,2,Affine> Affine2d</td> </tr> </table> </div><div class="memdoc"> </div> </div> <a class="anchor" id="gaaacee014168e2286115b9acd93de9c77"></a> <div class="memitem"> <div class="memproto"> <table class="memname"> <tr> <td class="memname">typedef Transform<float,2,Affine> Affine2f</td> </tr> </table> </div><div class="memdoc"> </div> </div> <a class="anchor" id="ga1cee3e45cad60050d3c94cf603e2ff45"></a> <div class="memitem"> <div class="memproto"> <table class="memname"> <tr> <td class="memname">typedef Transform<double,3,Affine> Affine3d</td> </tr> </table> </div><div class="memdoc"> </div> </div> <a class="anchor" id="gaca1e57a15de08fe20c934c5600f19314"></a> <div class="memitem"> <div class="memproto"> <table class="memname"> <tr> <td class="memname">typedef Transform<float,3,Affine> Affine3f</td> </tr> </table> </div><div class="memdoc"> </div> </div> <a class="anchor" id="gaa44476c8fb0e300b47e21111db7f8224"></a> <div class="memitem"> <div class="memproto"> <table class="memname"> <tr> <td class="memname">typedef Transform<double,2,AffineCompact> AffineCompact2d</td> </tr> </table> </div><div class="memdoc"> </div> </div> <a class="anchor" id="ga2ffd4a2121dd59a3955735299e5e6364"></a> <div class="memitem"> <div class="memproto"> <table class="memname"> <tr> <td class="memname">typedef Transform<float,2,AffineCompact> AffineCompact2f</td> </tr> </table> </div><div class="memdoc"> </div> </div> <a class="anchor" id="ga40d0ecae258b06c9633323da78162515"></a> <div class="memitem"> <div class="memproto"> <table class="memname"> <tr> <td class="memname">typedef Transform<double,3,AffineCompact> AffineCompact3d</td> </tr> </table> </div><div class="memdoc"> </div> </div> <a class="anchor" id="ga1c999eef555fe940f147f3e8d18be6de"></a> <div class="memitem"> <div class="memproto"> <table class="memname"> <tr> <td class="memname">typedef Transform<float,3,AffineCompact> AffineCompact3f</td> </tr> </table> </div><div class="memdoc"> </div> </div> <a class="anchor" id="gaae4796ca310c4b5b0d857333a3c99afc"></a> <div class="memitem"> <div class="memproto"> <table class="memname"> <tr> <td class="memname">typedef DiagonalMatrix<double,2> AlignedScaling2d</td> </tr> </table> </div><div class="memdoc"> </div> </div> <a class="anchor" id="ga2d228eea1fcc6870e75d608ef2aa9565"></a> <div class="memitem"> <div class="memproto"> <table class="memname"> <tr> <td class="memname">typedef DiagonalMatrix<float, 2> AlignedScaling2f</td> </tr> </table> </div><div class="memdoc"> </div> </div> <a class="anchor" id="gabe609cfffef979a1b2ab553486b9f48b"></a> <div class="memitem"> <div class="memproto"> <table class="memname"> <tr> <td class="memname">typedef DiagonalMatrix<double,3> AlignedScaling3d</td> </tr> </table> </div><div class="memdoc"> </div> </div> <a class="anchor" id="ga462a06cebdc958f79cf7b50df07461ff"></a> <div class="memitem"> <div class="memproto"> <table class="memname"> <tr> <td class="memname">typedef DiagonalMatrix<float, 3> AlignedScaling3f</td> </tr> </table> </div><div class="memdoc"> </div> </div> <a class="anchor" id="ga0db1cc067c51aaa6dedf5805ee0c53d7"></a> <div class="memitem"> <div class="memproto"> <table class="memname"> <tr> <td class="memname">typedef AngleAxis<double> AngleAxisd</td> </tr> </table> </div><div class="memdoc"> <p>double precision angle-axis type </p> </div> </div> <a class="anchor" id="ga811d6fdab2002723bc7a72f055ce8c1d"></a> <div class="memitem"> <div class="memproto"> <table class="memname"> <tr> <td class="memname">typedef AngleAxis<float> AngleAxisf</td> </tr> </table> </div><div class="memdoc"> <p>single precision angle-axis type </p> </div> </div> <a class="anchor" id="ga3802603a0c0017d3bc38e1c2483bc0c8"></a> <div class="memitem"> <div class="memproto"> <table class="memname"> <tr> <td class="memname">typedef Transform<double,2,Isometry> Isometry2d</td> </tr> </table> </div><div class="memdoc"> </div> </div> <a class="anchor" id="ga4ee356e527fdc5aa9b048542ff6d4278"></a> <div class="memitem"> <div class="memproto"> <table class="memname"> <tr> <td class="memname">typedef Transform<float,2,Isometry> Isometry2f</td> </tr> </table> </div><div class="memdoc"> </div> </div> <a class="anchor" id="gaf6a2d462b1cb12f96cb6baedc132a3f2"></a> <div class="memitem"> <div class="memproto"> <table class="memname"> <tr> <td class="memname">typedef Transform<double,3,Isometry> Isometry3d</td> </tr> </table> </div><div class="memdoc"> </div> </div> <a class="anchor" id="gaaab0779d16c81454c038cee3324585f2"></a> <div class="memitem"> <div class="memproto"> <table class="memname"> <tr> <td class="memname">typedef Transform<float,3,Isometry> Isometry3f</td> </tr> </table> </div><div class="memdoc"> </div> </div> <a class="anchor" id="gab2b63f5ecf343e58bee2eb482592bc49"></a> <div class="memitem"> <div class="memproto"> <table class="memname"> <tr> <td class="memname">typedef Transform<double,2,Projective> Projective2d</td> </tr> </table> </div><div class="memdoc"> </div> </div> <a class="anchor" id="ga0af228e7443e5447de55bc7ba10b4d8e"></a> <div class="memitem"> <div class="memproto"> <table class="memname"> <tr> <td class="memname">typedef Transform<float,2,Projective> Projective2f</td> </tr> </table> </div><div class="memdoc"> </div> </div> <a class="anchor" id="gad835f3f4d1783ee8b7ff85af90b9d13e"></a> <div class="memitem"> <div class="memproto"> <table class="memname"> <tr> <td class="memname">typedef Transform<double,3,Projective> Projective3d</td> </tr> </table> </div><div class="memdoc"> </div> </div> <a class="anchor" id="gaaf2985f9f12b96ba19c07a5ecd47056b"></a> <div class="memitem"> <div class="memproto"> <table class="memname"> <tr> <td class="memname">typedef Transform<float,3,Projective> Projective3f</td> </tr> </table> </div><div class="memdoc"> </div> </div> <a class="anchor" id="ga0d2bd45f1215359f8e7c0d7ab53c4acb"></a> <div class="memitem"> <div class="memproto"> <table class="memname"> <tr> <td class="memname">typedef Quaternion<double> Quaterniond</td> </tr> </table> </div><div class="memdoc"> <p>double precision quaternion type </p> </div> </div> <a class="anchor" id="gaf65cf6f803890e57488d7de750bef682"></a> <div class="memitem"> <div class="memproto"> <table class="memname"> <tr> <td class="memname">typedef Quaternion<float> Quaternionf</td> </tr> </table> </div><div class="memdoc"> <p>single precision quaternion type </p> </div> </div> <a class="anchor" id="gabfe832b0fc8c4c8b4dcf8b59ebe25603"></a> <div class="memitem"> <div class="memproto"> <table class="memname"> <tr> <td class="memname">typedef Map<Quaternion<double>, Aligned> QuaternionMapAlignedd</td> </tr> </table> </div><div class="memdoc"> <p><a class="el" href="classEigen_1_1Map.html" title="A matrix or vector expression mapping an existing array of data. ">Map</a> a 16-bits aligned array of double precision scalars as a quaternion </p> </div> </div> <a class="anchor" id="ga987bcae344d68a98892d0e9404138404"></a> <div class="memitem"> <div class="memproto"> <table class="memname"> <tr> <td class="memname">typedef Map<Quaternion<float>, Aligned> QuaternionMapAlignedf</td> </tr> </table> </div><div class="memdoc"> <p><a class="el" href="classEigen_1_1Map.html" title="A matrix or vector expression mapping an existing array of data. ">Map</a> a 16-bits aligned array of double precision scalars as a quaternion </p> </div> </div> <a class="anchor" id="gaca6aed0c662d8272c53663c0093aacaa"></a> <div class="memitem"> <div class="memproto"> <table class="memname"> <tr> <td class="memname">typedef Map<Quaternion<double>, 0> QuaternionMapd</td> </tr> </table> </div><div class="memdoc"> <p><a class="el" href="classEigen_1_1Map.html" title="A matrix or vector expression mapping an existing array of data. ">Map</a> an unaligned array of double precision scalar as a quaternion </p> </div> </div> <a class="anchor" id="ga8d125241c02e63c656b75b38b2382816"></a> <div class="memitem"> <div class="memproto"> <table class="memname"> <tr> <td class="memname">typedef Map<Quaternion<float>, 0> QuaternionMapf</td> </tr> </table> </div><div class="memdoc"> <p><a class="el" href="classEigen_1_1Map.html" title="A matrix or vector expression mapping an existing array of data. ">Map</a> an unaligned array of single precision scalar as a quaternion </p> </div> </div> <a class="anchor" id="ga197b898c965d08135ebfb21bf41b23e2"></a> <div class="memitem"> <div class="memproto"> <table class="memname"> <tr> <td class="memname">typedef Rotation2D<double> Rotation2Dd</td> </tr> </table> </div><div class="memdoc"> <p>double precision 2D rotation type </p> </div> </div> <a class="anchor" id="ga9a907f59280ce5650106a74904bab16d"></a> <div class="memitem"> <div class="memproto"> <table class="memname"> <tr> <td class="memname">typedef Rotation2D<float> Rotation2Df</td> </tr> </table> </div><div class="memdoc"> <p>single precision 2D rotation type </p> </div> </div> <h2 class="groupheader">Function Documentation</h2> <a class="anchor" id="gad118fececd448d7485ffea4858775e5a"></a> <div class="memitem"> <div class="memproto"> <table class="mlabels"> <tr> <td class="mlabels-left"> <table class="memname"> <tr> <td class="memname">Matrix< typename MatrixBase< Derived >::Scalar, 3, 1 > eulerAngles </td> <td>(</td> <td class="paramtype"><a class="el" href="classEigen_1_1DenseBase.html#a4d4873e91be950c079f067fa97fd5c40">Index</a> </td> <td class="paramname"><em>a0</em>, </td> </tr> <tr> <td class="paramkey"></td> <td></td> <td class="paramtype"><a class="el" href="classEigen_1_1DenseBase.html#a4d4873e91be950c079f067fa97fd5c40">Index</a> </td> <td class="paramname"><em>a1</em>, </td> </tr> <tr> <td class="paramkey"></td> <td></td> <td class="paramtype"><a class="el" href="classEigen_1_1DenseBase.html#a4d4873e91be950c079f067fa97fd5c40">Index</a> </td> <td class="paramname"><em>a2</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"> <p>This is defined in the Geometry module.</p> <div class="fragment"><div class="line"><span class="preprocessor">#include <Eigen/Geometry></span> </div> </div><!-- fragment --><dl class="section return"><dt>Returns</dt><dd>the Euler-angles of the rotation matrix <code>*this</code> using the convention defined by the triplet (<em>a0</em>,<em>a1</em>,<em>a2</em>)</dd></dl> <p>Each of the three parameters <em>a0</em>,<em>a1</em>,<em>a2</em> represents the respective rotation axis as an integer in {0,1,2}. For instance, in: </p> <div class="fragment"><div class="line"><a class="code" href="group__matrixtypedefs.html#ga5ac9fb0df0c4858477890cce1f998096">Vector3f</a> ea = mat.eulerAngles(2, 0, 2); </div> </div><!-- fragment --><p> "2" represents the z axis and "0" the x axis, etc. The returned angles are such that we have the following equality: </p> <div class="fragment"><div class="line">* mat == <a class="code" href="group__Geometry__Module.html#ga811d6fdab2002723bc7a72f055ce8c1d">AngleAxisf</a>(ea[0], <a class="code" href="classEigen_1_1MatrixBase.html#a135f14351a7213bf0b764272c14ca68c">Vector3f::UnitZ</a>())</div> <div class="line">* * <a class="code" href="group__Geometry__Module.html#ga811d6fdab2002723bc7a72f055ce8c1d">AngleAxisf</a>(ea[1], <a class="code" href="classEigen_1_1MatrixBase.html#a5e0e42c39330e4274d2d0479048ebc37">Vector3f::UnitX</a>())</div> <div class="line">* * <a class="code" href="group__Geometry__Module.html#ga811d6fdab2002723bc7a72f055ce8c1d">AngleAxisf</a>(ea[2], <a class="code" href="classEigen_1_1MatrixBase.html#a135f14351a7213bf0b764272c14ca68c">Vector3f::UnitZ</a>()); </div> </div><!-- fragment --><p> This corresponds to the right-multiply conventions (with right hand side frames).</p> <p>The returned angles are in the ranges [0:pi]x[0:pi]x[-pi:pi].</p> <dl class="section see"><dt>See Also</dt><dd>class <a class="el" href="classEigen_1_1AngleAxis.html" title="Represents a 3D rotation as a rotation angle around an arbitrary 3D axis. ">AngleAxis</a> </dd></dl> </div> </div> <a class="anchor" id="gab3f5a82a24490b936f8694cf8fef8e60"></a> <div class="memitem"> <div class="memproto"> <table class="memname"> <tr> <td class="memname">internal::umeyama_transform_matrix_type<Derived, OtherDerived>::type Eigen::umeyama </td> <td>(</td> <td class="paramtype">const MatrixBase< Derived > & </td> <td class="paramname"><em>src</em>, </td> </tr> <tr> <td class="paramkey"></td> <td></td> <td class="paramtype">const MatrixBase< OtherDerived > & </td> <td class="paramname"><em>dst</em>, </td> </tr> <tr> <td class="paramkey"></td> <td></td> <td class="paramtype">bool </td> <td class="paramname"><em>with_scaling</em> = <code>true</code> </td> </tr> <tr> <td></td> <td>)</td> <td></td><td></td> </tr> </table> </div><div class="memdoc"> <p>Returns the transformation between two point sets. </p> <p>This is defined in the Geometry module.</p> <div class="fragment"><div class="line"><span class="preprocessor">#include <Eigen/Geometry></span> </div> </div><!-- fragment --><p>The algorithm is based on: "Least-squares estimation of transformation parameters between two point patterns", Shinji Umeyama, PAMI 1991, DOI: 10.1109/34.88573</p> <p>It estimates parameters <img class="formulaInl" alt="$ c, \mathbf{R}, $" src="form_116.png"/> and <img class="formulaInl" alt="$ \mathbf{t} $" src="form_117.png"/> such that </p> <p class="formulaDsp"> <img class="formulaDsp" alt="\begin{align*} \frac{1}{n} \sum_{i=1}^n \vert\vert y_i - (c\mathbf{R}x_i + \mathbf{t}) \vert\vert_2^2 \end{align*}" src="form_118.png"/> </p> <p> is minimized.</p> <p>The algorithm is based on the analysis of the covariance matrix <img class="formulaInl" alt="$ \Sigma_{\mathbf{x}\mathbf{y}} \in \mathbb{R}^{d \times d} $" src="form_119.png"/> of the input point sets <img class="formulaInl" alt="$ \mathbf{x} $" src="form_120.png"/> and <img class="formulaInl" alt="$ \mathbf{y} $" src="form_121.png"/> where <img class="formulaInl" alt="$d$" src="form_122.png"/> is corresponding to the dimension (which is typically small). The analysis is involving the SVD having a complexity of <img class="formulaInl" alt="$O(d^3)$" src="form_123.png"/> though the actual computational effort lies in the covariance matrix computation which has an asymptotic lower bound of <img class="formulaInl" alt="$O(dm)$" src="form_124.png"/> when the input point sets have dimension <img class="formulaInl" alt="$d \times m$" src="form_125.png"/>.</p> <p>Currently the method is working only for floating point matrices.</p> <dl class="params"><dt>Parameters</dt><dd> <table class="params"> <tr><td class="paramname">src</td><td>Source points <img class="formulaInl" alt="$ \mathbf{x} = \left( x_1, \hdots, x_n \right) $" src="form_126.png"/>. </td></tr> <tr><td class="paramname">dst</td><td>Destination points <img class="formulaInl" alt="$ \mathbf{y} = \left( y_1, \hdots, y_n \right) $" src="form_127.png"/>. </td></tr> <tr><td class="paramname">with_scaling</td><td>Sets <img class="formulaInl" alt="$ c=1 $" src="form_128.png"/> when <code>false</code> is passed. </td></tr> </table> </dd> </dl> <dl class="section return"><dt>Returns</dt><dd>The homogeneous transformation <p class="formulaDsp"> <img class="formulaDsp" alt="\begin{align*} T = \begin{bmatrix} c\mathbf{R} & \mathbf{t} \\ \mathbf{0} & 1 \end{bmatrix} \end{align*}" src="form_129.png"/> </p> minimizing the resudiual above. This transformation is always returned as an <a class="el" href="classEigen_1_1Matrix.html" title="The matrix class, also used for vectors and row-vectors. ">Eigen::Matrix</a>. </dd></dl> <p>References <a class="el" href="classEigen_1_1DenseBase.html#a49a617f24129ca31a27fe8a67ec20370">DenseBase< Derived >::colwise()</a>, <a class="el" href="group__enums.html#gga2d78499b99ddc29b9494f7ea33864d52a75953ad8f5837a79d6fbd3c116e7d462">Eigen::ComputeFullU</a>, <a class="el" href="group__enums.html#gga2d78499b99ddc29b9494f7ea33864d52a1785ac1174dab733556ac572448984c7">Eigen::ComputeFullV</a>, <a class="el" href="classEigen_1_1JacobiSVD.html#a49e16a4adf4fe58a2d65a5e5a31e7654">JacobiSVD< MatrixType, QRPreconditioner >::matrixU()</a>, <a class="el" href="classEigen_1_1JacobiSVD.html#ae5158ab7ca44a705c2a3b56ec926b42a">JacobiSVD< MatrixType, QRPreconditioner >::matrixV()</a>, <a class="el" href="classEigen_1_1DenseBase.html#a3af2f03b1d2affcec24e0748edf892cd">DenseBase< Derived >::rowwise()</a>, and <a class="el" href="classEigen_1_1JacobiSVD.html#a48d4068b97dfbb83d62599e56e26797a">JacobiSVD< MatrixType, QRPreconditioner >::singularValues()</a>.</p> </div> </div> </div><!-- contents --> </div><!-- doc-content --> <!-- start footer part --> <div id="nav-path" class="navpath"><!-- id is needed for treeview function! --> <ul> <li class="footer">Generated on Mon Oct 28 2013 11:04:27 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>