changeset 21451:763e30285935

eliminate configuration macros from oct-fftw.h (bug #41027) * oct-fftw.h: Eliminate FFTW configuration macros. Don't include fftw.h. * oct-fftw.h, oct-fftw.cc (octave_fftw_planner, octave_float_fftw_planner): Use void pointer to represent fftw_plan and fftwf_plan and cast as needed. * fftw.cc: Include fftw.h
author John W. Eaton <jwe@octave.org>
date Tue, 15 Mar 2016 14:58:19 -0400
parents 2868abbc88eb
children 769f9a7c02ae
files libinterp/dldfcn/fftw.cc liboctave/numeric/oct-fftw.cc liboctave/numeric/oct-fftw.h
diffstat 3 files changed, 133 insertions(+), 119 deletions(-) [+]
line wrap: on
line diff
--- a/libinterp/dldfcn/fftw.cc	Mon Mar 14 23:22:32 2016 -0400
+++ b/libinterp/dldfcn/fftw.cc	Tue Mar 15 14:58:19 2016 -0400
@@ -26,6 +26,10 @@
 
 #include <algorithm>
 
+#if defined (HAVE_FFTW3_H)
+#  include <fftw3.h>
+#endif
+
 #include "oct-fftw.h"
 
 #include "defun-dld.h"
--- a/liboctave/numeric/oct-fftw.cc	Mon Mar 14 23:22:32 2016 -0400
+++ b/liboctave/numeric/oct-fftw.cc	Tue Mar 15 14:58:19 2016 -0400
@@ -29,6 +29,10 @@
 #include <iostream>
 #include <vector>
 
+#if defined (HAVE_FFTW3_H)
+#  include <fftw3.h>
+#endif
+
 #include "lo-error.h"
 #include "oct-fftw.h"
 #include "quit.h"
@@ -61,7 +65,7 @@
 
 octave_fftw_planner::octave_fftw_planner (void)
   : meth (ESTIMATE), rplan (0), rd (0), rs (0), rr (0), rh (0), rn (),
-    rsimd_align (false)
+    rsimd_align (false), nthreads (1)
 {
   plan[0] = plan[1] = 0;
   d[0] = d[1] = s[0] = s[1] = r[0] = r[1] = h[0] = h[1] = 0;
@@ -74,8 +78,8 @@
   if (! init_ret)
     (*current_liboctave_error_handler) ("Error initializing FFTW threads");
 
-  //Use number of processors available to the current process
-  //This can be later changed with fftw ("threads", nthreads)
+  // Use number of processors available to the current process
+  // This can be later changed with fftw ("threads", nthreads).
   nthreads = num_processors (NPROC_CURRENT);
   fftw_plan_with_nthreads (nthreads);
 #endif
@@ -88,15 +92,15 @@
 {
   fftw_plan *plan_p;
 
-  plan_p = &rplan;
+  plan_p = reinterpret_cast<fftw_plan *> (&rplan);
   if (*plan_p)
     fftw_destroy_plan (*plan_p);
 
-  plan_p = &plan[0];
+  plan_p = reinterpret_cast<fftw_plan *> (&plan[0]);
   if (*plan_p)
     fftw_destroy_plan (*plan_p);
 
-  plan_p = &plan[1];
+  plan_p = reinterpret_cast<fftw_plan *> (&plan[1]);
   if (*plan_p)
     fftw_destroy_plan (*plan_p);
 }
@@ -121,10 +125,27 @@
   return retval;
 }
 
-#define CHECK_SIMD_ALIGNMENT(x) \
+void
+octave_fftw_planner::threads (int nt)
+{
+#if defined (HAVE_FFTW3F_THREADS)
+  if (instance_ok () && nt != threads ())
+    {
+      instance->nthreads = nt;
+      fftw_plan_with_nthreads (nt);
+      // Clear the current plans.
+      instance->rplan = instance->plan[0] = instance->plan[1] = 0;
+    }
+#else
+  (*current_liboctave_warning_handler)
+    ("unable to change number of threads without FFTW thread support");
+#endif
+}
+
+#define CHECK_SIMD_ALIGNMENT(x)                         \
   (((reinterpret_cast<ptrdiff_t> (x)) & 0xF) == 0)
 
-fftw_plan
+void *
 octave_fftw_planner::do_create_plan (int dir, const int rank,
                                      const dim_vector dims,
                                      octave_idx_type howmany,
@@ -133,7 +154,7 @@
                                      const Complex *in, Complex *out)
 {
   int which = (dir == FFTW_FORWARD) ? 0 : 1;
-  fftw_plan *cur_plan_p = &plan[which];
+  fftw_plan *cur_plan_p = reinterpret_cast<fftw_plan *> (&plan[which]);
   bool create_new_plan = false;
   bool ioalign = CHECK_SIMD_ALIGNMENT (in) && CHECK_SIMD_ALIGNMENT (out);
   bool ioinplace = (in == out);
@@ -249,14 +270,14 @@
   return *cur_plan_p;
 }
 
-fftw_plan
+void *
 octave_fftw_planner::do_create_plan (const int rank, const dim_vector dims,
                                      octave_idx_type howmany,
                                      octave_idx_type stride,
                                      octave_idx_type dist,
                                      const double *in, Complex *out)
 {
-  fftw_plan *cur_plan_p = &rplan;
+  fftw_plan *cur_plan_p = reinterpret_cast<fftw_plan *> (&rplan);
   bool create_new_plan = false;
   bool ioalign = CHECK_SIMD_ALIGNMENT (in) && CHECK_SIMD_ALIGNMENT (out);
 
@@ -385,11 +406,11 @@
         {
           meth = _meth;
           if (rplan)
-            fftw_destroy_plan (rplan);
+            fftw_destroy_plan (reinterpret_cast<fftw_plan> (rplan));
           if (plan[0])
-            fftw_destroy_plan (plan[0]);
+            fftw_destroy_plan (reinterpret_cast<fftw_plan> (plan[0]));
           if (plan[1])
-            fftw_destroy_plan (plan[1]);
+            fftw_destroy_plan (reinterpret_cast<fftw_plan> (plan[1]));
           rplan = plan[0] = plan[1] = 0;
         }
     }
@@ -402,7 +423,7 @@
 
 octave_float_fftw_planner::octave_float_fftw_planner (void)
   : meth (ESTIMATE), rplan (0), rd (0), rs (0), rr (0), rh (0), rn (),
-    rsimd_align (false)
+    rsimd_align (false), nthreads (1)
 {
   plan[0] = plan[1] = 0;
   d[0] = d[1] = s[0] = s[1] = r[0] = r[1] = h[0] = h[1] = 0;
@@ -415,8 +436,8 @@
   if (! init_ret)
     (*current_liboctave_error_handler) ("Error initializing FFTW3F threads");
 
-  //Use number of processors available to the current process
-  //This can be later changed with fftw ("threads", nthreads)
+  // Use number of processors available to the current process
+  // This can be later changed with fftw ("threads", nthreads).
   nthreads = num_processors (NPROC_CURRENT);
   fftwf_plan_with_nthreads (nthreads);
 #endif
@@ -429,15 +450,15 @@
 {
   fftwf_plan *plan_p;
 
-  plan_p = &rplan;
+  plan_p = reinterpret_cast<fftwf_plan *> (&rplan);
   if (*plan_p)
     fftwf_destroy_plan (*plan_p);
 
-  plan_p = &plan[0];
+  plan_p = reinterpret_cast<fftwf_plan *> (&plan[0]);
   if (*plan_p)
     fftwf_destroy_plan (*plan_p);
 
-  plan_p = &plan[1];
+  plan_p = reinterpret_cast<fftwf_plan *> (&plan[1]);
   if (*plan_p)
     fftwf_destroy_plan (*plan_p);
 }
@@ -462,7 +483,24 @@
   return retval;
 }
 
-fftwf_plan
+void
+octave_float_fftw_planner::threads (int nt)
+{
+#if defined (HAVE_FFTW3F_THREADS)
+  if (instance_ok () && nt != threads ())
+    {
+      instance->nthreads = nt;
+      fftwf_plan_with_nthreads (nt);
+      // Clear the current plans.
+      instance->rplan = instance->plan[0] = instance->plan[1] = 0;
+    }
+#else
+  (*current_liboctave_warning_handler)
+    ("unable to change number of threads without FFTW thread support");
+#endif
+}
+
+void *
 octave_float_fftw_planner::do_create_plan (int dir, const int rank,
                                            const dim_vector dims,
                                            octave_idx_type howmany,
@@ -472,7 +510,7 @@
                                            FloatComplex *out)
 {
   int which = (dir == FFTW_FORWARD) ? 0 : 1;
-  fftwf_plan *cur_plan_p = &plan[which];
+  fftwf_plan *cur_plan_p = reinterpret_cast<fftwf_plan *> (&plan[which]);
   bool create_new_plan = false;
   bool ioalign = CHECK_SIMD_ALIGNMENT (in) && CHECK_SIMD_ALIGNMENT (out);
   bool ioinplace = (in == out);
@@ -588,7 +626,7 @@
   return *cur_plan_p;
 }
 
-fftwf_plan
+void *
 octave_float_fftw_planner::do_create_plan (const int rank,
                                            const dim_vector dims,
                                            octave_idx_type howmany,
@@ -596,7 +634,7 @@
                                            octave_idx_type dist,
                                            const float *in, FloatComplex *out)
 {
-  fftwf_plan *cur_plan_p = &rplan;
+  fftwf_plan *cur_plan_p = reinterpret_cast<fftwf_plan *> (&rplan);
   bool create_new_plan = false;
   bool ioalign = CHECK_SIMD_ALIGNMENT (in) && CHECK_SIMD_ALIGNMENT (out);
 
@@ -725,11 +763,11 @@
         {
           meth = _meth;
           if (rplan)
-            fftwf_destroy_plan (rplan);
+            fftwf_destroy_plan (reinterpret_cast<fftwf_plan> (rplan));
           if (plan[0])
-            fftwf_destroy_plan (plan[0]);
+            fftwf_destroy_plan (reinterpret_cast<fftwf_plan> (plan[0]));
           if (plan[1])
-            fftwf_destroy_plan (plan[1]);
+            fftwf_destroy_plan (reinterpret_cast<fftwf_plan> (plan[1]));
           rplan = plan[0] = plan[1] = 0;
         }
     }
@@ -824,8 +862,9 @@
   dist = (dist < 0 ? npts : dist);
 
   dim_vector dv (npts, 1);
-  fftw_plan plan = octave_fftw_planner::create_plan (1, dv, nsamples,
-                                                     stride, dist, in, out);
+  void *vplan = octave_fftw_planner::create_plan (1, dv, nsamples,
+                                                  stride, dist, in, out);
+  fftw_plan plan = reinterpret_cast<fftw_plan> (vplan);
 
   fftw_execute_dft_r2c (plan, (const_cast<double *>(in)),
                         reinterpret_cast<fftw_complex *> (out));
@@ -844,9 +883,10 @@
   dist = (dist < 0 ? npts : dist);
 
   dim_vector dv (npts, 1);
-  fftw_plan plan = octave_fftw_planner::create_plan (FFTW_FORWARD, 1, dv,
+  void *vplan = octave_fftw_planner::create_plan (FFTW_FORWARD, 1, dv,
                                                      nsamples, stride,
                                                      dist, in, out);
+  fftw_plan plan = reinterpret_cast<fftw_plan> (vplan);
 
   fftw_execute_dft (plan,
                     reinterpret_cast<fftw_complex *> (const_cast<Complex *>(in)),
@@ -863,9 +903,10 @@
   dist = (dist < 0 ? npts : dist);
 
   dim_vector dv (npts, 1);
-  fftw_plan plan = octave_fftw_planner::create_plan (FFTW_BACKWARD, 1, dv,
+  void *vplan = octave_fftw_planner::create_plan (FFTW_BACKWARD, 1, dv,
                                                      nsamples, stride,
                                                      dist, in, out);
+  fftw_plan plan = reinterpret_cast<fftw_plan> (vplan);
 
   fftw_execute_dft (plan,
                     reinterpret_cast<fftw_complex *> (const_cast<Complex *>(in)),
@@ -892,8 +933,9 @@
 
   octave_idx_type offset = (dv.numel () / dv(0)) * ((dv(0) - 1) / 2);
 
-  fftw_plan plan = octave_fftw_planner::create_plan (rank, dv, 1, 1, dist,
-                                                     in, out + offset);
+  void *vplan = octave_fftw_planner::create_plan (rank, dv, 1, 1, dist,
+                                                  in, out + offset);
+  fftw_plan plan = reinterpret_cast<fftw_plan> (vplan);
 
   fftw_execute_dft_r2c (plan, (const_cast<double *>(in)),
                         reinterpret_cast<fftw_complex *> (out+ offset));
@@ -913,8 +955,9 @@
   for (int i = 0; i < rank; i++)
     dist *= dv(i);
 
-  fftw_plan plan = octave_fftw_planner::create_plan (FFTW_FORWARD, rank,
-                                                     dv, 1, 1, dist, in, out);
+  void *vplan = octave_fftw_planner::create_plan (FFTW_FORWARD, rank,
+                                                  dv, 1, 1, dist, in, out);
+  fftw_plan plan = reinterpret_cast<fftw_plan> (vplan);
 
   fftw_execute_dft (plan,
                     reinterpret_cast<fftw_complex *> (const_cast<Complex *>(in)),
@@ -931,8 +974,9 @@
   for (int i = 0; i < rank; i++)
     dist *= dv(i);
 
-  fftw_plan plan = octave_fftw_planner::create_plan (FFTW_BACKWARD, rank,
-                                                     dv, 1, 1, dist, in, out);
+  void *vplan = octave_fftw_planner::create_plan (FFTW_BACKWARD, rank,
+                                                  dv, 1, 1, dist, in, out);
+  fftw_plan plan = reinterpret_cast<fftw_plan> (vplan);
 
   fftw_execute_dft (plan,
                     reinterpret_cast<fftw_complex *> (const_cast<Complex *>(in)),
@@ -953,9 +997,10 @@
   dist = (dist < 0 ? npts : dist);
 
   dim_vector dv (npts, 1);
-  fftwf_plan plan = octave_float_fftw_planner::create_plan (1, dv, nsamples,
-                                                            stride, dist,
-                                                            in, out);
+  void *vplan = octave_float_fftw_planner::create_plan (1, dv, nsamples,
+                                                        stride, dist,
+                                                        in, out);
+  fftwf_plan plan = reinterpret_cast<fftwf_plan> (vplan);
 
   fftwf_execute_dft_r2c (plan, (const_cast<float *>(in)),
                          reinterpret_cast<fftwf_complex *> (out));
@@ -974,10 +1019,11 @@
   dist = (dist < 0 ? npts : dist);
 
   dim_vector dv (npts, 1);
-  fftwf_plan plan = octave_float_fftw_planner::create_plan (FFTW_FORWARD, 1,
-                                                            dv, nsamples,
-                                                            stride, dist,
-                                                            in, out);
+  void *vplan = octave_float_fftw_planner::create_plan (FFTW_FORWARD, 1,
+                                                        dv, nsamples,
+                                                        stride, dist,
+                                                        in, out);
+  fftwf_plan plan = reinterpret_cast<fftwf_plan> (vplan);
 
   fftwf_execute_dft (plan,
                      reinterpret_cast<fftwf_complex *> (const_cast<FloatComplex *>(in)),
@@ -994,10 +1040,11 @@
   dist = (dist < 0 ? npts : dist);
 
   dim_vector dv (npts, 1);
-  fftwf_plan plan = octave_float_fftw_planner::create_plan (FFTW_BACKWARD, 1,
-                                                            dv, nsamples,
-                                                            stride, dist,
-                                                            in, out);
+  void *vplan = octave_float_fftw_planner::create_plan (FFTW_BACKWARD, 1,
+                                                        dv, nsamples,
+                                                        stride, dist,
+                                                        in, out);
+  fftwf_plan plan = reinterpret_cast<fftwf_plan> (vplan);
 
   fftwf_execute_dft (plan,
                      reinterpret_cast<fftwf_complex *> (const_cast<FloatComplex *>(in)),
@@ -1024,9 +1071,10 @@
 
   octave_idx_type offset = (dv.numel () / dv(0)) * ((dv(0) - 1) / 2);
 
-  fftwf_plan plan = octave_float_fftw_planner::create_plan (rank, dv, 1, 1,
-                                                            dist, in,
-                                                            out + offset);
+  void *vplan = octave_float_fftw_planner::create_plan (rank, dv, 1, 1,
+                                                        dist, in,
+                                                        out + offset);
+  fftwf_plan plan = reinterpret_cast<fftwf_plan> (vplan);
 
   fftwf_execute_dft_r2c (plan, (const_cast<float *>(in)),
                          reinterpret_cast<fftwf_complex *> (out+ offset));
@@ -1046,9 +1094,10 @@
   for (int i = 0; i < rank; i++)
     dist *= dv(i);
 
-  fftwf_plan plan = octave_float_fftw_planner::create_plan (FFTW_FORWARD,
-                                                            rank, dv, 1, 1,
-                                                            dist, in, out);
+  void *vplan = octave_float_fftw_planner::create_plan (FFTW_FORWARD,
+                                                        rank, dv, 1, 1,
+                                                        dist, in, out);
+  fftwf_plan plan = reinterpret_cast<fftwf_plan> (vplan);
 
   fftwf_execute_dft (plan,
                      reinterpret_cast<fftwf_complex *> (const_cast<FloatComplex *>(in)),
@@ -1065,9 +1114,10 @@
   for (int i = 0; i < rank; i++)
     dist *= dv(i);
 
-  fftwf_plan plan = octave_float_fftw_planner::create_plan (FFTW_BACKWARD,
+  void *vplan = octave_float_fftw_planner::create_plan (FFTW_BACKWARD,
                                                             rank, dv, 1, 1,
                                                             dist, in, out);
+  fftwf_plan plan = reinterpret_cast<fftwf_plan> (vplan);
 
   fftwf_execute_dft (plan,
                      reinterpret_cast<fftwf_complex *> (const_cast<FloatComplex *>(in)),
--- a/liboctave/numeric/oct-fftw.h	Mon Mar 14 23:22:32 2016 -0400
+++ b/liboctave/numeric/oct-fftw.h	Tue Mar 15 14:58:19 2016 -0400
@@ -27,15 +27,9 @@
 
 #include <cstddef>
 
-#if defined (HAVE_FFTW3_H)
-#  include <fftw3.h>
-#endif
-
 #include "oct-cmplx.h"
 #include "dim-vector.h"
 
-#if defined (HAVE_FFTW)
-
 class
 OCTAVE_API
 octave_fftw_planner
@@ -60,31 +54,27 @@
 
   static bool instance_ok (void);
 
-  static fftw_plan
+  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 fftw_plan dummy;
-
     return instance_ok ()
            ? instance->do_create_plan (dir, rank, dims, howmany, stride,
                                        dist, in, out)
-           : dummy;
+           : 0;
   }
 
-  static fftw_plan
+  static void *
   create_plan (const int rank, const dim_vector dims,
                octave_idx_type howmany, octave_idx_type stride,
                octave_idx_type dist, const double *in, Complex *out)
   {
-    static fftw_plan dummy;
-
     return instance_ok ()
            ? instance->do_create_plan (rank, dims, howmany, stride, dist,
                                        in, out)
-           : dummy;
+           : 0;
   }
 
   static FftwMethod method (void)
@@ -101,23 +91,12 @@
     return instance_ok () ? instance->do_method (_meth) : dummy;
   }
 
-#if defined (HAVE_FFTW3F_THREADS)
-  static void threads (int _nthreads)
-  {
-    if (instance_ok () && _nthreads != threads ())
-      {
-        instance->nthreads = _nthreads;
-        fftw_plan_with_nthreads (_nthreads);
-        //Clear the current plans
-        instance->rplan = instance->plan[0] = instance->plan[1] = 0;
-      }
-  }
+  static void threads (int nt);
 
-  static int threads ()
+  static int threads (void)
   {
     return instance_ok () ? instance->nthreads : 0;
   }
-#endif
 
 private:
 
@@ -131,13 +110,13 @@
 
   static void cleanup_instance (void) { delete instance; instance = 0; }
 
-  fftw_plan
+  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);
 
-  fftw_plan
+  void *
   do_create_plan (const int rank, const dim_vector dims,
                   octave_idx_type howmany, octave_idx_type stride,
                   octave_idx_type dist, const double *in, Complex *out);
@@ -151,7 +130,7 @@
   // FIXME: perhaps this should be split into two classes?
 
   // Plan for fft and ifft of complex values
-  fftw_plan plan[2];
+  void *plan[2];
 
   // dist
   octave_idx_type d[2];
@@ -172,7 +151,7 @@
   bool inplace[2];
 
   // Plan for fft of real values
-  fftw_plan rplan;
+  void *rplan;
 
   // dist
   octave_idx_type rd;
@@ -191,10 +170,9 @@
 
   bool rsimd_align;
 
-#if defined (HAVE_FFTW3_THREADS)
-  //number of threads when compiled with Multi-threading support
+  // number of threads.  Always 1 unless compiled with Multi-threading
+  // support.
   int nthreads;
-#endif
 };
 
 class
@@ -221,31 +199,27 @@
 
   static bool instance_ok (void);
 
-  static fftwf_plan
+  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)
   {
-    static fftwf_plan dummy;
-
     return instance_ok ()
            ? instance->do_create_plan (dir, rank, dims, howmany, stride,
                                        dist, in, out)
-           : dummy;
+           : 0;
   }
 
-  static fftwf_plan
+  static void *
   create_plan (const int rank, const dim_vector dims,
                octave_idx_type howmany, octave_idx_type stride,
                octave_idx_type dist, const float *in, FloatComplex *out)
   {
-    static fftwf_plan dummy;
-
     return instance_ok ()
            ? instance->do_create_plan (rank, dims, howmany, stride, dist,
                                        in, out)
-           : dummy;
+           : 0;
   }
 
   static FftwMethod method (void)
@@ -262,23 +236,12 @@
     return instance_ok () ? instance->do_method (_meth) : dummy;
   }
 
-#if defined (HAVE_FFTW3F_THREADS)
-  static void threads (int _nthreads)
-  {
-    if (instance_ok () && _nthreads != threads ())
-      {
-        instance->nthreads = _nthreads;
-        fftwf_plan_with_nthreads (_nthreads);
-        //Clear the current plans
-        instance->rplan = instance->plan[0] = instance->plan[1] = 0;
-      }
-  }
+  static void threads (int nt);
 
-  static int threads ()
+  static int threads (void)
   {
     return instance_ok () ? instance->nthreads : 0;
   }
-#endif
 
 private:
 
@@ -292,13 +255,13 @@
 
   static void cleanup_instance (void) { delete instance; instance = 0; }
 
-  fftwf_plan
+  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);
 
-  fftwf_plan
+  void *
   do_create_plan (const int rank, const dim_vector dims,
                   octave_idx_type howmany, octave_idx_type stride,
                   octave_idx_type dist, const float *in, FloatComplex *out);
@@ -312,7 +275,7 @@
   // FIXME: perhaps this should be split into two classes?
 
   // Plan for fft and ifft of complex values
-  fftwf_plan plan[2];
+  void *plan[2];
 
   // dist
   octave_idx_type d[2];
@@ -333,7 +296,7 @@
   bool inplace[2];
 
   // Plan for fft of real values
-  fftwf_plan rplan;
+  void *rplan;
 
   // dist
   octave_idx_type rd;
@@ -352,10 +315,9 @@
 
   bool rsimd_align;
 
-#if defined (HAVE_FFTW3F_THREADS)
-  //number of threads when compiled with Multi-threading support
+  // number of threads.  Always 1 unless compiled with Multi-threading
+  // support.
   int nthreads;
-#endif
 };
 
 class
@@ -403,5 +365,3 @@
 };
 
 #endif
-
-#endif