26#if defined (HAVE_CONFIG_H)
30#if defined (HAVE_FFTW3_H)
40#if defined (HAVE_FFTW3_THREADS) || defined (HAVE_FFTW3F_THREADS)
46#if defined (HAVE_FFTW)
69 : m_meth (ESTIMATE), m_rplan (nullptr), m_rd (0), m_rs (0), m_rr (0),
70 m_rh (0), m_rn (), m_rsimd_align (false), m_nthreads (1)
78#if defined (HAVE_FFTW3_THREADS)
79 int init_ret = fftw_init_threads ();
81 (*current_liboctave_error_handler) (
"Error initializing FFTW threads");
92 fftw_import_system_wisdom ();
99 plan_p =
reinterpret_cast<fftw_plan *
> (&
m_rplan);
101 fftw_destroy_plan (*plan_p);
103 plan_p =
reinterpret_cast<fftw_plan *
> (&
m_plan[0]);
105 fftw_destroy_plan (*plan_p);
107 plan_p =
reinterpret_cast<fftw_plan *
> (&
m_plan[1]);
109 fftw_destroy_plan (*plan_p);
129#if defined (HAVE_FFTW3_THREADS)
133 fftw_plan_with_nthreads (nt);
139 octave_unused_parameter (nt);
141 (*current_liboctave_warning_handler)
142 (
"unable to change number of threads without FFTW thread support");
146#define CHECK_SIMD_ALIGNMENT(x) \
147 (((reinterpret_cast<std::ptrdiff_t> (x)) & 0xF) == 0)
157 int which = (dir == FFTW_FORWARD) ? 0 : 1;
158 fftw_plan *cur_plan_p =
reinterpret_cast<fftw_plan *
> (&
m_plan[which]);
159 bool create_new_plan =
false;
161 bool ioinplace = (in == out);
167 if (
m_plan[which] ==
nullptr ||
m_d[which] != dist ||
m_s[which] != stride
168 ||
m_r[which] != rank ||
m_h[which] != howmany
170 || ((ioalign !=
m_simd_align[which]) ? ! ioalign :
false))
171 create_new_plan =
true;
176 for (
int i = 0; i < rank; i++)
177 if (dims(i) !=
m_n[which](i))
179 create_new_plan =
true;
189 m_h[which] = howmany;
198 for (
int i = 0, j = rank-1; i < rank; i++, j--)
205 bool plan_destroys_in =
true;
211 plan_flags |= FFTW_ESTIMATE;
212 plan_destroys_in =
false;
215 plan_flags |= FFTW_MEASURE;
218 plan_flags |= FFTW_PATIENT;
221 plan_flags |= FFTW_EXHAUSTIVE;
225 plan_flags |= FFTW_MEASURE;
228 plan_flags |= FFTW_ESTIMATE;
229 plan_destroys_in =
false;
235 plan_flags &= ~FFTW_UNALIGNED;
237 plan_flags |= FFTW_UNALIGNED;
240 fftw_destroy_plan (*cur_plan_p);
242 if (plan_destroys_in)
246 itmp =
reinterpret_cast<Complex *
>
247 (((
reinterpret_cast<std::ptrdiff_t
> (itmp) + 15) & ~ 0xF) +
248 ((
reinterpret_cast<std::ptrdiff_t
> (in)) & 0xF));
251 = fftw_plan_many_dft (rank, tmp, howmany,
252 reinterpret_cast<fftw_complex *
> (itmp),
253 nullptr, stride, dist,
254 reinterpret_cast<fftw_complex *
> (out),
255 nullptr, stride, dist, dir, plan_flags);
260 = fftw_plan_many_dft (rank, tmp, howmany,
261 reinterpret_cast<fftw_complex *
> (
const_cast<Complex *
> (in)),
262 nullptr, stride, dist,
263 reinterpret_cast<fftw_complex *
> (out),
264 nullptr, stride, dist, dir, plan_flags);
267 if (*cur_plan_p ==
nullptr)
268 (*current_liboctave_error_handler) (
"Error creating FFTW plan");
279 const double *in,
Complex *out)
281 fftw_plan *cur_plan_p =
reinterpret_cast<fftw_plan *
> (&
m_rplan);
282 bool create_new_plan =
false;
291 create_new_plan =
true;
296 for (
int i = 0; i < rank; i++)
297 if (dims(i) !=
m_rn(i))
299 create_new_plan =
true;
317 for (
int i = 0, j = rank-1; i < rank; i++, j--)
324 bool plan_destroys_in =
true;
330 plan_flags |= FFTW_ESTIMATE;
331 plan_destroys_in =
false;
334 plan_flags |= FFTW_MEASURE;
337 plan_flags |= FFTW_PATIENT;
340 plan_flags |= FFTW_EXHAUSTIVE;
344 plan_flags |= FFTW_MEASURE;
347 plan_flags |= FFTW_ESTIMATE;
348 plan_destroys_in =
false;
354 plan_flags &= ~FFTW_UNALIGNED;
356 plan_flags |= FFTW_UNALIGNED;
359 fftw_destroy_plan (*cur_plan_p);
361 if (plan_destroys_in)
365 itmp =
reinterpret_cast<double *
>
366 (((
reinterpret_cast<std::ptrdiff_t
> (itmp) + 15) & ~ 0xF) +
367 ((
reinterpret_cast<std::ptrdiff_t
> (in)) & 0xF));
370 = fftw_plan_many_dft_r2c (rank, tmp, howmany, itmp,
371 nullptr, stride, dist,
372 reinterpret_cast<fftw_complex *
> (out),
373 nullptr, stride, dist, plan_flags);
378 = fftw_plan_many_dft_r2c (rank, tmp, howmany,
379 (
const_cast<double *
> (in)),
380 nullptr, stride, dist,
381 reinterpret_cast<fftw_complex *
> (out),
382 nullptr, stride, dist, plan_flags);
385 if (*cur_plan_p ==
nullptr)
386 (*current_liboctave_error_handler) (
"Error creating FFTW plan");
410 fftw_destroy_plan (
reinterpret_cast<fftw_plan
> (
m_rplan));
412 fftw_destroy_plan (
reinterpret_cast<fftw_plan
> (
m_plan[0]));
414 fftw_destroy_plan (
reinterpret_cast<fftw_plan
> (
m_plan[1]));
426 : m_meth (ESTIMATE), m_rplan (nullptr), m_rd (0), m_rs (0), m_rr (0),
427 m_rh (0), m_rn (), m_rsimd_align (false), m_nthreads (1)
435#if defined (HAVE_FFTW3F_THREADS)
436 int init_ret = fftwf_init_threads ();
438 (*current_liboctave_error_handler) (
"Error initializing FFTW3F threads");
449 fftwf_import_system_wisdom ();
456 plan_p =
reinterpret_cast<fftwf_plan *
> (&
m_rplan);
458 fftwf_destroy_plan (*plan_p);
460 plan_p =
reinterpret_cast<fftwf_plan *
> (&
m_plan[0]);
462 fftwf_destroy_plan (*plan_p);
464 plan_p =
reinterpret_cast<fftwf_plan *
> (&
m_plan[1]);
466 fftwf_destroy_plan (*plan_p);
486#if defined (HAVE_FFTW3F_THREADS)
490 fftwf_plan_with_nthreads (nt);
496 octave_unused_parameter (nt);
498 (*current_liboctave_warning_handler)
499 (
"unable to change number of threads without FFTW thread support");
512 int which = (dir == FFTW_FORWARD) ? 0 : 1;
513 fftwf_plan *cur_plan_p =
reinterpret_cast<fftwf_plan *
> (&
m_plan[which]);
514 bool create_new_plan =
false;
516 bool ioinplace = (in == out);
522 if (
m_plan[which] ==
nullptr ||
m_d[which] != dist ||
m_s[which] != stride
523 ||
m_r[which] != rank ||
m_h[which] != howmany
525 || ((ioalign !=
m_simd_align[which]) ? ! ioalign :
false))
526 create_new_plan =
true;
531 for (
int i = 0; i < rank; i++)
532 if (dims(i) !=
m_n[which](i))
534 create_new_plan =
true;
544 m_h[which] = howmany;
553 for (
int i = 0, j = rank-1; i < rank; i++, j--)
560 bool plan_destroys_in =
true;
566 plan_flags |= FFTW_ESTIMATE;
567 plan_destroys_in =
false;
570 plan_flags |= FFTW_MEASURE;
573 plan_flags |= FFTW_PATIENT;
576 plan_flags |= FFTW_EXHAUSTIVE;
580 plan_flags |= FFTW_MEASURE;
583 plan_flags |= FFTW_ESTIMATE;
584 plan_destroys_in =
false;
590 plan_flags &= ~FFTW_UNALIGNED;
592 plan_flags |= FFTW_UNALIGNED;
595 fftwf_destroy_plan (*cur_plan_p);
597 if (plan_destroys_in)
602 (((
reinterpret_cast<std::ptrdiff_t
> (itmp) + 15) & ~ 0xF) +
603 ((
reinterpret_cast<std::ptrdiff_t
> (in)) & 0xF));
606 = fftwf_plan_many_dft (rank, tmp, howmany,
607 reinterpret_cast<fftwf_complex *
> (itmp),
608 nullptr, stride, dist,
609 reinterpret_cast<fftwf_complex *
> (out),
610 nullptr, stride, dist, dir, plan_flags);
615 = fftwf_plan_many_dft (rank, tmp, howmany,
616 reinterpret_cast<fftwf_complex *
> (
const_cast<FloatComplex *
> (in)),
617 nullptr, stride, dist,
618 reinterpret_cast<fftwf_complex *
> (out),
619 nullptr, stride, dist, dir, plan_flags);
622 if (*cur_plan_p ==
nullptr)
623 (*current_liboctave_error_handler) (
"Error creating FFTW plan");
636 fftwf_plan *cur_plan_p =
reinterpret_cast<fftwf_plan *
> (&
m_rplan);
637 bool create_new_plan =
false;
646 create_new_plan =
true;
651 for (
int i = 0; i < rank; i++)
652 if (dims(i) !=
m_rn(i))
654 create_new_plan =
true;
672 for (
int i = 0, j = rank-1; i < rank; i++, j--)
679 bool plan_destroys_in =
true;
685 plan_flags |= FFTW_ESTIMATE;
686 plan_destroys_in =
false;
689 plan_flags |= FFTW_MEASURE;
692 plan_flags |= FFTW_PATIENT;
695 plan_flags |= FFTW_EXHAUSTIVE;
699 plan_flags |= FFTW_MEASURE;
702 plan_flags |= FFTW_ESTIMATE;
703 plan_destroys_in =
false;
709 plan_flags &= ~FFTW_UNALIGNED;
711 plan_flags |= FFTW_UNALIGNED;
714 fftwf_destroy_plan (*cur_plan_p);
716 if (plan_destroys_in)
720 itmp =
reinterpret_cast<float *
>
721 (((
reinterpret_cast<std::ptrdiff_t
> (itmp) + 15) & ~ 0xF) +
722 ((
reinterpret_cast<std::ptrdiff_t
> (in)) & 0xF));
725 = fftwf_plan_many_dft_r2c (rank, tmp, howmany, itmp,
726 nullptr, stride, dist,
727 reinterpret_cast<fftwf_complex *
> (out),
728 nullptr, stride, dist, plan_flags);
733 = fftwf_plan_many_dft_r2c (rank, tmp, howmany,
734 (
const_cast<float *
> (in)),
735 nullptr, stride, dist,
736 reinterpret_cast<fftwf_complex *
> (out),
737 nullptr, stride, dist, plan_flags);
740 if (*cur_plan_p ==
nullptr)
741 (*current_liboctave_error_handler) (
"Error creating FFTW plan");
765 fftwf_destroy_plan (
reinterpret_cast<fftwf_plan
> (
m_rplan));
767 fftwf_destroy_plan (
reinterpret_cast<fftwf_plan
> (
m_plan[0]));
769 fftwf_destroy_plan (
reinterpret_cast<fftwf_plan
> (
m_plan[1]));
778 template <
typename T>
787 for (std::size_t i = 0; i < nr; i++)
788 for (std::size_t j = nc/2+1; j < nc; j++)
789 out[j*stride + i*dist] =
conj (out[(nc - j)*stride + i*dist]);
794 template <
typename T>
798 std::size_t nc = dv(0);
799 std::size_t nr = dv(1);
800 std::size_t np = (dv.
ndims () > 2 ? dv.
numel () / nc / nr : 1);
801 std::size_t nrp = nr * np;
808 for (std::size_t i = 0; i < nrp; i++)
810 ptr1 = out + i * (nc/2 + 1) + nrp*((nc-1)/2);
812 for (std::size_t j = 0; j < nc/2+1; j++)
820 for (std::size_t i = 0; i < np; i++)
822 for (std::size_t j = 1; j < nr; j++)
823 for (std::size_t k = nc/2+1; k < nc; k++)
824 out[k + (j + i*nr)*nc] =
conj (out[nc - k + ((i+1)*nr - j)*nc]);
826 for (std::size_t j = nc/2+1; j < nc; j++)
827 out[j + i*nr*nc] =
conj (out[(i*nr+1)*nc - j]);
834 std::size_t jstart = dv(0) * dv(1);
835 std::size_t kstep = dv(0);
836 std::size_t nel = dv.
numel ();
838 for (
int inner = 2; inner < dv.
ndims (); inner++)
840 std::size_t jmax = jstart * dv(inner);
841 for (std::size_t i = 0; i < nel; i+=jmax)
842 for (std::size_t j = jstart, jj = jmax-jstart; j < jj;
843 j+=jstart, jj-=jstart)
844 for (std::size_t k = 0; k < jstart; k+= kstep)
845 for (std::size_t l = nc/2+1; l < nc; l++)
847 T tmp = out[i+ j + k + l];
848 out[i + j + k + l] = out[i + jj + k + l];
849 out[i + jj + k + l] = tmp;
862 dist = (dist < 0 ? npts : dist);
866 stride, dist, in, out);
867 fftw_plan m_plan =
reinterpret_cast<fftw_plan
> (vplan);
869 fftw_execute_dft_r2c (m_plan, (
const_cast<double *
>(in)),
870 reinterpret_cast<fftw_complex *
> (out));
884 dist = (dist < 0 ? npts : dist);
890 fftw_plan m_plan =
reinterpret_cast<fftw_plan
> (vplan);
892 fftw_execute_dft (m_plan,
893 reinterpret_cast<fftw_complex *
> (
const_cast<Complex *
>(in)),
894 reinterpret_cast<fftw_complex *
> (out));
904 dist = (dist < 0 ? npts : dist);
908 stride, dist, in, out);
909 fftw_plan m_plan =
reinterpret_cast<fftw_plan
> (vplan);
911 fftw_execute_dft (m_plan,
912 reinterpret_cast<fftw_complex *
> (
const_cast<Complex *
>(in)),
913 reinterpret_cast<fftw_complex *
> (out));
916 for (std::size_t j = 0; j < nsamples; j++)
917 for (std::size_t i = 0; i < npts; i++)
918 out[i*stride + j*dist] /=
scale;
928 for (
int i = 0; i < rank; i++)
938 fftw_plan m_plan =
reinterpret_cast<fftw_plan
> (vplan);
940 fftw_execute_dft_r2c (m_plan, (
const_cast<double *
>(in)),
941 reinterpret_cast<fftw_complex *
> (out+ offset));
955 for (
int i = 0; i < rank; i++)
959 1, 1, dist, in, out);
960 fftw_plan m_plan =
reinterpret_cast<fftw_plan
> (vplan);
962 fftw_execute_dft (m_plan,
963 reinterpret_cast<fftw_complex *
> (
const_cast<Complex *
>(in)),
964 reinterpret_cast<fftw_complex *
> (out));
974 for (
int i = 0; i < rank; i++)
978 1, 1, dist, in, out);
979 fftw_plan m_plan =
reinterpret_cast<fftw_plan
> (vplan);
981 fftw_execute_dft (m_plan,
982 reinterpret_cast<fftw_complex *
> (
const_cast<Complex *
>(in)),
983 reinterpret_cast<fftw_complex *
> (out));
985 const std::size_t npts = dv.
numel ();
987 for (std::size_t i = 0; i < npts; i++)
998 dist = (dist < 0 ? npts : dist);
1003 fftwf_plan m_plan =
reinterpret_cast<fftwf_plan
> (vplan);
1005 fftwf_execute_dft_r2c (m_plan, (
const_cast<float *
>(in)),
1006 reinterpret_cast<fftwf_complex *
> (out));
1020 dist = (dist < 0 ? npts : dist);
1024 nsamples, stride, dist,
1026 fftwf_plan m_plan =
reinterpret_cast<fftwf_plan
> (vplan);
1028 fftwf_execute_dft (m_plan,
1029 reinterpret_cast<fftwf_complex *
> (
const_cast<FloatComplex *
>(in)),
1030 reinterpret_cast<fftwf_complex *
> (out));
1040 dist = (dist < 0 ? npts : dist);
1044 nsamples, stride, dist,
1046 fftwf_plan m_plan =
reinterpret_cast<fftwf_plan
> (vplan);
1048 fftwf_execute_dft (m_plan,
1049 reinterpret_cast<fftwf_complex *
> (
const_cast<FloatComplex *
>(in)),
1050 reinterpret_cast<fftwf_complex *
> (out));
1053 for (std::size_t j = 0; j < nsamples; j++)
1054 for (std::size_t i = 0; i < npts; i++)
1055 out[i*stride + j*dist] /=
scale;
1065 for (
int i = 0; i < rank; i++)
1075 fftwf_plan m_plan =
reinterpret_cast<fftwf_plan
> (vplan);
1077 fftwf_execute_dft_r2c (m_plan, (
const_cast<float *
>(in)),
1078 reinterpret_cast<fftwf_complex *
> (out+ offset));
1092 for (
int i = 0; i < rank; i++)
1096 1, 1, dist, in, out);
1097 fftwf_plan m_plan =
reinterpret_cast<fftwf_plan
> (vplan);
1099 fftwf_execute_dft (m_plan,
1100 reinterpret_cast<fftwf_complex *
> (
const_cast<FloatComplex *
>(in)),
1101 reinterpret_cast<fftwf_complex *
> (out));
1111 for (
int i = 0; i < rank; i++)
1115 1, 1, dist, in, out);
1116 fftwf_plan m_plan =
reinterpret_cast<fftwf_plan
> (vplan);
1118 fftwf_execute_dft (m_plan,
1119 reinterpret_cast<fftwf_complex *
> (
const_cast<FloatComplex *
>(in)),
1120 reinterpret_cast<fftwf_complex *
> (out));
1122 const std::size_t npts = dv.
numel ();
1124 for (std::size_t i = 0; i < npts; i++)
1135#if defined (HAVE_FFTW)
1145#if defined (HAVE_FFTW)
ComplexColumnVector conj(const ComplexColumnVector &a)
Vector representing the dimensions (size) of an Array.
octave_idx_type numel(int n=0) const
Number of elements that a matrix with this dimensions would have.
octave_idx_type ndims(void) const
Number of dimensions.
static fftw_planner * s_instance
FftwMethod do_method(void)
static void * create_plan(int dir, const int rank, const dim_vector &dims, octave_idx_type howmany, octave_idx_type stride, octave_idx_type dist, const Complex *in, Complex *out)
static bool instance_ok(void)
static void cleanup_instance(void)
void * do_create_plan(int dir, const int rank, const dim_vector &dims, octave_idx_type howmany, octave_idx_type stride, octave_idx_type dist, const Complex *in, Complex *out)
static int ifft(const Complex *in, Complex *out, std::size_t npts, std::size_t nsamples=1, octave_idx_type stride=1, octave_idx_type dist=-1)
static int fftNd(const double *, Complex *, const int, const dim_vector &)
static int ifftNd(const Complex *, Complex *, const int, const dim_vector &)
static int fft(const double *in, Complex *out, std::size_t npts, std::size_t nsamples=1, octave_idx_type stride=1, octave_idx_type dist=-1)
static void cleanup_instance(void)
void * do_create_plan(int dir, const int rank, const dim_vector &dims, octave_idx_type howmany, octave_idx_type stride, octave_idx_type dist, const FloatComplex *in, FloatComplex *out)
static bool instance_ok(void)
static void * create_plan(int dir, const int rank, const dim_vector &dims, octave_idx_type howmany, octave_idx_type stride, octave_idx_type dist, const FloatComplex *in, FloatComplex *out)
FftwMethod do_method(void)
static float_fftw_planner * s_instance
~float_fftw_planner(void)
void scale(Matrix &m, double x, double y, double z)
static void convert_packcomplex_Nd(T *out, const dim_vector &dv)
std::string fftw_version(void)
std::string fftwf_version(void)
static void convert_packcomplex_1d(T *out, std::size_t nr, std::size_t nc, octave_idx_type stride, octave_idx_type dist)
unsigned long int octave_num_processors_wrapper(enum octave_nproc_query octave_query)
@ OCTAVE_NPROC_CURRENT_OVERRIDABLE
std::complex< double > Complex
std::complex< float > FloatComplex
#define CHECK_SIMD_ALIGNMENT(x)
#define OCTAVE_LOCAL_BUFFER(T, buf, size)