changeset 32365:65e27afd86a2 stable

fft: Avoid segmentation fault with ND-arrays (bug #64729). * liboctave/numeric/act-fftw.cc (fftw_planner::do_create_plan, float_fftw_planner::do_create_plan): Allocate sufficiently large buffer for FFTW plan (as suggested by Markus Meisinger). * libinterp/corefcn/fft.cc: Add test cases.
author Markus Mützel <markus.muetzel@gmx.de>
date Sun, 01 Oct 2023 17:42:52 +0200
parents 4acb694518de
children e3c66ad99652 903c9100178b
files libinterp/corefcn/fft.cc liboctave/numeric/oct-fftw.cc
diffstat 2 files changed, 24 insertions(+), 2 deletions(-) [+]
line wrap: on
line diff
--- a/libinterp/corefcn/fft.cc	Wed Sep 20 12:45:45 2023 -0400
+++ b/libinterp/corefcn/fft.cc	Sun Oct 01 17:42:52 2023 +0200
@@ -312,6 +312,28 @@
 %! S(N-n+1) = N/2;
 %!
 %! assert (ifft (S), s, 4*N* eps ("single"));
+
+%!testif HAVE_FFTW <*64729>
+%! x = rand (5, 5, 5, 5);
+%! old_planner = fftw ('planner', 'measure');
+%! unwind_protect
+%!   y = fft (x);
+%!   y = fft (x);  # second invocation might crash Octave
+%!   assert (size (y), [5, 5, 5, 5]);
+%! unwind_protect_cleanup
+%!   fftw ('planner', old_planner);
+%! end_unwind_protect
+
+%!testif HAVE_FFTW <*64729>
+%! x = rand (5, 5, 5, 5, 'single');
+%! old_planner = fftw ('planner', 'measure');
+%! unwind_protect
+%!   y = fft (x);
+%!   y = fft (x);  # second invocation might crash Octave
+%!   assert (size (y), [5, 5, 5, 5]);
+%! unwind_protect_cleanup
+%!   fftw ('planner', old_planner);
+%! end_unwind_protect
 */
 
 OCTAVE_END_NAMESPACE(octave)
--- a/liboctave/numeric/oct-fftw.cc	Wed Sep 20 12:45:45 2023 -0400
+++ b/liboctave/numeric/oct-fftw.cc	Sun Oct 01 17:42:52 2023 +0200
@@ -366,7 +366,7 @@
       if (plan_destroys_in)
         {
           // Create matrix with the same size and 16-byte alignment as input
-          OCTAVE_LOCAL_BUFFER (double, itmp, nn + 32);
+          OCTAVE_LOCAL_BUFFER (double, itmp, nn * howmany + 32);
           itmp = reinterpret_cast<double *>
                  (((reinterpret_cast<std::ptrdiff_t> (itmp) + 15) & ~ 0xF) +
                   ((reinterpret_cast<std::ptrdiff_t> (in)) & 0xF));
@@ -721,7 +721,7 @@
       if (plan_destroys_in)
         {
           // Create matrix with the same size and 16-byte alignment as input
-          OCTAVE_LOCAL_BUFFER (float, itmp, nn + 32);
+          OCTAVE_LOCAL_BUFFER (float, itmp, nn * howmany + 32);
           itmp = reinterpret_cast<float *>
                  (((reinterpret_cast<std::ptrdiff_t> (itmp) + 15) & ~ 0xF) +
                   ((reinterpret_cast<std::ptrdiff_t> (in)) & 0xF));