<!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: OrthoMethods.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('OrthoMethods_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">OrthoMethods.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-2009 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) 2006-2008 Benoit Jacob <jacob.benoit.1@gmail.com></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_ORTHOMETHODS_H</span></div> <div class="line"><a name="l00012"></a><span class="lineno"> 12</span> <span class="preprocessor"></span><span class="preprocessor">#define EIGEN_ORTHOMETHODS_H</span></div> <div class="line"><a name="l00013"></a><span class="lineno"> 13</span> <span class="preprocessor"></span></div> <div class="line"><a name="l00014"></a><span class="lineno"> 14</span> <span class="keyword">namespace </span>Eigen { </div> <div class="line"><a name="l00015"></a><span class="lineno"> 15</span> </div> <div class="line"><a name="l00023"></a><span class="lineno"> 23</span> <span class="keyword">template</span><<span class="keyword">typename</span> Derived></div> <div class="line"><a name="l00024"></a><span class="lineno"> 24</span> <span class="keyword">template</span><<span class="keyword">typename</span> OtherDerived></div> <div class="line"><a name="l00025"></a><span class="lineno"> 25</span> <span class="keyword">inline</span> <span class="keyword">typename</span> MatrixBase<Derived>::template cross_product_return_type<OtherDerived>::type</div> <div class="line"><a name="l00026"></a><span class="lineno"><a class="line" href="classEigen_1_1MatrixBase.html#a8ba9f136bec0fc2e415fc1ebcd89ce50"> 26</a></span> <a class="code" href="classEigen_1_1MatrixBase.html">MatrixBase<Derived>::cross</a>(<span class="keyword">const</span> <a class="code" href="classEigen_1_1MatrixBase.html">MatrixBase<OtherDerived></a>& other)<span class="keyword"> const</span></div> <div class="line"><a name="l00027"></a><span class="lineno"> 27</span> <span class="keyword"></span>{</div> <div class="line"><a name="l00028"></a><span class="lineno"> 28</span>  EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(Derived,3)</div> <div class="line"><a name="l00029"></a><span class="lineno"> 29</span>  EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(OtherDerived,3)</div> <div class="line"><a name="l00030"></a><span class="lineno"> 30</span> </div> <div class="line"><a name="l00031"></a><span class="lineno"> 31</span>  <span class="comment">// Note that there is no need for an expression here since the compiler</span></div> <div class="line"><a name="l00032"></a><span class="lineno"> 32</span>  <span class="comment">// optimize such a small temporary very well (even within a complex expression)</span></div> <div class="line"><a name="l00033"></a><span class="lineno"> 33</span>  <span class="keyword">typename</span> internal::nested<Derived,2>::type lhs(derived());</div> <div class="line"><a name="l00034"></a><span class="lineno"> 34</span>  <span class="keyword">typename</span> internal::nested<OtherDerived,2>::type rhs(other.derived());</div> <div class="line"><a name="l00035"></a><span class="lineno"> 35</span>  <span class="keywordflow">return</span> <span class="keyword">typename</span> cross_product_return_type<OtherDerived>::type(</div> <div class="line"><a name="l00036"></a><span class="lineno"> 36</span>  numext::conj(lhs.coeff(1) * rhs.coeff(2) - lhs.coeff(2) * rhs.coeff(1)),</div> <div class="line"><a name="l00037"></a><span class="lineno"> 37</span>  numext::conj(lhs.coeff(2) * rhs.coeff(0) - lhs.coeff(0) * rhs.coeff(2)),</div> <div class="line"><a name="l00038"></a><span class="lineno"> 38</span>  numext::conj(lhs.coeff(0) * rhs.coeff(1) - lhs.coeff(1) * rhs.coeff(0))</div> <div class="line"><a name="l00039"></a><span class="lineno"> 39</span>  );</div> <div class="line"><a name="l00040"></a><span class="lineno"> 40</span> }</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">namespace </span>internal {</div> <div class="line"><a name="l00043"></a><span class="lineno"> 43</span> </div> <div class="line"><a name="l00044"></a><span class="lineno"> 44</span> <span class="keyword">template</span>< <span class="keywordtype">int</span> Arch,<span class="keyword">typename</span> VectorLhs,<span class="keyword">typename</span> VectorRhs,</div> <div class="line"><a name="l00045"></a><span class="lineno"> 45</span>  <span class="keyword">typename</span> Scalar = <span class="keyword">typename</span> VectorLhs::Scalar,</div> <div class="line"><a name="l00046"></a><span class="lineno"> 46</span>  <span class="keywordtype">bool</span> Vectorizable = bool((VectorLhs::Flags&VectorRhs::Flags)&<a class="code" href="group__flags.html#gaa780614dc11271c147db56e9c1524e76">PacketAccessBit</a>)></div> <div class="line"><a name="l00047"></a><span class="lineno"> 47</span> <span class="keyword">struct </span>cross3_impl {</div> <div class="line"><a name="l00048"></a><span class="lineno"> 48</span>  <span class="keyword">static</span> <span class="keyword">inline</span> <span class="keyword">typename</span> internal::plain_matrix_type<VectorLhs>::type</div> <div class="line"><a name="l00049"></a><span class="lineno"> 49</span>  run(<span class="keyword">const</span> VectorLhs& lhs, <span class="keyword">const</span> VectorRhs& rhs)</div> <div class="line"><a name="l00050"></a><span class="lineno"> 50</span>  {</div> <div class="line"><a name="l00051"></a><span class="lineno"> 51</span>  <span class="keywordflow">return</span> <span class="keyword">typename</span> internal::plain_matrix_type<VectorLhs>::type(</div> <div class="line"><a name="l00052"></a><span class="lineno"> 52</span>  numext::conj(lhs.coeff(1) * rhs.coeff(2) - lhs.coeff(2) * rhs.coeff(1)),</div> <div class="line"><a name="l00053"></a><span class="lineno"> 53</span>  numext::conj(lhs.coeff(2) * rhs.coeff(0) - lhs.coeff(0) * rhs.coeff(2)),</div> <div class="line"><a name="l00054"></a><span class="lineno"> 54</span>  numext::conj(lhs.coeff(0) * rhs.coeff(1) - lhs.coeff(1) * rhs.coeff(0)),</div> <div class="line"><a name="l00055"></a><span class="lineno"> 55</span>  0</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="l00059"></a><span class="lineno"> 59</span> </div> <div class="line"><a name="l00060"></a><span class="lineno"> 60</span> }</div> <div class="line"><a name="l00061"></a><span class="lineno"> 61</span> </div> <div class="line"><a name="l00071"></a><span class="lineno"> 71</span> <span class="keyword">template</span><<span class="keyword">typename</span> Derived></div> <div class="line"><a name="l00072"></a><span class="lineno"> 72</span> <span class="keyword">template</span><<span class="keyword">typename</span> OtherDerived></div> <div class="line"><a name="l00073"></a><span class="lineno"> 73</span> <span class="keyword">inline</span> <span class="keyword">typename</span> MatrixBase<Derived>::PlainObject</div> <div class="line"><a name="l00074"></a><span class="lineno"><a class="line" href="classEigen_1_1MatrixBase.html#aa74b2131193df7ff5590099b8db634af"> 74</a></span> <a class="code" href="classEigen_1_1MatrixBase.html">MatrixBase<Derived>::cross3</a>(<span class="keyword">const</span> <a class="code" href="classEigen_1_1MatrixBase.html">MatrixBase<OtherDerived></a>& other)<span class="keyword"> const</span></div> <div class="line"><a name="l00075"></a><span class="lineno"> 75</span> <span class="keyword"></span>{</div> <div class="line"><a name="l00076"></a><span class="lineno"> 76</span>  EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(Derived,4)</div> <div class="line"><a name="l00077"></a><span class="lineno"> 77</span>  EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(OtherDerived,4)</div> <div class="line"><a name="l00078"></a><span class="lineno"> 78</span> </div> <div class="line"><a name="l00079"></a><span class="lineno"> 79</span>  <span class="keyword">typedef</span> <span class="keyword">typename</span> internal::nested<Derived,2>::type DerivedNested;</div> <div class="line"><a name="l00080"></a><span class="lineno"> 80</span>  <span class="keyword">typedef</span> <span class="keyword">typename</span> internal::nested<OtherDerived,2>::type OtherDerivedNested;</div> <div class="line"><a name="l00081"></a><span class="lineno"> 81</span>  DerivedNested lhs(derived());</div> <div class="line"><a name="l00082"></a><span class="lineno"> 82</span>  OtherDerivedNested rhs(other.derived());</div> <div class="line"><a name="l00083"></a><span class="lineno"> 83</span> </div> <div class="line"><a name="l00084"></a><span class="lineno"> 84</span>  <span class="keywordflow">return</span> internal::cross3_impl<Architecture::Target,</div> <div class="line"><a name="l00085"></a><span class="lineno"> 85</span>  <span class="keyword">typename</span> internal::remove_all<DerivedNested>::type,</div> <div class="line"><a name="l00086"></a><span class="lineno"> 86</span>  <span class="keyword">typename</span> internal::remove_all<OtherDerivedNested>::type>::run(lhs,rhs);</div> <div class="line"><a name="l00087"></a><span class="lineno"> 87</span> }</div> <div class="line"><a name="l00088"></a><span class="lineno"> 88</span> </div> <div class="line"><a name="l00098"></a><span class="lineno"> 98</span> <span class="keyword">template</span><<span class="keyword">typename</span> ExpressionType, <span class="keywordtype">int</span> Direction></div> <div class="line"><a name="l00099"></a><span class="lineno"> 99</span> <span class="keyword">template</span><<span class="keyword">typename</span> OtherDerived></div> <div class="line"><a name="l00100"></a><span class="lineno"> 100</span> <span class="keyword">const</span> <span class="keyword">typename</span> VectorwiseOp<ExpressionType,Direction>::CrossReturnType</div> <div class="line"><a name="l00101"></a><span class="lineno"><a class="line" href="classEigen_1_1VectorwiseOp.html#a93a7a9d23fc28219f5dd82f3f90222ca"> 101</a></span> <a class="code" href="classEigen_1_1VectorwiseOp.html">VectorwiseOp<ExpressionType,Direction>::cross</a>(<span class="keyword">const</span> <a class="code" href="classEigen_1_1MatrixBase.html">MatrixBase<OtherDerived></a>& other)<span class="keyword"> const</span></div> <div class="line"><a name="l00102"></a><span class="lineno"> 102</span> <span class="keyword"></span>{</div> <div class="line"><a name="l00103"></a><span class="lineno"> 103</span>  EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(OtherDerived,3)</div> <div class="line"><a name="l00104"></a><span class="lineno"> 104</span>  EIGEN_STATIC_ASSERT((internal::is_same<Scalar, typename OtherDerived::Scalar>::value),</div> <div class="line"><a name="l00105"></a><span class="lineno"> 105</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="l00106"></a><span class="lineno"> 106</span> </div> <div class="line"><a name="l00107"></a><span class="lineno"> 107</span>  CrossReturnType res(_expression().rows(),_expression().cols());</div> <div class="line"><a name="l00108"></a><span class="lineno"> 108</span>  <span class="keywordflow">if</span>(Direction==<a class="code" href="group__enums.html#gga8ef30fa9c08e08c8706653571f9f5b81a2434cd8c1a594a4cdaa250f86639c600">Vertical</a>)</div> <div class="line"><a name="l00109"></a><span class="lineno"> 109</span>  {</div> <div class="line"><a name="l00110"></a><span class="lineno"> 110</span>  eigen_assert(CrossReturnType::RowsAtCompileTime==3 && <span class="stringliteral">"the matrix must have exactly 3 rows"</span>);</div> <div class="line"><a name="l00111"></a><span class="lineno"> 111</span>  res.row(0) = (_expression().row(1) * other.coeff(2) - _expression().<a class="code" href="classEigen_1_1DenseBase.html#aa8716d44f51321072ee5c88665c28813">row</a>(2) * other.coeff(1)).conjugate();</div> <div class="line"><a name="l00112"></a><span class="lineno"> 112</span>  res.row(1) = (_expression().row(2) * other.coeff(0) - _expression().<a class="code" href="classEigen_1_1DenseBase.html#aa8716d44f51321072ee5c88665c28813">row</a>(0) * other.coeff(2)).conjugate();</div> <div class="line"><a name="l00113"></a><span class="lineno"> 113</span>  res.row(2) = (_expression().row(0) * other.coeff(1) - _expression().<a class="code" href="classEigen_1_1DenseBase.html#aa8716d44f51321072ee5c88665c28813">row</a>(1) * other.coeff(0)).conjugate();</div> <div class="line"><a name="l00114"></a><span class="lineno"> 114</span>  }</div> <div class="line"><a name="l00115"></a><span class="lineno"> 115</span>  <span class="keywordflow">else</span></div> <div class="line"><a name="l00116"></a><span class="lineno"> 116</span>  {</div> <div class="line"><a name="l00117"></a><span class="lineno"> 117</span>  eigen_assert(CrossReturnType::ColsAtCompileTime==3 && <span class="stringliteral">"the matrix must have exactly 3 columns"</span>);</div> <div class="line"><a name="l00118"></a><span class="lineno"> 118</span>  res.col(0) = (_expression().col(1) * other.coeff(2) - _expression().<a class="code" href="classEigen_1_1DenseBase.html#a58c77695de3b33405f01f2fdf3dc389d">col</a>(2) * other.coeff(1)).conjugate();</div> <div class="line"><a name="l00119"></a><span class="lineno"> 119</span>  res.col(1) = (_expression().col(2) * other.coeff(0) - _expression().<a class="code" href="classEigen_1_1DenseBase.html#a58c77695de3b33405f01f2fdf3dc389d">col</a>(0) * other.coeff(2)).conjugate();</div> <div class="line"><a name="l00120"></a><span class="lineno"> 120</span>  res.col(2) = (_expression().col(0) * other.coeff(1) - _expression().<a class="code" href="classEigen_1_1DenseBase.html#a58c77695de3b33405f01f2fdf3dc389d">col</a>(1) * other.coeff(0)).conjugate();</div> <div class="line"><a name="l00121"></a><span class="lineno"> 121</span>  }</div> <div class="line"><a name="l00122"></a><span class="lineno"> 122</span>  <span class="keywordflow">return</span> res;</div> <div class="line"><a name="l00123"></a><span class="lineno"> 123</span> }</div> <div class="line"><a name="l00124"></a><span class="lineno"> 124</span> </div> <div class="line"><a name="l00125"></a><span class="lineno"> 125</span> <span class="keyword">namespace </span>internal {</div> <div class="line"><a name="l00126"></a><span class="lineno"> 126</span> </div> <div class="line"><a name="l00127"></a><span class="lineno"> 127</span> <span class="keyword">template</span><<span class="keyword">typename</span> Derived, <span class="keywordtype">int</span> Size = Derived::SizeAtCompileTime></div> <div class="line"><a name="l00128"></a><span class="lineno"> 128</span> <span class="keyword">struct </span>unitOrthogonal_selector</div> <div class="line"><a name="l00129"></a><span class="lineno"> 129</span> {</div> <div class="line"><a name="l00130"></a><span class="lineno"> 130</span>  <span class="keyword">typedef</span> <span class="keyword">typename</span> plain_matrix_type<Derived>::type VectorType;</div> <div class="line"><a name="l00131"></a><span class="lineno"> 131</span>  <span class="keyword">typedef</span> <span class="keyword">typename</span> traits<Derived>::Scalar Scalar;</div> <div class="line"><a name="l00132"></a><span class="lineno"> 132</span>  <span class="keyword">typedef</span> <span class="keyword">typename</span> <a class="code" href="structEigen_1_1NumTraits.html">NumTraits<Scalar>::Real</a> RealScalar;</div> <div class="line"><a name="l00133"></a><span class="lineno"> 133</span>  <span class="keyword">typedef</span> <span class="keyword">typename</span> Derived::Index Index;</div> <div class="line"><a name="l00134"></a><span class="lineno"> 134</span>  <span class="keyword">typedef</span> <a class="code" href="classEigen_1_1Matrix.html">Matrix<Scalar,2,1></a> Vector2;</div> <div class="line"><a name="l00135"></a><span class="lineno"> 135</span>  <span class="keyword">static</span> <span class="keyword">inline</span> VectorType run(<span class="keyword">const</span> Derived& src)</div> <div class="line"><a name="l00136"></a><span class="lineno"> 136</span>  {</div> <div class="line"><a name="l00137"></a><span class="lineno"> 137</span>  VectorType perp = VectorType::Zero(src.size());</div> <div class="line"><a name="l00138"></a><span class="lineno"> 138</span>  Index maxi = 0;</div> <div class="line"><a name="l00139"></a><span class="lineno"> 139</span>  Index sndi = 0;</div> <div class="line"><a name="l00140"></a><span class="lineno"> 140</span>  src.<a class="code" href="classEigen_1_1MatrixBase.html#ab4d8398a4497e4a888cfc11c51c14a81">cwiseAbs</a>().maxCoeff(&maxi);</div> <div class="line"><a name="l00141"></a><span class="lineno"> 141</span>  <span class="keywordflow">if</span> (maxi==0)</div> <div class="line"><a name="l00142"></a><span class="lineno"> 142</span>  sndi = 1;</div> <div class="line"><a name="l00143"></a><span class="lineno"> 143</span>  RealScalar invnm = RealScalar(1)/(Vector2() << src.coeff(sndi),src.coeff(maxi)).finished().norm();</div> <div class="line"><a name="l00144"></a><span class="lineno"> 144</span>  perp.coeffRef(maxi) = -numext::conj(src.coeff(sndi)) * invnm;</div> <div class="line"><a name="l00145"></a><span class="lineno"> 145</span>  perp.coeffRef(sndi) = numext::conj(src.coeff(maxi)) * invnm;</div> <div class="line"><a name="l00146"></a><span class="lineno"> 146</span> </div> <div class="line"><a name="l00147"></a><span class="lineno"> 147</span>  <span class="keywordflow">return</span> perp;</div> <div class="line"><a name="l00148"></a><span class="lineno"> 148</span>  }</div> <div class="line"><a name="l00149"></a><span class="lineno"> 149</span> };</div> <div class="line"><a name="l00150"></a><span class="lineno"> 150</span> </div> <div class="line"><a name="l00151"></a><span class="lineno"> 151</span> <span class="keyword">template</span><<span class="keyword">typename</span> Derived></div> <div class="line"><a name="l00152"></a><span class="lineno"> 152</span> <span class="keyword">struct </span>unitOrthogonal_selector<Derived,3></div> <div class="line"><a name="l00153"></a><span class="lineno"> 153</span> {</div> <div class="line"><a name="l00154"></a><span class="lineno"> 154</span>  <span class="keyword">typedef</span> <span class="keyword">typename</span> plain_matrix_type<Derived>::type VectorType;</div> <div class="line"><a name="l00155"></a><span class="lineno"> 155</span>  <span class="keyword">typedef</span> <span class="keyword">typename</span> traits<Derived>::Scalar Scalar;</div> <div class="line"><a name="l00156"></a><span class="lineno"> 156</span>  <span class="keyword">typedef</span> <span class="keyword">typename</span> NumTraits<Scalar>::Real RealScalar;</div> <div class="line"><a name="l00157"></a><span class="lineno"> 157</span>  <span class="keyword">static</span> <span class="keyword">inline</span> VectorType run(<span class="keyword">const</span> Derived& src)</div> <div class="line"><a name="l00158"></a><span class="lineno"> 158</span>  {</div> <div class="line"><a name="l00159"></a><span class="lineno"> 159</span>  VectorType perp;</div> <div class="line"><a name="l00160"></a><span class="lineno"> 160</span>  <span class="comment">/* Let us compute the crossed product of *this with a vector</span></div> <div class="line"><a name="l00161"></a><span class="lineno"> 161</span> <span class="comment"> * that is not too close to being colinear to *this.</span></div> <div class="line"><a name="l00162"></a><span class="lineno"> 162</span> <span class="comment"> */</span></div> <div class="line"><a name="l00163"></a><span class="lineno"> 163</span> </div> <div class="line"><a name="l00164"></a><span class="lineno"> 164</span>  <span class="comment">/* unless the x and y coords are both close to zero, we can</span></div> <div class="line"><a name="l00165"></a><span class="lineno"> 165</span> <span class="comment"> * simply take ( -y, x, 0 ) and normalize it.</span></div> <div class="line"><a name="l00166"></a><span class="lineno"> 166</span> <span class="comment"> */</span></div> <div class="line"><a name="l00167"></a><span class="lineno"> 167</span>  <span class="keywordflow">if</span>((!isMuchSmallerThan(src.x(), src.z()))</div> <div class="line"><a name="l00168"></a><span class="lineno"> 168</span>  || (!isMuchSmallerThan(src.y(), src.z())))</div> <div class="line"><a name="l00169"></a><span class="lineno"> 169</span>  {</div> <div class="line"><a name="l00170"></a><span class="lineno"> 170</span>  RealScalar invnm = RealScalar(1)/src.template head<2>().norm();</div> <div class="line"><a name="l00171"></a><span class="lineno"> 171</span>  perp.coeffRef(0) = -numext::conj(src.y())*invnm;</div> <div class="line"><a name="l00172"></a><span class="lineno"> 172</span>  perp.coeffRef(1) = numext::conj(src.x())*invnm;</div> <div class="line"><a name="l00173"></a><span class="lineno"> 173</span>  perp.coeffRef(2) = 0;</div> <div class="line"><a name="l00174"></a><span class="lineno"> 174</span>  }</div> <div class="line"><a name="l00175"></a><span class="lineno"> 175</span>  <span class="comment">/* if both x and y are close to zero, then the vector is close</span></div> <div class="line"><a name="l00176"></a><span class="lineno"> 176</span> <span class="comment"> * to the z-axis, so it's far from colinear to the x-axis for instance.</span></div> <div class="line"><a name="l00177"></a><span class="lineno"> 177</span> <span class="comment"> * So we take the crossed product with (1,0,0) and normalize it.</span></div> <div class="line"><a name="l00178"></a><span class="lineno"> 178</span> <span class="comment"> */</span></div> <div class="line"><a name="l00179"></a><span class="lineno"> 179</span>  <span class="keywordflow">else</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>  RealScalar invnm = RealScalar(1)/src.template tail<2>().norm();</div> <div class="line"><a name="l00182"></a><span class="lineno"> 182</span>  perp.coeffRef(0) = 0;</div> <div class="line"><a name="l00183"></a><span class="lineno"> 183</span>  perp.coeffRef(1) = -numext::conj(src.z())*invnm;</div> <div class="line"><a name="l00184"></a><span class="lineno"> 184</span>  perp.coeffRef(2) = numext::conj(src.y())*invnm;</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> </div> <div class="line"><a name="l00187"></a><span class="lineno"> 187</span>  <span class="keywordflow">return</span> perp;</div> <div class="line"><a name="l00188"></a><span class="lineno"> 188</span>  }</div> <div class="line"><a name="l00189"></a><span class="lineno"> 189</span> };</div> <div class="line"><a name="l00190"></a><span class="lineno"> 190</span> </div> <div class="line"><a name="l00191"></a><span class="lineno"> 191</span> <span class="keyword">template</span><<span class="keyword">typename</span> Derived></div> <div class="line"><a name="l00192"></a><span class="lineno"> 192</span> <span class="keyword">struct </span>unitOrthogonal_selector<Derived,2></div> <div class="line"><a name="l00193"></a><span class="lineno"> 193</span> {</div> <div class="line"><a name="l00194"></a><span class="lineno"> 194</span>  <span class="keyword">typedef</span> <span class="keyword">typename</span> plain_matrix_type<Derived>::type VectorType;</div> <div class="line"><a name="l00195"></a><span class="lineno"> 195</span>  <span class="keyword">static</span> <span class="keyword">inline</span> VectorType run(<span class="keyword">const</span> Derived& src)</div> <div class="line"><a name="l00196"></a><span class="lineno"> 196</span>  { <span class="keywordflow">return</span> VectorType(-numext::conj(src.y()), numext::conj(src.x())).normalized(); }</div> <div class="line"><a name="l00197"></a><span class="lineno"> 197</span> };</div> <div class="line"><a name="l00198"></a><span class="lineno"> 198</span> </div> <div class="line"><a name="l00199"></a><span class="lineno"> 199</span> } <span class="comment">// end namespace internal</span></div> <div class="line"><a name="l00200"></a><span class="lineno"> 200</span> </div> <div class="line"><a name="l00208"></a><span class="lineno"> 208</span> <span class="keyword">template</span><<span class="keyword">typename</span> Derived></div> <div class="line"><a name="l00209"></a><span class="lineno"> 209</span> <span class="keyword">typename</span> MatrixBase<Derived>::PlainObject</div> <div class="line"><a name="l00210"></a><span class="lineno"><a class="line" href="classEigen_1_1MatrixBase.html#a42fdca0c5854b85d9482853c8a085dc1"> 210</a></span> <a class="code" href="classEigen_1_1MatrixBase.html#a42fdca0c5854b85d9482853c8a085dc1">MatrixBase<Derived>::unitOrthogonal</a>()<span class="keyword"> const</span></div> <div class="line"><a name="l00211"></a><span class="lineno"> 211</span> <span class="keyword"></span>{</div> <div class="line"><a name="l00212"></a><span class="lineno"> 212</span>  EIGEN_STATIC_ASSERT_VECTOR_ONLY(Derived)</div> <div class="line"><a name="l00213"></a><span class="lineno"> 213</span>  <span class="keywordflow">return</span> internal::unitOrthogonal_selector<Derived>::run(derived());</div> <div class="line"><a name="l00214"></a><span class="lineno"> 214</span> }</div> <div class="line"><a name="l00215"></a><span class="lineno"> 215</span> </div> <div class="line"><a name="l00216"></a><span class="lineno"> 216</span> } <span class="comment">// end namespace Eigen</span></div> <div class="line"><a name="l00217"></a><span class="lineno"> 217</span> </div> <div class="line"><a name="l00218"></a><span class="lineno"> 218</span> <span class="preprocessor">#endif // EIGEN_ORTHOMETHODS_H</span></div> <div class="ttc" id="classEigen_1_1DenseBase_html_aa8716d44f51321072ee5c88665c28813"><div class="ttname"><a href="classEigen_1_1DenseBase.html#aa8716d44f51321072ee5c88665c28813">Eigen::DenseBase::row</a></div><div class="ttdeci">RowXpr row(Index i)</div><div class="ttdef"><b>Definition:</b> DenseBase.h:726</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_1VectorwiseOp_html"><div class="ttname"><a href="classEigen_1_1VectorwiseOp.html">Eigen::VectorwiseOp</a></div><div class="ttdoc">Pseudo expression providing partial reduction operations. </div><div class="ttdef"><b>Definition:</b> ForwardDeclarations.h:212</div></div> <div class="ttc" id="classEigen_1_1DenseBase_html_a58c77695de3b33405f01f2fdf3dc389d"><div class="ttname"><a href="classEigen_1_1DenseBase.html#a58c77695de3b33405f01f2fdf3dc389d">Eigen::DenseBase::col</a></div><div class="ttdeci">ColXpr col(Index i)</div><div class="ttdef"><b>Definition:</b> DenseBase.h:709</div></div> <div class="ttc" id="group__flags_html_gaa780614dc11271c147db56e9c1524e76"><div class="ttname"><a href="group__flags.html#gaa780614dc11271c147db56e9c1524e76">Eigen::PacketAccessBit</a></div><div class="ttdeci">const unsigned int PacketAccessBit</div><div class="ttdef"><b>Definition:</b> Constants.h:81</div></div> <div class="ttc" id="group__enums_html_gga8ef30fa9c08e08c8706653571f9f5b81a2434cd8c1a594a4cdaa250f86639c600"><div class="ttname"><a href="group__enums.html#gga8ef30fa9c08e08c8706653571f9f5b81a2434cd8c1a594a4cdaa250f86639c600">Eigen::Vertical</a></div><div class="ttdef"><b>Definition:</b> Constants.h:209</div></div> <div class="ttc" id="classEigen_1_1MatrixBase_html_a42fdca0c5854b85d9482853c8a085dc1"><div class="ttname"><a href="classEigen_1_1MatrixBase.html#a42fdca0c5854b85d9482853c8a085dc1">Eigen::MatrixBase::unitOrthogonal</a></div><div class="ttdeci">PlainObject unitOrthogonal(void) const </div><div class="ttdef"><b>Definition:</b> OrthoMethods.h:210</div></div> <div class="ttc" id="classEigen_1_1Matrix_html"><div class="ttname"><a href="classEigen_1_1Matrix.html">Eigen::Matrix</a></div><div class="ttdoc">The matrix class, also used for vectors and row-vectors. </div><div class="ttdef"><b>Definition:</b> Matrix.h:127</div></div> <div class="ttc" id="classEigen_1_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_1MatrixBase_html_ab4d8398a4497e4a888cfc11c51c14a81"><div class="ttname"><a href="classEigen_1_1MatrixBase.html#ab4d8398a4497e4a888cfc11c51c14a81">Eigen::MatrixBase::cwiseAbs</a></div><div class="ttdeci">const CwiseUnaryOp< internal::scalar_abs_op< Scalar >, const Derived > cwiseAbs() const </div><div class="ttdef"><b>Definition:</b> MatrixBase.h:22</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>OrthoMethods.h</b></li> <li class="footer">Generated on Mon Oct 28 2013 11:04:24 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>