<!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-unsupported: ei_kissfft_impl.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-unsupported</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('ei__kissfft__impl_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>Enumerator</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(7)"><span class="SelectionMark"> </span>Groups</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(8)"><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">ei_kissfft_impl.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) 2009 Mark Borgerding mark a borgerding net</span></div> <div class="line"><a name="l00005"></a><span class="lineno"> 5</span> <span class="comment">//</span></div> <div class="line"><a name="l00006"></a><span class="lineno"> 6</span> <span class="comment">// This Source Code Form is subject to the terms of the Mozilla</span></div> <div class="line"><a name="l00007"></a><span class="lineno"> 7</span> <span class="comment">// Public License v. 2.0. If a copy of the MPL was not distributed</span></div> <div class="line"><a name="l00008"></a><span class="lineno"> 8</span> <span class="comment">// with this file, You can obtain one at http://mozilla.org/MPL/2.0/.</span></div> <div class="line"><a name="l00009"></a><span class="lineno"> 9</span> </div> <div class="line"><a name="l00010"></a><span class="lineno"> 10</span> <span class="keyword">namespace </span>Eigen { </div> <div class="line"><a name="l00011"></a><span class="lineno"> 11</span> </div> <div class="line"><a name="l00012"></a><span class="lineno"> 12</span> <span class="keyword">namespace </span>internal {</div> <div class="line"><a name="l00013"></a><span class="lineno"> 13</span> </div> <div class="line"><a name="l00014"></a><span class="lineno"> 14</span>  <span class="comment">// This FFT implementation was derived from kissfft http:sourceforge.net/projects/kissfft</span></div> <div class="line"><a name="l00015"></a><span class="lineno"> 15</span>  <span class="comment">// Copyright 2003-2009 Mark Borgerding</span></div> <div class="line"><a name="l00016"></a><span class="lineno"> 16</span> </div> <div class="line"><a name="l00017"></a><span class="lineno"> 17</span> <span class="keyword">template</span> <<span class="keyword">typename</span> _Scalar></div> <div class="line"><a name="l00018"></a><span class="lineno"> 18</span> <span class="keyword">struct </span>kiss_cpx_fft</div> <div class="line"><a name="l00019"></a><span class="lineno"> 19</span> {</div> <div class="line"><a name="l00020"></a><span class="lineno"> 20</span>  <span class="keyword">typedef</span> _Scalar Scalar;</div> <div class="line"><a name="l00021"></a><span class="lineno"> 21</span>  <span class="keyword">typedef</span> std::complex<Scalar> Complex;</div> <div class="line"><a name="l00022"></a><span class="lineno"> 22</span>  std::vector<Complex> m_twiddles;</div> <div class="line"><a name="l00023"></a><span class="lineno"> 23</span>  std::vector<int> m_stageRadix;</div> <div class="line"><a name="l00024"></a><span class="lineno"> 24</span>  std::vector<int> m_stageRemainder;</div> <div class="line"><a name="l00025"></a><span class="lineno"> 25</span>  std::vector<Complex> m_scratchBuf;</div> <div class="line"><a name="l00026"></a><span class="lineno"> 26</span>  <span class="keywordtype">bool</span> m_inverse;</div> <div class="line"><a name="l00027"></a><span class="lineno"> 27</span> </div> <div class="line"><a name="l00028"></a><span class="lineno"> 28</span>  <span class="keyword">inline</span></div> <div class="line"><a name="l00029"></a><span class="lineno"> 29</span>  <span class="keywordtype">void</span> make_twiddles(<span class="keywordtype">int</span> nfft,<span class="keywordtype">bool</span> inverse)</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="keyword">using</span> std::acos;</div> <div class="line"><a name="l00032"></a><span class="lineno"> 32</span>  m_inverse = inverse;</div> <div class="line"><a name="l00033"></a><span class="lineno"> 33</span>  m_twiddles.resize(nfft);</div> <div class="line"><a name="l00034"></a><span class="lineno"> 34</span>  Scalar phinc = (inverse?2:-2)* acos( (Scalar) -1) / nfft;</div> <div class="line"><a name="l00035"></a><span class="lineno"> 35</span>  <span class="keywordflow">for</span> (<span class="keywordtype">int</span> i=0;i<nfft;++i)</div> <div class="line"><a name="l00036"></a><span class="lineno"> 36</span>  m_twiddles[i] = exp( Complex(0,i*phinc) );</div> <div class="line"><a name="l00037"></a><span class="lineno"> 37</span>  }</div> <div class="line"><a name="l00038"></a><span class="lineno"> 38</span> </div> <div class="line"><a name="l00039"></a><span class="lineno"> 39</span>  <span class="keywordtype">void</span> factorize(<span class="keywordtype">int</span> nfft)</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>  <span class="comment">//start factoring out 4's, then 2's, then 3,5,7,9,...</span></div> <div class="line"><a name="l00042"></a><span class="lineno"> 42</span>  <span class="keywordtype">int</span> n= nfft;</div> <div class="line"><a name="l00043"></a><span class="lineno"> 43</span>  <span class="keywordtype">int</span> p=4;</div> <div class="line"><a name="l00044"></a><span class="lineno"> 44</span>  <span class="keywordflow">do</span> {</div> <div class="line"><a name="l00045"></a><span class="lineno"> 45</span>  <span class="keywordflow">while</span> (n % p) {</div> <div class="line"><a name="l00046"></a><span class="lineno"> 46</span>  <span class="keywordflow">switch</span> (p) {</div> <div class="line"><a name="l00047"></a><span class="lineno"> 47</span>  <span class="keywordflow">case</span> 4: p = 2; <span class="keywordflow">break</span>;</div> <div class="line"><a name="l00048"></a><span class="lineno"> 48</span>  <span class="keywordflow">case</span> 2: p = 3; <span class="keywordflow">break</span>;</div> <div class="line"><a name="l00049"></a><span class="lineno"> 49</span>  <span class="keywordflow">default</span>: p += 2; <span class="keywordflow">break</span>;</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">if</span> (p*p>n)</div> <div class="line"><a name="l00052"></a><span class="lineno"> 52</span>  p=n;<span class="comment">// impossible to have a factor > sqrt(n)</span></div> <div class="line"><a name="l00053"></a><span class="lineno"> 53</span>  }</div> <div class="line"><a name="l00054"></a><span class="lineno"> 54</span>  n /= p;</div> <div class="line"><a name="l00055"></a><span class="lineno"> 55</span>  m_stageRadix.push_back(p);</div> <div class="line"><a name="l00056"></a><span class="lineno"> 56</span>  m_stageRemainder.push_back(n);</div> <div class="line"><a name="l00057"></a><span class="lineno"> 57</span>  <span class="keywordflow">if</span> ( p > 5 )</div> <div class="line"><a name="l00058"></a><span class="lineno"> 58</span>  m_scratchBuf.resize(p); <span class="comment">// scratchbuf will be needed in bfly_generic</span></div> <div class="line"><a name="l00059"></a><span class="lineno"> 59</span>  }<span class="keywordflow">while</span>(n>1);</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="l00062"></a><span class="lineno"> 62</span>  <span class="keyword">template</span> <<span class="keyword">typename</span> _Src></div> <div class="line"><a name="l00063"></a><span class="lineno"> 63</span>  <span class="keyword">inline</span></div> <div class="line"><a name="l00064"></a><span class="lineno"> 64</span>  <span class="keywordtype">void</span> work( <span class="keywordtype">int</span> stage,Complex * xout, <span class="keyword">const</span> _Src * xin, <span class="keywordtype">size_t</span> fstride,<span class="keywordtype">size_t</span> in_stride)</div> <div class="line"><a name="l00065"></a><span class="lineno"> 65</span>  {</div> <div class="line"><a name="l00066"></a><span class="lineno"> 66</span>  <span class="keywordtype">int</span> p = m_stageRadix[stage];</div> <div class="line"><a name="l00067"></a><span class="lineno"> 67</span>  <span class="keywordtype">int</span> m = m_stageRemainder[stage];</div> <div class="line"><a name="l00068"></a><span class="lineno"> 68</span>  Complex * Fout_beg = xout;</div> <div class="line"><a name="l00069"></a><span class="lineno"> 69</span>  Complex * Fout_end = xout + p*m;</div> <div class="line"><a name="l00070"></a><span class="lineno"> 70</span> </div> <div class="line"><a name="l00071"></a><span class="lineno"> 71</span>  <span class="keywordflow">if</span> (m>1) {</div> <div class="line"><a name="l00072"></a><span class="lineno"> 72</span>  <span class="keywordflow">do</span>{</div> <div class="line"><a name="l00073"></a><span class="lineno"> 73</span>  <span class="comment">// recursive call:</span></div> <div class="line"><a name="l00074"></a><span class="lineno"> 74</span>  <span class="comment">// DFT of size m*p performed by doing</span></div> <div class="line"><a name="l00075"></a><span class="lineno"> 75</span>  <span class="comment">// p instances of smaller DFTs of size m, </span></div> <div class="line"><a name="l00076"></a><span class="lineno"> 76</span>  <span class="comment">// each one takes a decimated version of the input</span></div> <div class="line"><a name="l00077"></a><span class="lineno"> 77</span>  work(stage+1, xout , xin, fstride*p,in_stride);</div> <div class="line"><a name="l00078"></a><span class="lineno"> 78</span>  xin += fstride*in_stride;</div> <div class="line"><a name="l00079"></a><span class="lineno"> 79</span>  }<span class="keywordflow">while</span>( (xout += m) != Fout_end );</div> <div class="line"><a name="l00080"></a><span class="lineno"> 80</span>  }<span class="keywordflow">else</span>{</div> <div class="line"><a name="l00081"></a><span class="lineno"> 81</span>  <span class="keywordflow">do</span>{</div> <div class="line"><a name="l00082"></a><span class="lineno"> 82</span>  *xout = *xin;</div> <div class="line"><a name="l00083"></a><span class="lineno"> 83</span>  xin += fstride*in_stride;</div> <div class="line"><a name="l00084"></a><span class="lineno"> 84</span>  }<span class="keywordflow">while</span>(++xout != Fout_end );</div> <div class="line"><a name="l00085"></a><span class="lineno"> 85</span>  }</div> <div class="line"><a name="l00086"></a><span class="lineno"> 86</span>  xout=Fout_beg;</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>  <span class="comment">// recombine the p smaller DFTs </span></div> <div class="line"><a name="l00089"></a><span class="lineno"> 89</span>  <span class="keywordflow">switch</span> (p) {</div> <div class="line"><a name="l00090"></a><span class="lineno"> 90</span>  <span class="keywordflow">case</span> 2: bfly2(xout,fstride,m); <span class="keywordflow">break</span>;</div> <div class="line"><a name="l00091"></a><span class="lineno"> 91</span>  <span class="keywordflow">case</span> 3: bfly3(xout,fstride,m); <span class="keywordflow">break</span>;</div> <div class="line"><a name="l00092"></a><span class="lineno"> 92</span>  <span class="keywordflow">case</span> 4: bfly4(xout,fstride,m); <span class="keywordflow">break</span>;</div> <div class="line"><a name="l00093"></a><span class="lineno"> 93</span>  <span class="keywordflow">case</span> 5: bfly5(xout,fstride,m); <span class="keywordflow">break</span>;</div> <div class="line"><a name="l00094"></a><span class="lineno"> 94</span>  <span class="keywordflow">default</span>: bfly_generic(xout,fstride,m,p); <span class="keywordflow">break</span>;</div> <div class="line"><a name="l00095"></a><span class="lineno"> 95</span>  }</div> <div class="line"><a name="l00096"></a><span class="lineno"> 96</span>  }</div> <div class="line"><a name="l00097"></a><span class="lineno"> 97</span> </div> <div class="line"><a name="l00098"></a><span class="lineno"> 98</span>  <span class="keyword">inline</span></div> <div class="line"><a name="l00099"></a><span class="lineno"> 99</span>  <span class="keywordtype">void</span> bfly2( Complex * Fout, <span class="keyword">const</span> <span class="keywordtype">size_t</span> fstride, <span class="keywordtype">int</span> m)</div> <div class="line"><a name="l00100"></a><span class="lineno"> 100</span>  {</div> <div class="line"><a name="l00101"></a><span class="lineno"> 101</span>  <span class="keywordflow">for</span> (<span class="keywordtype">int</span> k=0;k<m;++k) {</div> <div class="line"><a name="l00102"></a><span class="lineno"> 102</span>  Complex t = Fout[m+k] * m_twiddles[k*fstride];</div> <div class="line"><a name="l00103"></a><span class="lineno"> 103</span>  Fout[m+k] = Fout[k] - t;</div> <div class="line"><a name="l00104"></a><span class="lineno"> 104</span>  Fout[k] += t;</div> <div class="line"><a name="l00105"></a><span class="lineno"> 105</span>  }</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> </div> <div class="line"><a name="l00108"></a><span class="lineno"> 108</span>  <span class="keyword">inline</span></div> <div class="line"><a name="l00109"></a><span class="lineno"> 109</span>  <span class="keywordtype">void</span> bfly4( Complex * Fout, <span class="keyword">const</span> <span class="keywordtype">size_t</span> fstride, <span class="keyword">const</span> <span class="keywordtype">size_t</span> m)</div> <div class="line"><a name="l00110"></a><span class="lineno"> 110</span>  {</div> <div class="line"><a name="l00111"></a><span class="lineno"> 111</span>  Complex scratch[6];</div> <div class="line"><a name="l00112"></a><span class="lineno"> 112</span>  <span class="keywordtype">int</span> negative_if_inverse = m_inverse * -2 +1;</div> <div class="line"><a name="l00113"></a><span class="lineno"> 113</span>  <span class="keywordflow">for</span> (<span class="keywordtype">size_t</span> k=0;k<m;++k) {</div> <div class="line"><a name="l00114"></a><span class="lineno"> 114</span>  scratch[0] = Fout[k+m] * m_twiddles[k*fstride];</div> <div class="line"><a name="l00115"></a><span class="lineno"> 115</span>  scratch[1] = Fout[k+2*m] * m_twiddles[k*fstride*2];</div> <div class="line"><a name="l00116"></a><span class="lineno"> 116</span>  scratch[2] = Fout[k+3*m] * m_twiddles[k*fstride*3];</div> <div class="line"><a name="l00117"></a><span class="lineno"> 117</span>  scratch[5] = Fout[k] - scratch[1];</div> <div class="line"><a name="l00118"></a><span class="lineno"> 118</span> </div> <div class="line"><a name="l00119"></a><span class="lineno"> 119</span>  Fout[k] += scratch[1];</div> <div class="line"><a name="l00120"></a><span class="lineno"> 120</span>  scratch[3] = scratch[0] + scratch[2];</div> <div class="line"><a name="l00121"></a><span class="lineno"> 121</span>  scratch[4] = scratch[0] - scratch[2];</div> <div class="line"><a name="l00122"></a><span class="lineno"> 122</span>  scratch[4] = Complex( scratch[4].imag()*negative_if_inverse , -scratch[4].real()* negative_if_inverse );</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>  Fout[k+2*m] = Fout[k] - scratch[3];</div> <div class="line"><a name="l00125"></a><span class="lineno"> 125</span>  Fout[k] += scratch[3];</div> <div class="line"><a name="l00126"></a><span class="lineno"> 126</span>  Fout[k+m] = scratch[5] + scratch[4];</div> <div class="line"><a name="l00127"></a><span class="lineno"> 127</span>  Fout[k+3*m] = scratch[5] - scratch[4];</div> <div class="line"><a name="l00128"></a><span class="lineno"> 128</span>  }</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> </div> <div class="line"><a name="l00131"></a><span class="lineno"> 131</span>  <span class="keyword">inline</span></div> <div class="line"><a name="l00132"></a><span class="lineno"> 132</span>  <span class="keywordtype">void</span> bfly3( Complex * Fout, <span class="keyword">const</span> <span class="keywordtype">size_t</span> fstride, <span class="keyword">const</span> <span class="keywordtype">size_t</span> m)</div> <div class="line"><a name="l00133"></a><span class="lineno"> 133</span>  {</div> <div class="line"><a name="l00134"></a><span class="lineno"> 134</span>  <span class="keywordtype">size_t</span> k=m;</div> <div class="line"><a name="l00135"></a><span class="lineno"> 135</span>  <span class="keyword">const</span> <span class="keywordtype">size_t</span> m2 = 2*m;</div> <div class="line"><a name="l00136"></a><span class="lineno"> 136</span>  Complex *tw1,*tw2;</div> <div class="line"><a name="l00137"></a><span class="lineno"> 137</span>  Complex scratch[5];</div> <div class="line"><a name="l00138"></a><span class="lineno"> 138</span>  Complex epi3;</div> <div class="line"><a name="l00139"></a><span class="lineno"> 139</span>  epi3 = m_twiddles[fstride*m];</div> <div class="line"><a name="l00140"></a><span class="lineno"> 140</span> </div> <div class="line"><a name="l00141"></a><span class="lineno"> 141</span>  tw1=tw2=&m_twiddles[0];</div> <div class="line"><a name="l00142"></a><span class="lineno"> 142</span> </div> <div class="line"><a name="l00143"></a><span class="lineno"> 143</span>  <span class="keywordflow">do</span>{</div> <div class="line"><a name="l00144"></a><span class="lineno"> 144</span>  scratch[1]=Fout[m] * *tw1;</div> <div class="line"><a name="l00145"></a><span class="lineno"> 145</span>  scratch[2]=Fout[m2] * *tw2;</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>  scratch[3]=scratch[1]+scratch[2];</div> <div class="line"><a name="l00148"></a><span class="lineno"> 148</span>  scratch[0]=scratch[1]-scratch[2];</div> <div class="line"><a name="l00149"></a><span class="lineno"> 149</span>  tw1 += fstride;</div> <div class="line"><a name="l00150"></a><span class="lineno"> 150</span>  tw2 += fstride*2;</div> <div class="line"><a name="l00151"></a><span class="lineno"> 151</span>  Fout[m] = Complex( Fout->real() - Scalar(.5)*scratch[3].real() , Fout->imag() - Scalar(.5)*scratch[3].imag() );</div> <div class="line"><a name="l00152"></a><span class="lineno"> 152</span>  scratch[0] *= epi3.imag();</div> <div class="line"><a name="l00153"></a><span class="lineno"> 153</span>  *Fout += scratch[3];</div> <div class="line"><a name="l00154"></a><span class="lineno"> 154</span>  Fout[m2] = Complex( Fout[m].real() + scratch[0].imag() , Fout[m].imag() - scratch[0].real() );</div> <div class="line"><a name="l00155"></a><span class="lineno"> 155</span>  Fout[m] += Complex( -scratch[0].imag(),scratch[0].real() );</div> <div class="line"><a name="l00156"></a><span class="lineno"> 156</span>  ++Fout;</div> <div class="line"><a name="l00157"></a><span class="lineno"> 157</span>  }<span class="keywordflow">while</span>(--k);</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> </div> <div class="line"><a name="l00160"></a><span class="lineno"> 160</span>  <span class="keyword">inline</span></div> <div class="line"><a name="l00161"></a><span class="lineno"> 161</span>  <span class="keywordtype">void</span> bfly5( Complex * Fout, <span class="keyword">const</span> <span class="keywordtype">size_t</span> fstride, <span class="keyword">const</span> <span class="keywordtype">size_t</span> m)</div> <div class="line"><a name="l00162"></a><span class="lineno"> 162</span>  {</div> <div class="line"><a name="l00163"></a><span class="lineno"> 163</span>  Complex *Fout0,*Fout1,*Fout2,*Fout3,*Fout4;</div> <div class="line"><a name="l00164"></a><span class="lineno"> 164</span>  <span class="keywordtype">size_t</span> u;</div> <div class="line"><a name="l00165"></a><span class="lineno"> 165</span>  Complex scratch[13];</div> <div class="line"><a name="l00166"></a><span class="lineno"> 166</span>  Complex * twiddles = &m_twiddles[0];</div> <div class="line"><a name="l00167"></a><span class="lineno"> 167</span>  Complex *tw;</div> <div class="line"><a name="l00168"></a><span class="lineno"> 168</span>  Complex ya,yb;</div> <div class="line"><a name="l00169"></a><span class="lineno"> 169</span>  ya = twiddles[fstride*m];</div> <div class="line"><a name="l00170"></a><span class="lineno"> 170</span>  yb = twiddles[fstride*2*m];</div> <div class="line"><a name="l00171"></a><span class="lineno"> 171</span> </div> <div class="line"><a name="l00172"></a><span class="lineno"> 172</span>  Fout0=Fout;</div> <div class="line"><a name="l00173"></a><span class="lineno"> 173</span>  Fout1=Fout0+m;</div> <div class="line"><a name="l00174"></a><span class="lineno"> 174</span>  Fout2=Fout0+2*m;</div> <div class="line"><a name="l00175"></a><span class="lineno"> 175</span>  Fout3=Fout0+3*m;</div> <div class="line"><a name="l00176"></a><span class="lineno"> 176</span>  Fout4=Fout0+4*m;</div> <div class="line"><a name="l00177"></a><span class="lineno"> 177</span> </div> <div class="line"><a name="l00178"></a><span class="lineno"> 178</span>  tw=twiddles;</div> <div class="line"><a name="l00179"></a><span class="lineno"> 179</span>  <span class="keywordflow">for</span> ( u=0; u<m; ++u ) {</div> <div class="line"><a name="l00180"></a><span class="lineno"> 180</span>  scratch[0] = *Fout0;</div> <div class="line"><a name="l00181"></a><span class="lineno"> 181</span> </div> <div class="line"><a name="l00182"></a><span class="lineno"> 182</span>  scratch[1] = *Fout1 * tw[u*fstride];</div> <div class="line"><a name="l00183"></a><span class="lineno"> 183</span>  scratch[2] = *Fout2 * tw[2*u*fstride];</div> <div class="line"><a name="l00184"></a><span class="lineno"> 184</span>  scratch[3] = *Fout3 * tw[3*u*fstride];</div> <div class="line"><a name="l00185"></a><span class="lineno"> 185</span>  scratch[4] = *Fout4 * tw[4*u*fstride];</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>  scratch[7] = scratch[1] + scratch[4];</div> <div class="line"><a name="l00188"></a><span class="lineno"> 188</span>  scratch[10] = scratch[1] - scratch[4];</div> <div class="line"><a name="l00189"></a><span class="lineno"> 189</span>  scratch[8] = scratch[2] + scratch[3];</div> <div class="line"><a name="l00190"></a><span class="lineno"> 190</span>  scratch[9] = scratch[2] - scratch[3];</div> <div class="line"><a name="l00191"></a><span class="lineno"> 191</span> </div> <div class="line"><a name="l00192"></a><span class="lineno"> 192</span>  *Fout0 += scratch[7];</div> <div class="line"><a name="l00193"></a><span class="lineno"> 193</span>  *Fout0 += scratch[8];</div> <div class="line"><a name="l00194"></a><span class="lineno"> 194</span> </div> <div class="line"><a name="l00195"></a><span class="lineno"> 195</span>  scratch[5] = scratch[0] + Complex(</div> <div class="line"><a name="l00196"></a><span class="lineno"> 196</span>  (scratch[7].real()*ya.real() ) + (scratch[8].real() *yb.real() ),</div> <div class="line"><a name="l00197"></a><span class="lineno"> 197</span>  (scratch[7].imag()*ya.real()) + (scratch[8].imag()*yb.real())</div> <div class="line"><a name="l00198"></a><span class="lineno"> 198</span>  );</div> <div class="line"><a name="l00199"></a><span class="lineno"> 199</span> </div> <div class="line"><a name="l00200"></a><span class="lineno"> 200</span>  scratch[6] = Complex(</div> <div class="line"><a name="l00201"></a><span class="lineno"> 201</span>  (scratch[10].imag()*ya.imag()) + (scratch[9].imag()*yb.imag()),</div> <div class="line"><a name="l00202"></a><span class="lineno"> 202</span>  -(scratch[10].real()*ya.imag()) - (scratch[9].real()*yb.imag())</div> <div class="line"><a name="l00203"></a><span class="lineno"> 203</span>  );</div> <div class="line"><a name="l00204"></a><span class="lineno"> 204</span> </div> <div class="line"><a name="l00205"></a><span class="lineno"> 205</span>  *Fout1 = scratch[5] - scratch[6];</div> <div class="line"><a name="l00206"></a><span class="lineno"> 206</span>  *Fout4 = scratch[5] + scratch[6];</div> <div class="line"><a name="l00207"></a><span class="lineno"> 207</span> </div> <div class="line"><a name="l00208"></a><span class="lineno"> 208</span>  scratch[11] = scratch[0] +</div> <div class="line"><a name="l00209"></a><span class="lineno"> 209</span>  Complex(</div> <div class="line"><a name="l00210"></a><span class="lineno"> 210</span>  (scratch[7].real()*yb.real()) + (scratch[8].real()*ya.real()),</div> <div class="line"><a name="l00211"></a><span class="lineno"> 211</span>  (scratch[7].imag()*yb.real()) + (scratch[8].imag()*ya.real())</div> <div class="line"><a name="l00212"></a><span class="lineno"> 212</span>  );</div> <div class="line"><a name="l00213"></a><span class="lineno"> 213</span> </div> <div class="line"><a name="l00214"></a><span class="lineno"> 214</span>  scratch[12] = Complex(</div> <div class="line"><a name="l00215"></a><span class="lineno"> 215</span>  -(scratch[10].imag()*yb.imag()) + (scratch[9].imag()*ya.imag()),</div> <div class="line"><a name="l00216"></a><span class="lineno"> 216</span>  (scratch[10].real()*yb.imag()) - (scratch[9].real()*ya.imag())</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> </div> <div class="line"><a name="l00219"></a><span class="lineno"> 219</span>  *Fout2=scratch[11]+scratch[12];</div> <div class="line"><a name="l00220"></a><span class="lineno"> 220</span>  *Fout3=scratch[11]-scratch[12];</div> <div class="line"><a name="l00221"></a><span class="lineno"> 221</span> </div> <div class="line"><a name="l00222"></a><span class="lineno"> 222</span>  ++Fout0;++Fout1;++Fout2;++Fout3;++Fout4;</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>  <span class="comment">/* perform the butterfly for one stage of a mixed radix FFT */</span></div> <div class="line"><a name="l00227"></a><span class="lineno"> 227</span>  <span class="keyword">inline</span></div> <div class="line"><a name="l00228"></a><span class="lineno"> 228</span>  <span class="keywordtype">void</span> bfly_generic(</div> <div class="line"><a name="l00229"></a><span class="lineno"> 229</span>  Complex * Fout,</div> <div class="line"><a name="l00230"></a><span class="lineno"> 230</span>  <span class="keyword">const</span> <span class="keywordtype">size_t</span> fstride,</div> <div class="line"><a name="l00231"></a><span class="lineno"> 231</span>  <span class="keywordtype">int</span> m,</div> <div class="line"><a name="l00232"></a><span class="lineno"> 232</span>  <span class="keywordtype">int</span> p</div> <div class="line"><a name="l00233"></a><span class="lineno"> 233</span>  )</div> <div class="line"><a name="l00234"></a><span class="lineno"> 234</span>  {</div> <div class="line"><a name="l00235"></a><span class="lineno"> 235</span>  <span class="keywordtype">int</span> u,k,q1,q;</div> <div class="line"><a name="l00236"></a><span class="lineno"> 236</span>  Complex * twiddles = &m_twiddles[0];</div> <div class="line"><a name="l00237"></a><span class="lineno"> 237</span>  Complex t;</div> <div class="line"><a name="l00238"></a><span class="lineno"> 238</span>  <span class="keywordtype">int</span> Norig = <span class="keyword">static_cast<</span><span class="keywordtype">int</span><span class="keyword">></span>(m_twiddles.size());</div> <div class="line"><a name="l00239"></a><span class="lineno"> 239</span>  Complex * scratchbuf = &m_scratchBuf[0];</div> <div class="line"><a name="l00240"></a><span class="lineno"> 240</span> </div> <div class="line"><a name="l00241"></a><span class="lineno"> 241</span>  <span class="keywordflow">for</span> ( u=0; u<m; ++u ) {</div> <div class="line"><a name="l00242"></a><span class="lineno"> 242</span>  k=u;</div> <div class="line"><a name="l00243"></a><span class="lineno"> 243</span>  <span class="keywordflow">for</span> ( q1=0 ; q1<p ; ++q1 ) {</div> <div class="line"><a name="l00244"></a><span class="lineno"> 244</span>  scratchbuf[q1] = Fout[ k ];</div> <div class="line"><a name="l00245"></a><span class="lineno"> 245</span>  k += m;</div> <div class="line"><a name="l00246"></a><span class="lineno"> 246</span>  }</div> <div class="line"><a name="l00247"></a><span class="lineno"> 247</span> </div> <div class="line"><a name="l00248"></a><span class="lineno"> 248</span>  k=u;</div> <div class="line"><a name="l00249"></a><span class="lineno"> 249</span>  <span class="keywordflow">for</span> ( q1=0 ; q1<p ; ++q1 ) {</div> <div class="line"><a name="l00250"></a><span class="lineno"> 250</span>  <span class="keywordtype">int</span> twidx=0;</div> <div class="line"><a name="l00251"></a><span class="lineno"> 251</span>  Fout[ k ] = scratchbuf[0];</div> <div class="line"><a name="l00252"></a><span class="lineno"> 252</span>  <span class="keywordflow">for</span> (q=1;q<p;++q ) {</div> <div class="line"><a name="l00253"></a><span class="lineno"> 253</span>  twidx += <span class="keyword">static_cast<</span><span class="keywordtype">int</span><span class="keyword">></span>(fstride) * k;</div> <div class="line"><a name="l00254"></a><span class="lineno"> 254</span>  <span class="keywordflow">if</span> (twidx>=Norig) twidx-=Norig;</div> <div class="line"><a name="l00255"></a><span class="lineno"> 255</span>  t=scratchbuf[q] * twiddles[twidx];</div> <div class="line"><a name="l00256"></a><span class="lineno"> 256</span>  Fout[ k ] += t;</div> <div class="line"><a name="l00257"></a><span class="lineno"> 257</span>  }</div> <div class="line"><a name="l00258"></a><span class="lineno"> 258</span>  k += m;</div> <div class="line"><a name="l00259"></a><span class="lineno"> 259</span>  }</div> <div class="line"><a name="l00260"></a><span class="lineno"> 260</span>  }</div> <div class="line"><a name="l00261"></a><span class="lineno"> 261</span>  }</div> <div class="line"><a name="l00262"></a><span class="lineno"> 262</span> };</div> <div class="line"><a name="l00263"></a><span class="lineno"> 263</span> </div> <div class="line"><a name="l00264"></a><span class="lineno"> 264</span> <span class="keyword">template</span> <<span class="keyword">typename</span> _Scalar></div> <div class="line"><a name="l00265"></a><span class="lineno"> 265</span> <span class="keyword">struct </span>kissfft_impl</div> <div class="line"><a name="l00266"></a><span class="lineno"> 266</span> {</div> <div class="line"><a name="l00267"></a><span class="lineno"> 267</span>  <span class="keyword">typedef</span> _Scalar Scalar;</div> <div class="line"><a name="l00268"></a><span class="lineno"> 268</span>  <span class="keyword">typedef</span> std::complex<Scalar> Complex;</div> <div class="line"><a name="l00269"></a><span class="lineno"> 269</span> </div> <div class="line"><a name="l00270"></a><span class="lineno"> 270</span>  <span class="keywordtype">void</span> clear() </div> <div class="line"><a name="l00271"></a><span class="lineno"> 271</span>  {</div> <div class="line"><a name="l00272"></a><span class="lineno"> 272</span>  m_plans.clear();</div> <div class="line"><a name="l00273"></a><span class="lineno"> 273</span>  m_realTwiddles.clear();</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> </div> <div class="line"><a name="l00276"></a><span class="lineno"> 276</span>  <span class="keyword">inline</span></div> <div class="line"><a name="l00277"></a><span class="lineno"> 277</span>  <span class="keywordtype">void</span> fwd( Complex * dst,<span class="keyword">const</span> Complex *src,<span class="keywordtype">int</span> nfft)</div> <div class="line"><a name="l00278"></a><span class="lineno"> 278</span>  {</div> <div class="line"><a name="l00279"></a><span class="lineno"> 279</span>  get_plan(nfft,<span class="keyword">false</span>).work(0, dst, src, 1,1);</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> </div> <div class="line"><a name="l00282"></a><span class="lineno"> 282</span>  <span class="keyword">inline</span></div> <div class="line"><a name="l00283"></a><span class="lineno"> 283</span>  <span class="keywordtype">void</span> fwd2( Complex * dst,<span class="keyword">const</span> Complex *src,<span class="keywordtype">int</span> n0,<span class="keywordtype">int</span> n1)</div> <div class="line"><a name="l00284"></a><span class="lineno"> 284</span>  {</div> <div class="line"><a name="l00285"></a><span class="lineno"> 285</span>  EIGEN_UNUSED_VARIABLE(dst);</div> <div class="line"><a name="l00286"></a><span class="lineno"> 286</span>  EIGEN_UNUSED_VARIABLE(src);</div> <div class="line"><a name="l00287"></a><span class="lineno"> 287</span>  EIGEN_UNUSED_VARIABLE(n0);</div> <div class="line"><a name="l00288"></a><span class="lineno"> 288</span>  EIGEN_UNUSED_VARIABLE(n1);</div> <div class="line"><a name="l00289"></a><span class="lineno"> 289</span>  }</div> <div class="line"><a name="l00290"></a><span class="lineno"> 290</span> </div> <div class="line"><a name="l00291"></a><span class="lineno"> 291</span>  <span class="keyword">inline</span></div> <div class="line"><a name="l00292"></a><span class="lineno"> 292</span>  <span class="keywordtype">void</span> inv2( Complex * dst,<span class="keyword">const</span> Complex *src,<span class="keywordtype">int</span> n0,<span class="keywordtype">int</span> n1)</div> <div class="line"><a name="l00293"></a><span class="lineno"> 293</span>  {</div> <div class="line"><a name="l00294"></a><span class="lineno"> 294</span>  EIGEN_UNUSED_VARIABLE(dst);</div> <div class="line"><a name="l00295"></a><span class="lineno"> 295</span>  EIGEN_UNUSED_VARIABLE(src);</div> <div class="line"><a name="l00296"></a><span class="lineno"> 296</span>  EIGEN_UNUSED_VARIABLE(n0);</div> <div class="line"><a name="l00297"></a><span class="lineno"> 297</span>  EIGEN_UNUSED_VARIABLE(n1);</div> <div class="line"><a name="l00298"></a><span class="lineno"> 298</span>  }</div> <div class="line"><a name="l00299"></a><span class="lineno"> 299</span> </div> <div class="line"><a name="l00300"></a><span class="lineno"> 300</span>  <span class="comment">// real-to-complex forward FFT</span></div> <div class="line"><a name="l00301"></a><span class="lineno"> 301</span>  <span class="comment">// perform two FFTs of src even and src odd</span></div> <div class="line"><a name="l00302"></a><span class="lineno"> 302</span>  <span class="comment">// then twiddle to recombine them into the half-spectrum format</span></div> <div class="line"><a name="l00303"></a><span class="lineno"> 303</span>  <span class="comment">// then fill in the conjugate symmetric half</span></div> <div class="line"><a name="l00304"></a><span class="lineno"> 304</span>  <span class="keyword">inline</span></div> <div class="line"><a name="l00305"></a><span class="lineno"> 305</span>  <span class="keywordtype">void</span> fwd( Complex * dst,<span class="keyword">const</span> Scalar * src,<span class="keywordtype">int</span> nfft) </div> <div class="line"><a name="l00306"></a><span class="lineno"> 306</span>  {</div> <div class="line"><a name="l00307"></a><span class="lineno"> 307</span>  <span class="keywordflow">if</span> ( nfft&3 ) {</div> <div class="line"><a name="l00308"></a><span class="lineno"> 308</span>  <span class="comment">// use generic mode for odd</span></div> <div class="line"><a name="l00309"></a><span class="lineno"> 309</span>  m_tmpBuf1.resize(nfft);</div> <div class="line"><a name="l00310"></a><span class="lineno"> 310</span>  get_plan(nfft,<span class="keyword">false</span>).work(0, &m_tmpBuf1[0], src, 1,1);</div> <div class="line"><a name="l00311"></a><span class="lineno"> 311</span>  std::copy(m_tmpBuf1.begin(),m_tmpBuf1.begin()+(nfft>>1)+1,dst );</div> <div class="line"><a name="l00312"></a><span class="lineno"> 312</span>  }<span class="keywordflow">else</span>{</div> <div class="line"><a name="l00313"></a><span class="lineno"> 313</span>  <span class="keywordtype">int</span> ncfft = nfft>>1;</div> <div class="line"><a name="l00314"></a><span class="lineno"> 314</span>  <span class="keywordtype">int</span> ncfft2 = nfft>>2;</div> <div class="line"><a name="l00315"></a><span class="lineno"> 315</span>  Complex * rtw = real_twiddles(ncfft2);</div> <div class="line"><a name="l00316"></a><span class="lineno"> 316</span> </div> <div class="line"><a name="l00317"></a><span class="lineno"> 317</span>  <span class="comment">// use optimized mode for even real</span></div> <div class="line"><a name="l00318"></a><span class="lineno"> 318</span>  fwd( dst, reinterpret_cast<const Complex*> (src), ncfft);</div> <div class="line"><a name="l00319"></a><span class="lineno"> 319</span>  Complex dc = dst[0].real() + dst[0].imag();</div> <div class="line"><a name="l00320"></a><span class="lineno"> 320</span>  Complex nyquist = dst[0].real() - dst[0].imag();</div> <div class="line"><a name="l00321"></a><span class="lineno"> 321</span>  <span class="keywordtype">int</span> k;</div> <div class="line"><a name="l00322"></a><span class="lineno"> 322</span>  <span class="keywordflow">for</span> ( k=1;k <= ncfft2 ; ++k ) {</div> <div class="line"><a name="l00323"></a><span class="lineno"> 323</span>  Complex fpk = dst[k];</div> <div class="line"><a name="l00324"></a><span class="lineno"> 324</span>  Complex fpnk = conj(dst[ncfft-k]);</div> <div class="line"><a name="l00325"></a><span class="lineno"> 325</span>  Complex f1k = fpk + fpnk;</div> <div class="line"><a name="l00326"></a><span class="lineno"> 326</span>  Complex f2k = fpk - fpnk;</div> <div class="line"><a name="l00327"></a><span class="lineno"> 327</span>  Complex tw= f2k * rtw[k-1];</div> <div class="line"><a name="l00328"></a><span class="lineno"> 328</span>  dst[k] = (f1k + tw) * Scalar(.5);</div> <div class="line"><a name="l00329"></a><span class="lineno"> 329</span>  dst[ncfft-k] = conj(f1k -tw)*Scalar(.5);</div> <div class="line"><a name="l00330"></a><span class="lineno"> 330</span>  }</div> <div class="line"><a name="l00331"></a><span class="lineno"> 331</span>  dst[0] = dc;</div> <div class="line"><a name="l00332"></a><span class="lineno"> 332</span>  dst[ncfft] = nyquist;</div> <div class="line"><a name="l00333"></a><span class="lineno"> 333</span>  }</div> <div class="line"><a name="l00334"></a><span class="lineno"> 334</span>  }</div> <div class="line"><a name="l00335"></a><span class="lineno"> 335</span> </div> <div class="line"><a name="l00336"></a><span class="lineno"> 336</span>  <span class="comment">// inverse complex-to-complex</span></div> <div class="line"><a name="l00337"></a><span class="lineno"> 337</span>  <span class="keyword">inline</span></div> <div class="line"><a name="l00338"></a><span class="lineno"> 338</span>  <span class="keywordtype">void</span> inv(Complex * dst,<span class="keyword">const</span> Complex *src,<span class="keywordtype">int</span> nfft)</div> <div class="line"><a name="l00339"></a><span class="lineno"> 339</span>  {</div> <div class="line"><a name="l00340"></a><span class="lineno"> 340</span>  get_plan(nfft,<span class="keyword">true</span>).work(0, dst, src, 1,1);</div> <div class="line"><a name="l00341"></a><span class="lineno"> 341</span>  }</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="comment">// half-complex to scalar</span></div> <div class="line"><a name="l00344"></a><span class="lineno"> 344</span>  <span class="keyword">inline</span></div> <div class="line"><a name="l00345"></a><span class="lineno"> 345</span>  <span class="keywordtype">void</span> inv( Scalar * dst,<span class="keyword">const</span> Complex * src,<span class="keywordtype">int</span> nfft) </div> <div class="line"><a name="l00346"></a><span class="lineno"> 346</span>  {</div> <div class="line"><a name="l00347"></a><span class="lineno"> 347</span>  <span class="keywordflow">if</span> (nfft&3) {</div> <div class="line"><a name="l00348"></a><span class="lineno"> 348</span>  m_tmpBuf1.resize(nfft);</div> <div class="line"><a name="l00349"></a><span class="lineno"> 349</span>  m_tmpBuf2.resize(nfft);</div> <div class="line"><a name="l00350"></a><span class="lineno"> 350</span>  std::copy(src,src+(nfft>>1)+1,m_tmpBuf1.begin() );</div> <div class="line"><a name="l00351"></a><span class="lineno"> 351</span>  <span class="keywordflow">for</span> (<span class="keywordtype">int</span> k=1;k<(nfft>>1)+1;++k)</div> <div class="line"><a name="l00352"></a><span class="lineno"> 352</span>  m_tmpBuf1[nfft-k] = conj(m_tmpBuf1[k]);</div> <div class="line"><a name="l00353"></a><span class="lineno"> 353</span>  inv(&m_tmpBuf2[0],&m_tmpBuf1[0],nfft);</div> <div class="line"><a name="l00354"></a><span class="lineno"> 354</span>  <span class="keywordflow">for</span> (<span class="keywordtype">int</span> k=0;k<nfft;++k)</div> <div class="line"><a name="l00355"></a><span class="lineno"> 355</span>  dst[k] = m_tmpBuf2[k].real();</div> <div class="line"><a name="l00356"></a><span class="lineno"> 356</span>  }<span class="keywordflow">else</span>{</div> <div class="line"><a name="l00357"></a><span class="lineno"> 357</span>  <span class="comment">// optimized version for multiple of 4</span></div> <div class="line"><a name="l00358"></a><span class="lineno"> 358</span>  <span class="keywordtype">int</span> ncfft = nfft>>1;</div> <div class="line"><a name="l00359"></a><span class="lineno"> 359</span>  <span class="keywordtype">int</span> ncfft2 = nfft>>2;</div> <div class="line"><a name="l00360"></a><span class="lineno"> 360</span>  Complex * rtw = real_twiddles(ncfft2);</div> <div class="line"><a name="l00361"></a><span class="lineno"> 361</span>  m_tmpBuf1.resize(ncfft);</div> <div class="line"><a name="l00362"></a><span class="lineno"> 362</span>  m_tmpBuf1[0] = Complex( src[0].real() + src[ncfft].real(), src[0].real() - src[ncfft].real() );</div> <div class="line"><a name="l00363"></a><span class="lineno"> 363</span>  <span class="keywordflow">for</span> (<span class="keywordtype">int</span> k = 1; k <= ncfft / 2; ++k) {</div> <div class="line"><a name="l00364"></a><span class="lineno"> 364</span>  Complex fk = src[k];</div> <div class="line"><a name="l00365"></a><span class="lineno"> 365</span>  Complex fnkc = conj(src[ncfft-k]);</div> <div class="line"><a name="l00366"></a><span class="lineno"> 366</span>  Complex fek = fk + fnkc;</div> <div class="line"><a name="l00367"></a><span class="lineno"> 367</span>  Complex tmp = fk - fnkc;</div> <div class="line"><a name="l00368"></a><span class="lineno"> 368</span>  Complex fok = tmp * conj(rtw[k-1]);</div> <div class="line"><a name="l00369"></a><span class="lineno"> 369</span>  m_tmpBuf1[k] = fek + fok;</div> <div class="line"><a name="l00370"></a><span class="lineno"> 370</span>  m_tmpBuf1[ncfft-k] = conj(fek - fok);</div> <div class="line"><a name="l00371"></a><span class="lineno"> 371</span>  }</div> <div class="line"><a name="l00372"></a><span class="lineno"> 372</span>  get_plan(ncfft,<span class="keyword">true</span>).work(0, reinterpret_cast<Complex*>(dst), &m_tmpBuf1[0], 1,1);</div> <div class="line"><a name="l00373"></a><span class="lineno"> 373</span>  }</div> <div class="line"><a name="l00374"></a><span class="lineno"> 374</span>  }</div> <div class="line"><a name="l00375"></a><span class="lineno"> 375</span> </div> <div class="line"><a name="l00376"></a><span class="lineno"> 376</span>  <span class="keyword">protected</span>:</div> <div class="line"><a name="l00377"></a><span class="lineno"> 377</span>  <span class="keyword">typedef</span> kiss_cpx_fft<Scalar> PlanData;</div> <div class="line"><a name="l00378"></a><span class="lineno"> 378</span>  <span class="keyword">typedef</span> std::map<int,PlanData> PlanMap;</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>  PlanMap m_plans;</div> <div class="line"><a name="l00381"></a><span class="lineno"> 381</span>  std::map<int, std::vector<Complex> > m_realTwiddles;</div> <div class="line"><a name="l00382"></a><span class="lineno"> 382</span>  std::vector<Complex> m_tmpBuf1;</div> <div class="line"><a name="l00383"></a><span class="lineno"> 383</span>  std::vector<Complex> m_tmpBuf2;</div> <div class="line"><a name="l00384"></a><span class="lineno"> 384</span> </div> <div class="line"><a name="l00385"></a><span class="lineno"> 385</span>  <span class="keyword">inline</span></div> <div class="line"><a name="l00386"></a><span class="lineno"> 386</span>  <span class="keywordtype">int</span> PlanKey(<span class="keywordtype">int</span> nfft, <span class="keywordtype">bool</span> isinverse)<span class="keyword"> const </span>{ <span class="keywordflow">return</span> (nfft<<1) | int(isinverse); }</div> <div class="line"><a name="l00387"></a><span class="lineno"> 387</span> </div> <div class="line"><a name="l00388"></a><span class="lineno"> 388</span>  <span class="keyword">inline</span></div> <div class="line"><a name="l00389"></a><span class="lineno"> 389</span>  PlanData & get_plan(<span class="keywordtype">int</span> nfft, <span class="keywordtype">bool</span> inverse)</div> <div class="line"><a name="l00390"></a><span class="lineno"> 390</span>  {</div> <div class="line"><a name="l00391"></a><span class="lineno"> 391</span>  <span class="comment">// TODO look for PlanKey(nfft, ! inverse) and conjugate the twiddles</span></div> <div class="line"><a name="l00392"></a><span class="lineno"> 392</span>  PlanData & pd = m_plans[ PlanKey(nfft,inverse) ];</div> <div class="line"><a name="l00393"></a><span class="lineno"> 393</span>  <span class="keywordflow">if</span> ( pd.m_twiddles.size() == 0 ) {</div> <div class="line"><a name="l00394"></a><span class="lineno"> 394</span>  pd.make_twiddles(nfft,inverse);</div> <div class="line"><a name="l00395"></a><span class="lineno"> 395</span>  pd.factorize(nfft);</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="keywordflow">return</span> pd;</div> <div class="line"><a name="l00398"></a><span class="lineno"> 398</span>  }</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>  <span class="keyword">inline</span></div> <div class="line"><a name="l00401"></a><span class="lineno"> 401</span>  Complex * real_twiddles(<span class="keywordtype">int</span> ncfft2)</div> <div class="line"><a name="l00402"></a><span class="lineno"> 402</span>  {</div> <div class="line"><a name="l00403"></a><span class="lineno"> 403</span>  <span class="keyword">using</span> std::acos;</div> <div class="line"><a name="l00404"></a><span class="lineno"> 404</span>  std::vector<Complex> & twidref = m_realTwiddles[ncfft2];<span class="comment">// creates new if not there</span></div> <div class="line"><a name="l00405"></a><span class="lineno"> 405</span>  <span class="keywordflow">if</span> ( (<span class="keywordtype">int</span>)twidref.size() != ncfft2 ) {</div> <div class="line"><a name="l00406"></a><span class="lineno"> 406</span>  twidref.resize(ncfft2);</div> <div class="line"><a name="l00407"></a><span class="lineno"> 407</span>  <span class="keywordtype">int</span> ncfft= ncfft2<<1;</div> <div class="line"><a name="l00408"></a><span class="lineno"> 408</span>  Scalar pi = acos( Scalar(-1) );</div> <div class="line"><a name="l00409"></a><span class="lineno"> 409</span>  <span class="keywordflow">for</span> (<span class="keywordtype">int</span> k=1;k<=ncfft2;++k) </div> <div class="line"><a name="l00410"></a><span class="lineno"> 410</span>  twidref[k-1] = exp( Complex(0,-pi * (Scalar(k) / ncfft + Scalar(.5)) ) );</div> <div class="line"><a name="l00411"></a><span class="lineno"> 411</span>  }</div> <div class="line"><a name="l00412"></a><span class="lineno"> 412</span>  <span class="keywordflow">return</span> &twidref[0];</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> };</div> <div class="line"><a name="l00415"></a><span class="lineno"> 415</span> </div> <div class="line"><a name="l00416"></a><span class="lineno"> 416</span> } <span class="comment">// end namespace internal</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">// end namespace Eigen</span></div> <div class="line"><a name="l00419"></a><span class="lineno"> 419</span> </div> <div class="line"><a name="l00420"></a><span class="lineno"> 420</span> <span class="comment">/* vim: set filetype=cpp et sw=2 ts=2 ai: */</span></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_70b2be79c95c9d5bfaa4c2dafa46bf10.html">unsupported</a></li><li class="navelem"><a class="el" href="dir_f12b092121fb86d54df52b635b2d8129.html">Eigen</a></li><li class="navelem"><a class="el" href="dir_756fd3610c3abb5994ea9c814224d188.html">src</a></li><li class="navelem"><a class="el" href="dir_b0ccb2a6e9eea075eb5c733286d59031.html">FFT</a></li><li class="navelem"><b>ei_kissfft_impl.h</b></li> <li class="footer">Generated on Mon Oct 28 2013 11:05:27 for Eigen-unsupported 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>