27 #ifndef MLN_TRANSFORM_FFT_HH
28 # define MLN_TRANSFORM_FFT_HH
34 # include <mln/core/image/image2d.hh>
35 # include <mln/estim/min_max.hh>
36 # include <mln/opt/at.hh>
65 struct fft_trait< std::complex<T> >
120 bool ordered =
true)
const;
131 bool ordered =
true)
const;
168 bool ordered =
true)
const;
181 bool ordered =
true)
const;
210 typename fft_trait<T>::fftw_input*
in;
251 template <
typename D>
281 fft(
const image2d< std::complex<T> >& original_im);
302 # ifndef MLN_INCLUDE_ONLY
317 image2d< std::complex<T> >&
331 const unsigned nrows = trans_im.nrows();
332 const unsigned ncols = trans_im.ncols();
333 image2d<R> new_im(trans_im.domain());
335 for (
unsigned row = 0; row < new_im.nrows(); ++row)
336 for (
unsigned col = 0; col < new_im.ncols(); ++col)
339 (row + nrows / 2) % nrows,
340 (col + ncols / 2) % ncols));
343 mln_piter(image2d< std::complex<T> >) it(trans_im.domain());
345 new_im(it) = std::norm(trans_im(it));
353 fft<
T>::transformed_image_magn(
bool ordered)
const
355 return transformed_image_magn<T>(ordered);
367 double max = mln_min(
double);
368 mln_piter(image2d<T>) it(trans_im.domain());
370 if (std::norm(trans_im(it)) > max)
371 max = std::norm(trans_im(it));
373 const
unsigned nrows = trans_im.nrows();
374 const
unsigned ncols = trans_im.ncols();
375 image2d<R> new_im(trans_im.domain());
377 for (
unsigned row = 0; row < new_im.nrows(); ++row)
378 for (
unsigned col = 0; col < new_im.ncols(); ++col)
379 if (std::norm(opt::at(trans_im,
380 (row + nrows / 2) % nrows,
381 (col + ncols / 2) % ncols))
383 opt::at(new_im, row, col) = mln_max(R);
385 opt::at(new_im, row, col) =
386 (
double) mln_max(R) *
387 std::norm(opt::at(trans_im,
388 (row + nrows / 2) % nrows,
389 (col + ncols / 2) % ncols)) / (max * clip);
392 if (std::norm(trans_im(it)) >= max * clip)
393 new_im(it) = mln_max(R);
396 (
double) mln_max(R) * std::norm(trans_im(it)) / (max * clip);
403 fft<
T>::transformed_image_clipped_magn(
double clip,
bool ordered)
const
405 return transformed_image_clipped_magn<T>(clip, ordered);
414 return transformed_image_clipped_magn<R>(1, ordered);
422 return transformed_image_clipped_magn<T>(1, ordered);
438 image2d<R> new_im(trans_im.domain());
440 double max = mln_min(
double);
441 mln_piter(image2d<R>) it(trans_im.domain());
443 if (std::norm(trans_im(it)) > max)
444 max = std::norm(trans_im(it));
446 const
unsigned nrows = trans_im.nrows();
447 const
unsigned ncols = trans_im.ncols();
449 for (
unsigned row = 0; row < new_im.nrows(); ++row)
450 for (
unsigned col = 0; col < new_im.ncols(); ++col)
451 opt::at(new_im, row, col) =
452 log(a + b * std::norm(opt::at(trans_im,
453 (row + nrows / 2) % nrows,
454 (col + ncols / 2) % ncols)))
455 / log (a + b * max) * mln_max(R);
458 mln_piter(image2d< std::complex<T> >) it(trans_im.domain());
461 log(a + b * std::norm(trans_im(it)))
462 / log (a + b * max) * mln_max(R);
470 fft<
T>::transformed_image_log_magn(
double a,
double b,
bool ordered)
const
472 return transformed_image_log_magn<T>(a, b, ordered);
481 return transformed_image_log_magn<R>(1, 100, ordered);
489 return transformed_image_log_magn<T>(1, 100, ordered);
498 fftw_destroy_plan(
p);
499 fftw_destroy_plan(p_inv);
513 template <
typename D>
515 fft<T, internal::fft_real>::fft(
const image2d<D>& original_im)
517 const unsigned nrows = original_im.nrows();
518 const unsigned ncols = original_im.ncols();
520 this->in = (
T*) fftw_malloc(nrows * ncols *
sizeof(
T));
521 this->out = (std::complex<T>*) fftw_malloc(nrows
523 *
sizeof(std::complex<T>));
525 for (
unsigned row = 0; row <
nrows; ++row)
526 for (
unsigned col = 0; col <
ncols; ++col)
527 this->in[row * ncols + col] =
528 opt::at(original_im, row, col);
531 fftw_plan_dft_r2c_2d (nrows, ncols,
533 reinterpret_cast<fftw_complex*>(this->out),
536 fftw_plan_dft_c2r_2d (nrows, ncols,
537 reinterpret_cast<fftw_complex*>(this->out),
540 this->trans_im = image2d< std::complex<T> >(original_im.domain());
545 image2d< std::complex<T> >
548 fftw_execute(this->
p);
550 const unsigned nrows = this->trans_im.nrows();
551 const unsigned ncols = this->trans_im.ncols();
552 unsigned denom = nrows *
ncols;
554 for (
unsigned row = 0; row <
nrows; ++row)
555 for (
unsigned col = 0; col <= ncols / 2; ++col)
557 this->out[i] = std::complex<T> (this->out[i].real() / denom,
558 this->out[i].imag() / denom);
559 opt::at(this->trans_im, row, col) = this->out[i];
562 for (
unsigned row = 0; row <
nrows; ++row)
563 for (
unsigned col = ncols - 1; col > ncols / 2; --col)
564 opt::at(this->trans_im, row, col) =
565 opt::at(this->trans_im, nrows - row - 1, ncols - col - 1);
566 return this->trans_im;
573 fft<T, internal::fft_real>::transform_inv()
577 const unsigned nrows = this->trans_im.nrows();
578 const unsigned ncols = this->trans_im.ncols();
579 for (
unsigned row = 0; row <
nrows; ++row)
580 for (
unsigned col = 0; col <= ncols / 2; ++col)
581 this->out[row * (ncols / 2 + 1) + col] =
582 opt::at(this->trans_im, row, col);
584 fftw_execute(this->p_inv);
586 image2d<R> new_im(this->trans_im.domain());
588 for (
unsigned row = 0; row <
nrows; ++row)
589 for (
unsigned col = 0; col <
ncols; ++col)
591 opt::at(new_im, row, col) = (this->in[i] >= mln_min(R)
592 ? (this->in[i] <= mln_max(R)
604 fft<T, internal::fft_real>::transform_inv()
606 return transform_inv<T>();
616 fft<T, internal::fft_cplx>::fft(
const image2d< std::complex<T> >& original_im)
618 const unsigned nrows = original_im.nrows();
619 const unsigned ncols = original_im.ncols();
621 this->in = fftw_malloc(
sizeof(std::complex<T>) * nrows * ncols);
622 this->out = fftw_malloc(
sizeof(std::complex<T>) * nrows * ncols);
624 for (
unsigned row = 0; row <
nrows; ++row)
625 for (
unsigned col = 0; col <
ncols; ++col)
627 this->in[row * ncols + col].re = original_im(row, col).real();
628 this->in[row * ncols + col].im = original_im(row, col).imag();
631 this->
p = fftw_plan_dft_2d(nrows, ncols, this->in, this->out,
632 FFTW_FORWARD, FFTW_ESTIMATE);
633 this->p_inv = fftw_plan_dft_2d(nrows, ncols, this->out, this->in,
634 FFTW_BACKWARD, FFTW_ESTIMATE);
635 this->trans_im = image2d< std::complex<T> >(original_im.domain());
640 image2d< std::complex<T> >
643 fftw_execute(this->
p);
645 const unsigned nrows = this->trans_im.nrows();
646 const unsigned ncols = this->trans_im.ncols();
647 unsigned denom = nrows *
ncols;
649 for (
unsigned row = 0; row <
nrows; ++row)
650 for (
unsigned col = 0; col <
ncols; ++col)
652 this->out[i].re /= denom;
653 this->out[i].im /= denom;
654 this->trans_im(row, col) =
655 std::complex<double>(this->out[i].re, this->out[i].im);
658 return this->trans_im;
664 image2d< std::complex<R> >
665 fft<T, internal::fft_cplx>::transform_inv()
667 const unsigned nrows = this->trans_im.nrows();
668 const unsigned ncols = this->trans_im.ncols();
669 for (
unsigned row = 0; row <
nrows; ++row)
670 for (
unsigned col = 0; col <
ncols; ++col)
672 this->out[row * ncols + col].re =
673 this->trans_im(row, col).real();
674 this->out[row * ncols + col].im =
675 this->trans_im(row, col).imag();
678 fftw_execute(this->p_inv);
680 image2d< std::complex<R> > new_im(this->trans_im.size());
682 for (
unsigned row = 0; row <
nrows; ++row)
683 for (
unsigned col = 0; col <
ncols; ++col)
685 new_im(row, col) = this->in[i];
693 image2d< std::complex<T> >
694 fft<T, internal::fft_cplx>::transform_inv()
696 return transform_inv<T>();
701 # endif // ! MLN_INCLUDE_ONLY
705 #endif // ! MLN_TRANSFORM_FFT_HH