27 #if defined (HAVE_FFTW)
38 #if defined (HAVE_FFTW3_THREADS)
63 : meth (ESTIMATE), rplan (0), rd (0), rs (0), rr (0), rh (0), rn (),
67 d[0] =
d[1] =
s[0] =
s[1] =
r[0] =
r[1] =
h[0] =
h[1] = 0;
72 #if defined (HAVE_FFTW3_THREADS)
73 int init_ret = fftw_init_threads ();
75 (*current_liboctave_error_handler) (
"Error initializing FFTW threads");
78 nthreads = num_processors (NPROC_CURRENT);
79 fftw_plan_with_nthreads (nthreads);
83 fftw_import_system_wisdom ();
92 fftw_destroy_plan (*plan_p);
96 fftw_destroy_plan (*plan_p);
100 fftw_destroy_plan (*plan_p);
118 (*current_liboctave_error_handler)
119 (
"unable to create octave_fftw_planner object!");
127 #define CHECK_SIMD_ALIGNMENT(x) \
128 (((reinterpret_cast<ptrdiff_t> (x)) & 0xF) == 0)
138 int which = (dir == FFTW_FORWARD) ? 0 : 1;
139 fftw_plan *cur_plan_p = &
plan[which];
140 bool create_new_plan =
false;
142 bool ioinplace = (in == out);
148 if (
plan[which] == 0 ||
d[which] != dist ||
s[which] != stride
149 ||
r[which] != rank ||
h[which] != howmany
151 || ((ioalign !=
simd_align[which]) ? !ioalign :
false))
152 create_new_plan =
true;
157 for (
int i = 0; i < rank; i++)
158 if (dims(i) !=
n[which](i))
160 create_new_plan =
true;
179 for (
int i = 0, j = rank-1; i < rank; i++, j--)
186 bool plan_destroys_in =
true;
192 plan_flags |= FFTW_ESTIMATE;
193 plan_destroys_in =
false;
196 plan_flags |= FFTW_MEASURE;
199 plan_flags |= FFTW_PATIENT;
202 plan_flags |= FFTW_EXHAUSTIVE;
206 plan_flags |= FFTW_MEASURE;
209 plan_flags |= FFTW_ESTIMATE;
210 plan_destroys_in =
false;
216 plan_flags &= ~FFTW_UNALIGNED;
218 plan_flags |= FFTW_UNALIGNED;
221 fftw_destroy_plan (*cur_plan_p);
223 if (plan_destroys_in)
227 itmp =
reinterpret_cast<Complex *
>
228 (((
reinterpret_cast<ptrdiff_t
>(itmp) + 15) & ~ 0xF) +
229 ((reinterpret_cast<ptrdiff_t> (in)) & 0xF));
232 fftw_plan_many_dft (rank, tmp, howmany,
233 reinterpret_cast<fftw_complex *> (itmp),
235 reinterpret_cast<fftw_complex *> (out),
236 0, stride, dist, dir, plan_flags);
241 fftw_plan_many_dft (rank, tmp, howmany,
242 reinterpret_cast<fftw_complex *> (const_cast<Complex *> (in)),
244 reinterpret_cast<fftw_complex *> (out),
245 0, stride, dist, dir, plan_flags);
248 if (*cur_plan_p == 0)
249 (*current_liboctave_error_handler) (
"Error creating fftw plan");
260 const double *in,
Complex *out)
262 fftw_plan *cur_plan_p = &
rplan;
263 bool create_new_plan =
false;
270 if (
rplan == 0 ||
rd != dist ||
rs != stride ||
rr != rank
271 ||
rh != howmany || ((ioalign !=
rsimd_align) ? !ioalign :
false))
272 create_new_plan =
true;
277 for (
int i = 0; i < rank; i++)
278 if (dims(i) !=
rn(i))
280 create_new_plan =
true;
298 for (
int i = 0, j = rank-1; i < rank; i++, j--)
305 bool plan_destroys_in =
true;
311 plan_flags |= FFTW_ESTIMATE;
312 plan_destroys_in =
false;
315 plan_flags |= FFTW_MEASURE;
318 plan_flags |= FFTW_PATIENT;
321 plan_flags |= FFTW_EXHAUSTIVE;
325 plan_flags |= FFTW_MEASURE;
328 plan_flags |= FFTW_ESTIMATE;
329 plan_destroys_in =
false;
335 plan_flags &= ~FFTW_UNALIGNED;
337 plan_flags |= FFTW_UNALIGNED;
340 fftw_destroy_plan (*cur_plan_p);
342 if (plan_destroys_in)
346 itmp =
reinterpret_cast<double *
>
347 (((
reinterpret_cast<ptrdiff_t
>(itmp) + 15) & ~ 0xF) +
348 ((reinterpret_cast<ptrdiff_t> (in)) & 0xF));
351 fftw_plan_many_dft_r2c (rank, tmp, howmany, itmp,
353 reinterpret_cast<fftw_complex *> (out),
354 0, stride, dist, plan_flags);
359 fftw_plan_many_dft_r2c (rank, tmp, howmany,
360 (const_cast<double *> (in)),
362 reinterpret_cast<fftw_complex *> (out),
363 0, stride, dist, plan_flags);
366 if (*cur_plan_p == 0)
367 (*current_liboctave_error_handler) (
"Error creating fftw plan");
391 fftw_destroy_plan (
rplan);
393 fftw_destroy_plan (
plan[0]);
395 fftw_destroy_plan (
plan[1]);
407 : meth (ESTIMATE), rplan (0), rd (0), rs (0), rr (0), rh (0), rn (),
411 d[0] =
d[1] =
s[0] =
s[1] =
r[0] =
r[1] =
h[0] =
h[1] = 0;
416 #if defined (HAVE_FFTW3F_THREADS)
417 int init_ret = fftwf_init_threads ();
419 (*current_liboctave_error_handler) (
"Error initializing FFTW3F threads");
422 nthreads = num_processors (NPROC_CURRENT);
423 fftwf_plan_with_nthreads (nthreads);
427 fftwf_import_system_wisdom ();
436 fftwf_destroy_plan (*plan_p);
440 fftwf_destroy_plan (*plan_p);
444 fftwf_destroy_plan (*plan_p);
462 (*current_liboctave_error_handler)
463 (
"unable to create octave_fftw_planner object!");
480 int which = (dir == FFTW_FORWARD) ? 0 : 1;
481 fftwf_plan *cur_plan_p = &
plan[which];
482 bool create_new_plan =
false;
484 bool ioinplace = (in == out);
490 if (
plan[which] == 0 ||
d[which] != dist ||
s[which] != stride
491 ||
r[which] != rank ||
h[which] != howmany
493 || ((ioalign !=
simd_align[which]) ? !ioalign :
false))
494 create_new_plan =
true;
499 for (
int i = 0; i < rank; i++)
500 if (dims(i) !=
n[which](i))
502 create_new_plan =
true;
521 for (
int i = 0, j = rank-1; i < rank; i++, j--)
528 bool plan_destroys_in =
true;
534 plan_flags |= FFTW_ESTIMATE;
535 plan_destroys_in =
false;
538 plan_flags |= FFTW_MEASURE;
541 plan_flags |= FFTW_PATIENT;
544 plan_flags |= FFTW_EXHAUSTIVE;
548 plan_flags |= FFTW_MEASURE;
551 plan_flags |= FFTW_ESTIMATE;
552 plan_destroys_in =
false;
558 plan_flags &= ~FFTW_UNALIGNED;
560 plan_flags |= FFTW_UNALIGNED;
563 fftwf_destroy_plan (*cur_plan_p);
565 if (plan_destroys_in)
570 (((
reinterpret_cast<ptrdiff_t
>(itmp) + 15) & ~ 0xF) +
571 ((reinterpret_cast<ptrdiff_t> (in)) & 0xF));
574 fftwf_plan_many_dft (rank, tmp, howmany,
575 reinterpret_cast<fftwf_complex *> (itmp),
577 reinterpret_cast<fftwf_complex *> (out),
578 0, stride, dist, dir, plan_flags);
583 fftwf_plan_many_dft (rank, tmp, howmany,
584 reinterpret_cast<fftwf_complex *> (const_cast<FloatComplex *> (in)),
586 reinterpret_cast<fftwf_complex *> (out),
587 0, stride, dist, dir, plan_flags);
590 if (*cur_plan_p == 0)
591 (*current_liboctave_error_handler) (
"Error creating fftw plan");
605 fftwf_plan *cur_plan_p = &
rplan;
606 bool create_new_plan =
false;
613 if (
rplan == 0 ||
rd != dist ||
rs != stride ||
rr != rank
614 ||
rh != howmany || ((ioalign !=
rsimd_align) ? !ioalign :
false))
615 create_new_plan =
true;
620 for (
int i = 0; i < rank; i++)
621 if (dims(i) !=
rn(i))
623 create_new_plan =
true;
641 for (
int i = 0, j = rank-1; i < rank; i++, j--)
648 bool plan_destroys_in =
true;
654 plan_flags |= FFTW_ESTIMATE;
655 plan_destroys_in =
false;
658 plan_flags |= FFTW_MEASURE;
661 plan_flags |= FFTW_PATIENT;
664 plan_flags |= FFTW_EXHAUSTIVE;
668 plan_flags |= FFTW_MEASURE;
671 plan_flags |= FFTW_ESTIMATE;
672 plan_destroys_in =
false;
678 plan_flags &= ~FFTW_UNALIGNED;
680 plan_flags |= FFTW_UNALIGNED;
683 fftwf_destroy_plan (*cur_plan_p);
685 if (plan_destroys_in)
689 itmp =
reinterpret_cast<float *
>
690 (((
reinterpret_cast<ptrdiff_t
>(itmp) + 15) & ~ 0xF) +
691 ((reinterpret_cast<ptrdiff_t> (in)) & 0xF));
694 fftwf_plan_many_dft_r2c (rank, tmp, howmany, itmp,
696 reinterpret_cast<fftwf_complex *> (out),
697 0, stride, dist, plan_flags);
702 fftwf_plan_many_dft_r2c (rank, tmp, howmany,
703 (const_cast<float *> (in)),
705 reinterpret_cast<fftwf_complex *> (out),
706 0, stride, dist, plan_flags);
709 if (*cur_plan_p == 0)
710 (*current_liboctave_error_handler) (
"Error creating fftw plan");
734 fftwf_destroy_plan (
rplan);
736 fftwf_destroy_plan (
plan[0]);
738 fftwf_destroy_plan (
plan[1]);
756 for (
size_t i = 0; i < nr; i++)
757 for (
size_t j = nc/2+1; j < nc; j++)
758 out[j*stride + i*dist] =
conj (out[(nc - j)*stride + i*dist]);
769 size_t np = (dv.
length () > 2 ? dv.
numel () / nc / nr : 1);
770 size_t nrp = nr * np;
777 for (
size_t i = 0; i < nrp; i++)
779 ptr1 = out + i * (nc/2 + 1) + nrp*((nc-1)/2);
781 for (
size_t j = 0; j < nc/2+1; j++)
789 for (
size_t i = 0; i < np; i++)
791 for (
size_t j = 1; j < nr; j++)
792 for (
size_t k = nc/2+1; k < nc; k++)
793 out[k + (j + i*nr)*nc] =
conj (out[nc - k + ((i+1)*nr - j)*nc]);
795 for (
size_t j = nc/2+1; j < nc; j++)
796 out[j + i*nr*nc] =
conj (out[(i*nr+1)*nc - j]);
803 size_t jstart = dv(0) * dv(1);
804 size_t kstep = dv(0);
805 size_t nel = dv.
numel ();
807 for (
int inner = 2; inner < dv.
length (); inner++)
809 size_t jmax = jstart * dv(inner);
810 for (
size_t i = 0; i < nel; i+=jmax)
811 for (
size_t j = jstart, jj = jmax-jstart; j < jj;
812 j+=jstart, jj-=jstart)
813 for (
size_t k = 0; k < jstart; k+= kstep)
814 for (
size_t l = nc/2+1; l < nc; l++)
816 T tmp = out[i+ j + k + l];
817 out[i + j + k + l] = out[i + jj + k + l];
818 out[i + jj + k + l] = tmp;
830 dist = (dist < 0 ? npts : dist);
834 stride, dist, in, out);
836 fftw_execute_dft_r2c (plan, (const_cast<double *>(in)),
837 reinterpret_cast<fftw_complex *> (out));
850 dist = (dist < 0 ? npts : dist);
857 fftw_execute_dft (plan,
858 reinterpret_cast<fftw_complex *> (const_cast<Complex *>(in)),
859 reinterpret_cast<fftw_complex *> (out));
869 dist = (dist < 0 ? npts : dist);
876 fftw_execute_dft (plan,
877 reinterpret_cast<fftw_complex *> (const_cast<Complex *>(in)),
878 reinterpret_cast<fftw_complex *> (out));
881 for (
size_t j = 0; j < nsamples; j++)
882 for (
size_t i = 0; i < npts; i++)
883 out[i*stride + j*dist] /= scale;
893 for (
int i = 0; i < rank; i++)
904 fftw_execute_dft_r2c (plan, (const_cast<double *>(in)),
905 reinterpret_cast<fftw_complex *> (out+ offset));
919 for (
int i = 0; i < rank; i++)
923 dv, 1, 1, dist, in, out);
925 fftw_execute_dft (plan,
926 reinterpret_cast<fftw_complex *> (const_cast<Complex *>(in)),
927 reinterpret_cast<fftw_complex *> (out));
937 for (
int i = 0; i < rank; i++)
941 dv, 1, 1, dist, in, out);
943 fftw_execute_dft (plan,
944 reinterpret_cast<fftw_complex *> (const_cast<Complex *>(in)),
945 reinterpret_cast<fftw_complex *> (out));
947 const size_t npts = dv.
numel ();
949 for (
size_t i = 0; i < npts; i++)
959 dist = (dist < 0 ? npts : dist);
966 fftwf_execute_dft_r2c (plan, (const_cast<float *>(in)),
967 reinterpret_cast<fftwf_complex *> (out));
980 dist = (dist < 0 ? npts : dist);
988 fftwf_execute_dft (plan,
989 reinterpret_cast<fftwf_complex *> (const_cast<FloatComplex *>(in)),
990 reinterpret_cast<fftwf_complex *> (out));
1000 dist = (dist < 0 ? npts : dist);
1008 fftwf_execute_dft (plan,
1009 reinterpret_cast<fftwf_complex *> (const_cast<FloatComplex *>(in)),
1010 reinterpret_cast<fftwf_complex *> (out));
1013 for (
size_t j = 0; j < nsamples; j++)
1014 for (
size_t i = 0; i < npts; i++)
1015 out[i*stride + j*dist] /= scale;
1025 for (
int i = 0; i < rank; i++)
1037 fftwf_execute_dft_r2c (plan, (const_cast<float *>(in)),
1038 reinterpret_cast<fftwf_complex *> (out+ offset));
1052 for (
int i = 0; i < rank; i++)
1059 fftwf_execute_dft (plan,
1060 reinterpret_cast<fftwf_complex *> (const_cast<FloatComplex *>(in)),
1061 reinterpret_cast<fftwf_complex *> (out));
1071 for (
int i = 0; i < rank; i++)
1078 fftwf_execute_dft (plan,
1079 reinterpret_cast<fftwf_complex *> (const_cast<FloatComplex *>(in)),
1080 reinterpret_cast<fftwf_complex *> (out));
1082 const size_t npts = dv.
numel ();
1084 for (
size_t i = 0; i < npts; i++)