<!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: Quaternion.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('Quaternion_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">Quaternion.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-2010 Gael Guennebaud <gael.guennebaud@inria.fr></span></div> <div class="line"><a name="l00005"></a><span class="lineno"> 5</span> <span class="comment">// Copyright (C) 2009 Mathieu Gautier <mathieu.gautier@cea.fr></span></div> <div class="line"><a name="l00006"></a><span class="lineno"> 6</span> <span class="comment">//</span></div> <div class="line"><a name="l00007"></a><span class="lineno"> 7</span> <span class="comment">// This Source Code Form is subject to the terms of the Mozilla</span></div> <div class="line"><a name="l00008"></a><span class="lineno"> 8</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="l00009"></a><span class="lineno"> 9</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="l00010"></a><span class="lineno"> 10</span> </div> <div class="line"><a name="l00011"></a><span class="lineno"> 11</span> <span class="preprocessor">#ifndef EIGEN_QUATERNION_H</span></div> <div class="line"><a name="l00012"></a><span class="lineno"> 12</span> <span class="preprocessor"></span><span class="preprocessor">#define EIGEN_QUATERNION_H</span></div> <div class="line"><a name="l00013"></a><span class="lineno"> 13</span> <span class="preprocessor"></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="l00015"></a><span class="lineno"> 15</span> </div> <div class="line"><a name="l00016"></a><span class="lineno"> 16</span> <span class="comment">/***************************************************************************</span></div> <div class="line"><a name="l00017"></a><span class="lineno"> 17</span> <span class="comment">* Definition of QuaternionBase<Derived></span></div> <div class="line"><a name="l00018"></a><span class="lineno"> 18</span> <span class="comment">* The implementation is at the end of the file</span></div> <div class="line"><a name="l00019"></a><span class="lineno"> 19</span> <span class="comment">***************************************************************************/</span></div> <div class="line"><a name="l00020"></a><span class="lineno"> 20</span> </div> <div class="line"><a name="l00021"></a><span class="lineno"> 21</span> <span class="keyword">namespace </span>internal {</div> <div class="line"><a name="l00022"></a><span class="lineno"> 22</span> <span class="keyword">template</span><<span class="keyword">typename</span> Other,</div> <div class="line"><a name="l00023"></a><span class="lineno"> 23</span>  <span class="keywordtype">int</span> OtherRows=Other::RowsAtCompileTime,</div> <div class="line"><a name="l00024"></a><span class="lineno"> 24</span>  <span class="keywordtype">int</span> OtherCols=Other::ColsAtCompileTime></div> <div class="line"><a name="l00025"></a><span class="lineno"> 25</span> <span class="keyword">struct </span>quaternionbase_assign_impl;</div> <div class="line"><a name="l00026"></a><span class="lineno"> 26</span> }</div> <div class="line"><a name="l00027"></a><span class="lineno"> 27</span> </div> <div class="line"><a name="l00034"></a><span class="lineno"> 34</span> <span class="keyword">template</span><<span class="keyword">class</span> Derived></div> <div class="line"><a name="l00035"></a><span class="lineno"> 35</span> <span class="keyword">class </span>QuaternionBase : <span class="keyword">public</span> RotationBase<Derived, 3></div> <div class="line"><a name="l00036"></a><span class="lineno"> 36</span> {</div> <div class="line"><a name="l00037"></a><span class="lineno"> 37</span>  <span class="keyword">typedef</span> RotationBase<Derived, 3> Base;</div> <div class="line"><a name="l00038"></a><span class="lineno"> 38</span> <span class="keyword">public</span>:</div> <div class="line"><a name="l00039"></a><span class="lineno"> 39</span>  <span class="keyword">using</span> Base::operator*;</div> <div class="line"><a name="l00040"></a><span class="lineno"> 40</span>  <span class="keyword">using</span> Base::derived;</div> <div class="line"><a name="l00041"></a><span class="lineno"> 41</span> </div> <div class="line"><a name="l00042"></a><span class="lineno"> 42</span>  <span class="keyword">typedef</span> <span class="keyword">typename</span> internal::traits<Derived>::Scalar Scalar;</div> <div class="line"><a name="l00043"></a><span class="lineno"> 43</span>  <span class="keyword">typedef</span> <span class="keyword">typename</span> NumTraits<Scalar>::Real RealScalar;</div> <div class="line"><a name="l00044"></a><span class="lineno"> 44</span>  <span class="keyword">typedef</span> <span class="keyword">typename</span> internal::traits<Derived>::Coefficients Coefficients;</div> <div class="line"><a name="l00045"></a><span class="lineno"> 45</span>  <span class="keyword">enum</span> {</div> <div class="line"><a name="l00046"></a><span class="lineno"> 46</span>  Flags = Eigen::internal::traits<Derived>::Flags</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> </div> <div class="line"><a name="l00049"></a><span class="lineno"> 49</span>  <span class="comment">// typedef typename Matrix<Scalar,4,1> Coefficients;</span></div> <div class="line"><a name="l00051"></a><span class="lineno"><a class="line" href="classEigen_1_1QuaternionBase.html#ab05eeebdd63f50a52fecb6353018075f"> 51</a></span> <span class="comment"></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_1QuaternionBase.html#ab05eeebdd63f50a52fecb6353018075f">Vector3</a>;</div> <div class="line"><a name="l00053"></a><span class="lineno"><a class="line" href="classEigen_1_1QuaternionBase.html#abdb0743d0fe840b927780b2bd78fb29e"> 53</a></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_1QuaternionBase.html#abdb0743d0fe840b927780b2bd78fb29e">Matrix3</a>;</div> <div class="line"><a name="l00055"></a><span class="lineno"><a class="line" href="classEigen_1_1QuaternionBase.html#ac8cbb66bbccf5f2c9596004a14129c26"> 55</a></span>  <span class="keyword">typedef</span> <a class="code" href="classEigen_1_1AngleAxis.html">AngleAxis<Scalar></a> <a class="code" href="classEigen_1_1QuaternionBase.html#ac8cbb66bbccf5f2c9596004a14129c26">AngleAxisType</a>;</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> </div> <div class="line"><a name="l00058"></a><span class="lineno"> 58</span> </div> <div class="line"><a name="l00060"></a><span class="lineno"><a class="line" href="classEigen_1_1QuaternionBase.html#a5eb2596dc2a509b276d01578a9c3dd27"> 60</a></span>  <span class="keyword">inline</span> <a class="code" href="classEigen_1_1RotationBase.html#a844358c46408e878e60c4026c52eb1e9">Scalar</a> <a class="code" href="classEigen_1_1QuaternionBase.html#a5eb2596dc2a509b276d01578a9c3dd27">x</a>()<span class="keyword"> const </span>{ <span class="keywordflow">return</span> this->derived().coeffs().coeff(0); }</div> <div class="line"><a name="l00062"></a><span class="lineno"><a class="line" href="classEigen_1_1QuaternionBase.html#a7b19faa70d70b5f16a39f1dd4b69e7a0"> 62</a></span>  <span class="keyword">inline</span> <a class="code" href="classEigen_1_1RotationBase.html#a844358c46408e878e60c4026c52eb1e9">Scalar</a> <a class="code" href="classEigen_1_1QuaternionBase.html#a7b19faa70d70b5f16a39f1dd4b69e7a0">y</a>()<span class="keyword"> const </span>{ <span class="keywordflow">return</span> this->derived().coeffs().coeff(1); }</div> <div class="line"><a name="l00064"></a><span class="lineno"><a class="line" href="classEigen_1_1QuaternionBase.html#a94a075384f8f54dc97a69546e600b78e"> 64</a></span>  <span class="keyword">inline</span> <a class="code" href="classEigen_1_1RotationBase.html#a844358c46408e878e60c4026c52eb1e9">Scalar</a> <a class="code" href="classEigen_1_1QuaternionBase.html#a94a075384f8f54dc97a69546e600b78e">z</a>()<span class="keyword"> const </span>{ <span class="keywordflow">return</span> this->derived().coeffs().coeff(2); }</div> <div class="line"><a name="l00066"></a><span class="lineno"><a class="line" href="classEigen_1_1QuaternionBase.html#ae11dcbf26de9bc8c2afddf9a8b5c32a9"> 66</a></span>  <span class="keyword">inline</span> <a class="code" href="classEigen_1_1RotationBase.html#a844358c46408e878e60c4026c52eb1e9">Scalar</a> <a class="code" href="classEigen_1_1QuaternionBase.html#ae11dcbf26de9bc8c2afddf9a8b5c32a9">w</a>()<span class="keyword"> const </span>{ <span class="keywordflow">return</span> this->derived().coeffs().coeff(3); }</div> <div class="line"><a name="l00067"></a><span class="lineno"> 67</span> </div> <div class="line"><a name="l00069"></a><span class="lineno"><a class="line" href="classEigen_1_1QuaternionBase.html#a5e6fc074bced85f298c76edac1af6cd9"> 69</a></span>  <span class="keyword">inline</span> <a class="code" href="classEigen_1_1RotationBase.html#a844358c46408e878e60c4026c52eb1e9">Scalar</a>& <a class="code" href="classEigen_1_1QuaternionBase.html#a5e6fc074bced85f298c76edac1af6cd9">x</a>() { <span class="keywordflow">return</span> this->derived().coeffs().coeffRef(0); }</div> <div class="line"><a name="l00071"></a><span class="lineno"><a class="line" href="classEigen_1_1QuaternionBase.html#ad948bb3bd6866bdca15149880c34be62"> 71</a></span>  <span class="keyword">inline</span> <a class="code" href="classEigen_1_1RotationBase.html#a844358c46408e878e60c4026c52eb1e9">Scalar</a>& <a class="code" href="classEigen_1_1QuaternionBase.html#ad948bb3bd6866bdca15149880c34be62">y</a>() { <span class="keywordflow">return</span> this->derived().coeffs().coeffRef(1); }</div> <div class="line"><a name="l00073"></a><span class="lineno"><a class="line" href="classEigen_1_1QuaternionBase.html#a1abcf9e222c8e5c38cf63f6821cd6480"> 73</a></span>  <span class="keyword">inline</span> <a class="code" href="classEigen_1_1RotationBase.html#a844358c46408e878e60c4026c52eb1e9">Scalar</a>& <a class="code" href="classEigen_1_1QuaternionBase.html#a1abcf9e222c8e5c38cf63f6821cd6480">z</a>() { <span class="keywordflow">return</span> this->derived().coeffs().coeffRef(2); }</div> <div class="line"><a name="l00075"></a><span class="lineno"><a class="line" href="classEigen_1_1QuaternionBase.html#a81f3d09a835c62dbde6a0824f5db517b"> 75</a></span>  <span class="keyword">inline</span> <a class="code" href="classEigen_1_1RotationBase.html#a844358c46408e878e60c4026c52eb1e9">Scalar</a>& <a class="code" href="classEigen_1_1QuaternionBase.html#a81f3d09a835c62dbde6a0824f5db517b">w</a>() { <span class="keywordflow">return</span> this->derived().coeffs().coeffRef(3); }</div> <div class="line"><a name="l00076"></a><span class="lineno"> 76</span> </div> <div class="line"><a name="l00078"></a><span class="lineno"><a class="line" href="classEigen_1_1QuaternionBase.html#a154f40d7fec7be85112f7e7c12e4900b"> 78</a></span>  <span class="keyword">inline</span> <span class="keyword">const</span> <a class="code" href="classEigen_1_1VectorBlock.html">VectorBlock<const Coefficients,3></a> <a class="code" href="classEigen_1_1QuaternionBase.html#a154f40d7fec7be85112f7e7c12e4900b">vec</a>()<span class="keyword"> const </span>{ <span class="keywordflow">return</span> <a class="code" href="classEigen_1_1QuaternionBase.html#a316774ca4ac3e8d56a7e24e5e0ad8cc4">coeffs</a>().template head<3>(); }</div> <div class="line"><a name="l00079"></a><span class="lineno"> 79</span> </div> <div class="line"><a name="l00081"></a><span class="lineno"><a class="line" href="classEigen_1_1QuaternionBase.html#a5a8ee232e99cf184e74c6a8dfab5d43e"> 81</a></span>  <span class="keyword">inline</span> <a class="code" href="classEigen_1_1VectorBlock.html">VectorBlock<Coefficients,3></a> <a class="code" href="classEigen_1_1QuaternionBase.html#a5a8ee232e99cf184e74c6a8dfab5d43e">vec</a>() { <span class="keywordflow">return</span> <a class="code" href="classEigen_1_1QuaternionBase.html#a316774ca4ac3e8d56a7e24e5e0ad8cc4">coeffs</a>().template head<3>(); }</div> <div class="line"><a name="l00082"></a><span class="lineno"> 82</span> </div> <div class="line"><a name="l00084"></a><span class="lineno"><a class="line" href="classEigen_1_1QuaternionBase.html#a316774ca4ac3e8d56a7e24e5e0ad8cc4"> 84</a></span>  <span class="keyword">inline</span> <span class="keyword">const</span> <span class="keyword">typename</span> internal::traits<Derived>::Coefficients& <a class="code" href="classEigen_1_1QuaternionBase.html#a316774ca4ac3e8d56a7e24e5e0ad8cc4">coeffs</a>()<span class="keyword"> const </span>{ <span class="keywordflow">return</span> derived().coeffs(); }</div> <div class="line"><a name="l00085"></a><span class="lineno"> 85</span> </div> <div class="line"><a name="l00087"></a><span class="lineno"><a class="line" href="classEigen_1_1QuaternionBase.html#a33e390da91b3c4f42c2e7330c0a5338b"> 87</a></span>  <span class="keyword">inline</span> <span class="keyword">typename</span> internal::traits<Derived>::Coefficients& <a class="code" href="classEigen_1_1QuaternionBase.html#a33e390da91b3c4f42c2e7330c0a5338b">coeffs</a>() { <span class="keywordflow">return</span> derived().coeffs(); }</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>  EIGEN_STRONG_INLINE <a class="code" href="classEigen_1_1QuaternionBase.html">QuaternionBase<Derived></a>& operator=(<span class="keyword">const</span> <a class="code" href="classEigen_1_1QuaternionBase.html">QuaternionBase<Derived></a>& other);</div> <div class="line"><a name="l00090"></a><span class="lineno"> 90</span>  <span class="keyword">template</span><<span class="keyword">class</span> OtherDerived> EIGEN_STRONG_INLINE Derived& operator=(<span class="keyword">const</span> <a class="code" href="classEigen_1_1QuaternionBase.html">QuaternionBase<OtherDerived></a>& other);</div> <div class="line"><a name="l00091"></a><span class="lineno"> 91</span> </div> <div class="line"><a name="l00092"></a><span class="lineno"> 92</span> <span class="comment">// disabled this copy operator as it is giving very strange compilation errors when compiling</span></div> <div class="line"><a name="l00093"></a><span class="lineno"> 93</span> <span class="comment">// test_stdvector with GCC 4.4.2. This looks like a GCC bug though, so feel free to re-enable it if it's</span></div> <div class="line"><a name="l00094"></a><span class="lineno"> 94</span> <span class="comment">// useful; however notice that we already have the templated operator= above and e.g. in MatrixBase</span></div> <div class="line"><a name="l00095"></a><span class="lineno"> 95</span> <span class="comment">// we didn't have to add, in addition to templated operator=, such a non-templated copy operator.</span></div> <div class="line"><a name="l00096"></a><span class="lineno"> 96</span> <span class="comment">// Derived& operator=(const QuaternionBase& other)</span></div> <div class="line"><a name="l00097"></a><span class="lineno"> 97</span> <span class="comment">// { return operator=<Derived>(other); }</span></div> <div class="line"><a name="l00098"></a><span class="lineno"> 98</span> </div> <div class="line"><a name="l00099"></a><span class="lineno"> 99</span>  Derived& operator=(<span class="keyword">const</span> <a class="code" href="classEigen_1_1QuaternionBase.html#ac8cbb66bbccf5f2c9596004a14129c26">AngleAxisType</a>& aa);</div> <div class="line"><a name="l00100"></a><span class="lineno"> 100</span>  <span class="keyword">template</span><<span class="keyword">class</span> OtherDerived> Derived& operator=(<span class="keyword">const</span> <a class="code" href="classEigen_1_1MatrixBase.html">MatrixBase<OtherDerived></a>& m);</div> <div class="line"><a name="l00101"></a><span class="lineno"> 101</span> </div> <div class="line"><a name="l00105"></a><span class="lineno"><a class="line" href="classEigen_1_1QuaternionBase.html#a480a9376ddbb4d075331a3cefca97bc7"> 105</a></span>  <span class="keyword">static</span> <span class="keyword">inline</span> <a class="code" href="classEigen_1_1Quaternion.html">Quaternion<Scalar></a> <a class="code" href="classEigen_1_1QuaternionBase.html#a480a9376ddbb4d075331a3cefca97bc7">Identity</a>() { <span class="keywordflow">return</span> <a class="code" href="classEigen_1_1Quaternion.html">Quaternion<Scalar></a>(1, 0, 0, 0); }</div> <div class="line"><a name="l00106"></a><span class="lineno"> 106</span> </div> <div class="line"><a name="l00109"></a><span class="lineno"><a class="line" href="classEigen_1_1QuaternionBase.html#aa07f38d31b32610f5296f81d438b60fb"> 109</a></span>  <span class="keyword">inline</span> <a class="code" href="classEigen_1_1QuaternionBase.html">QuaternionBase</a>& <a class="code" href="classEigen_1_1QuaternionBase.html#aa07f38d31b32610f5296f81d438b60fb">setIdentity</a>() { <a class="code" href="classEigen_1_1QuaternionBase.html#a316774ca4ac3e8d56a7e24e5e0ad8cc4">coeffs</a>() << 0, 0, 0, 1; <span class="keywordflow">return</span> *<span class="keyword">this</span>; }</div> <div class="line"><a name="l00110"></a><span class="lineno"> 110</span> </div> <div class="line"><a name="l00114"></a><span class="lineno"><a class="line" href="classEigen_1_1QuaternionBase.html#a8699d72c996ca6cb4673e810fe3a616c"> 114</a></span>  <span class="keyword">inline</span> <a class="code" href="classEigen_1_1RotationBase.html#a844358c46408e878e60c4026c52eb1e9">Scalar</a> <a class="code" href="classEigen_1_1QuaternionBase.html#a8699d72c996ca6cb4673e810fe3a616c">squaredNorm</a>()<span class="keyword"> const </span>{ <span class="keywordflow">return</span> <a class="code" href="classEigen_1_1QuaternionBase.html#a316774ca4ac3e8d56a7e24e5e0ad8cc4">coeffs</a>().squaredNorm(); }</div> <div class="line"><a name="l00115"></a><span class="lineno"> 115</span> </div> <div class="line"><a name="l00119"></a><span class="lineno"><a class="line" href="classEigen_1_1QuaternionBase.html#aa2bb33b2ed836350b02a6b9083bf8b5b"> 119</a></span>  <span class="keyword">inline</span> <a class="code" href="classEigen_1_1RotationBase.html#a844358c46408e878e60c4026c52eb1e9">Scalar</a> <a class="code" href="classEigen_1_1QuaternionBase.html#aa2bb33b2ed836350b02a6b9083bf8b5b">norm</a>()<span class="keyword"> const </span>{ <span class="keywordflow">return</span> <a class="code" href="classEigen_1_1QuaternionBase.html#a316774ca4ac3e8d56a7e24e5e0ad8cc4">coeffs</a>().norm(); }</div> <div class="line"><a name="l00120"></a><span class="lineno"> 120</span> </div> <div class="line"><a name="l00123"></a><span class="lineno"><a class="line" href="classEigen_1_1QuaternionBase.html#acd0de676568888d848beb97dcc53ae47"> 123</a></span>  <span class="keyword">inline</span> <span class="keywordtype">void</span> <a class="code" href="classEigen_1_1QuaternionBase.html#acd0de676568888d848beb97dcc53ae47">normalize</a>() { <a class="code" href="classEigen_1_1QuaternionBase.html#a316774ca4ac3e8d56a7e24e5e0ad8cc4">coeffs</a>().normalize(); }</div> <div class="line"><a name="l00126"></a><span class="lineno"><a class="line" href="classEigen_1_1QuaternionBase.html#a4b12380cc06df343c220d0745b50759e"> 126</a></span>  <span class="keyword">inline</span> <a class="code" href="classEigen_1_1Quaternion.html">Quaternion<Scalar></a> <a class="code" href="classEigen_1_1QuaternionBase.html#a4b12380cc06df343c220d0745b50759e">normalized</a>()<span class="keyword"> const </span>{ <span class="keywordflow">return</span> <a class="code" href="classEigen_1_1Quaternion.html">Quaternion<Scalar></a>(<a class="code" href="classEigen_1_1QuaternionBase.html#a316774ca4ac3e8d56a7e24e5e0ad8cc4">coeffs</a>().normalized()); }</div> <div class="line"><a name="l00127"></a><span class="lineno"> 127</span> </div> <div class="line"><a name="l00133"></a><span class="lineno"><a class="line" href="classEigen_1_1QuaternionBase.html#a20b9b376e143961250055470585ac1d0"> 133</a></span>  <span class="keyword">template</span><<span class="keyword">class</span> OtherDerived> <span class="keyword">inline</span> <a class="code" href="classEigen_1_1RotationBase.html#a844358c46408e878e60c4026c52eb1e9">Scalar</a> <a class="code" href="classEigen_1_1QuaternionBase.html#a20b9b376e143961250055470585ac1d0">dot</a>(<span class="keyword">const</span> <a class="code" href="classEigen_1_1QuaternionBase.html">QuaternionBase<OtherDerived></a>& other)<span class="keyword"> const </span>{ <span class="keywordflow">return</span> <a class="code" href="classEigen_1_1QuaternionBase.html#a316774ca4ac3e8d56a7e24e5e0ad8cc4">coeffs</a>().dot(other.coeffs()); }</div> <div class="line"><a name="l00134"></a><span class="lineno"> 134</span> </div> <div class="line"><a name="l00135"></a><span class="lineno"> 135</span>  <span class="keyword">template</span><<span class="keyword">class</span> OtherDerived> Scalar angularDistance(<span class="keyword">const</span> <a class="code" href="classEigen_1_1QuaternionBase.html">QuaternionBase<OtherDerived></a>& other) <span class="keyword">const</span>;</div> <div class="line"><a name="l00136"></a><span class="lineno"> 136</span> </div> <div class="line"><a name="l00138"></a><span class="lineno"> 138</span>  <a class="code" href="classEigen_1_1QuaternionBase.html#abdb0743d0fe840b927780b2bd78fb29e">Matrix3</a> <a class="code" href="classEigen_1_1QuaternionBase.html#a7bdc6d3cce92c44d32281aa2f85b56f7">toRotationMatrix</a>() <span class="keyword">const</span>;</div> <div class="line"><a name="l00139"></a><span class="lineno"> 139</span> </div> <div class="line"><a name="l00141"></a><span class="lineno"> 141</span>  <span class="keyword">template</span><<span class="keyword">typename</span> Derived1, <span class="keyword">typename</span> Derived2></div> <div class="line"><a name="l00142"></a><span class="lineno"> 142</span>  Derived& <a class="code" href="classEigen_1_1QuaternionBase.html#ac35460294d855096e9b687cadf821452">setFromTwoVectors</a>(<span class="keyword">const</span> <a class="code" href="classEigen_1_1MatrixBase.html">MatrixBase<Derived1></a>& a, <span class="keyword">const</span> <a class="code" href="classEigen_1_1MatrixBase.html">MatrixBase<Derived2></a>& b);</div> <div class="line"><a name="l00143"></a><span class="lineno"> 143</span> </div> <div class="line"><a name="l00144"></a><span class="lineno"> 144</span>  <span class="keyword">template</span><<span class="keyword">class</span> OtherDerived> EIGEN_STRONG_INLINE <a class="code" href="classEigen_1_1Quaternion.html">Quaternion<Scalar></a> operator* (<span class="keyword">const</span> <a class="code" href="classEigen_1_1QuaternionBase.html">QuaternionBase<OtherDerived></a>& q) <span class="keyword">const</span>;</div> <div class="line"><a name="l00145"></a><span class="lineno"> 145</span>  <span class="keyword">template</span><<span class="keyword">class</span> OtherDerived> EIGEN_STRONG_INLINE Derived& <a class="code" href="classEigen_1_1QuaternionBase.html#aa76d91cad6945c0d8f57d11289b3e8c2">operator*= </a>(<span class="keyword">const</span> <a class="code" href="classEigen_1_1QuaternionBase.html">QuaternionBase<OtherDerived></a>& q);</div> <div class="line"><a name="l00146"></a><span class="lineno"> 146</span> </div> <div class="line"><a name="l00148"></a><span class="lineno"> 148</span>  <a class="code" href="classEigen_1_1Quaternion.html">Quaternion<Scalar></a> <a class="code" href="classEigen_1_1QuaternionBase.html#a666e60218567ecf0ec8e2cd5040e1d87">inverse</a>() <span class="keyword">const</span>;</div> <div class="line"><a name="l00149"></a><span class="lineno"> 149</span> </div> <div class="line"><a name="l00151"></a><span class="lineno"> 151</span>  <a class="code" href="classEigen_1_1Quaternion.html">Quaternion<Scalar></a> <a class="code" href="classEigen_1_1QuaternionBase.html#a6761d13accfe34e38d907505221d081c">conjugate</a>() <span class="keyword">const</span>;</div> <div class="line"><a name="l00152"></a><span class="lineno"> 152</span> </div> <div class="line"><a name="l00157"></a><span class="lineno"> 157</span>  <span class="keyword">template</span><<span class="keyword">class</span> OtherDerived> <a class="code" href="classEigen_1_1Quaternion.html">Quaternion<Scalar></a> <a class="code" href="classEigen_1_1QuaternionBase.html#a8cf5be4be8788856ca8dfe97b3cde1c4">slerp</a>(<span class="keyword">const</span> Scalar& t, <span class="keyword">const</span> <a class="code" href="classEigen_1_1QuaternionBase.html">QuaternionBase<OtherDerived></a>& other) <span class="keyword">const</span>;</div> <div class="line"><a name="l00158"></a><span class="lineno"> 158</span> </div> <div class="line"><a name="l00163"></a><span class="lineno"> 163</span>  <span class="keyword">template</span><<span class="keyword">class</span> OtherDerived></div> <div class="line"><a name="l00164"></a><span class="lineno"><a class="line" href="classEigen_1_1QuaternionBase.html#a3956cb2fddca9905ac6ed71a42ee2d2a"> 164</a></span>  <span class="keywordtype">bool</span> <a class="code" href="classEigen_1_1QuaternionBase.html#a3956cb2fddca9905ac6ed71a42ee2d2a">isApprox</a>(<span class="keyword">const</span> <a class="code" href="classEigen_1_1QuaternionBase.html">QuaternionBase<OtherDerived></a>& other, <span class="keyword">const</span> RealScalar& 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="l00165"></a><span class="lineno"> 165</span> <span class="keyword"> </span>{ <span class="keywordflow">return</span> <a class="code" href="classEigen_1_1QuaternionBase.html#a316774ca4ac3e8d56a7e24e5e0ad8cc4">coeffs</a>().isApprox(other.coeffs(), prec); }</div> <div class="line"><a name="l00166"></a><span class="lineno"> 166</span> </div> <div class="line"><a name="l00168"></a><span class="lineno"> 168</span>  EIGEN_STRONG_INLINE <a class="code" href="classEigen_1_1QuaternionBase.html#ab05eeebdd63f50a52fecb6353018075f">Vector3</a> <a class="code" href="classEigen_1_1QuaternionBase.html#ab45aef6f20921b0997f4e8d75ef03300">_transformVector</a>(<a class="code" href="classEigen_1_1QuaternionBase.html#ab05eeebdd63f50a52fecb6353018075f">Vector3</a> v) <span class="keyword">const</span>;</div> <div class="line"><a name="l00169"></a><span class="lineno"> 169</span> </div> <div class="line"><a name="l00175"></a><span class="lineno"> 175</span>  <span class="keyword">template</span><<span class="keyword">typename</span> NewScalarType></div> <div class="line"><a name="l00176"></a><span class="lineno"><a class="line" href="classEigen_1_1QuaternionBase.html#a65f6d2422debd1a64a2452641497e2f4"> 176</a></span>  <span class="keyword">inline</span> <span class="keyword">typename</span> internal::cast_return_type<Derived,Quaternion<NewScalarType> >::type <a class="code" href="classEigen_1_1QuaternionBase.html#a65f6d2422debd1a64a2452641497e2f4">cast</a>()<span class="keyword"> const</span></div> <div class="line"><a name="l00177"></a><span class="lineno"> 177</span> <span class="keyword"> </span>{</div> <div class="line"><a name="l00178"></a><span class="lineno"> 178</span>  <span class="keywordflow">return</span> <span class="keyword">typename</span> internal::cast_return_type<Derived,Quaternion<NewScalarType> >::type(derived());</div> <div class="line"><a name="l00179"></a><span class="lineno"> 179</span>  }</div> <div class="line"><a name="l00180"></a><span class="lineno"> 180</span> </div> <div class="line"><a name="l00181"></a><span class="lineno"> 181</span> <span class="preprocessor">#ifdef EIGEN_QUATERNIONBASE_PLUGIN</span></div> <div class="line"><a name="l00182"></a><span class="lineno"> 182</span> <span class="preprocessor"></span><span class="preprocessor"># include EIGEN_QUATERNIONBASE_PLUGIN</span></div> <div class="line"><a name="l00183"></a><span class="lineno"> 183</span> <span class="preprocessor"></span><span class="preprocessor">#endif</span></div> <div class="line"><a name="l00184"></a><span class="lineno"> 184</span> <span class="preprocessor"></span>};</div> <div class="line"><a name="l00185"></a><span class="lineno"> 185</span> </div> <div class="line"><a name="l00186"></a><span class="lineno"> 186</span> <span class="comment">/***************************************************************************</span></div> <div class="line"><a name="l00187"></a><span class="lineno"> 187</span> <span class="comment">* Definition/implementation of Quaternion<Scalar></span></div> <div class="line"><a name="l00188"></a><span class="lineno"> 188</span> <span class="comment">***************************************************************************/</span></div> <div class="line"><a name="l00189"></a><span class="lineno"> 189</span> </div> <div class="line"><a name="l00213"></a><span class="lineno"> 213</span> <span class="keyword">namespace </span>internal {</div> <div class="line"><a name="l00214"></a><span class="lineno"> 214</span> <span class="keyword">template</span><<span class="keyword">typename</span> _Scalar,<span class="keywordtype">int</span> _Options></div> <div class="line"><a name="l00215"></a><span class="lineno"> 215</span> <span class="keyword">struct </span>traits<Quaternion<_Scalar,_Options> ></div> <div class="line"><a name="l00216"></a><span class="lineno"> 216</span> {</div> <div class="line"><a name="l00217"></a><span class="lineno"> 217</span>  <span class="keyword">typedef</span> Quaternion<_Scalar,_Options> PlainObject;</div> <div class="line"><a name="l00218"></a><span class="lineno"> 218</span>  <span class="keyword">typedef</span> _Scalar Scalar;</div> <div class="line"><a name="l00219"></a><span class="lineno"> 219</span>  <span class="keyword">typedef</span> Matrix<_Scalar,4,1,_Options> Coefficients;</div> <div class="line"><a name="l00220"></a><span class="lineno"> 220</span>  <span class="keyword">enum</span>{</div> <div class="line"><a name="l00221"></a><span class="lineno"> 221</span>  IsAligned = internal::traits<Coefficients>::Flags & <a class="code" href="group__flags.html#ga972a2dcb6603215fa53e0b9e82051426">AlignedBit</a>,</div> <div class="line"><a name="l00222"></a><span class="lineno"> 222</span>  Flags = IsAligned ? (AlignedBit | <a class="code" href="group__flags.html#ga64e21b7543bdedce27f013512a4403a3">LvalueBit</a>) : <a class="code" href="group__flags.html#ga64e21b7543bdedce27f013512a4403a3">LvalueBit</a></div> <div class="line"><a name="l00223"></a><span class="lineno"> 223</span>  };</div> <div class="line"><a name="l00224"></a><span class="lineno"> 224</span> };</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> </div> <div class="line"><a name="l00227"></a><span class="lineno"> 227</span> <span class="keyword">template</span><<span class="keyword">typename</span> _Scalar, <span class="keywordtype">int</span> _Options></div> <div class="line"><a name="l00228"></a><span class="lineno"> 228</span> <span class="keyword">class </span>Quaternion : <span class="keyword">public</span> QuaternionBase<Quaternion<_Scalar,_Options> ></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>  <span class="keyword">typedef</span> QuaternionBase<Quaternion<_Scalar,_Options> > Base;</div> <div class="line"><a name="l00231"></a><span class="lineno"> 231</span>  <span class="keyword">enum</span> { IsAligned = internal::traits<Quaternion>::IsAligned };</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="keyword">public</span>:</div> <div class="line"><a name="l00234"></a><span class="lineno"> 234</span>  <span class="keyword">typedef</span> _Scalar Scalar;</div> <div class="line"><a name="l00235"></a><span class="lineno"> 235</span> </div> <div class="line"><a name="l00236"></a><span class="lineno"> 236</span>  EIGEN_INHERIT_ASSIGNMENT_EQUAL_OPERATOR(Quaternion)</div> <div class="line"><a name="l00237"></a><span class="lineno"> 237</span>  using Base::operator*=;</div> <div class="line"><a name="l00238"></a><span class="lineno"> 238</span> </div> <div class="line"><a name="l00239"></a><span class="lineno"> 239</span>  typedef typename internal::traits<Quaternion>::Coefficients Coefficients;</div> <div class="line"><a name="l00240"></a><span class="lineno"> 240</span>  typedef typename Base::AngleAxisType AngleAxisType;</div> <div class="line"><a name="l00241"></a><span class="lineno"> 241</span> </div> <div class="line"><a name="l00243"></a><span class="lineno"><a class="line" href="classEigen_1_1Quaternion.html#a65ed15cc19af958b5933b5c522f10e66"> 243</a></span>  inline <a class="code" href="classEigen_1_1Quaternion.html">Quaternion</a>() {}</div> <div class="line"><a name="l00244"></a><span class="lineno"> 244</span> </div> <div class="line"><a name="l00252"></a><span class="lineno"><a class="line" href="classEigen_1_1Quaternion.html#a91b6ea2cac13ab2d33b6e74818ee1490"> 252</a></span>  <span class="keyword">inline</span> <a class="code" href="classEigen_1_1Quaternion.html#a91b6ea2cac13ab2d33b6e74818ee1490">Quaternion</a>(<span class="keyword">const</span> Scalar& w, <span class="keyword">const</span> Scalar& x, <span class="keyword">const</span> Scalar& y, <span class="keyword">const</span> Scalar& z) : m_coeffs(x, y, z, w){}</div> <div class="line"><a name="l00253"></a><span class="lineno"> 253</span> </div> <div class="line"><a name="l00255"></a><span class="lineno"><a class="line" href="classEigen_1_1Quaternion.html#afb78170e1d8b745832839d4585becf85"> 255</a></span>  <span class="keyword">inline</span> <a class="code" href="classEigen_1_1Quaternion.html#afb78170e1d8b745832839d4585becf85">Quaternion</a>(<span class="keyword">const</span> Scalar* data) : m_coeffs(data) {}</div> <div class="line"><a name="l00256"></a><span class="lineno"> 256</span> </div> <div class="line"><a name="l00258"></a><span class="lineno"><a class="line" href="classEigen_1_1Quaternion.html#a5f7c4a058b9ca7bf6fb5d40c7c047a4a"> 258</a></span>  <span class="keyword">template</span><<span class="keyword">class</span> Derived> EIGEN_STRONG_INLINE <a class="code" href="classEigen_1_1Quaternion.html#a5f7c4a058b9ca7bf6fb5d40c7c047a4a">Quaternion</a>(<span class="keyword">const</span> <a class="code" href="classEigen_1_1QuaternionBase.html">QuaternionBase<Derived></a>& other) { this->Base::operator=(other); }</div> <div class="line"><a name="l00259"></a><span class="lineno"> 259</span> </div> <div class="line"><a name="l00261"></a><span class="lineno"><a class="line" href="classEigen_1_1Quaternion.html#a344a7c038d0b89e798baf2baa3b4a592"> 261</a></span>  <span class="keyword">explicit</span> <span class="keyword">inline</span> <a class="code" href="classEigen_1_1Quaternion.html#a344a7c038d0b89e798baf2baa3b4a592">Quaternion</a>(<span class="keyword">const</span> <a class="code" href="classEigen_1_1AngleAxis.html">AngleAxisType</a>& aa) { *<span class="keyword">this</span> = aa; }</div> <div class="line"><a name="l00262"></a><span class="lineno"> 262</span> </div> <div class="line"><a name="l00267"></a><span class="lineno"> 267</span>  <span class="keyword">template</span><<span class="keyword">typename</span> Derived></div> <div class="line"><a name="l00268"></a><span class="lineno"><a class="line" href="classEigen_1_1Quaternion.html#ad90ae48f7378bb94dfbc6436e3a66aa2"> 268</a></span>  <span class="keyword">explicit</span> <span class="keyword">inline</span> <a class="code" href="classEigen_1_1Quaternion.html#ad90ae48f7378bb94dfbc6436e3a66aa2">Quaternion</a>(<span class="keyword">const</span> <a class="code" href="classEigen_1_1MatrixBase.html">MatrixBase<Derived></a>& other) { *<span class="keyword">this</span> = other; }</div> <div class="line"><a name="l00269"></a><span class="lineno"> 269</span> </div> <div class="line"><a name="l00271"></a><span class="lineno"> 271</span>  <span class="keyword">template</span><<span class="keyword">typename</span> OtherScalar, <span class="keywordtype">int</span> OtherOptions></div> <div class="line"><a name="l00272"></a><span class="lineno"><a class="line" href="classEigen_1_1Quaternion.html#a801ef1cd3194d3f9e3067c35d883ba4b"> 272</a></span>  <span class="keyword">explicit</span> <span class="keyword">inline</span> <a class="code" href="classEigen_1_1Quaternion.html#a801ef1cd3194d3f9e3067c35d883ba4b">Quaternion</a>(<span class="keyword">const</span> <a class="code" href="classEigen_1_1Quaternion.html">Quaternion<OtherScalar, OtherOptions></a>& other)</div> <div class="line"><a name="l00273"></a><span class="lineno"> 273</span>  { m_coeffs = other.coeffs().template cast<Scalar>(); }</div> <div class="line"><a name="l00274"></a><span class="lineno"> 274</span> </div> <div class="line"><a name="l00275"></a><span class="lineno"> 275</span>  <span class="keyword">template</span><<span class="keyword">typename</span> Derived1, <span class="keyword">typename</span> Derived2></div> <div class="line"><a name="l00276"></a><span class="lineno"> 276</span>  <span class="keyword">static</span> <a class="code" href="classEigen_1_1Quaternion.html">Quaternion</a> FromTwoVectors(<span class="keyword">const</span> <a class="code" href="classEigen_1_1MatrixBase.html">MatrixBase<Derived1></a>& a, <span class="keyword">const</span> <a class="code" href="classEigen_1_1MatrixBase.html">MatrixBase<Derived2></a>& b);</div> <div class="line"><a name="l00277"></a><span class="lineno"> 277</span> </div> <div class="line"><a name="l00278"></a><span class="lineno"> 278</span>  <span class="keyword">inline</span> Coefficients& coeffs() { <span class="keywordflow">return</span> m_coeffs;}</div> <div class="line"><a name="l00279"></a><span class="lineno"> 279</span>  <span class="keyword">inline</span> <span class="keyword">const</span> Coefficients& coeffs()<span class="keyword"> const </span>{ <span class="keywordflow">return</span> m_coeffs;}</div> <div class="line"><a name="l00280"></a><span class="lineno"> 280</span> </div> <div class="line"><a name="l00281"></a><span class="lineno"> 281</span>  EIGEN_MAKE_ALIGNED_OPERATOR_NEW_IF(IsAligned)</div> <div class="line"><a name="l00282"></a><span class="lineno"> 282</span> </div> <div class="line"><a name="l00283"></a><span class="lineno"> 283</span> protected:</div> <div class="line"><a name="l00284"></a><span class="lineno"> 284</span>  Coefficients m_coeffs;</div> <div class="line"><a name="l00285"></a><span class="lineno"> 285</span>  </div> <div class="line"><a name="l00286"></a><span class="lineno"> 286</span> <span class="preprocessor">#ifndef EIGEN_PARSED_BY_DOXYGEN</span></div> <div class="line"><a name="l00287"></a><span class="lineno"> 287</span> <span class="preprocessor"></span> <span class="keyword">static</span> EIGEN_STRONG_INLINE <span class="keywordtype">void</span> _check_template_params()</div> <div class="line"><a name="l00288"></a><span class="lineno"> 288</span>  {</div> <div class="line"><a name="l00289"></a><span class="lineno"> 289</span>  EIGEN_STATIC_ASSERT( (_Options & <a class="code" href="group__enums.html#gga0c5bde183ecefe103f70b49ad9740bcda761c0504a518c7450ed6dfe9eaeab8a6">DontAlign</a>) == _Options,</div> <div class="line"><a name="l00290"></a><span class="lineno"> 290</span>  INVALID_MATRIX_TEMPLATE_PARAMETERS)</div> <div class="line"><a name="l00291"></a><span class="lineno"> 291</span>  }</div> <div class="line"><a name="l00292"></a><span class="lineno"> 292</span> <span class="preprocessor">#endif</span></div> <div class="line"><a name="l00293"></a><span class="lineno"> 293</span> <span class="preprocessor"></span>};</div> <div class="line"><a name="l00294"></a><span class="lineno"> 294</span> </div> <div class="line"><a name="l00297"></a><span class="lineno"><a class="line" href="group__Geometry__Module.html#gaf65cf6f803890e57488d7de750bef682"> 297</a></span> <span class="keyword">typedef</span> <a class="code" href="classEigen_1_1Quaternion.html">Quaternion<float></a> <a class="code" href="group__Geometry__Module.html#gaf65cf6f803890e57488d7de750bef682">Quaternionf</a>;</div> <div class="line"><a name="l00300"></a><span class="lineno"><a class="line" href="group__Geometry__Module.html#ga0d2bd45f1215359f8e7c0d7ab53c4acb"> 300</a></span> <span class="keyword">typedef</span> <a class="code" href="classEigen_1_1Quaternion.html">Quaternion<double></a> <a class="code" href="group__Geometry__Module.html#ga0d2bd45f1215359f8e7c0d7ab53c4acb">Quaterniond</a>;</div> <div class="line"><a name="l00301"></a><span class="lineno"> 301</span> </div> <div class="line"><a name="l00302"></a><span class="lineno"> 302</span> <span class="comment">/***************************************************************************</span></div> <div class="line"><a name="l00303"></a><span class="lineno"> 303</span> <span class="comment">* Specialization of Map<Quaternion<Scalar>></span></div> <div class="line"><a name="l00304"></a><span class="lineno"> 304</span> <span class="comment">***************************************************************************/</span></div> <div class="line"><a name="l00305"></a><span class="lineno"> 305</span> </div> <div class="line"><a name="l00306"></a><span class="lineno"> 306</span> <span class="keyword">namespace </span>internal {</div> <div class="line"><a name="l00307"></a><span class="lineno"> 307</span>  <span class="keyword">template</span><<span class="keyword">typename</span> _Scalar, <span class="keywordtype">int</span> _Options></div> <div class="line"><a name="l00308"></a><span class="lineno"> 308</span>  <span class="keyword">struct </span>traits<Map<Quaternion<_Scalar>, _Options> > : traits<Quaternion<_Scalar, (int(_Options)&Aligned)==Aligned ? AutoAlign : DontAlign> ></div> <div class="line"><a name="l00309"></a><span class="lineno"> 309</span>  {</div> <div class="line"><a name="l00310"></a><span class="lineno"> 310</span>  <span class="keyword">typedef</span> <a class="code" href="classEigen_1_1Map.html">Map<Matrix<_Scalar,4,1></a>, _Options> Coefficients;</div> <div class="line"><a name="l00311"></a><span class="lineno"> 311</span>  };</div> <div class="line"><a name="l00312"></a><span class="lineno"> 312</span> }</div> <div class="line"><a name="l00313"></a><span class="lineno"> 313</span> </div> <div class="line"><a name="l00314"></a><span class="lineno"> 314</span> <span class="keyword">namespace </span>internal {</div> <div class="line"><a name="l00315"></a><span class="lineno"> 315</span>  <span class="keyword">template</span><<span class="keyword">typename</span> _Scalar, <span class="keywordtype">int</span> _Options></div> <div class="line"><a name="l00316"></a><span class="lineno"> 316</span>  <span class="keyword">struct </span>traits<Map<const Quaternion<_Scalar>, _Options> > : traits<Quaternion<_Scalar, (int(_Options)&Aligned)==Aligned ? AutoAlign : DontAlign> ></div> <div class="line"><a name="l00317"></a><span class="lineno"> 317</span>  {</div> <div class="line"><a name="l00318"></a><span class="lineno"> 318</span>  <span class="keyword">typedef</span> Map<const Matrix<_Scalar,4,1>, _Options> Coefficients;</div> <div class="line"><a name="l00319"></a><span class="lineno"> 319</span>  <span class="keyword">typedef</span> traits<Quaternion<_Scalar, (int(_Options)&Aligned)==Aligned ? AutoAlign : DontAlign> > TraitsBase;</div> <div class="line"><a name="l00320"></a><span class="lineno"> 320</span>  <span class="keyword">enum</span> {</div> <div class="line"><a name="l00321"></a><span class="lineno"> 321</span>  Flags = TraitsBase::Flags & ~<a class="code" href="group__flags.html#ga64e21b7543bdedce27f013512a4403a3">LvalueBit</a></div> <div class="line"><a name="l00322"></a><span class="lineno"> 322</span>  };</div> <div class="line"><a name="l00323"></a><span class="lineno"> 323</span>  };</div> <div class="line"><a name="l00324"></a><span class="lineno"> 324</span> }</div> <div class="line"><a name="l00325"></a><span class="lineno"> 325</span> </div> <div class="line"><a name="l00337"></a><span class="lineno"> 337</span> <span class="keyword">template</span><<span class="keyword">typename</span> _Scalar, <span class="keywordtype">int</span> _Options></div> <div class="line"><a name="l00338"></a><span class="lineno"><a class="line" href="classEigen_1_1Map_3_01const_01Quaternion_3_01__Scalar_01_4_00_01__Options_01_4.html"> 338</a></span> <span class="keyword">class </span><a class="code" href="classEigen_1_1Map.html">Map</a><const <a class="code" href="classEigen_1_1Quaternion.html">Quaternion</a><_Scalar>, _Options ></div> <div class="line"><a name="l00339"></a><span class="lineno"> 339</span>  : <span class="keyword">public</span> <a class="code" href="classEigen_1_1QuaternionBase.html">QuaternionBase</a><Map<const Quaternion<_Scalar>, _Options> ></div> <div class="line"><a name="l00340"></a><span class="lineno"> 340</span> {</div> <div class="line"><a name="l00341"></a><span class="lineno"> 341</span>  <span class="keyword">typedef</span> <a class="code" href="classEigen_1_1QuaternionBase.html">QuaternionBase<Map<const Quaternion<_Scalar></a>, _Options> > <a class="code" href="classEigen_1_1QuaternionBase.html">Base</a>;</div> <div class="line"><a name="l00342"></a><span class="lineno"> 342</span> </div> <div class="line"><a name="l00343"></a><span class="lineno"> 343</span>  <span class="keyword">public</span>:</div> <div class="line"><a name="l00344"></a><span class="lineno"> 344</span>  <span class="keyword">typedef</span> _Scalar Scalar;</div> <div class="line"><a name="l00345"></a><span class="lineno"> 345</span>  <span class="keyword">typedef</span> <span class="keyword">typename</span> internal::traits<Map>::Coefficients Coefficients;</div> <div class="line"><a name="l00346"></a><span class="lineno"> 346</span>  EIGEN_INHERIT_ASSIGNMENT_EQUAL_OPERATOR(<a class="code" href="classEigen_1_1Map.html">Map</a>)</div> <div class="line"><a name="l00347"></a><span class="lineno"> 347</span>  <span class="keyword">using</span> Base::operator*=;</div> <div class="line"><a name="l00348"></a><span class="lineno"> 348</span> </div> <div class="line"><a name="l00355"></a><span class="lineno"><a class="line" href="classEigen_1_1Map_3_01const_01Quaternion_3_01__Scalar_01_4_00_01__Options_01_4.html#ac040b01854bd9a68cf897936e8ba1efb"> 355</a></span>  EIGEN_STRONG_INLINE <a class="code" href="classEigen_1_1Map_3_01const_01Quaternion_3_01__Scalar_01_4_00_01__Options_01_4.html#ac040b01854bd9a68cf897936e8ba1efb">Map</a>(<span class="keyword">const</span> Scalar* coeffs) : m_coeffs(coeffs) {}</div> <div class="line"><a name="l00356"></a><span class="lineno"> 356</span> </div> <div class="line"><a name="l00357"></a><span class="lineno"> 357</span>  <span class="keyword">inline</span> <span class="keyword">const</span> Coefficients& coeffs()<span class="keyword"> const </span>{ <span class="keywordflow">return</span> m_coeffs;}</div> <div class="line"><a name="l00358"></a><span class="lineno"> 358</span> </div> <div class="line"><a name="l00359"></a><span class="lineno"> 359</span>  <span class="keyword">protected</span>:</div> <div class="line"><a name="l00360"></a><span class="lineno"> 360</span>  <span class="keyword">const</span> Coefficients m_coeffs;</div> <div class="line"><a name="l00361"></a><span class="lineno"> 361</span> };</div> <div class="line"><a name="l00362"></a><span class="lineno"> 362</span> </div> <div class="line"><a name="l00374"></a><span class="lineno"> 374</span> <span class="keyword">template</span><<span class="keyword">typename</span> _Scalar, <span class="keywordtype">int</span> _Options></div> <div class="line"><a name="l00375"></a><span class="lineno"><a class="line" href="classEigen_1_1Map_3_01Quaternion_3_01__Scalar_01_4_00_01__Options_01_4.html"> 375</a></span> <span class="keyword">class </span><a class="code" href="classEigen_1_1Map.html">Map</a><<a class="code" href="classEigen_1_1Quaternion.html">Quaternion</a><_Scalar>, _Options ></div> <div class="line"><a name="l00376"></a><span class="lineno"> 376</span>  : <span class="keyword">public</span> <a class="code" href="classEigen_1_1QuaternionBase.html">QuaternionBase</a><Map<Quaternion<_Scalar>, _Options> ></div> <div class="line"><a name="l00377"></a><span class="lineno"> 377</span> {</div> <div class="line"><a name="l00378"></a><span class="lineno"> 378</span>  <span class="keyword">typedef</span> <a class="code" href="classEigen_1_1QuaternionBase.html">QuaternionBase<Map<Quaternion<_Scalar></a>, _Options> > <a class="code" href="classEigen_1_1QuaternionBase.html">Base</a>;</div> <div class="line"><a name="l00379"></a><span class="lineno"> 379</span> </div> <div class="line"><a name="l00380"></a><span class="lineno"> 380</span>  <span class="keyword">public</span>:</div> <div class="line"><a name="l00381"></a><span class="lineno"> 381</span>  <span class="keyword">typedef</span> _Scalar Scalar;</div> <div class="line"><a name="l00382"></a><span class="lineno"> 382</span>  <span class="keyword">typedef</span> <span class="keyword">typename</span> internal::traits<Map>::Coefficients Coefficients;</div> <div class="line"><a name="l00383"></a><span class="lineno"> 383</span>  EIGEN_INHERIT_ASSIGNMENT_EQUAL_OPERATOR(<a class="code" href="classEigen_1_1Map.html">Map</a>)</div> <div class="line"><a name="l00384"></a><span class="lineno"> 384</span>  <span class="keyword">using</span> Base::operator*=;</div> <div class="line"><a name="l00385"></a><span class="lineno"> 385</span> </div> <div class="line"><a name="l00392"></a><span class="lineno"><a class="line" href="classEigen_1_1Map_3_01Quaternion_3_01__Scalar_01_4_00_01__Options_01_4.html#aa5742a3a943c3adca0ea1f74f3339cb7"> 392</a></span>  EIGEN_STRONG_INLINE <a class="code" href="classEigen_1_1Map_3_01Quaternion_3_01__Scalar_01_4_00_01__Options_01_4.html#aa5742a3a943c3adca0ea1f74f3339cb7">Map</a>(Scalar* coeffs) : m_coeffs(coeffs) {}</div> <div class="line"><a name="l00393"></a><span class="lineno"> 393</span> </div> <div class="line"><a name="l00394"></a><span class="lineno"> 394</span>  <span class="keyword">inline</span> Coefficients& coeffs() { <span class="keywordflow">return</span> m_coeffs; }</div> <div class="line"><a name="l00395"></a><span class="lineno"> 395</span>  <span class="keyword">inline</span> <span class="keyword">const</span> Coefficients& coeffs()<span class="keyword"> const </span>{ <span class="keywordflow">return</span> m_coeffs; }</div> <div class="line"><a name="l00396"></a><span class="lineno"> 396</span> </div> <div class="line"><a name="l00397"></a><span class="lineno"> 397</span>  <span class="keyword">protected</span>:</div> <div class="line"><a name="l00398"></a><span class="lineno"> 398</span>  Coefficients m_coeffs;</div> <div class="line"><a name="l00399"></a><span class="lineno"> 399</span> };</div> <div class="line"><a name="l00400"></a><span class="lineno"> 400</span> </div> <div class="line"><a name="l00403"></a><span class="lineno"><a class="line" href="group__Geometry__Module.html#ga8d125241c02e63c656b75b38b2382816"> 403</a></span> <span class="keyword">typedef</span> <a class="code" href="classEigen_1_1Map.html">Map<Quaternion<float></a>, 0> <a class="code" href="group__Geometry__Module.html#ga8d125241c02e63c656b75b38b2382816">QuaternionMapf</a>;</div> <div class="line"><a name="l00406"></a><span class="lineno"><a class="line" href="group__Geometry__Module.html#gaca6aed0c662d8272c53663c0093aacaa"> 406</a></span> <span class="keyword">typedef</span> <a class="code" href="classEigen_1_1Map.html">Map<Quaternion<double></a>, 0> <a class="code" href="group__Geometry__Module.html#gaca6aed0c662d8272c53663c0093aacaa">QuaternionMapd</a>;</div> <div class="line"><a name="l00409"></a><span class="lineno"><a class="line" href="group__Geometry__Module.html#ga987bcae344d68a98892d0e9404138404"> 409</a></span> <span class="keyword">typedef</span> <a class="code" href="classEigen_1_1Map.html">Map<Quaternion<float></a>, <a class="code" href="group__enums.html#gga456ac33d49271d3e2c371351cd1d6371ad5380ca00f3d74b38593adf8a0d06d3e">Aligned</a>> <a class="code" href="group__Geometry__Module.html#ga987bcae344d68a98892d0e9404138404">QuaternionMapAlignedf</a>;</div> <div class="line"><a name="l00412"></a><span class="lineno"><a class="line" href="group__Geometry__Module.html#gabfe832b0fc8c4c8b4dcf8b59ebe25603"> 412</a></span> <span class="keyword">typedef</span> <a class="code" href="classEigen_1_1Map.html">Map<Quaternion<double></a>, <a class="code" href="group__enums.html#gga456ac33d49271d3e2c371351cd1d6371ad5380ca00f3d74b38593adf8a0d06d3e">Aligned</a>> <a class="code" href="group__Geometry__Module.html#gabfe832b0fc8c4c8b4dcf8b59ebe25603">QuaternionMapAlignedd</a>;</div> <div class="line"><a name="l00413"></a><span class="lineno"> 413</span> </div> <div class="line"><a name="l00414"></a><span class="lineno"> 414</span> <span class="comment">/***************************************************************************</span></div> <div class="line"><a name="l00415"></a><span class="lineno"> 415</span> <span class="comment">* Implementation of QuaternionBase methods</span></div> <div class="line"><a name="l00416"></a><span class="lineno"> 416</span> <span class="comment">***************************************************************************/</span></div> <div class="line"><a name="l00417"></a><span class="lineno"> 417</span> </div> <div class="line"><a name="l00418"></a><span class="lineno"> 418</span> <span class="comment">// Generic Quaternion * Quaternion product</span></div> <div class="line"><a name="l00419"></a><span class="lineno"> 419</span> <span class="comment">// This product can be specialized for a given architecture via the Arch template argument.</span></div> <div class="line"><a name="l00420"></a><span class="lineno"> 420</span> <span class="keyword">namespace </span>internal {</div> <div class="line"><a name="l00421"></a><span class="lineno"> 421</span> <span class="keyword">template</span><<span class="keywordtype">int</span> Arch, <span class="keyword">class</span> Derived1, <span class="keyword">class</span> Derived2, <span class="keyword">typename</span> Scalar, <span class="keywordtype">int</span> _Options> <span class="keyword">struct </span>quat_product</div> <div class="line"><a name="l00422"></a><span class="lineno"> 422</span> {</div> <div class="line"><a name="l00423"></a><span class="lineno"> 423</span>  <span class="keyword">static</span> EIGEN_STRONG_INLINE <a class="code" href="classEigen_1_1Quaternion.html">Quaternion<Scalar></a> run(<span class="keyword">const</span> <a class="code" href="classEigen_1_1QuaternionBase.html">QuaternionBase<Derived1></a>& a, <span class="keyword">const</span> <a class="code" href="classEigen_1_1QuaternionBase.html">QuaternionBase<Derived2></a>& b){</div> <div class="line"><a name="l00424"></a><span class="lineno"> 424</span>  <span class="keywordflow">return</span> <a class="code" href="classEigen_1_1Quaternion.html">Quaternion<Scalar></a></div> <div class="line"><a name="l00425"></a><span class="lineno"> 425</span>  (</div> <div class="line"><a name="l00426"></a><span class="lineno"> 426</span>  a.w() * b.w() - a.x() * b.x() - a.y() * b.y() - a.z() * b.z(),</div> <div class="line"><a name="l00427"></a><span class="lineno"> 427</span>  a.w() * b.x() + a.x() * b.w() + a.y() * b.z() - a.z() * b.y(),</div> <div class="line"><a name="l00428"></a><span class="lineno"> 428</span>  a.w() * b.y() + a.y() * b.w() + a.z() * b.x() - a.x() * b.z(),</div> <div class="line"><a name="l00429"></a><span class="lineno"> 429</span>  a.w() * b.z() + a.z() * b.w() + a.x() * b.y() - a.y() * b.x()</div> <div class="line"><a name="l00430"></a><span class="lineno"> 430</span>  );</div> <div class="line"><a name="l00431"></a><span class="lineno"> 431</span>  }</div> <div class="line"><a name="l00432"></a><span class="lineno"> 432</span> };</div> <div class="line"><a name="l00433"></a><span class="lineno"> 433</span> }</div> <div class="line"><a name="l00434"></a><span class="lineno"> 434</span> </div> <div class="line"><a name="l00436"></a><span class="lineno"> 436</span> <span class="keyword">template</span> <<span class="keyword">class</span> Derived></div> <div class="line"><a name="l00437"></a><span class="lineno"> 437</span> <span class="keyword">template</span> <<span class="keyword">class</span> OtherDerived></div> <div class="line"><a name="l00438"></a><span class="lineno"> 438</span> EIGEN_STRONG_INLINE Quaternion<typename internal::traits<Derived>::Scalar></div> <div class="line"><a name="l00439"></a><span class="lineno"><a class="line" href="classEigen_1_1QuaternionBase.html#abcbdbe1d235761acf941c7f1d29704d2"> 439</a></span> <a class="code" href="classEigen_1_1QuaternionBase.html">QuaternionBase<Derived>::operator* </a>(<span class="keyword">const</span> <a class="code" href="classEigen_1_1QuaternionBase.html">QuaternionBase<OtherDerived></a>& other)<span class="keyword"> const</span></div> <div class="line"><a name="l00440"></a><span class="lineno"> 440</span> <span class="keyword"></span>{</div> <div class="line"><a name="l00441"></a><span class="lineno"> 441</span>  EIGEN_STATIC_ASSERT((internal::is_same<typename Derived::Scalar, typename OtherDerived::Scalar>::value),</div> <div class="line"><a name="l00442"></a><span class="lineno"> 442</span>  YOU_MIXED_DIFFERENT_NUMERIC_TYPES__YOU_NEED_TO_USE_THE_CAST_METHOD_OF_MATRIXBASE_TO_CAST_NUMERIC_TYPES_EXPLICITLY)</div> <div class="line"><a name="l00443"></a><span class="lineno"> 443</span>  <span class="keywordflow">return</span> internal::quat_product<Architecture::Target, Derived, OtherDerived,</div> <div class="line"><a name="l00444"></a><span class="lineno"> 444</span>  <span class="keyword">typename</span> internal::traits<Derived>::Scalar,</div> <div class="line"><a name="l00445"></a><span class="lineno"> 445</span>  internal::traits<Derived>::IsAligned && internal::traits<OtherDerived>::IsAligned>::run(*<span class="keyword">this</span>, other);</div> <div class="line"><a name="l00446"></a><span class="lineno"> 446</span> }</div> <div class="line"><a name="l00447"></a><span class="lineno"> 447</span> </div> <div class="line"><a name="l00449"></a><span class="lineno"> 449</span> <span class="keyword">template</span> <<span class="keyword">class</span> Derived></div> <div class="line"><a name="l00450"></a><span class="lineno"> 450</span> <span class="keyword">template</span> <<span class="keyword">class</span> OtherDerived></div> <div class="line"><a name="l00451"></a><span class="lineno"><a class="line" href="classEigen_1_1QuaternionBase.html#aa76d91cad6945c0d8f57d11289b3e8c2"> 451</a></span> EIGEN_STRONG_INLINE Derived& <a class="code" href="classEigen_1_1QuaternionBase.html#aa76d91cad6945c0d8f57d11289b3e8c2">QuaternionBase<Derived>::operator*= </a>(<span class="keyword">const</span> <a class="code" href="classEigen_1_1QuaternionBase.html">QuaternionBase<OtherDerived></a>& other)</div> <div class="line"><a name="l00452"></a><span class="lineno"> 452</span> {</div> <div class="line"><a name="l00453"></a><span class="lineno"> 453</span>  derived() = derived() * other.derived();</div> <div class="line"><a name="l00454"></a><span class="lineno"> 454</span>  <span class="keywordflow">return</span> derived();</div> <div class="line"><a name="l00455"></a><span class="lineno"> 455</span> }</div> <div class="line"><a name="l00456"></a><span class="lineno"> 456</span> </div> <div class="line"><a name="l00464"></a><span class="lineno"> 464</span> <span class="keyword">template</span> <<span class="keyword">class</span> Derived></div> <div class="line"><a name="l00465"></a><span class="lineno"> 465</span> EIGEN_STRONG_INLINE <span class="keyword">typename</span> <a class="code" href="classEigen_1_1Matrix.html">QuaternionBase<Derived>::Vector3</a></div> <div class="line"><a name="l00466"></a><span class="lineno"><a class="line" href="classEigen_1_1QuaternionBase.html#ab45aef6f20921b0997f4e8d75ef03300"> 466</a></span> <a class="code" href="classEigen_1_1QuaternionBase.html#ab45aef6f20921b0997f4e8d75ef03300">QuaternionBase<Derived>::_transformVector</a>(<a class="code" href="classEigen_1_1Matrix.html">Vector3</a> v)<span class="keyword"> const</span></div> <div class="line"><a name="l00467"></a><span class="lineno"> 467</span> <span class="keyword"></span>{</div> <div class="line"><a name="l00468"></a><span class="lineno"> 468</span>  <span class="comment">// Note that this algorithm comes from the optimization by hand</span></div> <div class="line"><a name="l00469"></a><span class="lineno"> 469</span>  <span class="comment">// of the conversion to a Matrix followed by a Matrix/Vector product.</span></div> <div class="line"><a name="l00470"></a><span class="lineno"> 470</span>  <span class="comment">// It appears to be much faster than the common algorithm found</span></div> <div class="line"><a name="l00471"></a><span class="lineno"> 471</span>  <span class="comment">// in the litterature (30 versus 39 flops). It also requires two</span></div> <div class="line"><a name="l00472"></a><span class="lineno"> 472</span>  <span class="comment">// Vector3 as temporaries.</span></div> <div class="line"><a name="l00473"></a><span class="lineno"> 473</span>  <a class="code" href="classEigen_1_1Matrix.html">Vector3</a> uv = this->vec().cross(v);</div> <div class="line"><a name="l00474"></a><span class="lineno"> 474</span>  uv += uv;</div> <div class="line"><a name="l00475"></a><span class="lineno"> 475</span>  <span class="keywordflow">return</span> v + this->w() * uv + this->vec().cross(uv);</div> <div class="line"><a name="l00476"></a><span class="lineno"> 476</span> }</div> <div class="line"><a name="l00477"></a><span class="lineno"> 477</span> </div> <div class="line"><a name="l00478"></a><span class="lineno"> 478</span> <span class="keyword">template</span><<span class="keyword">class</span> Derived></div> <div class="line"><a name="l00479"></a><span class="lineno"> 479</span> EIGEN_STRONG_INLINE <a class="code" href="classEigen_1_1QuaternionBase.html">QuaternionBase<Derived></a>& <a class="code" href="classEigen_1_1QuaternionBase.html">QuaternionBase<Derived>::operator=</a>(<span class="keyword">const</span> <a class="code" href="classEigen_1_1QuaternionBase.html">QuaternionBase<Derived></a>& other)</div> <div class="line"><a name="l00480"></a><span class="lineno"> 480</span> {</div> <div class="line"><a name="l00481"></a><span class="lineno"> 481</span>  coeffs() = other.coeffs();</div> <div class="line"><a name="l00482"></a><span class="lineno"> 482</span>  <span class="keywordflow">return</span> derived();</div> <div class="line"><a name="l00483"></a><span class="lineno"> 483</span> }</div> <div class="line"><a name="l00484"></a><span class="lineno"> 484</span> </div> <div class="line"><a name="l00485"></a><span class="lineno"> 485</span> <span class="keyword">template</span><<span class="keyword">class</span> Derived></div> <div class="line"><a name="l00486"></a><span class="lineno"> 486</span> <span class="keyword">template</span><<span class="keyword">class</span> OtherDerived></div> <div class="line"><a name="l00487"></a><span class="lineno"> 487</span> EIGEN_STRONG_INLINE Derived& QuaternionBase<Derived>::operator=(<span class="keyword">const</span> QuaternionBase<OtherDerived>& other)</div> <div class="line"><a name="l00488"></a><span class="lineno"> 488</span> {</div> <div class="line"><a name="l00489"></a><span class="lineno"> 489</span>  coeffs() = other.coeffs();</div> <div class="line"><a name="l00490"></a><span class="lineno"> 490</span>  <span class="keywordflow">return</span> derived();</div> <div class="line"><a name="l00491"></a><span class="lineno"> 491</span> }</div> <div class="line"><a name="l00492"></a><span class="lineno"> 492</span> </div> <div class="line"><a name="l00495"></a><span class="lineno"> 495</span> <span class="keyword">template</span><<span class="keyword">class</span> Derived></div> <div class="line"><a name="l00496"></a><span class="lineno"><a class="line" href="classEigen_1_1QuaternionBase.html#a8ffbc66c68cc438bbf5cb0caa50598f5"> 496</a></span> EIGEN_STRONG_INLINE Derived& <a class="code" href="classEigen_1_1QuaternionBase.html">QuaternionBase<Derived>::operator=</a>(<span class="keyword">const</span> <a class="code" href="classEigen_1_1AngleAxis.html">AngleAxisType</a>& aa)</div> <div class="line"><a name="l00497"></a><span class="lineno"> 497</span> {</div> <div class="line"><a name="l00498"></a><span class="lineno"> 498</span>  <span class="keyword">using</span> std::cos;</div> <div class="line"><a name="l00499"></a><span class="lineno"> 499</span>  <span class="keyword">using</span> std::sin;</div> <div class="line"><a name="l00500"></a><span class="lineno"> 500</span>  Scalar ha = Scalar(0.5)*aa.angle(); <span class="comment">// Scalar(0.5) to suppress precision loss warnings</span></div> <div class="line"><a name="l00501"></a><span class="lineno"> 501</span>  this->w() = cos(ha);</div> <div class="line"><a name="l00502"></a><span class="lineno"> 502</span>  this->vec() = sin(ha) * aa.axis();</div> <div class="line"><a name="l00503"></a><span class="lineno"> 503</span>  <span class="keywordflow">return</span> derived();</div> <div class="line"><a name="l00504"></a><span class="lineno"> 504</span> }</div> <div class="line"><a name="l00505"></a><span class="lineno"> 505</span> </div> <div class="line"><a name="l00512"></a><span class="lineno"> 512</span> <span class="keyword">template</span><<span class="keyword">class</span> Derived></div> <div class="line"><a name="l00513"></a><span class="lineno"> 513</span> <span class="keyword">template</span><<span class="keyword">class</span> MatrixDerived></div> <div class="line"><a name="l00514"></a><span class="lineno"><a class="line" href="classEigen_1_1QuaternionBase.html#a320e5057b6f2d1dcdd4d6438c0c2036c"> 514</a></span> <span class="keyword">inline</span> Derived& <a class="code" href="classEigen_1_1QuaternionBase.html">QuaternionBase<Derived>::operator=</a>(<span class="keyword">const</span> <a class="code" href="classEigen_1_1MatrixBase.html">MatrixBase<MatrixDerived></a>& xpr)</div> <div class="line"><a name="l00515"></a><span class="lineno"> 515</span> {</div> <div class="line"><a name="l00516"></a><span class="lineno"> 516</span>  EIGEN_STATIC_ASSERT((internal::is_same<typename Derived::Scalar, typename MatrixDerived::Scalar>::value),</div> <div class="line"><a name="l00517"></a><span class="lineno"> 517</span>  YOU_MIXED_DIFFERENT_NUMERIC_TYPES__YOU_NEED_TO_USE_THE_CAST_METHOD_OF_MATRIXBASE_TO_CAST_NUMERIC_TYPES_EXPLICITLY)</div> <div class="line"><a name="l00518"></a><span class="lineno"> 518</span>  internal::quaternionbase_assign_impl<MatrixDerived>::run(*<span class="keyword">this</span>, xpr.derived());</div> <div class="line"><a name="l00519"></a><span class="lineno"> 519</span>  <span class="keywordflow">return</span> derived();</div> <div class="line"><a name="l00520"></a><span class="lineno"> 520</span> }</div> <div class="line"><a name="l00521"></a><span class="lineno"> 521</span> </div> <div class="line"><a name="l00525"></a><span class="lineno"> 525</span> <span class="keyword">template</span><<span class="keyword">class</span> Derived></div> <div class="line"><a name="l00526"></a><span class="lineno"> 526</span> <span class="keyword">inline</span> <span class="keyword">typename</span> <a class="code" href="classEigen_1_1Matrix.html">QuaternionBase<Derived>::Matrix3</a></div> <div class="line"><a name="l00527"></a><span class="lineno"><a class="line" href="classEigen_1_1QuaternionBase.html#a7bdc6d3cce92c44d32281aa2f85b56f7"> 527</a></span> <a class="code" href="classEigen_1_1QuaternionBase.html#a7bdc6d3cce92c44d32281aa2f85b56f7">QuaternionBase<Derived>::toRotationMatrix</a>(<span class="keywordtype">void</span>)<span class="keyword"> const</span></div> <div class="line"><a name="l00528"></a><span class="lineno"> 528</span> <span class="keyword"></span>{</div> <div class="line"><a name="l00529"></a><span class="lineno"> 529</span>  <span class="comment">// NOTE if inlined, then gcc 4.2 and 4.4 get rid of the temporary (not gcc 4.3 !!)</span></div> <div class="line"><a name="l00530"></a><span class="lineno"> 530</span>  <span class="comment">// if not inlined then the cost of the return by value is huge ~ +35%,</span></div> <div class="line"><a name="l00531"></a><span class="lineno"> 531</span>  <span class="comment">// however, not inlining this function is an order of magnitude slower, so</span></div> <div class="line"><a name="l00532"></a><span class="lineno"> 532</span>  <span class="comment">// it has to be inlined, and so the return by value is not an issue</span></div> <div class="line"><a name="l00533"></a><span class="lineno"> 533</span>  <a class="code" href="classEigen_1_1Matrix.html">Matrix3</a> res;</div> <div class="line"><a name="l00534"></a><span class="lineno"> 534</span> </div> <div class="line"><a name="l00535"></a><span class="lineno"> 535</span>  <span class="keyword">const</span> Scalar tx = Scalar(2)*this->x();</div> <div class="line"><a name="l00536"></a><span class="lineno"> 536</span>  <span class="keyword">const</span> Scalar ty = Scalar(2)*this->y();</div> <div class="line"><a name="l00537"></a><span class="lineno"> 537</span>  <span class="keyword">const</span> Scalar tz = Scalar(2)*this->z();</div> <div class="line"><a name="l00538"></a><span class="lineno"> 538</span>  <span class="keyword">const</span> Scalar twx = tx*this->w();</div> <div class="line"><a name="l00539"></a><span class="lineno"> 539</span>  <span class="keyword">const</span> Scalar twy = ty*this->w();</div> <div class="line"><a name="l00540"></a><span class="lineno"> 540</span>  <span class="keyword">const</span> Scalar twz = tz*this->w();</div> <div class="line"><a name="l00541"></a><span class="lineno"> 541</span>  <span class="keyword">const</span> Scalar txx = tx*this->x();</div> <div class="line"><a name="l00542"></a><span class="lineno"> 542</span>  <span class="keyword">const</span> Scalar txy = ty*this->x();</div> <div class="line"><a name="l00543"></a><span class="lineno"> 543</span>  <span class="keyword">const</span> Scalar txz = tz*this->x();</div> <div class="line"><a name="l00544"></a><span class="lineno"> 544</span>  <span class="keyword">const</span> Scalar tyy = ty*this->y();</div> <div class="line"><a name="l00545"></a><span class="lineno"> 545</span>  <span class="keyword">const</span> Scalar tyz = tz*this->y();</div> <div class="line"><a name="l00546"></a><span class="lineno"> 546</span>  <span class="keyword">const</span> Scalar tzz = tz*this->z();</div> <div class="line"><a name="l00547"></a><span class="lineno"> 547</span> </div> <div class="line"><a name="l00548"></a><span class="lineno"> 548</span>  res.coeffRef(0,0) = Scalar(1)-(tyy+tzz);</div> <div class="line"><a name="l00549"></a><span class="lineno"> 549</span>  res.coeffRef(0,1) = txy-twz;</div> <div class="line"><a name="l00550"></a><span class="lineno"> 550</span>  res.coeffRef(0,2) = txz+twy;</div> <div class="line"><a name="l00551"></a><span class="lineno"> 551</span>  res.coeffRef(1,0) = txy+twz;</div> <div class="line"><a name="l00552"></a><span class="lineno"> 552</span>  res.coeffRef(1,1) = Scalar(1)-(txx+tzz);</div> <div class="line"><a name="l00553"></a><span class="lineno"> 553</span>  res.coeffRef(1,2) = tyz-twx;</div> <div class="line"><a name="l00554"></a><span class="lineno"> 554</span>  res.coeffRef(2,0) = txz-twy;</div> <div class="line"><a name="l00555"></a><span class="lineno"> 555</span>  res.coeffRef(2,1) = tyz+twx;</div> <div class="line"><a name="l00556"></a><span class="lineno"> 556</span>  res.coeffRef(2,2) = Scalar(1)-(txx+tyy);</div> <div class="line"><a name="l00557"></a><span class="lineno"> 557</span> </div> <div class="line"><a name="l00558"></a><span class="lineno"> 558</span>  <span class="keywordflow">return</span> res;</div> <div class="line"><a name="l00559"></a><span class="lineno"> 559</span> }</div> <div class="line"><a name="l00560"></a><span class="lineno"> 560</span> </div> <div class="line"><a name="l00571"></a><span class="lineno"> 571</span> <span class="keyword">template</span><<span class="keyword">class</span> Derived></div> <div class="line"><a name="l00572"></a><span class="lineno"> 572</span> <span class="keyword">template</span><<span class="keyword">typename</span> Derived1, <span class="keyword">typename</span> Derived2></div> <div class="line"><a name="l00573"></a><span class="lineno"><a class="line" href="classEigen_1_1QuaternionBase.html#ac35460294d855096e9b687cadf821452"> 573</a></span> <span class="keyword">inline</span> Derived& <a class="code" href="classEigen_1_1QuaternionBase.html#ac35460294d855096e9b687cadf821452">QuaternionBase<Derived>::setFromTwoVectors</a>(<span class="keyword">const</span> <a class="code" href="classEigen_1_1MatrixBase.html">MatrixBase<Derived1></a>& a, <span class="keyword">const</span> <a class="code" href="classEigen_1_1MatrixBase.html">MatrixBase<Derived2></a>& b)</div> <div class="line"><a name="l00574"></a><span class="lineno"> 574</span> {</div> <div class="line"><a name="l00575"></a><span class="lineno"> 575</span>  <span class="keyword">using</span> std::max;</div> <div class="line"><a name="l00576"></a><span class="lineno"> 576</span>  <span class="keyword">using</span> std::sqrt;</div> <div class="line"><a name="l00577"></a><span class="lineno"> 577</span>  <a class="code" href="classEigen_1_1Matrix.html">Vector3</a> v0 = a.<a class="code" href="classEigen_1_1MatrixBase.html#a8ed1fb2e792b1079639a74e3581fbc74">normalized</a>();</div> <div class="line"><a name="l00578"></a><span class="lineno"> 578</span>  <a class="code" href="classEigen_1_1Matrix.html">Vector3</a> v1 = b.<a class="code" href="classEigen_1_1MatrixBase.html#a8ed1fb2e792b1079639a74e3581fbc74">normalized</a>();</div> <div class="line"><a name="l00579"></a><span class="lineno"> 579</span>  Scalar c = v1.dot(v0);</div> <div class="line"><a name="l00580"></a><span class="lineno"> 580</span> </div> <div class="line"><a name="l00581"></a><span class="lineno"> 581</span>  <span class="comment">// if dot == -1, vectors are nearly opposites</span></div> <div class="line"><a name="l00582"></a><span class="lineno"> 582</span>  <span class="comment">// => accuraletly compute the rotation axis by computing the</span></div> <div class="line"><a name="l00583"></a><span class="lineno"> 583</span>  <span class="comment">// intersection of the two planes. This is done by solving:</span></div> <div class="line"><a name="l00584"></a><span class="lineno"> 584</span>  <span class="comment">// x^T v0 = 0</span></div> <div class="line"><a name="l00585"></a><span class="lineno"> 585</span>  <span class="comment">// x^T v1 = 0</span></div> <div class="line"><a name="l00586"></a><span class="lineno"> 586</span>  <span class="comment">// under the constraint:</span></div> <div class="line"><a name="l00587"></a><span class="lineno"> 587</span>  <span class="comment">// ||x|| = 1</span></div> <div class="line"><a name="l00588"></a><span class="lineno"> 588</span>  <span class="comment">// which yields a singular value problem</span></div> <div class="line"><a name="l00589"></a><span class="lineno"> 589</span>  <span class="keywordflow">if</span> (c < Scalar(-1)+<a class="code" href="structEigen_1_1NumTraits.html">NumTraits<Scalar>::dummy_precision</a>())</div> <div class="line"><a name="l00590"></a><span class="lineno"> 590</span>  {</div> <div class="line"><a name="l00591"></a><span class="lineno"> 591</span>  c = max<Scalar>(c,-1);</div> <div class="line"><a name="l00592"></a><span class="lineno"> 592</span>  <a class="code" href="classEigen_1_1Matrix.html">Matrix<Scalar,2,3></a> m; m << v0.transpose(), v1.transpose();</div> <div class="line"><a name="l00593"></a><span class="lineno"> 593</span>  <a class="code" href="classEigen_1_1JacobiSVD.html">JacobiSVD<Matrix<Scalar,2,3></a> > svd(m, <a class="code" href="group__enums.html#gga2d78499b99ddc29b9494f7ea33864d52a1785ac1174dab733556ac572448984c7">ComputeFullV</a>);</div> <div class="line"><a name="l00594"></a><span class="lineno"> 594</span>  <a class="code" href="classEigen_1_1Matrix.html">Vector3</a> axis = svd.<a class="code" href="classEigen_1_1JacobiSVD.html#ae5158ab7ca44a705c2a3b56ec926b42a">matrixV</a>().col(2);</div> <div class="line"><a name="l00595"></a><span class="lineno"> 595</span> </div> <div class="line"><a name="l00596"></a><span class="lineno"> 596</span>  Scalar w2 = (Scalar(1)+c)*Scalar(0.5);</div> <div class="line"><a name="l00597"></a><span class="lineno"> 597</span>  this->w() = sqrt(w2);</div> <div class="line"><a name="l00598"></a><span class="lineno"> 598</span>  this->vec() = axis * sqrt(Scalar(1) - w2);</div> <div class="line"><a name="l00599"></a><span class="lineno"> 599</span>  <span class="keywordflow">return</span> derived();</div> <div class="line"><a name="l00600"></a><span class="lineno"> 600</span>  }</div> <div class="line"><a name="l00601"></a><span class="lineno"> 601</span>  <a class="code" href="classEigen_1_1Matrix.html">Vector3</a> axis = v0.cross(v1);</div> <div class="line"><a name="l00602"></a><span class="lineno"> 602</span>  Scalar s = sqrt((Scalar(1)+c)*Scalar(2));</div> <div class="line"><a name="l00603"></a><span class="lineno"> 603</span>  Scalar invs = Scalar(1)/s;</div> <div class="line"><a name="l00604"></a><span class="lineno"> 604</span>  this->vec() = axis * invs;</div> <div class="line"><a name="l00605"></a><span class="lineno"> 605</span>  this->w() = s * Scalar(0.5);</div> <div class="line"><a name="l00606"></a><span class="lineno"> 606</span> </div> <div class="line"><a name="l00607"></a><span class="lineno"> 607</span>  <span class="keywordflow">return</span> derived();</div> <div class="line"><a name="l00608"></a><span class="lineno"> 608</span> }</div> <div class="line"><a name="l00609"></a><span class="lineno"> 609</span> </div> <div class="line"><a name="l00610"></a><span class="lineno"> 610</span> </div> <div class="line"><a name="l00621"></a><span class="lineno"> 621</span> <span class="keyword">template</span><<span class="keyword">typename</span> Scalar, <span class="keywordtype">int</span> Options></div> <div class="line"><a name="l00622"></a><span class="lineno"> 622</span> <span class="keyword">template</span><<span class="keyword">typename</span> Derived1, <span class="keyword">typename</span> Derived2></div> <div class="line"><a name="l00623"></a><span class="lineno"><a class="line" href="classEigen_1_1Quaternion.html#a4a0f4bb36de68ca238c88bb263d19c20"> 623</a></span> <a class="code" href="classEigen_1_1Quaternion.html">Quaternion<Scalar,Options></a> <a class="code" href="classEigen_1_1Quaternion.html">Quaternion<Scalar,Options>::FromTwoVectors</a>(<span class="keyword">const</span> <a class="code" href="classEigen_1_1MatrixBase.html">MatrixBase<Derived1></a>& a, <span class="keyword">const</span> <a class="code" href="classEigen_1_1MatrixBase.html">MatrixBase<Derived2></a>& b)</div> <div class="line"><a name="l00624"></a><span class="lineno"> 624</span> {</div> <div class="line"><a name="l00625"></a><span class="lineno"> 625</span>  <a class="code" href="classEigen_1_1Quaternion.html">Quaternion</a> quat;</div> <div class="line"><a name="l00626"></a><span class="lineno"> 626</span>  quat.setFromTwoVectors(a, b);</div> <div class="line"><a name="l00627"></a><span class="lineno"> 627</span>  <span class="keywordflow">return</span> quat;</div> <div class="line"><a name="l00628"></a><span class="lineno"> 628</span> }</div> <div class="line"><a name="l00629"></a><span class="lineno"> 629</span> </div> <div class="line"><a name="l00630"></a><span class="lineno"> 630</span> </div> <div class="line"><a name="l00637"></a><span class="lineno"> 637</span> <span class="keyword">template</span> <<span class="keyword">class</span> Derived></div> <div class="line"><a name="l00638"></a><span class="lineno"><a class="line" href="classEigen_1_1QuaternionBase.html#a666e60218567ecf0ec8e2cd5040e1d87"> 638</a></span> <span class="keyword">inline</span> <a class="code" href="classEigen_1_1Quaternion.html">Quaternion<typename internal::traits<Derived>::Scalar</a>> <a class="code" href="classEigen_1_1QuaternionBase.html#a666e60218567ecf0ec8e2cd5040e1d87">QuaternionBase<Derived>::inverse</a>()<span class="keyword"> const</span></div> <div class="line"><a name="l00639"></a><span class="lineno"> 639</span> <span class="keyword"></span>{</div> <div class="line"><a name="l00640"></a><span class="lineno"> 640</span>  <span class="comment">// FIXME should this function be called multiplicativeInverse and conjugate() be called inverse() or opposite() ??</span></div> <div class="line"><a name="l00641"></a><span class="lineno"> 641</span>  Scalar n2 = this->squaredNorm();</div> <div class="line"><a name="l00642"></a><span class="lineno"> 642</span>  <span class="keywordflow">if</span> (n2 > 0)</div> <div class="line"><a name="l00643"></a><span class="lineno"> 643</span>  <span class="keywordflow">return</span> <a class="code" href="classEigen_1_1Quaternion.html">Quaternion<Scalar></a>(conjugate().coeffs() / n2);</div> <div class="line"><a name="l00644"></a><span class="lineno"> 644</span>  <span class="keywordflow">else</span></div> <div class="line"><a name="l00645"></a><span class="lineno"> 645</span>  {</div> <div class="line"><a name="l00646"></a><span class="lineno"> 646</span>  <span class="comment">// return an invalid result to flag the error</span></div> <div class="line"><a name="l00647"></a><span class="lineno"> 647</span>  <span class="keywordflow">return</span> <a class="code" href="classEigen_1_1Quaternion.html">Quaternion<Scalar></a>(Coefficients::Zero());</div> <div class="line"><a name="l00648"></a><span class="lineno"> 648</span>  }</div> <div class="line"><a name="l00649"></a><span class="lineno"> 649</span> }</div> <div class="line"><a name="l00650"></a><span class="lineno"> 650</span> </div> <div class="line"><a name="l00657"></a><span class="lineno"> 657</span> <span class="keyword">template</span> <<span class="keyword">class</span> Derived></div> <div class="line"><a name="l00658"></a><span class="lineno"> 658</span> <span class="keyword">inline</span> <a class="code" href="classEigen_1_1Quaternion.html">Quaternion<typename internal::traits<Derived>::Scalar</a>></div> <div class="line"><a name="l00659"></a><span class="lineno"><a class="line" href="classEigen_1_1QuaternionBase.html#a6761d13accfe34e38d907505221d081c"> 659</a></span> <a class="code" href="classEigen_1_1QuaternionBase.html#a6761d13accfe34e38d907505221d081c">QuaternionBase<Derived>::conjugate</a>()<span class="keyword"> const</span></div> <div class="line"><a name="l00660"></a><span class="lineno"> 660</span> <span class="keyword"></span>{</div> <div class="line"><a name="l00661"></a><span class="lineno"> 661</span>  <span class="keywordflow">return</span> <a class="code" href="classEigen_1_1Quaternion.html">Quaternion<Scalar></a>(this->w(),-this->x(),-this->y(),-this->z());</div> <div class="line"><a name="l00662"></a><span class="lineno"> 662</span> }</div> <div class="line"><a name="l00663"></a><span class="lineno"> 663</span> </div> <div class="line"><a name="l00667"></a><span class="lineno"> 667</span> <span class="keyword">template</span> <<span class="keyword">class</span> Derived></div> <div class="line"><a name="l00668"></a><span class="lineno"> 668</span> <span class="keyword">template</span> <<span class="keyword">class</span> OtherDerived></div> <div class="line"><a name="l00669"></a><span class="lineno"> 669</span> <span class="keyword">inline</span> <span class="keyword">typename</span> internal::traits<Derived>::Scalar</div> <div class="line"><a name="l00670"></a><span class="lineno"><a class="line" href="classEigen_1_1QuaternionBase.html#ab769d24022d1b0045cc8071f7b55861e"> 670</a></span> <a class="code" href="classEigen_1_1QuaternionBase.html">QuaternionBase<Derived>::angularDistance</a>(<span class="keyword">const</span> <a class="code" href="classEigen_1_1QuaternionBase.html">QuaternionBase<OtherDerived></a>& other)<span class="keyword"> const</span></div> <div class="line"><a name="l00671"></a><span class="lineno"> 671</span> <span class="keyword"></span>{</div> <div class="line"><a name="l00672"></a><span class="lineno"> 672</span>  <span class="keyword">using</span> std::acos;</div> <div class="line"><a name="l00673"></a><span class="lineno"> 673</span>  <span class="keyword">using</span> std::abs;</div> <div class="line"><a name="l00674"></a><span class="lineno"> 674</span>  <span class="keywordtype">double</span> d = abs(this->dot(other));</div> <div class="line"><a name="l00675"></a><span class="lineno"> 675</span>  <span class="keywordflow">if</span> (d>=1.0)</div> <div class="line"><a name="l00676"></a><span class="lineno"> 676</span>  <span class="keywordflow">return</span> Scalar(0);</div> <div class="line"><a name="l00677"></a><span class="lineno"> 677</span>  <span class="keywordflow">return</span> <span class="keyword">static_cast<</span>Scalar<span class="keyword">></span>(2 * acos(d));</div> <div class="line"><a name="l00678"></a><span class="lineno"> 678</span> }</div> <div class="line"><a name="l00679"></a><span class="lineno"> 679</span> </div> <div class="line"><a name="l00683"></a><span class="lineno"> 683</span> <span class="keyword">template</span> <<span class="keyword">class</span> Derived></div> <div class="line"><a name="l00684"></a><span class="lineno"> 684</span> <span class="keyword">template</span> <<span class="keyword">class</span> OtherDerived></div> <div class="line"><a name="l00685"></a><span class="lineno"> 685</span> <a class="code" href="classEigen_1_1Quaternion.html">Quaternion<typename internal::traits<Derived>::Scalar</a>></div> <div class="line"><a name="l00686"></a><span class="lineno"><a class="line" href="classEigen_1_1QuaternionBase.html#aa29a81b780c3012d0fd126a4525781c2"> 686</a></span> <a class="code" href="classEigen_1_1QuaternionBase.html#a8cf5be4be8788856ca8dfe97b3cde1c4">QuaternionBase<Derived>::slerp</a>(<span class="keyword">const</span> Scalar& t, <span class="keyword">const</span> <a class="code" href="classEigen_1_1QuaternionBase.html">QuaternionBase<OtherDerived></a>& other)<span class="keyword"> const</span></div> <div class="line"><a name="l00687"></a><span class="lineno"> 687</span> <span class="keyword"></span>{</div> <div class="line"><a name="l00688"></a><span class="lineno"> 688</span>  <span class="keyword">using</span> std::acos;</div> <div class="line"><a name="l00689"></a><span class="lineno"> 689</span>  <span class="keyword">using</span> std::sin;</div> <div class="line"><a name="l00690"></a><span class="lineno"> 690</span>  <span class="keyword">using</span> std::abs;</div> <div class="line"><a name="l00691"></a><span class="lineno"> 691</span>  <span class="keyword">static</span> <span class="keyword">const</span> Scalar one = Scalar(1) - <a class="code" href="structEigen_1_1NumTraits.html">NumTraits<Scalar>::epsilon</a>();</div> <div class="line"><a name="l00692"></a><span class="lineno"> 692</span>  Scalar d = this->dot(other);</div> <div class="line"><a name="l00693"></a><span class="lineno"> 693</span>  Scalar absD = abs(d);</div> <div class="line"><a name="l00694"></a><span class="lineno"> 694</span> </div> <div class="line"><a name="l00695"></a><span class="lineno"> 695</span>  Scalar scale0;</div> <div class="line"><a name="l00696"></a><span class="lineno"> 696</span>  Scalar scale1;</div> <div class="line"><a name="l00697"></a><span class="lineno"> 697</span> </div> <div class="line"><a name="l00698"></a><span class="lineno"> 698</span>  <span class="keywordflow">if</span>(absD>=one)</div> <div class="line"><a name="l00699"></a><span class="lineno"> 699</span>  {</div> <div class="line"><a name="l00700"></a><span class="lineno"> 700</span>  scale0 = Scalar(1) - t;</div> <div class="line"><a name="l00701"></a><span class="lineno"> 701</span>  scale1 = t;</div> <div class="line"><a name="l00702"></a><span class="lineno"> 702</span>  }</div> <div class="line"><a name="l00703"></a><span class="lineno"> 703</span>  <span class="keywordflow">else</span></div> <div class="line"><a name="l00704"></a><span class="lineno"> 704</span>  {</div> <div class="line"><a name="l00705"></a><span class="lineno"> 705</span>  <span class="comment">// theta is the angle between the 2 quaternions</span></div> <div class="line"><a name="l00706"></a><span class="lineno"> 706</span>  Scalar theta = acos(absD);</div> <div class="line"><a name="l00707"></a><span class="lineno"> 707</span>  Scalar sinTheta = sin(theta);</div> <div class="line"><a name="l00708"></a><span class="lineno"> 708</span> </div> <div class="line"><a name="l00709"></a><span class="lineno"> 709</span>  scale0 = sin( ( Scalar(1) - t ) * theta) / sinTheta;</div> <div class="line"><a name="l00710"></a><span class="lineno"> 710</span>  scale1 = sin( ( t * theta) ) / sinTheta;</div> <div class="line"><a name="l00711"></a><span class="lineno"> 711</span>  }</div> <div class="line"><a name="l00712"></a><span class="lineno"> 712</span>  <span class="keywordflow">if</span>(d<0) scale1 = -scale1;</div> <div class="line"><a name="l00713"></a><span class="lineno"> 713</span> </div> <div class="line"><a name="l00714"></a><span class="lineno"> 714</span>  <span class="keywordflow">return</span> <a class="code" href="classEigen_1_1Quaternion.html">Quaternion<Scalar></a>(scale0 * coeffs() + scale1 * other.coeffs());</div> <div class="line"><a name="l00715"></a><span class="lineno"> 715</span> }</div> <div class="line"><a name="l00716"></a><span class="lineno"> 716</span> </div> <div class="line"><a name="l00717"></a><span class="lineno"> 717</span> <span class="keyword">namespace </span>internal {</div> <div class="line"><a name="l00718"></a><span class="lineno"> 718</span> </div> <div class="line"><a name="l00719"></a><span class="lineno"> 719</span> <span class="comment">// set from a rotation matrix</span></div> <div class="line"><a name="l00720"></a><span class="lineno"> 720</span> <span class="keyword">template</span><<span class="keyword">typename</span> Other></div> <div class="line"><a name="l00721"></a><span class="lineno"> 721</span> <span class="keyword">struct </span>quaternionbase_assign_impl<Other,3,3></div> <div class="line"><a name="l00722"></a><span class="lineno"> 722</span> {</div> <div class="line"><a name="l00723"></a><span class="lineno"> 723</span>  <span class="keyword">typedef</span> <span class="keyword">typename</span> Other::Scalar Scalar;</div> <div class="line"><a name="l00724"></a><span class="lineno"> 724</span>  <span class="keyword">typedef</span> DenseIndex Index;</div> <div class="line"><a name="l00725"></a><span class="lineno"> 725</span>  <span class="keyword">template</span><<span class="keyword">class</span> Derived> <span class="keyword">static</span> <span class="keyword">inline</span> <span class="keywordtype">void</span> run(<a class="code" href="classEigen_1_1QuaternionBase.html">QuaternionBase<Derived></a>& q, <span class="keyword">const</span> Other& mat)</div> <div class="line"><a name="l00726"></a><span class="lineno"> 726</span>  {</div> <div class="line"><a name="l00727"></a><span class="lineno"> 727</span>  <span class="keyword">using</span> std::sqrt;</div> <div class="line"><a name="l00728"></a><span class="lineno"> 728</span>  <span class="comment">// This algorithm comes from "Quaternion Calculus and Fast Animation",</span></div> <div class="line"><a name="l00729"></a><span class="lineno"> 729</span>  <span class="comment">// Ken Shoemake, 1987 SIGGRAPH course notes</span></div> <div class="line"><a name="l00730"></a><span class="lineno"> 730</span>  Scalar t = mat.trace();</div> <div class="line"><a name="l00731"></a><span class="lineno"> 731</span>  <span class="keywordflow">if</span> (t > Scalar(0))</div> <div class="line"><a name="l00732"></a><span class="lineno"> 732</span>  {</div> <div class="line"><a name="l00733"></a><span class="lineno"> 733</span>  t = sqrt(t + Scalar(1.0));</div> <div class="line"><a name="l00734"></a><span class="lineno"> 734</span>  q.w() = Scalar(0.5)*t;</div> <div class="line"><a name="l00735"></a><span class="lineno"> 735</span>  t = Scalar(0.5)/t;</div> <div class="line"><a name="l00736"></a><span class="lineno"> 736</span>  q.x() = (mat.coeff(2,1) - mat.coeff(1,2)) * t;</div> <div class="line"><a name="l00737"></a><span class="lineno"> 737</span>  q.y() = (mat.coeff(0,2) - mat.coeff(2,0)) * t;</div> <div class="line"><a name="l00738"></a><span class="lineno"> 738</span>  q.z() = (mat.coeff(1,0) - mat.coeff(0,1)) * t;</div> <div class="line"><a name="l00739"></a><span class="lineno"> 739</span>  }</div> <div class="line"><a name="l00740"></a><span class="lineno"> 740</span>  <span class="keywordflow">else</span></div> <div class="line"><a name="l00741"></a><span class="lineno"> 741</span>  {</div> <div class="line"><a name="l00742"></a><span class="lineno"> 742</span>  DenseIndex i = 0;</div> <div class="line"><a name="l00743"></a><span class="lineno"> 743</span>  <span class="keywordflow">if</span> (mat.coeff(1,1) > mat.coeff(0,0))</div> <div class="line"><a name="l00744"></a><span class="lineno"> 744</span>  i = 1;</div> <div class="line"><a name="l00745"></a><span class="lineno"> 745</span>  <span class="keywordflow">if</span> (mat.coeff(2,2) > mat.coeff(i,i))</div> <div class="line"><a name="l00746"></a><span class="lineno"> 746</span>  i = 2;</div> <div class="line"><a name="l00747"></a><span class="lineno"> 747</span>  DenseIndex j = (i+1)%3;</div> <div class="line"><a name="l00748"></a><span class="lineno"> 748</span>  DenseIndex k = (j+1)%3;</div> <div class="line"><a name="l00749"></a><span class="lineno"> 749</span> </div> <div class="line"><a name="l00750"></a><span class="lineno"> 750</span>  t = sqrt(mat.coeff(i,i)-mat.coeff(j,j)-mat.coeff(k,k) + Scalar(1.0));</div> <div class="line"><a name="l00751"></a><span class="lineno"> 751</span>  q.coeffs().coeffRef(i) = Scalar(0.5) * t;</div> <div class="line"><a name="l00752"></a><span class="lineno"> 752</span>  t = Scalar(0.5)/t;</div> <div class="line"><a name="l00753"></a><span class="lineno"> 753</span>  q.w() = (mat.coeff(k,j)-mat.coeff(j,k))*t;</div> <div class="line"><a name="l00754"></a><span class="lineno"> 754</span>  q.coeffs().coeffRef(j) = (mat.coeff(j,i)+mat.coeff(i,j))*t;</div> <div class="line"><a name="l00755"></a><span class="lineno"> 755</span>  q.coeffs().coeffRef(k) = (mat.coeff(k,i)+mat.coeff(i,k))*t;</div> <div class="line"><a name="l00756"></a><span class="lineno"> 756</span>  }</div> <div class="line"><a name="l00757"></a><span class="lineno"> 757</span>  }</div> <div class="line"><a name="l00758"></a><span class="lineno"> 758</span> };</div> <div class="line"><a name="l00759"></a><span class="lineno"> 759</span> </div> <div class="line"><a name="l00760"></a><span class="lineno"> 760</span> <span class="comment">// set from a vector of coefficients assumed to be a quaternion</span></div> <div class="line"><a name="l00761"></a><span class="lineno"> 761</span> <span class="keyword">template</span><<span class="keyword">typename</span> Other></div> <div class="line"><a name="l00762"></a><span class="lineno"> 762</span> <span class="keyword">struct </span>quaternionbase_assign_impl<Other,4,1></div> <div class="line"><a name="l00763"></a><span class="lineno"> 763</span> {</div> <div class="line"><a name="l00764"></a><span class="lineno"> 764</span>  <span class="keyword">typedef</span> <span class="keyword">typename</span> Other::Scalar Scalar;</div> <div class="line"><a name="l00765"></a><span class="lineno"> 765</span>  <span class="keyword">template</span><<span class="keyword">class</span> Derived> <span class="keyword">static</span> <span class="keyword">inline</span> <span class="keywordtype">void</span> run(QuaternionBase<Derived>& q, <span class="keyword">const</span> Other& vec)</div> <div class="line"><a name="l00766"></a><span class="lineno"> 766</span>  {</div> <div class="line"><a name="l00767"></a><span class="lineno"> 767</span>  q.coeffs() = vec;</div> <div class="line"><a name="l00768"></a><span class="lineno"> 768</span>  }</div> <div class="line"><a name="l00769"></a><span class="lineno"> 769</span> };</div> <div class="line"><a name="l00770"></a><span class="lineno"> 770</span> </div> <div class="line"><a name="l00771"></a><span class="lineno"> 771</span> } <span class="comment">// end namespace internal</span></div> <div class="line"><a name="l00772"></a><span class="lineno"> 772</span> </div> <div class="line"><a name="l00773"></a><span class="lineno"> 773</span> } <span class="comment">// end namespace Eigen</span></div> <div class="line"><a name="l00774"></a><span class="lineno"> 774</span> </div> <div class="line"><a name="l00775"></a><span class="lineno"> 775</span> <span class="preprocessor">#endif // EIGEN_QUATERNION_H</span></div> <div class="ttc" id="classEigen_1_1Map_3_01const_01Quaternion_3_01__Scalar_01_4_00_01__Options_01_4_html_ac040b01854bd9a68cf897936e8ba1efb"><div class="ttname"><a href="classEigen_1_1Map_3_01const_01Quaternion_3_01__Scalar_01_4_00_01__Options_01_4.html#ac040b01854bd9a68cf897936e8ba1efb">Eigen::Map< const Quaternion< _Scalar >, _Options >::Map</a></div><div class="ttdeci">Map(const Scalar *coeffs)</div><div class="ttdef"><b>Definition:</b> Quaternion.h:355</div></div> <div class="ttc" id="classEigen_1_1QuaternionBase_html_a5eb2596dc2a509b276d01578a9c3dd27"><div class="ttname"><a href="classEigen_1_1QuaternionBase.html#a5eb2596dc2a509b276d01578a9c3dd27">Eigen::QuaternionBase::x</a></div><div class="ttdeci">Scalar x() const </div><div class="ttdef"><b>Definition:</b> Quaternion.h:60</div></div> <div class="ttc" id="group__Geometry__Module_html_gaf65cf6f803890e57488d7de750bef682"><div class="ttname"><a href="group__Geometry__Module.html#gaf65cf6f803890e57488d7de750bef682">Eigen::Quaternionf</a></div><div class="ttdeci">Quaternion< float > Quaternionf</div><div class="ttdef"><b>Definition:</b> Quaternion.h:297</div></div> <div class="ttc" id="classEigen_1_1QuaternionBase_html_aa2bb33b2ed836350b02a6b9083bf8b5b"><div class="ttname"><a href="classEigen_1_1QuaternionBase.html#aa2bb33b2ed836350b02a6b9083bf8b5b">Eigen::QuaternionBase::norm</a></div><div class="ttdeci">Scalar norm() const </div><div class="ttdef"><b>Definition:</b> Quaternion.h:119</div></div> <div class="ttc" id="classEigen_1_1QuaternionBase_html_a20b9b376e143961250055470585ac1d0"><div class="ttname"><a href="classEigen_1_1QuaternionBase.html#a20b9b376e143961250055470585ac1d0">Eigen::QuaternionBase::dot</a></div><div class="ttdeci">Scalar dot(const QuaternionBase< OtherDerived > &other) const </div><div class="ttdef"><b>Definition:</b> Quaternion.h:133</div></div> <div class="ttc" id="classEigen_1_1Map_html"><div class="ttname"><a href="classEigen_1_1Map.html">Eigen::Map</a></div><div class="ttdoc">A matrix or vector expression mapping an existing array of data. </div><div class="ttdef"><b>Definition:</b> Map.h:104</div></div> <div class="ttc" id="classEigen_1_1QuaternionBase_html_a8699d72c996ca6cb4673e810fe3a616c"><div class="ttname"><a href="classEigen_1_1QuaternionBase.html#a8699d72c996ca6cb4673e810fe3a616c">Eigen::QuaternionBase::squaredNorm</a></div><div class="ttdeci">Scalar squaredNorm() const </div><div class="ttdef"><b>Definition:</b> Quaternion.h:114</div></div> <div class="ttc" id="classEigen_1_1RotationBase_html_a844358c46408e878e60c4026c52eb1e9"><div class="ttname"><a href="classEigen_1_1RotationBase.html#a844358c46408e878e60c4026c52eb1e9">Eigen::RotationBase< Quaternion< _Scalar, _Options >, 3 >::Scalar</a></div><div class="ttdeci">internal::traits< Quaternion< _Scalar, _Options > >::Scalar Scalar</div><div class="ttdef"><b>Definition:</b> RotationBase.h:34</div></div> <div class="ttc" id="classEigen_1_1QuaternionBase_html_a316774ca4ac3e8d56a7e24e5e0ad8cc4"><div class="ttname"><a href="classEigen_1_1QuaternionBase.html#a316774ca4ac3e8d56a7e24e5e0ad8cc4">Eigen::QuaternionBase::coeffs</a></div><div class="ttdeci">const internal::traits< Derived >::Coefficients & coeffs() const </div><div class="ttdef"><b>Definition:</b> Quaternion.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_1MatrixBase_html_a8ed1fb2e792b1079639a74e3581fbc74"><div class="ttname"><a href="classEigen_1_1MatrixBase.html#a8ed1fb2e792b1079639a74e3581fbc74">Eigen::MatrixBase::normalized</a></div><div class="ttdeci">const PlainObject normalized() const </div><div class="ttdef"><b>Definition:</b> Dot.h:139</div></div> <div class="ttc" id="classEigen_1_1QuaternionBase_html_a81f3d09a835c62dbde6a0824f5db517b"><div class="ttname"><a href="classEigen_1_1QuaternionBase.html#a81f3d09a835c62dbde6a0824f5db517b">Eigen::QuaternionBase::w</a></div><div class="ttdeci">Scalar & w()</div><div class="ttdef"><b>Definition:</b> Quaternion.h:75</div></div> <div class="ttc" id="group__enums_html_gga2d78499b99ddc29b9494f7ea33864d52a1785ac1174dab733556ac572448984c7"><div class="ttname"><a href="group__enums.html#gga2d78499b99ddc29b9494f7ea33864d52a1785ac1174dab733556ac572448984c7">Eigen::ComputeFullV</a></div><div class="ttdef"><b>Definition:</b> Constants.h:331</div></div> <div class="ttc" id="classEigen_1_1QuaternionBase_html_a666e60218567ecf0ec8e2cd5040e1d87"><div class="ttname"><a href="classEigen_1_1QuaternionBase.html#a666e60218567ecf0ec8e2cd5040e1d87">Eigen::QuaternionBase::inverse</a></div><div class="ttdeci">Quaternion< Scalar > inverse() const </div><div class="ttdef"><b>Definition:</b> Quaternion.h:638</div></div> <div class="ttc" id="classEigen_1_1QuaternionBase_html_a5a8ee232e99cf184e74c6a8dfab5d43e"><div class="ttname"><a href="classEigen_1_1QuaternionBase.html#a5a8ee232e99cf184e74c6a8dfab5d43e">Eigen::QuaternionBase::vec</a></div><div class="ttdeci">VectorBlock< Coefficients, 3 > vec()</div><div class="ttdef"><b>Definition:</b> Quaternion.h:81</div></div> <div class="ttc" id="classEigen_1_1Map_3_01Quaternion_3_01__Scalar_01_4_00_01__Options_01_4_html_aa5742a3a943c3adca0ea1f74f3339cb7"><div class="ttname"><a href="classEigen_1_1Map_3_01Quaternion_3_01__Scalar_01_4_00_01__Options_01_4.html#aa5742a3a943c3adca0ea1f74f3339cb7">Eigen::Map< Quaternion< _Scalar >, _Options >::Map</a></div><div class="ttdeci">Map(Scalar *coeffs)</div><div class="ttdef"><b>Definition:</b> Quaternion.h:392</div></div> <div class="ttc" id="classEigen_1_1VectorBlock_html"><div class="ttname"><a href="classEigen_1_1VectorBlock.html">Eigen::VectorBlock</a></div><div class="ttdoc">Expression of a fixed-size or dynamic-size sub-vector. </div><div class="ttdef"><b>Definition:</b> ForwardDeclarations.h:83</div></div> <div class="ttc" id="classEigen_1_1Quaternion_html_a344a7c038d0b89e798baf2baa3b4a592"><div class="ttname"><a href="classEigen_1_1Quaternion.html#a344a7c038d0b89e798baf2baa3b4a592">Eigen::Quaternion::Quaternion</a></div><div class="ttdeci">Quaternion(const AngleAxisType &aa)</div><div class="ttdef"><b>Definition:</b> Quaternion.h:261</div></div> <div class="ttc" id="classEigen_1_1QuaternionBase_html_a7b19faa70d70b5f16a39f1dd4b69e7a0"><div class="ttname"><a href="classEigen_1_1QuaternionBase.html#a7b19faa70d70b5f16a39f1dd4b69e7a0">Eigen::QuaternionBase::y</a></div><div class="ttdeci">Scalar y() const </div><div class="ttdef"><b>Definition:</b> Quaternion.h:62</div></div> <div class="ttc" id="classEigen_1_1QuaternionBase_html_a1abcf9e222c8e5c38cf63f6821cd6480"><div class="ttname"><a href="classEigen_1_1QuaternionBase.html#a1abcf9e222c8e5c38cf63f6821cd6480">Eigen::QuaternionBase::z</a></div><div class="ttdeci">Scalar & z()</div><div class="ttdef"><b>Definition:</b> Quaternion.h:73</div></div> <div class="ttc" id="classEigen_1_1QuaternionBase_html_a33e390da91b3c4f42c2e7330c0a5338b"><div class="ttname"><a href="classEigen_1_1QuaternionBase.html#a33e390da91b3c4f42c2e7330c0a5338b">Eigen::QuaternionBase::coeffs</a></div><div class="ttdeci">internal::traits< Derived >::Coefficients & coeffs()</div><div class="ttdef"><b>Definition:</b> Quaternion.h:87</div></div> <div class="ttc" id="group__flags_html_ga64e21b7543bdedce27f013512a4403a3"><div class="ttname"><a href="group__flags.html#ga64e21b7543bdedce27f013512a4403a3">Eigen::LvalueBit</a></div><div class="ttdeci">const unsigned int LvalueBit</div><div class="ttdef"><b>Definition:</b> Constants.h:131</div></div> <div class="ttc" id="classEigen_1_1Quaternion_html_afb78170e1d8b745832839d4585becf85"><div class="ttname"><a href="classEigen_1_1Quaternion.html#afb78170e1d8b745832839d4585becf85">Eigen::Quaternion::Quaternion</a></div><div class="ttdeci">Quaternion(const Scalar *data)</div><div class="ttdef"><b>Definition:</b> Quaternion.h:255</div></div> <div class="ttc" id="classEigen_1_1QuaternionBase_html_a6761d13accfe34e38d907505221d081c"><div class="ttname"><a href="classEigen_1_1QuaternionBase.html#a6761d13accfe34e38d907505221d081c">Eigen::QuaternionBase::conjugate</a></div><div class="ttdeci">Quaternion< Scalar > conjugate() const </div><div class="ttdef"><b>Definition:</b> Quaternion.h:659</div></div> <div class="ttc" id="group__enums_html_gga0c5bde183ecefe103f70b49ad9740bcda761c0504a518c7450ed6dfe9eaeab8a6"><div class="ttname"><a href="group__enums.html#gga0c5bde183ecefe103f70b49ad9740bcda761c0504a518c7450ed6dfe9eaeab8a6">Eigen::DontAlign</a></div><div class="ttdef"><b>Definition:</b> Constants.h:270</div></div> <div class="ttc" id="group__Geometry__Module_html_ga0d2bd45f1215359f8e7c0d7ab53c4acb"><div class="ttname"><a href="group__Geometry__Module.html#ga0d2bd45f1215359f8e7c0d7ab53c4acb">Eigen::Quaterniond</a></div><div class="ttdeci">Quaternion< double > Quaterniond</div><div class="ttdef"><b>Definition:</b> Quaternion.h:300</div></div> <div class="ttc" id="classEigen_1_1QuaternionBase_html_ac35460294d855096e9b687cadf821452"><div class="ttname"><a href="classEigen_1_1QuaternionBase.html#ac35460294d855096e9b687cadf821452">Eigen::QuaternionBase::setFromTwoVectors</a></div><div class="ttdeci">Derived & setFromTwoVectors(const MatrixBase< Derived1 > &a, const MatrixBase< Derived2 > &b)</div><div class="ttdef"><b>Definition:</b> Quaternion.h:573</div></div> <div class="ttc" id="classEigen_1_1QuaternionBase_html_a4b12380cc06df343c220d0745b50759e"><div class="ttname"><a href="classEigen_1_1QuaternionBase.html#a4b12380cc06df343c220d0745b50759e">Eigen::QuaternionBase::normalized</a></div><div class="ttdeci">Quaternion< Scalar > normalized() const </div><div class="ttdef"><b>Definition:</b> Quaternion.h:126</div></div> <div class="ttc" id="classEigen_1_1QuaternionBase_html_a65f6d2422debd1a64a2452641497e2f4"><div class="ttname"><a href="classEigen_1_1QuaternionBase.html#a65f6d2422debd1a64a2452641497e2f4">Eigen::QuaternionBase::cast</a></div><div class="ttdeci">internal::cast_return_type< Derived, Quaternion< NewScalarType > >::type cast() const </div><div class="ttdef"><b>Definition:</b> Quaternion.h:176</div></div> <div class="ttc" id="classEigen_1_1Quaternion_html_a801ef1cd3194d3f9e3067c35d883ba4b"><div class="ttname"><a href="classEigen_1_1Quaternion.html#a801ef1cd3194d3f9e3067c35d883ba4b">Eigen::Quaternion::Quaternion</a></div><div class="ttdeci">Quaternion(const Quaternion< OtherScalar, OtherOptions > &other)</div><div class="ttdef"><b>Definition:</b> Quaternion.h:272</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_1QuaternionBase_html_ab05eeebdd63f50a52fecb6353018075f"><div class="ttname"><a href="classEigen_1_1QuaternionBase.html#ab05eeebdd63f50a52fecb6353018075f">Eigen::QuaternionBase::Vector3</a></div><div class="ttdeci">Matrix< Scalar, 3, 1 > Vector3</div><div class="ttdef"><b>Definition:</b> Quaternion.h:51</div></div> <div class="ttc" id="classEigen_1_1QuaternionBase_html_a8cf5be4be8788856ca8dfe97b3cde1c4"><div class="ttname"><a href="classEigen_1_1QuaternionBase.html#a8cf5be4be8788856ca8dfe97b3cde1c4">Eigen::QuaternionBase::slerp</a></div><div class="ttdeci">Quaternion< Scalar > slerp(const Scalar &t, const QuaternionBase< OtherDerived > &other) const </div></div> <div class="ttc" id="classEigen_1_1Quaternion_html_a91b6ea2cac13ab2d33b6e74818ee1490"><div class="ttname"><a href="classEigen_1_1Quaternion.html#a91b6ea2cac13ab2d33b6e74818ee1490">Eigen::Quaternion::Quaternion</a></div><div class="ttdeci">Quaternion(const Scalar &w, const Scalar &x, const Scalar &y, const Scalar &z)</div><div class="ttdef"><b>Definition:</b> Quaternion.h:252</div></div> <div class="ttc" id="classEigen_1_1QuaternionBase_html_ac8cbb66bbccf5f2c9596004a14129c26"><div class="ttname"><a href="classEigen_1_1QuaternionBase.html#ac8cbb66bbccf5f2c9596004a14129c26">Eigen::QuaternionBase::AngleAxisType</a></div><div class="ttdeci">AngleAxis< Scalar > AngleAxisType</div><div class="ttdef"><b>Definition:</b> Quaternion.h:55</div></div> <div class="ttc" id="classEigen_1_1QuaternionBase_html_a5e6fc074bced85f298c76edac1af6cd9"><div class="ttname"><a href="classEigen_1_1QuaternionBase.html#a5e6fc074bced85f298c76edac1af6cd9">Eigen::QuaternionBase::x</a></div><div class="ttdeci">Scalar & x()</div><div class="ttdef"><b>Definition:</b> Quaternion.h:69</div></div> <div class="ttc" id="classEigen_1_1QuaternionBase_html_a7bdc6d3cce92c44d32281aa2f85b56f7"><div class="ttname"><a href="classEigen_1_1QuaternionBase.html#a7bdc6d3cce92c44d32281aa2f85b56f7">Eigen::QuaternionBase::toRotationMatrix</a></div><div class="ttdeci">Matrix3 toRotationMatrix() const </div><div class="ttdef"><b>Definition:</b> Quaternion.h:527</div></div> <div class="ttc" id="classEigen_1_1JacobiSVD_html_ae5158ab7ca44a705c2a3b56ec926b42a"><div class="ttname"><a href="classEigen_1_1JacobiSVD.html#ae5158ab7ca44a705c2a3b56ec926b42a">Eigen::JacobiSVD::matrixV</a></div><div class="ttdeci">const MatrixVType & matrixV() const </div><div class="ttdef"><b>Definition:</b> JacobiSVD.h:618</div></div> <div class="ttc" id="classEigen_1_1QuaternionBase_html_a94a075384f8f54dc97a69546e600b78e"><div class="ttname"><a href="classEigen_1_1QuaternionBase.html#a94a075384f8f54dc97a69546e600b78e">Eigen::QuaternionBase::z</a></div><div class="ttdeci">Scalar z() const </div><div class="ttdef"><b>Definition:</b> Quaternion.h:64</div></div> <div class="ttc" id="group__Geometry__Module_html_gabfe832b0fc8c4c8b4dcf8b59ebe25603"><div class="ttname"><a href="group__Geometry__Module.html#gabfe832b0fc8c4c8b4dcf8b59ebe25603">Eigen::QuaternionMapAlignedd</a></div><div class="ttdeci">Map< Quaternion< double >, Aligned > QuaternionMapAlignedd</div><div class="ttdef"><b>Definition:</b> Quaternion.h:412</div></div> <div class="ttc" id="classEigen_1_1QuaternionBase_html_a3956cb2fddca9905ac6ed71a42ee2d2a"><div class="ttname"><a href="classEigen_1_1QuaternionBase.html#a3956cb2fddca9905ac6ed71a42ee2d2a">Eigen::QuaternionBase::isApprox</a></div><div class="ttdeci">bool isApprox(const QuaternionBase< OtherDerived > &other, const RealScalar &prec=NumTraits< Scalar >::dummy_precision()) const </div><div class="ttdef"><b>Definition:</b> Quaternion.h:164</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_1QuaternionBase_html_ab45aef6f20921b0997f4e8d75ef03300"><div class="ttname"><a href="classEigen_1_1QuaternionBase.html#ab45aef6f20921b0997f4e8d75ef03300">Eigen::QuaternionBase::_transformVector</a></div><div class="ttdeci">Vector3 _transformVector(Vector3 v) const </div><div class="ttdef"><b>Definition:</b> Quaternion.h:466</div></div> <div class="ttc" id="group__enums_html_gga456ac33d49271d3e2c371351cd1d6371ad5380ca00f3d74b38593adf8a0d06d3e"><div class="ttname"><a href="group__enums.html#gga456ac33d49271d3e2c371351cd1d6371ad5380ca00f3d74b38593adf8a0d06d3e">Eigen::Aligned</a></div><div class="ttdef"><b>Definition:</b> Constants.h:194</div></div> <div class="ttc" id="classEigen_1_1Quaternion_html_ad90ae48f7378bb94dfbc6436e3a66aa2"><div class="ttname"><a href="classEigen_1_1Quaternion.html#ad90ae48f7378bb94dfbc6436e3a66aa2">Eigen::Quaternion::Quaternion</a></div><div class="ttdeci">Quaternion(const MatrixBase< Derived > &other)</div><div class="ttdef"><b>Definition:</b> Quaternion.h:268</div></div> <div class="ttc" id="classEigen_1_1JacobiSVD_html"><div class="ttname"><a href="classEigen_1_1JacobiSVD.html">Eigen::JacobiSVD</a></div><div class="ttdoc">Two-sided Jacobi SVD decomposition of a rectangular matrix. </div><div class="ttdef"><b>Definition:</b> ForwardDeclarations.h:224</div></div> <div class="ttc" id="group__Geometry__Module_html_ga987bcae344d68a98892d0e9404138404"><div class="ttname"><a href="group__Geometry__Module.html#ga987bcae344d68a98892d0e9404138404">Eigen::QuaternionMapAlignedf</a></div><div class="ttdeci">Map< Quaternion< float >, Aligned > QuaternionMapAlignedf</div><div class="ttdef"><b>Definition:</b> Quaternion.h:409</div></div> <div class="ttc" id="classEigen_1_1QuaternionBase_html_ad948bb3bd6866bdca15149880c34be62"><div class="ttname"><a href="classEigen_1_1QuaternionBase.html#ad948bb3bd6866bdca15149880c34be62">Eigen::QuaternionBase::y</a></div><div class="ttdeci">Scalar & y()</div><div class="ttdef"><b>Definition:</b> Quaternion.h:71</div></div> <div class="ttc" id="classEigen_1_1QuaternionBase_html_acd0de676568888d848beb97dcc53ae47"><div class="ttname"><a href="classEigen_1_1QuaternionBase.html#acd0de676568888d848beb97dcc53ae47">Eigen::QuaternionBase::normalize</a></div><div class="ttdeci">void normalize()</div><div class="ttdef"><b>Definition:</b> Quaternion.h:123</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< Scalar, 3, 1 ></a></div></div> <div class="ttc" id="classEigen_1_1QuaternionBase_html_abdb0743d0fe840b927780b2bd78fb29e"><div class="ttname"><a href="classEigen_1_1QuaternionBase.html#abdb0743d0fe840b927780b2bd78fb29e">Eigen::QuaternionBase::Matrix3</a></div><div class="ttdeci">Matrix< Scalar, 3, 3 > Matrix3</div><div class="ttdef"><b>Definition:</b> Quaternion.h:53</div></div> <div class="ttc" id="group__Geometry__Module_html_ga8d125241c02e63c656b75b38b2382816"><div class="ttname"><a href="group__Geometry__Module.html#ga8d125241c02e63c656b75b38b2382816">Eigen::QuaternionMapf</a></div><div class="ttdeci">Map< Quaternion< float >, 0 > QuaternionMapf</div><div class="ttdef"><b>Definition:</b> Quaternion.h:403</div></div> <div class="ttc" id="group__Geometry__Module_html_gaca6aed0c662d8272c53663c0093aacaa"><div class="ttname"><a href="group__Geometry__Module.html#gaca6aed0c662d8272c53663c0093aacaa">Eigen::QuaternionMapd</a></div><div class="ttdeci">Map< Quaternion< double >, 0 > QuaternionMapd</div><div class="ttdef"><b>Definition:</b> Quaternion.h:406</div></div> <div class="ttc" id="group__flags_html_ga972a2dcb6603215fa53e0b9e82051426"><div class="ttname"><a href="group__flags.html#ga972a2dcb6603215fa53e0b9e82051426">Eigen::AlignedBit</a></div><div class="ttdeci">const unsigned int AlignedBit</div><div class="ttdef"><b>Definition:</b> Constants.h:147</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_1Quaternion_html_a5f7c4a058b9ca7bf6fb5d40c7c047a4a"><div class="ttname"><a href="classEigen_1_1Quaternion.html#a5f7c4a058b9ca7bf6fb5d40c7c047a4a">Eigen::Quaternion::Quaternion</a></div><div class="ttdeci">Quaternion(const QuaternionBase< Derived > &other)</div><div class="ttdef"><b>Definition:</b> Quaternion.h:258</div></div> <div class="ttc" id="classEigen_1_1QuaternionBase_html_aa76d91cad6945c0d8f57d11289b3e8c2"><div class="ttname"><a href="classEigen_1_1QuaternionBase.html#aa76d91cad6945c0d8f57d11289b3e8c2">Eigen::QuaternionBase::operator*=</a></div><div class="ttdeci">Derived & operator*=(const QuaternionBase< OtherDerived > &q)</div><div class="ttdef"><b>Definition:</b> Quaternion.h:451</div></div> <div class="ttc" id="classEigen_1_1QuaternionBase_html_a480a9376ddbb4d075331a3cefca97bc7"><div class="ttname"><a href="classEigen_1_1QuaternionBase.html#a480a9376ddbb4d075331a3cefca97bc7">Eigen::QuaternionBase::Identity</a></div><div class="ttdeci">static Quaternion< Scalar > Identity()</div><div class="ttdef"><b>Definition:</b> Quaternion.h:105</div></div> <div class="ttc" id="classEigen_1_1QuaternionBase_html_aa07f38d31b32610f5296f81d438b60fb"><div class="ttname"><a href="classEigen_1_1QuaternionBase.html#aa07f38d31b32610f5296f81d438b60fb">Eigen::QuaternionBase::setIdentity</a></div><div class="ttdeci">QuaternionBase & setIdentity()</div><div class="ttdef"><b>Definition:</b> Quaternion.h:109</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>Quaternion.h</b></li> <li class="footer">Generated on Mon Oct 28 2013 11:04:25 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>