120 for(
int i = 0; i < log2n; i++)
135 unsigned int inline log2(
unsigned int x)
164 template<
class ForwardIterator>
167 ForwardIterator last)
169 ForwardIterator middle = first + ((slip::cardinal<std::size_t>(first,last) + 1) / 2);
170 std::rotate(first,middle,last);
189 template<
class ForwardIterator2D>
192 ForwardIterator2D bottom_right)
195 typename ForwardIterator2D::difference_type size2d = bottom_right - upper_left;
197 std::size_t nb_lig = size2d[0];
198 std::size_t nb_col = size2d[1];
202 for(
size_t i = 0; i < nb_lig; ++i)
206 for(
size_t j = 0; j < nb_col; ++j)
228 template<
class ForwardIterator3D>
231 ForwardIterator3D back_bottom_right)
234 typename ForwardIterator3D::difference_type size3d = back_bottom_right - front_upper_left;
236 std::size_t nb_sli = size3d[0];
237 std::size_t nb_lig = size3d[1];
238 std::size_t nb_col = size3d[2];
240 for(
size_t k = 0; k < nb_sli; ++k)
241 for(
size_t i = 0; i < nb_lig; ++i)
243 slip::fftshift(front_upper_left.row_begin(k,i),front_upper_left.row_end(k,i));
245 for(
size_t k = 0; k < nb_sli; ++k)
246 for(
size_t j = 0; j < nb_col; ++j)
248 slip::fftshift(front_upper_left.col_begin(k,j),front_upper_left.col_end(k,j));
250 for(
size_t i = 0; i < nb_lig; ++i)
251 for(
size_t j = 0; j < nb_col; ++j)
253 slip::fftshift(front_upper_left.slice_begin(i,j),front_upper_left.slice_end(i,j));
290 template<
class InputIter_T,
class Iter_T>
294 typedef typename std::iterator_traits<Iter_T>::value_type complex;
295 typedef typename complex::value_type
real;
299 const complex J(0,1);
300 unsigned int n = 1 << log2n;
301 for(
unsigned int i = 0; i < n; ++i)
305 for(
int s = 1; s <= log2n; ++s)
307 unsigned int m = 1 << s;
308 unsigned int m2 = m >> 1;
310 complex wm =
exp(-J*(pi/m2));
311 for(
unsigned int j = 0; j < m2; ++j)
313 for(
unsigned int k = j ; k < n; k += m)
315 complex t = w * b[k + m2];
349 template<
class InputIter_T,
class Iter_T>
353 typedef typename std::iterator_traits<InputIter_T>::value_type complex;
354 typedef typename complex::value_type
real;
356 const complex J(0,1);
357 unsigned int n = 1 << log2n;
358 real coef = 1.0 / (
real)n;
359 for(
unsigned int i = 0; i < n; ++i)
363 for(
int s = 1; s <= log2n; ++s)
365 unsigned int m = 1 << s;
366 unsigned int m2 = m >> 1;
368 complex wm =
exp(J*(pi/m2));
369 for(
unsigned int j = 0; j < m2; ++j)
371 for(
unsigned int k = j ; k < n; k += m)
373 complex t = w * b[k + m2];
381 for(
unsigned int i = 0; i < n; ++i)
405 template<
class InputIter_T,
class Iter_T>
448 template<
typename InputIterator1d,
typename OuputIterator1d>
450 void radix2_fft1d2d(InputIterator1d begin1,OuputIterator1d begin2,std::size_t nb_lig, std::size_t nb_col)
454 typename std::iterator_traits<InputIterator1d>::difference_type step_horiz1 = 1;
455 typename std::iterator_traits<OuputIterator1d>::difference_type step_horiz2 = 1;
456 typename std::iterator_traits<OuputIterator1d>::difference_type step_vert2 = nb_col;
459 for(
size_t i = 0; i < nb_lig; ++i)
466 for(
size_t j = 0; j < nb_col; ++j)
505 template<
typename InputIterator1d,
typename OuputIterator1d>
507 void radix2_ifft1d2d(InputIterator1d begin1,OuputIterator1d begin2,std::size_t nb_lig, std::size_t nb_col)
511 typename std::iterator_traits<InputIterator1d>::difference_type step_horiz1 = 1;
512 typename std::iterator_traits<OuputIterator1d>::difference_type step_horiz2 = 1;
513 typename std::iterator_traits<OuputIterator1d>::difference_type step_vert2 = nb_col;
515 for(
size_t i = 0; i < nb_lig; ++i)
522 for(
size_t j = 0; j < nb_col; ++j)
556 template<
typename InputBidirectionalIterator2d,
557 typename OutputBidirectionalIterator2d>
560 InputBidirectionalIterator2d in_bottom_right,
561 OutputBidirectionalIterator2d out_upper_left)
564 typename InputBidirectionalIterator2d::difference_type size2d = in_bottom_right - in_upper_left;
566 std::size_t nb_lig = size2d[0];
567 std::size_t nb_col = size2d[1];
575 for(
size_t i = 0; i < nb_lig; ++i)
577 slip::radix2_fft(in_upper_left.row_begin(i),out_upper_left.row_begin(i),n2);
581 for(
size_t j = 0; j < nb_col; ++j)
617 template<
typename InputBidirectionalIterator2d,
618 typename OutputBidirectionalIterator2d>
621 InputBidirectionalIterator2d in_bottom_right,
622 OutputBidirectionalIterator2d out_upper_left)
625 typename InputBidirectionalIterator2d::difference_type size2d = in_bottom_right - in_upper_left;
627 std::size_t nb_lig = size2d[0];
628 std::size_t nb_col = size2d[1];
639 for(
size_t i = 0; i < nb_lig; ++i)
645 for(
size_t j = 0; j < nb_col; ++j)
684 template<
typename Matrix1,
typename Matrix2>
688 assert(datain.cols() == dataout.cols());
689 assert(datain.rows() == dataout.rows());
723 template<
typename Matrix1,
typename Matrix2>
727 assert(datain.cols() == dataout.cols());
728 assert(datain.rows() == dataout.rows());
772 template<
class InputIter,
class OutputIter>
776 typedef typename std::iterator_traits<InputIter>::value_type
real;
777 typedef typename std::iterator_traits<OutputIter>::value_type complex;
778 const complex J(0,1);
779 size_t length = in_end - in_begin;
784 std::copy(in_begin,in_end,itst_in_begin);
788 size_t double_length = 2*length;
789 slip::rffti(&double_length,r.
begin(),ifac.
begin());
795 out_begin[0] = in[0];
796 for(
size_t i = 1; i < length; i++)
798 out_begin[i] = in[2*i-1] + J * in[2*i];
800 if(in[0] != in[2*length - 1])
801 out_begin[0]+= J * in[2*length - 1];
826 template<
class InputIter,
class OutputIter>
830 typedef typename std::iterator_traits<InputIter>::value_type complex;
831 typedef typename complex::value_type
real;
832 const complex J(0,1);
833 size_t length = in_end - in_begin;
852 for(
size_t i = 0; i < length; i++)
854 out_begin[i] = in[2*i] + J * in[2*i+1];
878 template<
class InputIter,
class OutputIter>
882 typedef typename std::iterator_traits<InputIter>::value_type complex;
883 typedef typename complex::value_type
real;
884 const complex J(0,1);
885 size_t length = in_end - in_begin;
904 for(
size_t i = 0; i < length; i++)
906 out_begin[i] = (in[2*i] + J * in[2*i+1]) / (real)length;
941 template<
typename InputBidirectionalIterator2d,
942 typename OutputBidirectionalIterator2d>
945 InputBidirectionalIterator2d in_bottom_right,
946 OutputBidirectionalIterator2d out_upper_left)
949 typename InputBidirectionalIterator2d::difference_type size2d = in_bottom_right - in_upper_left;
951 std::size_t nb_lig = size2d[0];
952 std::size_t nb_col = size2d[1];
955 for(
size_t i = 0; i < nb_lig; ++i)
961 for(
size_t j = 0; j < nb_col; ++j)
981 template<
typename InputBidirectionalIterator2d,
982 typename OutputBidirectionalIterator2d>
985 InputBidirectionalIterator2d in_bottom_right,
986 OutputBidirectionalIterator2d out_upper_left)
989 typename InputBidirectionalIterator2d::difference_type size2d = in_bottom_right - in_upper_left;
991 std::size_t nb_lig = size2d[0];
992 std::size_t nb_col = size2d[1];
995 for(
size_t i = 0; i < nb_lig; ++i)
1001 for(
size_t j = 0; j < nb_col; ++j)
1036 template<
typename InputBidirectionalIterator2d,
1037 typename OutputBidirectionalIterator2d>
1040 InputBidirectionalIterator2d in_bottom_right,
1041 OutputBidirectionalIterator2d out_upper_left)
1044 typename InputBidirectionalIterator2d::difference_type size2d = in_bottom_right - in_upper_left;
1046 std::size_t nb_lig = size2d[0];
1047 std::size_t nb_col = size2d[1];
1054 for(
size_t i = 0; i < nb_lig; ++i)
1060 for(
size_t j = 0; j < nb_col; ++j)
1099 template<
typename InputMatrix,
typename OutputMatrix>
1102 assert(datain.cols() == dataout.cols());
1103 assert(datain.rows() == dataout.rows());
1136 template<
typename InputMatrix,
typename OutputMatrix>
1140 assert(datain.cols() == dataout.cols());
1141 assert(datain.rows() == dataout.rows());
1177 template<
typename InputBidirectionalIterator3d,
1178 typename OutputBidirectionalIterator3d>
1181 InputBidirectionalIterator3d in_back_bottom_right,
1182 OutputBidirectionalIterator3d out_front_upper_left)
1185 typename InputBidirectionalIterator3d::difference_type size3d = in_back_bottom_right - in_front_upper_left;
1186 typedef typename std::iterator_traits<InputBidirectionalIterator3d>::value_type
real;
1187 typedef typename std::iterator_traits<OutputBidirectionalIterator3d>::value_type complex;
1189 std::size_t nb_sli = size3d[0];
1190 std::size_t nb_lig = size3d[1];
1191 std::size_t nb_col = size3d[2];
1193 for(
size_t k = 0; k < nb_sli; ++k)
1194 for(
size_t i = 0; i < nb_lig; ++i)
1196 slip::real_split_radix_fft(in_front_upper_left.row_begin(k,i),in_front_upper_left.row_end(k,i),out_front_upper_left.row_begin(k,i));
1200 for(
size_t k = 0; k < nb_sli; ++k)
1201 for(
size_t j = 0; j < nb_col; ++j)
1208 for(
size_t i = 0; i < nb_lig; ++i)
1209 for(
size_t j = 0; j < nb_col; ++j)
1248 template<
typename InputBidirectionalIterator3d,
1249 typename OutputBidirectionalIterator3d>
1252 InputBidirectionalIterator3d in_back_bottom_right,
1253 OutputBidirectionalIterator3d out_front_upper_left)
1256 typename InputBidirectionalIterator3d::difference_type size3d = in_back_bottom_right - in_front_upper_left;
1258 std::size_t nb_sli = size3d[0];
1259 std::size_t nb_lig = size3d[1];
1260 std::size_t nb_col = size3d[2];
1263 for(
size_t k = 0; k < nb_sli; ++k)
1264 for(
size_t i = 0; i < nb_lig; ++i)
1266 slip::split_radix_ifft(in_front_upper_left.row_begin(k,i),in_front_upper_left.row_end(k,i),out_front_upper_left.row_begin(k,i));
1270 for(
size_t k = 0; k < nb_sli; ++k)
1271 for(
size_t j = 0; j < nb_col; ++j)
1278 for(
size_t i = 0; i < nb_lig; ++i)
1279 for(
size_t j = 0; j < nb_col; ++j)
1300 template<
typename InputBidirectionalIterator3d,
1301 typename OutputBidirectionalIterator3d>
1304 InputBidirectionalIterator3d in_back_bottom_right,
1305 OutputBidirectionalIterator3d out_front_upper_left)
1308 typename InputBidirectionalIterator3d::difference_type size3d = in_back_bottom_right - in_front_upper_left;
1309 typedef typename std::iterator_traits<InputBidirectionalIterator3d>::value_type complex;
1310 typedef typename complex::value_type
real;
1313 std::size_t nb_sli = size3d[0];
1314 std::size_t nb_lig = size3d[1];
1315 std::size_t nb_col = size3d[2];
1317 for(
size_t k = 0; k < nb_sli; ++k)
1318 for(
size_t i = 0; i < nb_lig; ++i)
1321 in_front_upper_left.row_end(k,i),out_front_upper_left.row_begin(k,i));
1325 for(
size_t k = 0; k < nb_sli; ++k)
1326 for(
size_t j = 0; j < nb_col; ++j)
1333 for(
size_t i = 0; i < nb_lig; ++i)
1334 for(
size_t j = 0; j < nb_col; ++j)
1356 template<
typename InputMatrix3d,
1357 typename OutputMatrix3d>
1361 assert(datain.cols() == dataout.cols());
1362 assert(datain.rows() == dataout.rows());
1363 assert(datain.slices() == dataout.slices());
1401 template<
typename InputMatrix3d,
1402 typename OutputMatrix3d>
1406 assert(datain.cols() == dataout.cols());
1407 assert(datain.rows() == dataout.rows());
1408 assert(datain.slices() == dataout.slices());
1452 template<
typename InputMatrix3d,
1453 typename OutputMatrix3d>
1457 assert(datain.cols() == dataout.cols());
1458 assert(datain.rows() == dataout.rows());
1459 assert(datain.slices() == dataout.slices());
1490 template<
class InputIter,
class OutputIter>
1494 typedef typename std::iterator_traits<InputIter>::value_type
real;
1495 size_t length = in_end - in_begin;
1507 slip::cosqf(&length,out_begin,r.
begin(),ifac.
begin());
1508 std::transform(out_begin, out_begin+length, out_begin,std::bind2nd(std::divides<real>(),
std::sqrt(2*length)));
1534 template<
class InputIter,
class OutputIter>
1538 typedef typename std::iterator_traits<InputIter>::value_type
real;
1539 size_t length = in_end - in_begin;
1551 slip::cosqb(&length,out_begin,r.
begin(),ifac.
begin());
1552 std::transform(out_begin, out_begin+length, out_begin,std::bind2nd(std::divides<real>(),2*
std::sqrt(2*length)));
1580 template<
class ForwardIterator1,
class ForwardIterator2>
1583 ForwardIterator2 end2)
1585 typedef typename std::iterator_traits<ForwardIterator1>::value_type complex;
1590 while (begin2 < end2)
1592 *end2 = std::real(*begin1) - J*std::imag(*begin1);
1611 template<
class ForwardIterator2D>
1614 ForwardIterator2D bottom_right)
1616 typename ForwardIterator2D::difference_type size2d = bottom_right - upper_left;
1617 std::size_t nb_lig = size2d[0];
1618 slip::fftduplicate(upper_left.row_begin(0),upper_left.row_begin(0),upper_left.row_end(0));
1619 for(
size_t i = 1; i < nb_lig; ++i)
1621 slip::fftduplicate(upper_left.row_begin(i),upper_left.row_begin(nb_lig - i),upper_left.row_end(nb_lig - i));
1639 template<
class InputIter,
class OutputIter>
1641 void real_fftw(InputIter in_begin, InputIter in_end, OutputIter out_begin)
1644 typedef typename std::iterator_traits<OutputIter>::value_type complex;
1645 typedef typename complex::value_type
real;
1647 size_t length = in_end - in_begin;
1648 size_t new_length = length/2 + 1;
1657 plan = fftw_plan_dft_r2c_1d(
int(length),(
double*)&(*in.
begin()),(fftw_complex*)&(*out.
begin()),FFTW_ESTIMATE);
1661 OutputIter itout = out_begin;
1662 while(it < out.
end())
1664 *itout = (
real)std::real(*it) + J * (
real)std::imag(*it);
1671 fftw_destroy_plan(plan);
1688 template<
class InputIter,
class OutputIter>
1690 void complex_fftw(InputIter in_begin, InputIter in_end, OutputIter out_begin)
1692 typedef typename std::iterator_traits<OutputIter>::value_type complex;
1693 typedef typename complex::value_type
real;
1696 size_t length = in_end - in_begin;
1704 plan = fftw_plan_dft_1d(
int(length),(fftw_complex*)&(*in.
begin()),(fftw_complex*)&(*out.
begin()),-1,FFTW_ESTIMATE);
1707 OutputIter itout = out_begin;
1708 while(it < out.
end())
1710 *itout = ((
real)std::real(*it) + J * (
real)std::imag(*it));
1715 fftw_destroy_plan(plan);
1731 template<
class InputIter,
class OutputIter>
1733 void ifftw(InputIter in_begin, InputIter in_end, OutputIter out_begin)
1735 typedef typename std::iterator_traits<OutputIter>::value_type complex;
1736 typedef typename complex::value_type
real;
1738 size_t length = in_end - in_begin;
1746 plan = fftw_plan_dft_1d(
int(length),(fftw_complex*)&(*in.
begin()),(fftw_complex*)&(*out.
begin()),1,FFTW_ESTIMATE);
1749 OutputIter itout = out_begin;
1750 while(it < out.
end())
1752 *itout = ((
real)std::real(*it) + J * (
real)std::imag(*it)) / (real)length;
1757 fftw_destroy_plan(plan);
1775 template<
typename InputBidirectionalIterator2d,
1776 typename OutputBidirectionalIterator2d>
1779 InputBidirectionalIterator2d in_bottom_right,
1780 OutputBidirectionalIterator2d out_upper_left)
1782 typedef typename std::iterator_traits<OutputBidirectionalIterator2d>::value_type complex;
1783 typedef typename complex::value_type
real;
1786 typename InputBidirectionalIterator2d::difference_type size2d = in_bottom_right - in_upper_left;
1787 std::size_t nb_lig = size2d[0];
1788 std::size_t nb_col = size2d[1];
1789 int new_nb_col = (nb_col/2)+1;
1790 assert((nb_lig > 1)&&(nb_col>1));
1799 plan = fftw_plan_dft_r2c_2d(
int(nb_lig),
int(nb_col),(
double*)&(*in.
upper_left()),(fftw_complex*)&(*out.
upper_left()),FFTW_ESTIMATE);
1802 for(
int i=0; i<int(nb_lig);i++)
1803 for(
int j=0; j<new_nb_col;j++)
1806 *(out_upper_left + dp) = (real)std::real(out[i][j]) + J * (
real)std::imag(out[i][j]);
1810 fftw_destroy_plan(plan);
1828 template<
typename InputBidirectionalIterator2d,
1829 typename OutputBidirectionalIterator2d>
1832 InputBidirectionalIterator2d in_bottom_right,
1833 OutputBidirectionalIterator2d out_upper_left)
1835 typedef typename std::iterator_traits<OutputBidirectionalIterator2d>::value_type complex;
1836 typedef typename complex::value_type
real;
1839 typename InputBidirectionalIterator2d::difference_type size2d = in_bottom_right - in_upper_left;
1840 std::size_t nb_lig = size2d[0];
1841 std::size_t nb_col = size2d[1];
1842 assert((nb_lig > 1)&&(nb_col>1));
1849 plan = fftw_plan_dft_2d(
int(nb_lig),
int(nb_col),(fftw_complex*)&(*in.
upper_left()),(fftw_complex*)&(*out.
upper_left()),-1,FFTW_ESTIMATE);
1852 OutputBidirectionalIterator2d itout = out_upper_left;
1853 while(it < out.
end())
1855 *itout = ((
real)std::real(*it) + J * (
real)std::imag(*it));
1860 fftw_destroy_plan(plan);
1878 template<
typename InputBidirectionalIterator2d,
1879 typename OutputBidirectionalIterator2d>
1881 void ifftw2d(InputBidirectionalIterator2d in_upper_left,
1882 InputBidirectionalIterator2d in_bottom_right,
1883 OutputBidirectionalIterator2d out_upper_left)
1885 typedef typename std::iterator_traits<OutputBidirectionalIterator2d>::value_type complex;
1886 typedef typename complex::value_type
real;
1889 typename InputBidirectionalIterator2d::difference_type size2d = in_bottom_right - in_upper_left;
1890 std::size_t nb_lig = size2d[0];
1891 std::size_t nb_col = size2d[1];
1892 int length(nb_lig*nb_col);
1893 assert((nb_lig > 1)&&(nb_col>1));
1900 plan = fftw_plan_dft_2d(
int(nb_lig),
int(nb_col),(fftw_complex*)&(*in.
upper_left()),(fftw_complex*)&(*out.
upper_left()),1,FFTW_ESTIMATE);
1904 OutputBidirectionalIterator2d itout = out_upper_left;
1905 while(it < out.
end())
1907 *itout = ((
real)std::real(*it) + J * (
real)std::imag(*it)) / (real)length;
1912 fftw_destroy_plan(plan);
1932 template<
typename InputMatrix,
typename OutputMatrix>
1935 assert(datain.cols() == dataout.cols());
1936 assert(datain.rows() == dataout.rows());
1937 slip::real_fftw2d(datain.upper_left(),datain.bottom_right(),dataout.upper_left());
1955 template<
typename InputMatrix,
typename OutputMatrix>
1957 void ifftw2d(InputMatrix &datain, OutputMatrix &dataout)
1959 assert(datain.cols() == dataout.cols());
1960 assert(datain.rows() == dataout.rows());
1961 slip::ifftw2d(datain.upper_left(),datain.bottom_right(),dataout.upper_left());
1979 template<
typename InputBidirectionalIterator3d,
1980 typename OutputBidirectionalIterator3d>
1983 InputBidirectionalIterator3d in_back_bottom_right,
1984 OutputBidirectionalIterator3d out_front_upper_left)
1987 typename InputBidirectionalIterator3d::difference_type size3d = in_back_bottom_right - in_front_upper_left;
1988 typedef typename std::iterator_traits<InputBidirectionalIterator3d>::value_type
real;
1989 typedef typename std::iterator_traits<OutputBidirectionalIterator3d>::value_type complex;
1991 std::size_t nb_sli = size3d[0];
1992 std::size_t nb_lig = size3d[1];
1993 std::size_t nb_col = size3d[2];
1994 assert((nb_sli > 1)&&(nb_lig > 1)&&(nb_col>1));
1996 for(
size_t k = 0; k < nb_sli; ++k)
1997 for(
size_t i = 0; i < nb_lig; ++i)
1999 slip::real_fftw(in_front_upper_left.row_begin(k,i),in_front_upper_left.row_end(k,i),out_front_upper_left.row_begin(k,i));
2003 for(
size_t k = 0; k < nb_sli; ++k)
2004 for(
size_t j = 0; j < nb_col; ++j)
2011 for(
size_t i = 0; i < nb_lig; ++i)
2012 for(
size_t j = 0; j < nb_col; ++j)
2034 template<
typename InputBidirectionalIterator3d,
2035 typename OutputBidirectionalIterator3d>
2038 InputBidirectionalIterator3d in_back_bottom_right,
2039 OutputBidirectionalIterator3d out_front_upper_left)
2042 typename InputBidirectionalIterator3d::difference_type size3d = in_back_bottom_right - in_front_upper_left;
2043 typedef typename std::iterator_traits<InputBidirectionalIterator3d>::value_type complex;
2044 typedef typename complex::value_type
real;
2046 std::size_t nb_sli = size3d[0];
2047 std::size_t nb_lig = size3d[1];
2048 std::size_t nb_col = size3d[2];
2049 assert((nb_sli > 1)&&(nb_lig > 1)&&(nb_col>1));
2052 for(
size_t k = 0; k < nb_sli; ++k)
2053 for(
size_t i = 0; i < nb_lig; ++i)
2055 slip::complex_fftw(in_front_upper_left.row_begin(k,i),in_front_upper_left.row_end(k,i),out_front_upper_left.row_begin(k,i));
2059 for(
size_t k = 0; k < nb_sli; ++k)
2060 for(
size_t j = 0; j < nb_col; ++j)
2067 for(
size_t i = 0; i < nb_lig; ++i)
2068 for(
size_t j = 0; j < nb_col; ++j)
2110 template<
typename InputBidirectionalIterator3d,
2111 typename OutputBidirectionalIterator3d>
2114 InputBidirectionalIterator3d in_back_bottom_right,
2115 OutputBidirectionalIterator3d out_front_upper_left)
2162 template<
typename Volume1,
typename Volume2>
2166 assert(datain.cols() == dataout.cols());
2167 assert(datain.rows() == dataout.rows());
2168 assert(datain.slices() == dataout.slices());
2170 slip::complex_fftw3d(datain.front_upper_left(),datain.back_bottom_right(),dataout.front_upper_left());
2191 template<
typename InputBidirectionalIterator3d,
2192 typename OutputBidirectionalIterator3d>
2194 void ifftw3d(InputBidirectionalIterator3d in_front_upper_left,
2195 InputBidirectionalIterator3d in_back_bottom_right,
2196 OutputBidirectionalIterator3d out_front_upper_left)
2199 typename InputBidirectionalIterator3d::difference_type size3d = in_back_bottom_right - in_front_upper_left;
2201 std::size_t nb_sli = size3d[0];
2202 std::size_t nb_lig = size3d[1];
2203 std::size_t nb_col = size3d[2];
2206 for(
size_t k = 0; k < nb_sli; ++k)
2207 for(
size_t i = 0; i < nb_lig; ++i)
2209 slip::ifftw(in_front_upper_left.row_begin(k,i),in_front_upper_left.row_end(k,i),out_front_upper_left.row_begin(k,i));
2213 for(
size_t k = 0; k < nb_sli; ++k)
2214 for(
size_t j = 0; j < nb_col; ++j)
2216 slip::ifftw(out_front_upper_left.col_begin(k,j),out_front_upper_left.col_end(k,j),Vtmp.
begin());
2221 for(
size_t i = 0; i < nb_lig; ++i)
2222 for(
size_t j = 0; j < nb_col; ++j)
2224 slip::ifftw(out_front_upper_left.slice_begin(i,j),out_front_upper_left.slice_end(i,j),Vtmp.
begin());
2245 template<
typename InputMatrix3d,
2246 typename OutputMatrix3d>
2250 assert(datain.cols() == dataout.cols());
2251 assert(datain.rows() == dataout.rows());
2252 assert(datain.slices() == dataout.slices());
2253 slip::real_fftw3d(datain.front_upper_left(),datain.back_bottom_right(),dataout.front_upper_left());
2272 template<
typename InputMatrix3d,
2273 typename OutputMatrix3d>
2277 assert(datain.cols() == dataout.cols());
2278 assert(datain.rows() == dataout.rows());
2279 assert(datain.slices() == dataout.slices());
2280 slip::complex_fftw3d(datain.front_upper_left(),datain.back_bottom_right(),dataout.front_upper_left());
2301 template<
typename InputMatrix3d,
2302 typename OutputMatrix3d>
2304 void ifftw3d(InputMatrix3d & datain,OutputMatrix3d & dataout)
2306 assert(datain.cols() == dataout.cols());
2307 assert(datain.rows() == dataout.rows());
2308 assert(datain.slices() == dataout.slices());
2309 slip::ifftw3d(datain.front_upper_left(),datain.back_bottom_right(),dataout.front_upper_left());
2325 template<
class InputIter,
class OutputIter>
2327 void fftw_dct(InputIter in_begin, InputIter in_end, OutputIter out_begin)
2329 size_t length = in_end - in_begin;
2332 double *in = (
double*) fftw_malloc(
sizeof(
double) * int(length));
2333 double *out = (
double*) fftw_malloc(
sizeof(
double) * int(length));
2335 fftw_r2r_kind kind = FFTW_REDFT10;
2338 plan = fftw_plan_r2r_1d(
int(length),in,out,kind,FFTW_ESTIMATE);
2340 std::transform(out,out+length,out_begin,std::bind2nd(std::divides<double>(),
std::sqrt(2*length)));
2342 fftw_destroy_plan(plan);
2360 template<
class InputIter,
class OutputIter>
2362 void fftw_idct(InputIter in_begin, InputIter in_end, OutputIter out_begin)
2364 size_t length = in_end - in_begin;
2367 double *in = (
double*) fftw_malloc(
sizeof(
double) * int(length));
2368 double *out = (
double*) fftw_malloc(
sizeof(
double) * int(length));
2370 fftw_r2r_kind kind = FFTW_REDFT01;
2373 plan = fftw_plan_r2r_1d(
int(length),in,out,kind,FFTW_ESTIMATE);
2375 std::transform(out,out+length,out_begin,std::bind2nd(std::divides<double>(),
std::sqrt(2*length)));
2377 fftw_destroy_plan(plan);
2428 template<
class InputIter,
class OutputIter>
2430 void real_fft(InputIter in_begin, InputIter in_end, OutputIter out_begin)
2468 template<
class InputIter,
class OutputIter>
2470 void complex_fft(InputIter in_begin, InputIter in_end, OutputIter out_begin)
2505 template<
class InputIter,
class OutputIter>
2507 void ifft(InputIter in_begin, InputIter in_end, OutputIter out_begin)
2550 template<
typename InputBidirectionalIterator2d,
2551 typename OutputBidirectionalIterator2d>
2554 InputBidirectionalIterator2d in_bottom_right,
2555 OutputBidirectionalIterator2d out_upper_left)
2582 template<
typename InputBidirectionalIterator2d,
2583 typename OutputBidirectionalIterator2d>
2586 InputBidirectionalIterator2d in_bottom_right,
2587 OutputBidirectionalIterator2d out_upper_left)
2626 template<
typename InputBidirectionalIterator2d,
2627 typename OutputBidirectionalIterator2d>
2629 void ifft2d(InputBidirectionalIterator2d in_upper_left,
2630 InputBidirectionalIterator2d in_bottom_right,
2631 OutputBidirectionalIterator2d out_upper_left)
2634 slip::ifftw2d(in_upper_left,in_bottom_right,out_upper_left);
2674 template<
typename Matrix1,
typename Matrix2>
2678 assert(datain.cols() == dataout.cols());
2679 assert(datain.rows() == dataout.rows());
2681 slip::real_fftw2d(datain.upper_left(),datain.bottom_right(),dataout.upper_left());
2705 template<
typename Matrix1,
typename Matrix2>
2709 assert(datain.cols() == dataout.cols());
2710 assert(datain.rows() == dataout.rows());
2751 template<
typename Matrix1,
typename Matrix2>
2753 void ifft2d(Matrix1 &datain, Matrix2 &dataout)
2755 assert(datain.cols() == dataout.cols());
2756 assert(datain.rows() == dataout.rows());
2758 slip::ifftw2d(datain.upper_left(),datain.bottom_right(),dataout.upper_left());
2801 template<
typename InputBidirectionalIterator3d,
2802 typename OutputBidirectionalIterator3d>
2804 void real_fft3d(InputBidirectionalIterator3d in_front_upper_left,
2805 InputBidirectionalIterator3d in_back_bottom_right,
2806 OutputBidirectionalIterator3d out_front_upper_left)
2809 slip::real_fftw3d(in_front_upper_left,in_back_bottom_right,out_front_upper_left);
2849 template<
typename InputBidirectionalIterator3d,
2850 typename OutputBidirectionalIterator3d>
2852 void ifft3d(InputBidirectionalIterator3d in_front_upper_left,
2853 InputBidirectionalIterator3d in_back_bottom_right,
2854 OutputBidirectionalIterator3d out_front_upper_left)
2857 slip::ifftw3d(in_front_upper_left,in_back_bottom_right,out_front_upper_left);
2901 template<
typename Volume1,
typename Volume2>
2905 assert(datain.cols() == dataout.cols());
2906 assert(datain.rows() == dataout.rows());
2907 assert(datain.slices() == dataout.slices());
2909 slip::real_fftw3d(datain.front_upper_left(),datain.back_bottom_right(),dataout.front_upper_left());
2953 template<
typename Volume1,
typename Volume2>
2955 void ifft3d(Volume1 &datain, Volume2 &dataout)
2957 assert(datain.cols() == dataout.cols());
2958 assert(datain.rows() == dataout.rows());
2959 assert(datain.slices() == dataout.slices());
2961 slip::ifftw3d(datain.front_upper_left(),datain.back_bottom_right(),dataout.front_upper_left());
2991 template<
class InputIter,
class OutputIter>
2993 void dct(InputIter in_begin, InputIter in_end, OutputIter out_begin)
3027 template<
class InputIter,
class OutputIter>
3029 void idct(InputIter in_begin, InputIter in_end, OutputIter out_begin)
3056 template<
class ForwardIterator4D>
3059 ForwardIterator4D last_back_bottom_right)
3061 typename ForwardIterator4D::difference_type size4d = last_back_bottom_right - first_front_upper_left;
3063 std::size_t nb_sla = size4d[0];
3064 std::size_t nb_sli = size4d[1];
3065 std::size_t nb_lig = size4d[2];
3066 std::size_t nb_col = size4d[3];
3068 for(
size_t t = 0; t < nb_sla; ++t)
3069 for(
size_t k = 0; k < nb_sli; ++k)
3070 for(
size_t i = 0; i < nb_lig; ++i)
3072 slip::fftshift(first_front_upper_left.row_begin(t,k,i),first_front_upper_left.row_end(t,k,i));
3075 for(
size_t t = 0; t < nb_sla; ++t)
3076 for(
size_t k = 0; k < nb_sli; ++k)
3077 for(
size_t j = 0; j < nb_col; ++j)
3079 slip::fftshift(first_front_upper_left.col_begin(t,k,j),first_front_upper_left.col_end(t,k,j));
3081 for(
size_t t = 0; t < nb_sla; ++t)
3082 for(
size_t i = 0; i < nb_lig; ++i)
3083 for(
size_t j = 0; j < nb_col; ++j)
3085 slip::fftshift(first_front_upper_left.slice_begin(t,i,j),first_front_upper_left.slice_end(t,i,j));
3087 for(
size_t k = 0; k < nb_sli; ++k)
3088 for(
size_t i = 0; i < nb_lig; ++i)
3089 for(
size_t j = 0; j < nb_col; ++j)
3091 slip::fftshift(first_front_upper_left.slab_begin(k,i,j),first_front_upper_left.slab_end(k,i,j));
3134 template<
typename InputBidirectionalIterator4d,
3135 typename OutputBidirectionalIterator4d>
3138 InputBidirectionalIterator4d in_last_back_bottom_right,
3139 OutputBidirectionalIterator4d out_first_front_upper_left)
3142 typename InputBidirectionalIterator4d::difference_type size4d = in_last_back_bottom_right - in_first_front_upper_left;
3143 typedef typename std::iterator_traits<InputBidirectionalIterator4d>::value_type
real;
3144 typedef typename std::iterator_traits<OutputBidirectionalIterator4d>::value_type complex;
3146 std::size_t nb_sla = size4d[0];
3147 std::size_t nb_sli = size4d[1];
3148 std::size_t nb_lig = size4d[2];
3149 std::size_t nb_col = size4d[3];
3151 for(
size_t t = 0; t < nb_sla; ++t)
3152 for(
size_t k = 0; k < nb_sli; ++k)
3153 for(
size_t i = 0; i < nb_lig; ++i)
3156 in_first_front_upper_left.row_end(t,k,i),out_first_front_upper_left.row_begin(t,k,i));
3160 for(
size_t t = 0; t < nb_sla; ++t)
3161 for(
size_t k = 0; k < nb_sli; ++k)
3162 for(
size_t j = 0; j < nb_col; ++j)
3165 out_first_front_upper_left.col_end(t,k,j),Vtmp.
begin());
3170 for(
size_t t = 0; t < nb_sla; ++t)
3171 for(
size_t i = 0; i < nb_lig; ++i)
3172 for(
size_t j = 0; j < nb_col; ++j)
3175 out_first_front_upper_left.slice_end(t,i,j),Vtmp.
begin());
3176 std::copy(Vtmp.
begin(),Vtmp.
end(),out_first_front_upper_left.slice_begin(t,i,j));
3180 for(
size_t k = 0; k < nb_sli; ++k)
3181 for(
size_t i = 0; i < nb_lig; ++i)
3182 for(
size_t j = 0; j < nb_col; ++j)
3185 out_first_front_upper_left.slab_end(k,i,j),Vtmp.
begin());
3186 std::copy(Vtmp.
begin(),Vtmp.
end(),out_first_front_upper_left.slab_begin(k,i,j));
3203 template<
typename InputBidirectionalIterator4d,
3204 typename OutputBidirectionalIterator4d>
3207 InputBidirectionalIterator4d in_last_back_bottom_right,
3208 OutputBidirectionalIterator4d out_first_front_upper_left)
3211 typename InputBidirectionalIterator4d::difference_type size4d = in_last_back_bottom_right - in_first_front_upper_left;
3212 typedef typename std::iterator_traits<InputBidirectionalIterator4d>::value_type complex;
3214 std::size_t nb_sla = size4d[0];
3215 std::size_t nb_sli = size4d[1];
3216 std::size_t nb_lig = size4d[2];
3217 std::size_t nb_col = size4d[3];
3219 for(
size_t t = 0; t < nb_sla; ++t)
3220 for(
size_t k = 0; k < nb_sli; ++k)
3221 for(
size_t i = 0; i < nb_lig; ++i)
3224 in_first_front_upper_left.row_end(t,k,i),out_first_front_upper_left.row_begin(t,k,i));
3228 for(
size_t t = 0; t < nb_sla; ++t)
3229 for(
size_t k = 0; k < nb_sli; ++k)
3230 for(
size_t j = 0; j < nb_col; ++j)
3233 out_first_front_upper_left.col_end(t,k,j),Vtmp.
begin());
3238 for(
size_t t = 0; t < nb_sla; ++t)
3239 for(
size_t i = 0; i < nb_lig; ++i)
3240 for(
size_t j = 0; j < nb_col; ++j)
3243 out_first_front_upper_left.slice_end(t,i,j),Vtmp.
begin());
3244 std::copy(Vtmp.
begin(),Vtmp.
end(),out_first_front_upper_left.slice_begin(t,i,j));
3248 for(
size_t k = 0; k < nb_sli; ++k)
3249 for(
size_t i = 0; i < nb_lig; ++i)
3250 for(
size_t j = 0; j < nb_col; ++j)
3253 out_first_front_upper_left.slab_end(k,i,j),Vtmp.
begin());
3254 std::copy(Vtmp.
begin(),Vtmp.
end(),out_first_front_upper_left.slab_begin(k,i,j));
3290 template<
typename InputBidirectionalIterator4d,
3291 typename OutputBidirectionalIterator4d>
3294 InputBidirectionalIterator4d in_last_back_bottom_right,
3295 OutputBidirectionalIterator4d out_first_front_upper_left)
3298 typename InputBidirectionalIterator4d::difference_type size4d = in_last_back_bottom_right - in_first_front_upper_left;
3300 std::size_t nb_sla = size4d[0];
3301 std::size_t nb_sli = size4d[1];
3302 std::size_t nb_lig = size4d[2];
3303 std::size_t nb_col = size4d[3];
3306 for(
size_t t = 0; t < nb_sla; ++t)
3307 for(
size_t k = 0; k < nb_sli; ++k)
3308 for(
size_t i = 0; i < nb_lig; ++i)
3311 in_first_front_upper_left.row_end(t,k,i),out_first_front_upper_left.row_begin(t,k,i));
3315 for(
size_t t = 0; t < nb_sla; ++t)
3316 for(
size_t k = 0; k < nb_sli; ++k)
3317 for(
size_t j = 0; j < nb_col; ++j)
3320 out_first_front_upper_left.col_end(t,k,j),Vtmp.
begin());
3325 for(
size_t t = 0; t < nb_sla; ++t)
3326 for(
size_t i = 0; i < nb_lig; ++i)
3327 for(
size_t j = 0; j < nb_col; ++j)
3330 out_first_front_upper_left.slice_end(t,i,j),Vtmp.
begin());
3331 std::copy(Vtmp.
begin(),Vtmp.
end(),out_first_front_upper_left.slice_begin(t,i,j));
3335 for(
size_t k = 0; k < nb_sli; ++k)
3336 for(
size_t i = 0; i < nb_lig; ++i)
3337 for(
size_t j = 0; j < nb_col; ++j)
3340 out_first_front_upper_left.slab_end(k,i,j),Vtmp.
begin());
3341 std::copy(Vtmp.
begin(),Vtmp.
end(),out_first_front_upper_left.slab_begin(k,i,j));
3380 template<
typename InputMatrix4d,
3381 typename OutputMatrix4d>
3385 assert(datain.cols() == dataout.cols());
3386 assert(datain.rows() == dataout.rows());
3387 assert(datain.slices() == dataout.slices());
3389 dataout.first_front_upper_left());
3406 template<
typename InputMatrix4d,
3407 typename OutputMatrix4d>
3411 assert(datain.cols() == dataout.cols());
3412 assert(datain.rows() == dataout.rows());
3413 assert(datain.slices() == dataout.slices());
3415 dataout.first_front_upper_left());
3453 template<
typename InputMatrix4d,
3454 typename OutputMatrix4d>
3458 assert(datain.cols() == dataout.cols());
3459 assert(datain.rows() == dataout.rows());
3460 assert(datain.slices() == dataout.slices());
3461 split_radix_ifft4d(datain.first_front_upper_left(),datain.last_back_bottom_right(),dataout.first_front_upper_left());
3488 template<
typename InputBidirectionalIterator4d,
3489 typename OutputBidirectionalIterator4d>
3491 void real_fftw4d(InputBidirectionalIterator4d in_first_front_upper_left,
3492 InputBidirectionalIterator4d in_last_back_bottom_right,
3493 OutputBidirectionalIterator4d out_first_front_upper_left)
3496 typename InputBidirectionalIterator4d::difference_type size4d = in_last_back_bottom_right - in_first_front_upper_left;
3497 typedef typename std::iterator_traits<InputBidirectionalIterator4d>::value_type
real;
3498 typedef typename std::iterator_traits<OutputBidirectionalIterator4d>::value_type complex;
3500 std::size_t nb_sla = size4d[0];
3501 std::size_t nb_sli = size4d[1];
3502 std::size_t nb_lig = size4d[2];
3503 std::size_t nb_col = size4d[3];
3504 assert((nb_sla > 1)&&(nb_sli > 1)&&(nb_lig > 1)&&(nb_col>1));
3506 for(
size_t t = 0; t < nb_sla; ++t)
3507 for(
size_t k = 0; k < nb_sli; ++k)
3508 for(
size_t i = 0; i < nb_lig; ++i)
3511 in_first_front_upper_left.row_end(t,k,i),out_first_front_upper_left.row_begin(t,k,i));
3515 for(
size_t t = 0; t < nb_sla; ++t)
3516 for(
size_t k = 0; k < nb_sli; ++k)
3517 for(
size_t j = 0; j < nb_col; ++j)
3520 out_first_front_upper_left.col_end(t,k,j),Vtmp.
begin());
3525 for(
size_t t = 0; t < nb_sla; ++t)
3526 for(
size_t i = 0; i < nb_lig; ++i)
3527 for(
size_t j = 0; j < nb_col; ++j)
3530 out_first_front_upper_left.slice_end(t,i,j),Vtmp.
begin());
3531 std::copy(Vtmp.
begin(),Vtmp.
end(),out_first_front_upper_left.slice_begin(t,i,j));
3535 for(
size_t k = 0; k < nb_sli; ++k)
3536 for(
size_t i = 0; i < nb_lig; ++i)
3537 for(
size_t j = 0; j < nb_col; ++j)
3540 out_first_front_upper_left.slab_end(k,i,j),Vtmp.
begin());
3541 std::copy(Vtmp.
begin(),Vtmp.
end(),out_first_front_upper_left.slab_begin(k,i,j));
3560 template<
typename InputBidirectionalIterator4d,
3561 typename OutputBidirectionalIterator4d>
3564 InputBidirectionalIterator4d in_last_back_bottom_right,
3565 OutputBidirectionalIterator4d out_first_front_upper_left)
3568 typename InputBidirectionalIterator4d::difference_type size4d = in_last_back_bottom_right - in_first_front_upper_left;
3569 typedef typename std::iterator_traits<InputBidirectionalIterator4d>::value_type complex;
3570 typedef typename complex::value_type
real;
3572 std::size_t nb_sla = size4d[0];
3573 std::size_t nb_sli = size4d[1];
3574 std::size_t nb_lig = size4d[2];
3575 std::size_t nb_col = size4d[3];
3576 assert((nb_sla > 1)&&(nb_sli > 1)&&(nb_lig > 1)&&(nb_col>1));
3579 for(
size_t t = 0; t < nb_sla; ++t)
3580 for(
size_t k = 0; k < nb_sli; ++k)
3581 for(
size_t i = 0; i < nb_lig; ++i)
3584 in_first_front_upper_left.row_end(t,k,i),out_first_front_upper_left.row_begin(t,k,i));
3588 for(
size_t t = 0; t < nb_sla; ++t)
3589 for(
size_t k = 0; k < nb_sli; ++k)
3590 for(
size_t j = 0; j < nb_col; ++j)
3593 out_first_front_upper_left.col_end(t,k,j),Vtmp.
begin());
3598 for(
size_t t = 0; t < nb_sla; ++t)
3599 for(
size_t i = 0; i < nb_lig; ++i)
3600 for(
size_t j = 0; j < nb_col; ++j)
3603 out_first_front_upper_left.slice_end(t,i,j),Vtmp.
begin());
3604 std::copy(Vtmp.
begin(),Vtmp.
end(),out_first_front_upper_left.slice_begin(t,i,j));
3608 for(
size_t k = 0; k < nb_sli; ++k)
3609 for(
size_t i = 0; i < nb_lig; ++i)
3610 for(
size_t j = 0; j < nb_col; ++j)
3613 out_first_front_upper_left.slab_end(k,i,j),Vtmp.
begin());
3614 std::copy(Vtmp.
begin(),Vtmp.
end(),out_first_front_upper_left.slab_begin(k,i,j));
3633 template<
typename InputBidirectionalIterator4d,
3634 typename OutputBidirectionalIterator4d>
3636 void ifftw4d(InputBidirectionalIterator4d in_first_front_upper_left,
3637 InputBidirectionalIterator4d in_last_back_bottom_right,
3638 OutputBidirectionalIterator4d out_first_front_upper_left)
3641 typename InputBidirectionalIterator4d::difference_type size4d = in_last_back_bottom_right - in_first_front_upper_left;
3643 std::size_t nb_sla = size4d[0];
3644 std::size_t nb_sli = size4d[1];
3645 std::size_t nb_lig = size4d[2];
3646 std::size_t nb_col = size4d[3];
3649 for(
size_t t = 0; t < nb_sla; ++t)
3650 for(
size_t k = 0; k < nb_sli; ++k)
3651 for(
size_t i = 0; i < nb_lig; ++i)
3653 slip::ifftw(in_first_front_upper_left.row_begin(t,k,i),
3654 in_first_front_upper_left.row_end(t,k,i),out_first_front_upper_left.row_begin(t,k,i));
3658 for(
size_t t = 0; t < nb_sla; ++t)
3659 for(
size_t k = 0; k < nb_sli; ++k)
3660 for(
size_t j = 0; j < nb_col; ++j)
3662 slip::ifftw(out_first_front_upper_left.col_begin(t,k,j),
3663 out_first_front_upper_left.col_end(t,k,j),Vtmp.
begin());
3668 for(
size_t t = 0; t < nb_sla; ++t)
3669 for(
size_t i = 0; i < nb_lig; ++i)
3670 for(
size_t j = 0; j < nb_col; ++j)
3672 slip::ifftw(out_first_front_upper_left.slice_begin(t,i,j),
3673 out_first_front_upper_left.slice_end(t,i,j),Vtmp.
begin());
3674 std::copy(Vtmp.
begin(),Vtmp.
end(),out_first_front_upper_left.slice_begin(t,i,j));
3678 for(
size_t k = 0; k < nb_sli; ++k)
3679 for(
size_t i = 0; i < nb_lig; ++i)
3680 for(
size_t j = 0; j < nb_col; ++j)
3682 slip::ifftw(out_first_front_upper_left.slab_begin(k,i,j),
3683 out_first_front_upper_left.slab_end(k,i,j),Vtmp.
begin());
3684 std::copy(Vtmp.
begin(),Vtmp.
end(),out_first_front_upper_left.slab_begin(k,i,j));
3703 template<
typename InputMatrix4d,
3704 typename OutputMatrix4d>
3706 void real_fftw4d(
const InputMatrix4d & datain,OutputMatrix4d & dataout)
3708 assert(datain.cols() == dataout.cols());
3709 assert(datain.rows() == dataout.rows());
3710 assert(datain.slices() == dataout.slices());
3711 assert(datain.slabs() == dataout.slabs());
3712 real_fftw4d(datain.first_front_upper_left(),datain.last_back_bottom_right(),dataout.first_front_upper_left());
3731 template<
typename InputMatrix4d,
3732 typename OutputMatrix4d>
3736 assert(datain.cols() == dataout.cols());
3737 assert(datain.rows() == dataout.rows());
3738 assert(datain.slices() == dataout.slices());
3739 assert(datain.slabs() == dataout.slabs());
3740 complex_fftw4d(datain.first_front_upper_left(),datain.last_back_bottom_right(),dataout.first_front_upper_left());
3762 template<
typename InputMatrix4d,
3763 typename OutputMatrix4d>
3765 void ifftw4d(
const InputMatrix4d & datain,OutputMatrix4d & dataout)
3767 assert(datain.cols() == dataout.cols());
3768 assert(datain.rows() == dataout.rows());
3769 assert(datain.slices() == dataout.slices());
3770 assert(datain.slabs() == dataout.slabs());
3771 ifftw4d(datain.first_front_upper_left(),datain.last_back_bottom_right(),dataout.first_front_upper_left());
3820 template<
typename InputBidirectionalIterator4d,
3821 typename OutputBidirectionalIterator4d>
3823 void real_fft4d(InputBidirectionalIterator4d in_first_front_upper_left,
3824 InputBidirectionalIterator4d in_last_back_bottom_right,
3825 OutputBidirectionalIterator4d out_first_front_upper_left)
3828 real_fftw4d(in_first_front_upper_left,in_last_back_bottom_right,out_first_front_upper_left);
3870 template<
typename InputBidirectionalIterator4d,
3871 typename OutputBidirectionalIterator4d>
3874 InputBidirectionalIterator4d in_last_back_bottom_right,
3875 OutputBidirectionalIterator4d out_first_front_upper_left)
3878 complex_fftw4d(in_first_front_upper_left,in_last_back_bottom_right,out_first_front_upper_left);
3918 template<
typename InputBidirectionalIterator4d,
3919 typename OutputBidirectionalIterator4d>
3921 void ifft4d(InputBidirectionalIterator4d in_first_front_upper_left,
3922 InputBidirectionalIterator4d in_last_back_bottom_right,
3923 OutputBidirectionalIterator4d out_first_front_upper_left)
3926 ifftw4d(in_first_front_upper_left,in_last_back_bottom_right,out_first_front_upper_left);
3928 split_radix_ifft4d(in_first_front_upper_left,in_last_back_bottom_right,out_first_front_upper_left);
3970 template<
typename HyperVolume1,
typename HyperVolume2>
3972 void real_fft4d(
const HyperVolume1 &datain, HyperVolume2 &dataout)
3974 assert(datain.cols() == dataout.cols());
3975 assert(datain.rows() == dataout.rows());
3976 assert(datain.slices() == dataout.slices());
3978 real_fftw4d(datain.first_front_upper_left(),datain.last_back_bottom_right(),dataout.first_front_upper_left());
3980 real_split_radix_fft4d(datain.first_front_upper_left(),datain.last_back_bottom_right(),dataout.first_front_upper_left());
4021 template<
typename HyperVolume1,
typename HyperVolume2>
4025 assert(datain.cols() == dataout.cols());
4026 assert(datain.rows() == dataout.rows());
4027 assert(datain.slices() == dataout.slices());
4029 complex_fftw4d(datain.first_front_upper_left(),datain.last_back_bottom_right(),dataout.first_front_upper_left());
4031 complex_split_radix_fft4d(datain.first_front_upper_left(),datain.last_back_bottom_right(),dataout.first_front_upper_left());
4073 template<
typename HyperVolume1,
typename HyperVolume2>
4075 void ifft4d(
const HyperVolume1 &datain, HyperVolume2 &dataout)
4077 assert(datain.cols() == dataout.cols());
4078 assert(datain.rows() == dataout.rows());
4079 assert(datain.slices() == dataout.slices());
4081 ifftw4d(datain.first_front_upper_left(),datain.last_back_bottom_right(),dataout.first_front_upper_left());
4083 split_radix_ifft4d(datain.first_front_upper_left(),datain.last_back_bottom_right(),dataout.first_front_upper_left());
void complex_split_radix_fft4d(InputBidirectionalIterator4d in_first_front_upper_left, InputBidirectionalIterator4d in_last_back_bottom_right, OutputBidirectionalIterator4d out_first_front_upper_left)
Computes the complex split-radix fft4d of a container with 4d iterators.
void split_radix_ifft2d(InputBidirectionalIterator2d in_upper_left, InputBidirectionalIterator2d in_bottom_right, OutputBidirectionalIterator2d out_upper_left)
Computes the split-radix ifft2d of a container with 2d iterators.
void ifftw(InputIter in_begin, InputIter in_end, OutputIter out_begin)
Computes the fftw backward of a container (from fftw3.h).
iterator2d upper_left()
Returns a read/write iterator2d that points to the first element of the Array2d. It points to the upp...
void complex_fft4d(InputBidirectionalIterator4d in_first_front_upper_left, InputBidirectionalIterator4d in_last_back_bottom_right, OutputBidirectionalIterator4d out_first_front_upper_left)
Computes the complex fft4d of a container with 4d iterators.
void fftshift(ForwardIterator first, ForwardIterator last)
Performs a shift of a container, for use with the fft and ifft functions, in order to move the freque...
void fftshift2d(ForwardIterator2D upper_left, ForwardIterator2D bottom_right)
Performs a shift of a 2d container, for use with the fft2d and ifft2d functions, in order to move the...
void complex_fftw2d(InputBidirectionalIterator2d in_upper_left, InputBidirectionalIterator2d in_bottom_right, OutputBidirectionalIterator2d out_upper_left)
Computes the complex fftw2d of a container with 2d iterators.
iterator end()
Returns a read/write iterator that points one past the last element in the Vector. Iteration is done in ordinary element order.
void real_split_radix_fft(InputIter in_begin, InputIter in_end, OutputIter out_begin)
Computes the real split radix fft of a container (from FFTPACK.hpp).
iterator begin()
Returns a read/write iterator that points to the first element in the Array. Iteration is done in ord...
void radix2_fft(InputIter_T a, Iter_T b, const int log2n)
Computes the Cooley-Tukey (radix2) fft of a container.
void fftshift4d(ForwardIterator4D first_front_upper_left, ForwardIterator4D last_back_bottom_right)
Performs a shift of a 4d container, for use with the fft4d and ifft4d functions, in order to move the...
void ifftw2d(InputBidirectionalIterator2d in_upper_left, InputBidirectionalIterator2d in_bottom_right, OutputBidirectionalIterator2d out_upper_left)
Computes the ifftw2d of a container with 2d iterators.
void real_fftw4d(InputBidirectionalIterator4d in_first_front_upper_left, InputBidirectionalIterator4d in_last_back_bottom_right, OutputBidirectionalIterator4d out_first_front_upper_left)
Computes the real fftw4d of a container with 4d iterators.
Real functor. Return the real part of x.
void complex_fft2d(InputBidirectionalIterator2d in_upper_left, InputBidirectionalIterator2d in_bottom_right, OutputBidirectionalIterator2d out_upper_left)
Computes the complex fft2d of a container with 2d iterators.
void ifft4d(InputBidirectionalIterator4d in_first_front_upper_left, InputBidirectionalIterator4d in_last_back_bottom_right, OutputBidirectionalIterator4d out_first_front_upper_left)
Computes the ifft4d of a container with 4d iterators.
HyperVolume< T > exp(const HyperVolume< T > &M)
void real_fft2d(InputBidirectionalIterator2d in_upper_left, InputBidirectionalIterator2d in_bottom_right, OutputBidirectionalIterator2d out_upper_left)
Computes the real fft2d of a container with 2d iterators.
void radix2_fft1d2d(InputIterator1d begin1, OuputIterator1d begin2, std::size_t nb_lig, std::size_t nb_col)
Computes the radix 2 fft2d of a container with 1d iterator (compatible with simple pointers)...
void ifftw4d(InputBidirectionalIterator4d in_first_front_upper_left, InputBidirectionalIterator4d in_last_back_bottom_right, OutputBidirectionalIterator4d out_first_front_upper_left)
Computes the ifftw4d of a container with 4d iterators.
void real(const Matrix1 &C, Matrix2 &R)
Extract the real Matrix of a complex Matrix. .
void complex_split_radix_fft2d(InputBidirectionalIterator2d in_upper_left, InputBidirectionalIterator2d in_bottom_right, OutputBidirectionalIterator2d out_upper_left)
Computes the complex split-radix fft2d of a container with 2d iterators.
static _Tp pi()
Constant .
void fftduplicate(ForwardIterator1 begin1, ForwardIterator2 begin2, ForwardIterator2 end2)
Calculates the conjugates of the begin std::complex elements of a first container. Writes them at the end of a second container. This function is used to calculate the negative frequencies after a real fftw.
Provides a class to manipulate 1d dynamic and generic arrays.
void split_radix_ifft3d(InputBidirectionalIterator3d in_front_upper_left, InputBidirectionalIterator3d in_back_bottom_right, OutputBidirectionalIterator3d out_front_upper_left)
Computes the split-radix ifft3d of a container with 3d iterators.
void ifftw3d(InputBidirectionalIterator3d in_front_upper_left, InputBidirectionalIterator3d in_back_bottom_right, OutputBidirectionalIterator3d out_front_upper_left)
Computes the ifftw3d of a container with 3d iterators.
void ifft3d(InputBidirectionalIterator3d in_front_upper_left, InputBidirectionalIterator3d in_back_bottom_right, OutputBidirectionalIterator3d out_front_upper_left)
Computes the ifft3d of a container with 3d iterators.
void complex_split_radix_fft(InputIter in_begin, InputIter in_end, OutputIter out_begin)
Computes the complex split radix fft of a container (from FFTPACK.hpp).
void complex_split_radix_fft3d(InputBidirectionalIterator3d in_front_upper_left, InputBidirectionalIterator3d in_back_bottom_right, OutputBidirectionalIterator3d out_front_upper_left)
Computes the complex split-radix fft3d of a container with 3d iterators.
void fftduplicate2d(ForwardIterator2D upper_left, ForwardIterator2D bottom_right)
Performs a duplicate of a 2D container, for use with the fft2d and ifft2d functions, in order to calculate the negative frequencies with the real fftw2d.
void complex_fftw(InputIter in_begin, InputIter in_end, OutputIter out_begin)
complex fftw of a container (from fftw3.h).
void dct(InputIter in_begin, InputIter in_end, OutputIter out_begin)
Computes the Discrete Cosinus Transform II (forward) of a container.
unsigned int log2(unsigned int x)
Calculates the log2 of an integer.
void ifft(InputIter in_begin, InputIter in_end, OutputIter out_begin)
Computes the inverse fft of a container.
void idct(InputIter in_begin, InputIter in_end, OutputIter out_begin)
Computes the Discrete Cosinus Transform III (backward) of a container.
void real_fftw3d(InputBidirectionalIterator3d in_front_upper_left, InputBidirectionalIterator3d in_back_bottom_right, OutputBidirectionalIterator3d out_front_upper_left)
Computes the real fftw3d of a container with 3d iterators.
void split_radix_ifft(InputIter in_begin, InputIter in_end, OutputIter out_begin)
Computes the split radix fft backward of a container (from FFTPACK.hpp).
iterator end()
Returns a read/write iterator that points one past the last element in the Array. Iteration is done i...
Provides some statistics algorithms.
void radix2_fft2d(InputBidirectionalIterator2d in_upper_left, InputBidirectionalIterator2d in_bottom_right, OutputBidirectionalIterator2d out_upper_left)
Computes the radix2 (Cooley-Tukey) fft2d of a container with 2d iterators.
void real_fftw(InputIter in_begin, InputIter in_end, OutputIter out_begin)
Computes the real fftw of a container (from fftw3.h).
void real_split_radix_fft4d(InputBidirectionalIterator4d in_first_front_upper_left, InputBidirectionalIterator4d in_last_back_bottom_right, OutputBidirectionalIterator4d out_first_front_upper_left)
Computes the real split-radix fft4d of a container with 4d iterators.
void copy(_II first, _II last, _OI output_first)
Copy algorithm optimized for slip iterators.
Difference of Point2D class, specialization of DPoint<CoordType,DIM> with DIM = 2.
void split_radix_idct(InputIter in_begin, InputIter in_end, OutputIter out_begin)
Computes the split radix Discrete Cosinus Transform III (i.e. idct) of a container (from FFTPACK...
Numerical Vector class. This container statisfies the RandomAccessContainer concepts of the Standard ...
Imag functor. Return the imaginary part of x.
unsigned int bitReverse(unsigned int x, const int log2n)
bitReverse function used in the Cooley-Tukey FFT algorithm
Provides some mathematical functors and constants.
iterator end()
Returns a read/write iterator that points one past the last element in the Array2d. Iteration is done in ordinary element order.
HyperVolume< T > sqrt(const HyperVolume< T > &M)
void complex_fft(InputIter in_begin, InputIter in_end, OutputIter out_begin)
Computes the complex fft of a container.
void resize(const size_type new_n, const T &val=T())
Resizes a Vector.
void real_split_radix_fft2d(InputBidirectionalIterator2d in_upper_left, InputBidirectionalIterator2d in_bottom_right, OutputBidirectionalIterator2d out_upper_left)
Computes the real split-radix fft2d of a container with 2d iterators.
void fftw_idct(InputIter in_begin, InputIter in_end, OutputIter out_begin)
Computes the Discrete Cosinus Transform III (backward) of a container (from fftw).
void split_radix_dct(InputIter in_begin, InputIter in_end, OutputIter out_begin)
Computes the split radix Discrete Cosinus Transform II of a container (from FFTPACK.hpp).
void real_fft4d(InputBidirectionalIterator4d in_first_front_upper_left, InputBidirectionalIterator4d in_last_back_bottom_right, OutputBidirectionalIterator4d out_first_front_upper_left)
Computes the real fft4d of a container with 4d iterators.
void complex_fftw3d(InputBidirectionalIterator3d in_front_upper_left, InputBidirectionalIterator3d in_back_bottom_right, OutputBidirectionalIterator3d out_front_upper_left)
Computes the complex fftw3d of a container with 3d iterators.
void ifft2d(InputBidirectionalIterator2d in_upper_left, InputBidirectionalIterator2d in_bottom_right, OutputBidirectionalIterator2d out_upper_left)
Computes the ifft2d of a container with 2d iterators.
void complex_fft3d(InputBidirectionalIterator3d in_front_upper_left, InputBidirectionalIterator3d in_back_bottom_right, OutputBidirectionalIterator3d out_front_upper_left)
Computes the complex fft3d of a container with 3d iterators.
void fftshift3d(ForwardIterator3D front_upper_left, ForwardIterator3D back_bottom_right)
Performs a shift of a 2d container, for use with the fft2d and ifft2d functions, in order to move the...
void real_fftw2d(InputBidirectionalIterator2d in_upper_left, InputBidirectionalIterator2d in_bottom_right, OutputBidirectionalIterator2d out_upper_left)
Computes the real fftw2d of a container with 2d iterators.
void radix2_ifft1d2d(InputIterator1d begin1, OuputIterator1d begin2, std::size_t nb_lig, std::size_t nb_col)
Computes the inverse radix 2 fft2d of a container with 1d iterators (compatible with simple pointers)...
void real_fft(InputIter in_begin, InputIter in_end, OutputIter out_begin)
Computes the real fft of a container.
Provides a class to manipulate Numerical Matrix.
void real_radix2_fft(InputIter_T a, Iter_T b, const int log2n)
Computes the Cooley-Tukey (radix2) real_fft of a container.
iterator begin()
Returns a read/write iterator that points to the first element in the Array2d. Iteration is done in o...
void radix2_ifft2d(InputBidirectionalIterator2d in_upper_left, InputBidirectionalIterator2d in_bottom_right, OutputBidirectionalIterator2d out_upper_left)
Computes the radix2 ifft2d of a container with 2d iterators.
void complex_fftw4d(InputBidirectionalIterator4d in_first_front_upper_left, InputBidirectionalIterator4d in_last_back_bottom_right, OutputBidirectionalIterator4d out_first_front_upper_left)
Computes the complex fftw4d of a container with 4d iterators.
This is a linear (one-dimensional) dynamic template container. This container statisfies the RandomAc...
void split_radix_ifft4d(InputBidirectionalIterator4d in_first_front_upper_left, InputBidirectionalIterator4d in_last_back_bottom_right, OutputBidirectionalIterator4d out_first_front_upper_left)
Computes the split-radix ifft4d of a container with 4d iterators.
void fftw_dct(InputIter in_begin, InputIter in_end, OutputIter out_begin)
Computes the Discrete Cosinus Transform II (forward) of a container (from fftw).
void real_fft3d(InputBidirectionalIterator3d in_front_upper_left, InputBidirectionalIterator3d in_back_bottom_right, OutputBidirectionalIterator3d out_front_upper_left)
Computes the real fft3d of a container with 3d iterators.
void radix2_ifft(InputIter_T a, Iter_T b, const int log2n)
Computes the inverse Cooley-Tukey (radix2) fft of a container.
void real_split_radix_fft3d(InputBidirectionalIterator3d in_front_upper_left, InputBidirectionalIterator3d in_back_bottom_right, OutputBidirectionalIterator3d out_front_upper_left)
Computes the real split-radix fft3d of a container with 3d iterators.
Provides a class to manipulate numerical vectors.
This is a two-dimensional dynamic and generic container. This container statisfies the Bidirectionnal...
iterator begin()
Returns a read/write iterator that points to the first element in the Vector. Iteration is done in or...