<!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: AngleAxis.h Source File</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('AngleAxis_8h_source.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="headertitle"> <div class="title">AngleAxis.h</div> </div> </div><!--header--> <div class="contents"> <div class="fragment"><div class="line"><a name="l00001"></a><span class="lineno"> 1</span> <span class="comment">// This file is part of Eigen, a lightweight C++ template library</span></div> <div class="line"><a name="l00002"></a><span class="lineno"> 2</span> <span class="comment">// for linear algebra.</span></div> <div class="line"><a name="l00003"></a><span class="lineno"> 3</span> <span class="comment">//</span></div> <div class="line"><a name="l00004"></a><span class="lineno"> 4</span> <span class="comment">// Copyright (C) 2008 Gael Guennebaud <gael.guennebaud@inria.fr></span></div> <div class="line"><a name="l00005"></a><span class="lineno"> 5</span> <span class="comment">//</span></div> <div class="line"><a name="l00006"></a><span class="lineno"> 6</span> <span class="comment">// This Source Code Form is subject to the terms of the Mozilla</span></div> <div class="line"><a name="l00007"></a><span class="lineno"> 7</span> <span class="comment">// Public License v. 2.0. If a copy of the MPL was not distributed</span></div> <div class="line"><a name="l00008"></a><span class="lineno"> 8</span> <span class="comment">// with this file, You can obtain one at http://mozilla.org/MPL/2.0/.</span></div> <div class="line"><a name="l00009"></a><span class="lineno"> 9</span> </div> <div class="line"><a name="l00010"></a><span class="lineno"> 10</span> <span class="preprocessor">#ifndef EIGEN_ANGLEAXIS_H</span></div> <div class="line"><a name="l00011"></a><span class="lineno"> 11</span> <span class="preprocessor"></span><span class="preprocessor">#define EIGEN_ANGLEAXIS_H</span></div> <div class="line"><a name="l00012"></a><span class="lineno"> 12</span> <span class="preprocessor"></span></div> <div class="line"><a name="l00013"></a><span class="lineno"> 13</span> <span class="keyword">namespace </span>Eigen { </div> <div class="line"><a name="l00014"></a><span class="lineno"> 14</span> </div> <div class="line"><a name="l00041"></a><span class="lineno"> 41</span> <span class="keyword">namespace </span>internal {</div> <div class="line"><a name="l00042"></a><span class="lineno"> 42</span> <span class="keyword">template</span><<span class="keyword">typename</span> _Scalar> <span class="keyword">struct </span>traits<AngleAxis<_Scalar> ></div> <div class="line"><a name="l00043"></a><span class="lineno"> 43</span> {</div> <div class="line"><a name="l00044"></a><span class="lineno"> 44</span>  <span class="keyword">typedef</span> _Scalar Scalar;</div> <div class="line"><a name="l00045"></a><span class="lineno"> 45</span> };</div> <div class="line"><a name="l00046"></a><span class="lineno"> 46</span> }</div> <div class="line"><a name="l00047"></a><span class="lineno"> 47</span> </div> <div class="line"><a name="l00048"></a><span class="lineno"> 48</span> <span class="keyword">template</span><<span class="keyword">typename</span> _Scalar></div> <div class="line"><a name="l00049"></a><span class="lineno"> 49</span> <span class="keyword">class </span>AngleAxis : <span class="keyword">public</span> RotationBase<AngleAxis<_Scalar>,3></div> <div class="line"><a name="l00050"></a><span class="lineno"> 50</span> {</div> <div class="line"><a name="l00051"></a><span class="lineno"> 51</span>  <span class="keyword">typedef</span> RotationBase<AngleAxis<_Scalar>,3> Base;</div> <div class="line"><a name="l00052"></a><span class="lineno"> 52</span> </div> <div class="line"><a name="l00053"></a><span class="lineno"> 53</span> <span class="keyword">public</span>:</div> <div class="line"><a name="l00054"></a><span class="lineno"> 54</span> </div> <div class="line"><a name="l00055"></a><span class="lineno"> 55</span>  <span class="keyword">using</span> Base::operator*;</div> <div class="line"><a name="l00056"></a><span class="lineno"> 56</span> </div> <div class="line"><a name="l00057"></a><span class="lineno"> 57</span>  <span class="keyword">enum</span> { Dim = 3 };</div> <div class="line"><a name="l00059"></a><span class="lineno"><a class="line" href="classEigen_1_1AngleAxis.html#aea29b63b6c32046ae9a471d82c5cf869"> 59</a></span>  <span class="keyword">typedef</span> _Scalar <a class="code" href="classEigen_1_1AngleAxis.html#aea29b63b6c32046ae9a471d82c5cf869">Scalar</a>;</div> <div class="line"><a name="l00060"></a><span class="lineno"> 60</span>  <span class="keyword">typedef</span> <a class="code" href="classEigen_1_1Matrix.html">Matrix<Scalar,3,3></a> <a class="code" href="classEigen_1_1Matrix.html">Matrix3</a>;</div> <div class="line"><a name="l00061"></a><span class="lineno"> 61</span>  <span class="keyword">typedef</span> <a class="code" href="classEigen_1_1Matrix.html">Matrix<Scalar,3,1></a> <a class="code" href="classEigen_1_1Matrix.html">Vector3</a>;</div> <div class="line"><a name="l00062"></a><span class="lineno"> 62</span>  <span class="keyword">typedef</span> <a class="code" href="classEigen_1_1Quaternion.html">Quaternion<Scalar></a> <a class="code" href="classEigen_1_1Quaternion.html">QuaternionType</a>;</div> <div class="line"><a name="l00063"></a><span class="lineno"> 63</span> </div> <div class="line"><a name="l00064"></a><span class="lineno"> 64</span> <span class="keyword">protected</span>:</div> <div class="line"><a name="l00065"></a><span class="lineno"> 65</span> </div> <div class="line"><a name="l00066"></a><span class="lineno"> 66</span>  <a class="code" href="classEigen_1_1Matrix.html">Vector3</a> m_axis;</div> <div class="line"><a name="l00067"></a><span class="lineno"> 67</span>  <a class="code" href="classEigen_1_1AngleAxis.html#aea29b63b6c32046ae9a471d82c5cf869">Scalar</a> m_angle;</div> <div class="line"><a name="l00068"></a><span class="lineno"> 68</span> </div> <div class="line"><a name="l00069"></a><span class="lineno"> 69</span> <span class="keyword">public</span>:</div> <div class="line"><a name="l00070"></a><span class="lineno"> 70</span> </div> <div class="line"><a name="l00072"></a><span class="lineno"><a class="line" href="classEigen_1_1AngleAxis.html#a38772216b8d784ddbe99d2d3df961f32"> 72</a></span>  <a class="code" href="classEigen_1_1AngleAxis.html#a38772216b8d784ddbe99d2d3df961f32">AngleAxis</a>() {}</div> <div class="line"><a name="l00078"></a><span class="lineno"> 78</span>  <span class="keyword">template</span><<span class="keyword">typename</span> Derived></div> <div class="line"><a name="l00079"></a><span class="lineno"><a class="line" href="classEigen_1_1AngleAxis.html#aff009fcfc0377acd2db48fb754ae1f62"> 79</a></span>  <span class="keyword">inline</span> <a class="code" href="classEigen_1_1AngleAxis.html#aff009fcfc0377acd2db48fb754ae1f62">AngleAxis</a>(<span class="keyword">const</span> <a class="code" href="classEigen_1_1AngleAxis.html#aea29b63b6c32046ae9a471d82c5cf869">Scalar</a>& angle, <span class="keyword">const</span> <a class="code" href="classEigen_1_1MatrixBase.html">MatrixBase<Derived></a>& axis) : m_axis(axis), m_angle(angle) {}</div> <div class="line"><a name="l00081"></a><span class="lineno"><a class="line" href="classEigen_1_1AngleAxis.html#aae1f1d7628b3f3d29e26dc480f4d037d"> 81</a></span>  <span class="keyword">template</span><<span class="keyword">typename</span> QuatDerived> <span class="keyword">inline</span> <span class="keyword">explicit</span> <a class="code" href="classEigen_1_1AngleAxis.html#aae1f1d7628b3f3d29e26dc480f4d037d">AngleAxis</a>(<span class="keyword">const</span> <a class="code" href="classEigen_1_1QuaternionBase.html">QuaternionBase<QuatDerived></a>& q) { *<span class="keyword">this</span> = q; }</div> <div class="line"><a name="l00083"></a><span class="lineno"> 83</span>  <span class="keyword">template</span><<span class="keyword">typename</span> Derived></div> <div class="line"><a name="l00084"></a><span class="lineno"><a class="line" href="classEigen_1_1AngleAxis.html#a45ad10816a418d3d86130608fd6ca996"> 84</a></span>  <span class="keyword">inline</span> <span class="keyword">explicit</span> <a class="code" href="classEigen_1_1AngleAxis.html#a45ad10816a418d3d86130608fd6ca996">AngleAxis</a>(<span class="keyword">const</span> <a class="code" href="classEigen_1_1MatrixBase.html">MatrixBase<Derived></a>& m) { *<span class="keyword">this</span> = m; }</div> <div class="line"><a name="l00085"></a><span class="lineno"> 85</span> </div> <div class="line"><a name="l00086"></a><span class="lineno"> 86</span>  <a class="code" href="classEigen_1_1AngleAxis.html#aea29b63b6c32046ae9a471d82c5cf869">Scalar</a> angle()<span class="keyword"> const </span>{ <span class="keywordflow">return</span> m_angle; }</div> <div class="line"><a name="l00087"></a><span class="lineno"> 87</span>  <a class="code" href="classEigen_1_1AngleAxis.html#aea29b63b6c32046ae9a471d82c5cf869">Scalar</a>& angle() { <span class="keywordflow">return</span> m_angle; }</div> <div class="line"><a name="l00088"></a><span class="lineno"> 88</span> </div> <div class="line"><a name="l00089"></a><span class="lineno"> 89</span>  <span class="keyword">const</span> Vector3& axis()<span class="keyword"> const </span>{ <span class="keywordflow">return</span> m_axis; }</div> <div class="line"><a name="l00090"></a><span class="lineno"> 90</span>  Vector3& axis() { <span class="keywordflow">return</span> m_axis; }</div> <div class="line"><a name="l00091"></a><span class="lineno"> 91</span> </div> <div class="line"><a name="l00093"></a><span class="lineno"><a class="line" href="classEigen_1_1AngleAxis.html#afbaeb13e13184d381a9e443166fc66b8"> 93</a></span>  <span class="keyword">inline</span> <a class="code" href="classEigen_1_1Quaternion.html">QuaternionType</a> <a class="code" href="classEigen_1_1AngleAxis.html#afbaeb13e13184d381a9e443166fc66b8">operator* </a>(<span class="keyword">const</span> <a class="code" href="classEigen_1_1AngleAxis.html">AngleAxis</a>& other)<span class="keyword"> const</span></div> <div class="line"><a name="l00094"></a><span class="lineno"> 94</span> <span class="keyword"> </span>{ <span class="keywordflow">return</span> <a class="code" href="classEigen_1_1Quaternion.html">QuaternionType</a>(*<span class="keyword">this</span>) * <a class="code" href="classEigen_1_1Quaternion.html">QuaternionType</a>(other); }</div> <div class="line"><a name="l00095"></a><span class="lineno"> 95</span> </div> <div class="line"><a name="l00097"></a><span class="lineno"><a class="line" href="classEigen_1_1AngleAxis.html#a83533759bf90b1514c756cd538101736"> 97</a></span>  <span class="keyword">inline</span> <a class="code" href="classEigen_1_1Quaternion.html">QuaternionType</a> <a class="code" href="classEigen_1_1AngleAxis.html#afbaeb13e13184d381a9e443166fc66b8">operator* </a>(<span class="keyword">const</span> <a class="code" href="classEigen_1_1Quaternion.html">QuaternionType</a>& other)<span class="keyword"> const</span></div> <div class="line"><a name="l00098"></a><span class="lineno"> 98</span> <span class="keyword"> </span>{ <span class="keywordflow">return</span> <a class="code" href="classEigen_1_1Quaternion.html">QuaternionType</a>(*<span class="keyword">this</span>) * other; }</div> <div class="line"><a name="l00099"></a><span class="lineno"> 99</span> </div> <div class="line"><a name="l00101"></a><span class="lineno"><a class="line" href="classEigen_1_1AngleAxis.html#acdec7200dd68b70b615bcda9659fed26"> 101</a></span>  <span class="keyword">friend</span> <span class="keyword">inline</span> <a class="code" href="classEigen_1_1Quaternion.html">QuaternionType</a> <a class="code" href="classEigen_1_1AngleAxis.html#afbaeb13e13184d381a9e443166fc66b8">operator* </a>(<span class="keyword">const</span> <a class="code" href="classEigen_1_1Quaternion.html">QuaternionType</a>& a, <span class="keyword">const</span> <a class="code" href="classEigen_1_1AngleAxis.html">AngleAxis</a>& b)</div> <div class="line"><a name="l00102"></a><span class="lineno"> 102</span>  { <span class="keywordflow">return</span> a * <a class="code" href="classEigen_1_1Quaternion.html">QuaternionType</a>(b); }</div> <div class="line"><a name="l00103"></a><span class="lineno"> 103</span> </div> <div class="line"><a name="l00105"></a><span class="lineno"><a class="line" href="classEigen_1_1AngleAxis.html#a68dbca4080079a7376bc8d6cfd8fe098"> 105</a></span>  <a class="code" href="classEigen_1_1AngleAxis.html">AngleAxis</a> <a class="code" href="classEigen_1_1AngleAxis.html#a68dbca4080079a7376bc8d6cfd8fe098">inverse</a>()<span class="keyword"> const</span></div> <div class="line"><a name="l00106"></a><span class="lineno"> 106</span> <span class="keyword"> </span>{ <span class="keywordflow">return</span> <a class="code" href="classEigen_1_1AngleAxis.html#a38772216b8d784ddbe99d2d3df961f32">AngleAxis</a>(-m_angle, m_axis); }</div> <div class="line"><a name="l00107"></a><span class="lineno"> 107</span> </div> <div class="line"><a name="l00108"></a><span class="lineno"> 108</span>  <span class="keyword">template</span><<span class="keyword">class</span> QuatDerived></div> <div class="line"><a name="l00109"></a><span class="lineno"> 109</span>  <a class="code" href="classEigen_1_1AngleAxis.html">AngleAxis</a>& operator=(<span class="keyword">const</span> <a class="code" href="classEigen_1_1QuaternionBase.html">QuaternionBase<QuatDerived></a>& q);</div> <div class="line"><a name="l00110"></a><span class="lineno"> 110</span>  <span class="keyword">template</span><<span class="keyword">typename</span> Derived></div> <div class="line"><a name="l00111"></a><span class="lineno"> 111</span>  <a class="code" href="classEigen_1_1AngleAxis.html">AngleAxis</a>& operator=(<span class="keyword">const</span> <a class="code" href="classEigen_1_1MatrixBase.html">MatrixBase<Derived></a>& m);</div> <div class="line"><a name="l00112"></a><span class="lineno"> 112</span> </div> <div class="line"><a name="l00113"></a><span class="lineno"> 113</span>  <span class="keyword">template</span><<span class="keyword">typename</span> Derived></div> <div class="line"><a name="l00114"></a><span class="lineno"> 114</span>  <a class="code" href="classEigen_1_1AngleAxis.html">AngleAxis</a>& fromRotationMatrix(<span class="keyword">const</span> <a class="code" href="classEigen_1_1MatrixBase.html">MatrixBase<Derived></a>& m);</div> <div class="line"><a name="l00115"></a><span class="lineno"> 115</span>  Matrix3 <a class="code" href="classEigen_1_1AngleAxis.html#a34adf4403bbc909c67309fbc7b6dbc16">toRotationMatrix</a>(<span class="keywordtype">void</span>) <span class="keyword">const</span>;</div> <div class="line"><a name="l00116"></a><span class="lineno"> 116</span> </div> <div class="line"><a name="l00122"></a><span class="lineno"> 122</span>  <span class="keyword">template</span><<span class="keyword">typename</span> NewScalarType></div> <div class="line"><a name="l00123"></a><span class="lineno"><a class="line" href="classEigen_1_1AngleAxis.html#a0469e863c2bf6b5c2d5655c8f4956d65"> 123</a></span>  <span class="keyword">inline</span> <span class="keyword">typename</span> internal::cast_return_type<AngleAxis,AngleAxis<NewScalarType> >::type <a class="code" href="classEigen_1_1AngleAxis.html#a0469e863c2bf6b5c2d5655c8f4956d65">cast</a>()<span class="keyword"> const</span></div> <div class="line"><a name="l00124"></a><span class="lineno"> 124</span> <span class="keyword"> </span>{ <span class="keywordflow">return</span> <span class="keyword">typename</span> internal::cast_return_type<AngleAxis,AngleAxis<NewScalarType> >::type(*<span class="keyword">this</span>); }</div> <div class="line"><a name="l00125"></a><span class="lineno"> 125</span> </div> <div class="line"><a name="l00127"></a><span class="lineno"> 127</span>  <span class="keyword">template</span><<span class="keyword">typename</span> OtherScalarType></div> <div class="line"><a name="l00128"></a><span class="lineno"><a class="line" href="classEigen_1_1AngleAxis.html#ad71c1c9e6a581420cba5205f837af1b5"> 128</a></span>  <span class="keyword">inline</span> <span class="keyword">explicit</span> <a class="code" href="classEigen_1_1AngleAxis.html#ad71c1c9e6a581420cba5205f837af1b5">AngleAxis</a>(<span class="keyword">const</span> <a class="code" href="classEigen_1_1AngleAxis.html">AngleAxis<OtherScalarType></a>& other)</div> <div class="line"><a name="l00129"></a><span class="lineno"> 129</span>  {</div> <div class="line"><a name="l00130"></a><span class="lineno"> 130</span>  m_axis = other.axis().template cast<Scalar>();</div> <div class="line"><a name="l00131"></a><span class="lineno"> 131</span>  m_angle = <a class="code" href="classEigen_1_1AngleAxis.html#aea29b63b6c32046ae9a471d82c5cf869">Scalar</a>(other.angle());</div> <div class="line"><a name="l00132"></a><span class="lineno"> 132</span>  }</div> <div class="line"><a name="l00133"></a><span class="lineno"> 133</span> </div> <div class="line"><a name="l00134"></a><span class="lineno"> 134</span>  <span class="keyword">static</span> <span class="keyword">inline</span> <span class="keyword">const</span> <a class="code" href="classEigen_1_1AngleAxis.html">AngleAxis</a> Identity() { <span class="keywordflow">return</span> <a class="code" href="classEigen_1_1AngleAxis.html#a38772216b8d784ddbe99d2d3df961f32">AngleAxis</a>(0, Vector3::UnitX()); }</div> <div class="line"><a name="l00135"></a><span class="lineno"> 135</span> </div> <div class="line"><a name="l00140"></a><span class="lineno"><a class="line" href="classEigen_1_1AngleAxis.html#a729107ea8555246d34f872febb8b14ce"> 140</a></span>  <span class="keywordtype">bool</span> <a class="code" href="classEigen_1_1AngleAxis.html#a729107ea8555246d34f872febb8b14ce">isApprox</a>(<span class="keyword">const</span> <a class="code" href="classEigen_1_1AngleAxis.html">AngleAxis</a>& other, <span class="keyword">const</span> <span class="keyword">typename</span> <a class="code" href="structEigen_1_1NumTraits.html">NumTraits<Scalar>::Real</a>& prec = <a class="code" href="structEigen_1_1NumTraits.html">NumTraits<Scalar>::dummy_precision</a>())<span class="keyword"> const</span></div> <div class="line"><a name="l00141"></a><span class="lineno"> 141</span> <span class="keyword"> </span>{ <span class="keywordflow">return</span> m_axis.isApprox(other.m_axis, prec) && internal::isApprox(m_angle,other.m_angle, prec); }</div> <div class="line"><a name="l00142"></a><span class="lineno"> 142</span> };</div> <div class="line"><a name="l00143"></a><span class="lineno"> 143</span> </div> <div class="line"><a name="l00146"></a><span class="lineno"><a class="line" href="group__Geometry__Module.html#ga811d6fdab2002723bc7a72f055ce8c1d"> 146</a></span> <span class="keyword">typedef</span> <a class="code" href="classEigen_1_1AngleAxis.html">AngleAxis<float></a> <a class="code" href="group__Geometry__Module.html#ga811d6fdab2002723bc7a72f055ce8c1d">AngleAxisf</a>;</div> <div class="line"><a name="l00149"></a><span class="lineno"><a class="line" href="group__Geometry__Module.html#ga0db1cc067c51aaa6dedf5805ee0c53d7"> 149</a></span> <span class="keyword">typedef</span> <a class="code" href="classEigen_1_1AngleAxis.html">AngleAxis<double></a> <a class="code" href="group__Geometry__Module.html#ga0db1cc067c51aaa6dedf5805ee0c53d7">AngleAxisd</a>;</div> <div class="line"><a name="l00150"></a><span class="lineno"> 150</span> </div> <div class="line"><a name="l00157"></a><span class="lineno"> 157</span> <span class="keyword">template</span><<span class="keyword">typename</span> Scalar></div> <div class="line"><a name="l00158"></a><span class="lineno"> 158</span> <span class="keyword">template</span><<span class="keyword">typename</span> QuatDerived></div> <div class="line"><a name="l00159"></a><span class="lineno"><a class="line" href="classEigen_1_1AngleAxis.html#a9694b8b0ff05dc4673d64af0e566d9d9"> 159</a></span> <a class="code" href="classEigen_1_1AngleAxis.html">AngleAxis<Scalar></a>& <a class="code" href="classEigen_1_1AngleAxis.html">AngleAxis<Scalar>::operator=</a>(<span class="keyword">const</span> <a class="code" href="classEigen_1_1QuaternionBase.html">QuaternionBase<QuatDerived></a>& q)</div> <div class="line"><a name="l00160"></a><span class="lineno"> 160</span> {</div> <div class="line"><a name="l00161"></a><span class="lineno"> 161</span>  <span class="keyword">using</span> std::acos;</div> <div class="line"><a name="l00162"></a><span class="lineno"> 162</span>  <span class="keyword">using</span> std::min;</div> <div class="line"><a name="l00163"></a><span class="lineno"> 163</span>  <span class="keyword">using</span> std::max;</div> <div class="line"><a name="l00164"></a><span class="lineno"> 164</span>  <span class="keyword">using</span> std::sqrt;</div> <div class="line"><a name="l00165"></a><span class="lineno"> 165</span>  <a class="code" href="classEigen_1_1AngleAxis.html#aea29b63b6c32046ae9a471d82c5cf869">Scalar</a> n2 = q.<a class="code" href="classEigen_1_1QuaternionBase.html#a154f40d7fec7be85112f7e7c12e4900b">vec</a>().squaredNorm();</div> <div class="line"><a name="l00166"></a><span class="lineno"> 166</span>  <span class="keywordflow">if</span> (n2 < <a class="code" href="structEigen_1_1NumTraits.html">NumTraits<Scalar>::dummy_precision</a>()*<a class="code" href="structEigen_1_1NumTraits.html">NumTraits<Scalar>::dummy_precision</a>())</div> <div class="line"><a name="l00167"></a><span class="lineno"> 167</span>  {</div> <div class="line"><a name="l00168"></a><span class="lineno"> 168</span>  m_angle = 0;</div> <div class="line"><a name="l00169"></a><span class="lineno"> 169</span>  m_axis << 1, 0, 0;</div> <div class="line"><a name="l00170"></a><span class="lineno"> 170</span>  }</div> <div class="line"><a name="l00171"></a><span class="lineno"> 171</span>  <span class="keywordflow">else</span></div> <div class="line"><a name="l00172"></a><span class="lineno"> 172</span>  {</div> <div class="line"><a name="l00173"></a><span class="lineno"> 173</span>  m_angle = <a class="code" href="classEigen_1_1AngleAxis.html#aea29b63b6c32046ae9a471d82c5cf869">Scalar</a>(2)*acos((min)((max)(<a class="code" href="classEigen_1_1AngleAxis.html#aea29b63b6c32046ae9a471d82c5cf869">Scalar</a>(-1),q.<a class="code" href="classEigen_1_1QuaternionBase.html#ae11dcbf26de9bc8c2afddf9a8b5c32a9">w</a>()),<a class="code" href="classEigen_1_1AngleAxis.html#aea29b63b6c32046ae9a471d82c5cf869">Scalar</a>(1)));</div> <div class="line"><a name="l00174"></a><span class="lineno"> 174</span>  m_axis = q.<a class="code" href="classEigen_1_1QuaternionBase.html#a154f40d7fec7be85112f7e7c12e4900b">vec</a>() / sqrt(n2);</div> <div class="line"><a name="l00175"></a><span class="lineno"> 175</span>  }</div> <div class="line"><a name="l00176"></a><span class="lineno"> 176</span>  <span class="keywordflow">return</span> *<span class="keyword">this</span>;</div> <div class="line"><a name="l00177"></a><span class="lineno"> 177</span> }</div> <div class="line"><a name="l00178"></a><span class="lineno"> 178</span> </div> <div class="line"><a name="l00181"></a><span class="lineno"> 181</span> <span class="keyword">template</span><<span class="keyword">typename</span> Scalar></div> <div class="line"><a name="l00182"></a><span class="lineno"> 182</span> <span class="keyword">template</span><<span class="keyword">typename</span> Derived></div> <div class="line"><a name="l00183"></a><span class="lineno"><a class="line" href="classEigen_1_1AngleAxis.html#aa472ded3ba078c186ac5cbac1ad4ad14"> 183</a></span> <a class="code" href="classEigen_1_1AngleAxis.html">AngleAxis<Scalar></a>& <a class="code" href="classEigen_1_1AngleAxis.html">AngleAxis<Scalar>::operator=</a>(<span class="keyword">const</span> <a class="code" href="classEigen_1_1MatrixBase.html">MatrixBase<Derived></a>& mat)</div> <div class="line"><a name="l00184"></a><span class="lineno"> 184</span> {</div> <div class="line"><a name="l00185"></a><span class="lineno"> 185</span>  <span class="comment">// Since a direct conversion would not be really faster,</span></div> <div class="line"><a name="l00186"></a><span class="lineno"> 186</span>  <span class="comment">// let's use the robust Quaternion implementation:</span></div> <div class="line"><a name="l00187"></a><span class="lineno"> 187</span>  <span class="keywordflow">return</span> *<span class="keyword">this</span> = <a class="code" href="classEigen_1_1Quaternion.html">QuaternionType</a>(mat);</div> <div class="line"><a name="l00188"></a><span class="lineno"> 188</span> }</div> <div class="line"><a name="l00189"></a><span class="lineno"> 189</span> </div> <div class="line"><a name="l00193"></a><span class="lineno"> 193</span> <span class="keyword">template</span><<span class="keyword">typename</span> Scalar></div> <div class="line"><a name="l00194"></a><span class="lineno"> 194</span> <span class="keyword">template</span><<span class="keyword">typename</span> Derived></div> <div class="line"><a name="l00195"></a><span class="lineno"><a class="line" href="classEigen_1_1AngleAxis.html#a45563d717551a8d18438d2481349536a"> 195</a></span> <a class="code" href="classEigen_1_1AngleAxis.html">AngleAxis<Scalar></a>& <a class="code" href="classEigen_1_1AngleAxis.html">AngleAxis<Scalar>::fromRotationMatrix</a>(<span class="keyword">const</span> <a class="code" href="classEigen_1_1MatrixBase.html">MatrixBase<Derived></a>& mat)</div> <div class="line"><a name="l00196"></a><span class="lineno"> 196</span> {</div> <div class="line"><a name="l00197"></a><span class="lineno"> 197</span>  <span class="keywordflow">return</span> *<span class="keyword">this</span> = <a class="code" href="classEigen_1_1Quaternion.html">QuaternionType</a>(mat);</div> <div class="line"><a name="l00198"></a><span class="lineno"> 198</span> }</div> <div class="line"><a name="l00199"></a><span class="lineno"> 199</span> </div> <div class="line"><a name="l00202"></a><span class="lineno"> 202</span> <span class="keyword">template</span><<span class="keyword">typename</span> Scalar></div> <div class="line"><a name="l00203"></a><span class="lineno"> 203</span> <span class="keyword">typename</span> <a class="code" href="classEigen_1_1Matrix.html">AngleAxis<Scalar>::Matrix3</a></div> <div class="line"><a name="l00204"></a><span class="lineno"><a class="line" href="classEigen_1_1AngleAxis.html#a34adf4403bbc909c67309fbc7b6dbc16"> 204</a></span> <a class="code" href="classEigen_1_1AngleAxis.html#a34adf4403bbc909c67309fbc7b6dbc16">AngleAxis<Scalar>::toRotationMatrix</a>(<span class="keywordtype">void</span>)<span class="keyword"> const</span></div> <div class="line"><a name="l00205"></a><span class="lineno"> 205</span> <span class="keyword"></span>{</div> <div class="line"><a name="l00206"></a><span class="lineno"> 206</span>  <span class="keyword">using</span> std::sin;</div> <div class="line"><a name="l00207"></a><span class="lineno"> 207</span>  <span class="keyword">using</span> std::cos;</div> <div class="line"><a name="l00208"></a><span class="lineno"> 208</span>  <a class="code" href="classEigen_1_1Matrix.html">Matrix3</a> res;</div> <div class="line"><a name="l00209"></a><span class="lineno"> 209</span>  <a class="code" href="classEigen_1_1Matrix.html">Vector3</a> sin_axis = sin(m_angle) * m_axis;</div> <div class="line"><a name="l00210"></a><span class="lineno"> 210</span>  <a class="code" href="classEigen_1_1AngleAxis.html#aea29b63b6c32046ae9a471d82c5cf869">Scalar</a> c = cos(m_angle);</div> <div class="line"><a name="l00211"></a><span class="lineno"> 211</span>  <a class="code" href="classEigen_1_1Matrix.html">Vector3</a> cos1_axis = (<a class="code" href="classEigen_1_1AngleAxis.html#aea29b63b6c32046ae9a471d82c5cf869">Scalar</a>(1)-c) * m_axis;</div> <div class="line"><a name="l00212"></a><span class="lineno"> 212</span> </div> <div class="line"><a name="l00213"></a><span class="lineno"> 213</span>  <a class="code" href="classEigen_1_1AngleAxis.html#aea29b63b6c32046ae9a471d82c5cf869">Scalar</a> tmp;</div> <div class="line"><a name="l00214"></a><span class="lineno"> 214</span>  tmp = cos1_axis.x() * m_axis.y();</div> <div class="line"><a name="l00215"></a><span class="lineno"> 215</span>  res.coeffRef(0,1) = tmp - sin_axis.z();</div> <div class="line"><a name="l00216"></a><span class="lineno"> 216</span>  res.coeffRef(1,0) = tmp + sin_axis.z();</div> <div class="line"><a name="l00217"></a><span class="lineno"> 217</span> </div> <div class="line"><a name="l00218"></a><span class="lineno"> 218</span>  tmp = cos1_axis.x() * m_axis.z();</div> <div class="line"><a name="l00219"></a><span class="lineno"> 219</span>  res.coeffRef(0,2) = tmp + sin_axis.y();</div> <div class="line"><a name="l00220"></a><span class="lineno"> 220</span>  res.coeffRef(2,0) = tmp - sin_axis.y();</div> <div class="line"><a name="l00221"></a><span class="lineno"> 221</span> </div> <div class="line"><a name="l00222"></a><span class="lineno"> 222</span>  tmp = cos1_axis.y() * m_axis.z();</div> <div class="line"><a name="l00223"></a><span class="lineno"> 223</span>  res.coeffRef(1,2) = tmp - sin_axis.x();</div> <div class="line"><a name="l00224"></a><span class="lineno"> 224</span>  res.coeffRef(2,1) = tmp + sin_axis.x();</div> <div class="line"><a name="l00225"></a><span class="lineno"> 225</span> </div> <div class="line"><a name="l00226"></a><span class="lineno"> 226</span>  res.<a class="code" href="classEigen_1_1MatrixBase.html#a0a45dd0ed5a44ec3f8f43239f2e4ac25">diagonal</a>() = (cos1_axis.cwiseProduct(m_axis)).array() + c;</div> <div class="line"><a name="l00227"></a><span class="lineno"> 227</span> </div> <div class="line"><a name="l00228"></a><span class="lineno"> 228</span>  <span class="keywordflow">return</span> res;</div> <div class="line"><a name="l00229"></a><span class="lineno"> 229</span> }</div> <div class="line"><a name="l00230"></a><span class="lineno"> 230</span> </div> <div class="line"><a name="l00231"></a><span class="lineno"> 231</span> } <span class="comment">// end namespace Eigen</span></div> <div class="line"><a name="l00232"></a><span class="lineno"> 232</span> </div> <div class="line"><a name="l00233"></a><span class="lineno"> 233</span> <span class="preprocessor">#endif // EIGEN_ANGLEAXIS_H</span></div> <div class="ttc" id="classEigen_1_1AngleAxis_html_a38772216b8d784ddbe99d2d3df961f32"><div class="ttname"><a href="classEigen_1_1AngleAxis.html#a38772216b8d784ddbe99d2d3df961f32">Eigen::AngleAxis::AngleAxis</a></div><div class="ttdeci">AngleAxis()</div><div class="ttdef"><b>Definition:</b> AngleAxis.h:72</div></div> <div class="ttc" id="classEigen_1_1AngleAxis_html_a68dbca4080079a7376bc8d6cfd8fe098"><div class="ttname"><a href="classEigen_1_1AngleAxis.html#a68dbca4080079a7376bc8d6cfd8fe098">Eigen::AngleAxis::inverse</a></div><div class="ttdeci">AngleAxis inverse() const </div><div class="ttdef"><b>Definition:</b> AngleAxis.h:105</div></div> <div class="ttc" id="classEigen_1_1AngleAxis_html_a45ad10816a418d3d86130608fd6ca996"><div class="ttname"><a href="classEigen_1_1AngleAxis.html#a45ad10816a418d3d86130608fd6ca996">Eigen::AngleAxis::AngleAxis</a></div><div class="ttdeci">AngleAxis(const MatrixBase< Derived > &m)</div><div class="ttdef"><b>Definition:</b> AngleAxis.h:84</div></div> <div class="ttc" id="structEigen_1_1NumTraits_html"><div class="ttname"><a href="structEigen_1_1NumTraits.html">Eigen::NumTraits</a></div><div class="ttdoc">Holds information about the various numeric (i.e. scalar) types allowed by Eigen. ...</div><div class="ttdef"><b>Definition:</b> NumTraits.h:88</div></div> <div class="ttc" id="classEigen_1_1AngleAxis_html_a34adf4403bbc909c67309fbc7b6dbc16"><div class="ttname"><a href="classEigen_1_1AngleAxis.html#a34adf4403bbc909c67309fbc7b6dbc16">Eigen::AngleAxis::toRotationMatrix</a></div><div class="ttdeci">Matrix3 toRotationMatrix(void) const </div><div class="ttdef"><b>Definition:</b> AngleAxis.h:204</div></div> <div class="ttc" id="classEigen_1_1AngleAxis_html_aff009fcfc0377acd2db48fb754ae1f62"><div class="ttname"><a href="classEigen_1_1AngleAxis.html#aff009fcfc0377acd2db48fb754ae1f62">Eigen::AngleAxis::AngleAxis</a></div><div class="ttdeci">AngleAxis(const Scalar &angle, const MatrixBase< Derived > &axis)</div><div class="ttdef"><b>Definition:</b> AngleAxis.h:79</div></div> <div class="ttc" id="group__Geometry__Module_html_ga0db1cc067c51aaa6dedf5805ee0c53d7"><div class="ttname"><a href="group__Geometry__Module.html#ga0db1cc067c51aaa6dedf5805ee0c53d7">Eigen::AngleAxisd</a></div><div class="ttdeci">AngleAxis< double > AngleAxisd</div><div class="ttdef"><b>Definition:</b> AngleAxis.h:149</div></div> <div class="ttc" id="group__Geometry__Module_html_ga811d6fdab2002723bc7a72f055ce8c1d"><div class="ttname"><a href="group__Geometry__Module.html#ga811d6fdab2002723bc7a72f055ce8c1d">Eigen::AngleAxisf</a></div><div class="ttdeci">AngleAxis< float > AngleAxisf</div><div class="ttdef"><b>Definition:</b> AngleAxis.h:146</div></div> <div class="ttc" id="classEigen_1_1AngleAxis_html_aae1f1d7628b3f3d29e26dc480f4d037d"><div class="ttname"><a href="classEigen_1_1AngleAxis.html#aae1f1d7628b3f3d29e26dc480f4d037d">Eigen::AngleAxis::AngleAxis</a></div><div class="ttdeci">AngleAxis(const QuaternionBase< QuatDerived > &q)</div><div class="ttdef"><b>Definition:</b> AngleAxis.h:81</div></div> <div class="ttc" id="classEigen_1_1MatrixBase_html_a0a45dd0ed5a44ec3f8f43239f2e4ac25"><div class="ttname"><a href="classEigen_1_1MatrixBase.html#a0a45dd0ed5a44ec3f8f43239f2e4ac25">Eigen::MatrixBase::diagonal</a></div><div class="ttdeci">DiagonalReturnType diagonal()</div><div class="ttdef"><b>Definition:</b> Diagonal.h:168</div></div> <div class="ttc" id="classEigen_1_1AngleAxis_html_afbaeb13e13184d381a9e443166fc66b8"><div class="ttname"><a href="classEigen_1_1AngleAxis.html#afbaeb13e13184d381a9e443166fc66b8">Eigen::AngleAxis::operator*</a></div><div class="ttdeci">QuaternionType operator*(const AngleAxis &other) const </div><div class="ttdef"><b>Definition:</b> AngleAxis.h:93</div></div> <div class="ttc" id="classEigen_1_1QuaternionBase_html"><div class="ttname"><a href="classEigen_1_1QuaternionBase.html">Eigen::QuaternionBase</a></div><div class="ttdoc">Base class for quaternion expressions. </div><div class="ttdef"><b>Definition:</b> ForwardDeclarations.h:233</div></div> <div class="ttc" id="classEigen_1_1AngleAxis_html_a729107ea8555246d34f872febb8b14ce"><div class="ttname"><a href="classEigen_1_1AngleAxis.html#a729107ea8555246d34f872febb8b14ce">Eigen::AngleAxis::isApprox</a></div><div class="ttdeci">bool isApprox(const AngleAxis &other, const typename NumTraits< Scalar >::Real &prec=NumTraits< Scalar >::dummy_precision()) const </div><div class="ttdef"><b>Definition:</b> AngleAxis.h:140</div></div> <div class="ttc" id="classEigen_1_1AngleAxis_html_a0469e863c2bf6b5c2d5655c8f4956d65"><div class="ttname"><a href="classEigen_1_1AngleAxis.html#a0469e863c2bf6b5c2d5655c8f4956d65">Eigen::AngleAxis::cast</a></div><div class="ttdeci">internal::cast_return_type< AngleAxis, AngleAxis< NewScalarType > >::type cast() const </div><div class="ttdef"><b>Definition:</b> AngleAxis.h:123</div></div> <div class="ttc" id="classEigen_1_1Quaternion_html"><div class="ttname"><a href="classEigen_1_1Quaternion.html">Eigen::Quaternion</a></div><div class="ttdoc">The quaternion class used to represent 3D orientations and rotations. </div><div class="ttdef"><b>Definition:</b> ForwardDeclarations.h:260</div></div> <div class="ttc" id="classEigen_1_1AngleAxis_html_aea29b63b6c32046ae9a471d82c5cf869"><div class="ttname"><a href="classEigen_1_1AngleAxis.html#aea29b63b6c32046ae9a471d82c5cf869">Eigen::AngleAxis::Scalar</a></div><div class="ttdeci">_Scalar Scalar</div><div class="ttdef"><b>Definition:</b> AngleAxis.h:59</div></div> <div class="ttc" id="classEigen_1_1QuaternionBase_html_ae11dcbf26de9bc8c2afddf9a8b5c32a9"><div class="ttname"><a href="classEigen_1_1QuaternionBase.html#ae11dcbf26de9bc8c2afddf9a8b5c32a9">Eigen::QuaternionBase::w</a></div><div class="ttdeci">Scalar w() const </div><div class="ttdef"><b>Definition:</b> Quaternion.h:66</div></div> <div class="ttc" id="classEigen_1_1Matrix_html"><div class="ttname"><a href="classEigen_1_1Matrix.html">Eigen::Matrix</a></div><div class="ttdoc">The matrix class, also used for vectors and row-vectors. </div><div class="ttdef"><b>Definition:</b> Matrix.h:127</div></div> <div class="ttc" id="classEigen_1_1AngleAxis_html_ad71c1c9e6a581420cba5205f837af1b5"><div class="ttname"><a href="classEigen_1_1AngleAxis.html#ad71c1c9e6a581420cba5205f837af1b5">Eigen::AngleAxis::AngleAxis</a></div><div class="ttdeci">AngleAxis(const AngleAxis< OtherScalarType > &other)</div><div class="ttdef"><b>Definition:</b> AngleAxis.h:128</div></div> <div class="ttc" id="classEigen_1_1MatrixBase_html"><div class="ttname"><a href="classEigen_1_1MatrixBase.html">Eigen::MatrixBase</a></div><div class="ttdoc">Base class for all dense matrices, vectors, and expressions. </div><div class="ttdef"><b>Definition:</b> MatrixBase.h:48</div></div> <div class="ttc" id="classEigen_1_1AngleAxis_html"><div class="ttname"><a href="classEigen_1_1AngleAxis.html">Eigen::AngleAxis</a></div><div class="ttdoc">Represents a 3D rotation as a rotation angle around an arbitrary 3D axis. </div><div class="ttdef"><b>Definition:</b> ForwardDeclarations.h:235</div></div> <div class="ttc" id="classEigen_1_1QuaternionBase_html_a154f40d7fec7be85112f7e7c12e4900b"><div class="ttname"><a href="classEigen_1_1QuaternionBase.html#a154f40d7fec7be85112f7e7c12e4900b">Eigen::QuaternionBase::vec</a></div><div class="ttdeci">const VectorBlock< const Coefficients, 3 > vec() const </div><div class="ttdef"><b>Definition:</b> Quaternion.h:78</div></div> </div><!-- fragment --></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="dir_e49d68e3078f12dfcf157021597ad168.html">Eigen</a></li><li class="navelem"><a class="el" href="dir_64b228556dc7f9fe757d43bb57fbfc24.html">src</a></li><li class="navelem"><a class="el" href="dir_d3830aef5e02ff83b0a743a58bf9d99e.html">Geometry</a></li><li class="navelem"><b>AngleAxis.h</b></li> <li class="footer">Generated on Mon Oct 28 2013 11:04:22 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>