octave-bug-tracker
[Top][All Lists]
Advanced

[Date Prev][Date Next][Thread Prev][Thread Next][Date Index][Thread Index]

[Octave-bug-tracker] [bug #43742] ifft, ifft2, ifftn should accept 'nons


From: Hg200
Subject: [Octave-bug-tracker] [bug #43742] ifft, ifft2, ifftn should accept 'nonsymmetric' and 'symmetric' options
Date: Sat, 25 Dec 2021 10:44:46 -0500 (EST)
User-agent: Mozilla/5.0 (X11; Fedora; Linux x86_64; rv:92.0) Gecko/20100101 Firefox/92.0

Follow-up Comment #7, bug #43742 (project octave):

1.) The attached patch "patch_ifft_symmetric.diff" implements the
"real(ifftx(..., 'symmetric'))" solution described in comment #6 and comment
#5. The changes made are simple and limited to the "fft.cc", "fft2.cc",
"fftn.cc" files. Because of the simplicity and because we don't break any
architecture, my recommendation is to push the file
"patch_ifft_symmetric.diff" to the mercurial repository now. If we decide to
use "FFTW c2r" transformations later on, these changes are still necessary. So
we don't make anything worse with this patch, but still fix the bug (see also
comment #6 and comment #5).

2.) The alternative is to use the FFTW routines "fftw_plan_many_dft_c2r" and
"fftw_execute_dft_c2r" according [1],[2] instead of "real(ifft())". I
implemented the "c2r" routines as an example for the "ComplexNDArray" data
type and attached a patch as "test_runtime_fftw_c2r_functions.diff". With this
patch we can do runtime measurements and we can see where the difficulties in
the implementation are: A lot of new code is added. In particular, a new
method must be introduced in CNDArray.cc to transport the "symflag"
information from the top software layer (fft.c) to the bottom software layer
(oct-fft.cc):


NDArray
ComplexNDArray::ifouriersym (int dim) const


In oct-fft.cc there are two new methods. A new (third) FFTW planner is called
to control the "c2r" transformation. This is in addition to the existing
planners used for the "r2c" or "c2c" transformations.


fftw::ifftsym (const Complex *in, double *out, std::size_t npts,
               std::size_t nsamples, octave_idx_type stride,
               octave_idx_type dist)



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 complex *in, double *out)


Regarding the runtime, we have to be careful what is compared to what. The
IFFT routine is very fast and the difference is hardly noticeable when running
"normal code". Are we comparing runtime when the planner is already set up, or
are we including the runtime of the planner? Are we comparing runtime between
processing large or small data sets? Are we comparing to "c2c" with or without
taking the real part? I would say, depending on the use case, we have no
advantage at worst and 30% runtime improvement at best if we compare "c2r"
with "c2c" without taking the real part. If we include taking the real part,
then its more.

Conclusion:
ifouriersym() needs to be also introduced in dNDArray.cc, fNDArray.cc and
fCNDArray.cc. Additionally, there must be also methods for ifouriersym2() and
ifouriersymn(). Where else does the code need to be changed?  Are we happy
adding those new methods to the existing datatypes for just supporting a
"flag"? IMO this is quite a big change.

[1]: https://www.fftw.org/fftw3_doc/Advanced-Real_002ddata-DFTs.html
[2]: http://fftw.org/fftw3_doc/One_002dDimensional-DFTs-of-Real-Data.html



    _______________________________________________________

Reply to this item at:

  <https://savannah.gnu.org/bugs/?43742>

_______________________________________________
  Message sent via Savannah
  https://savannah.gnu.org/




reply via email to

[Prev in Thread] Current Thread [Next in Thread]