00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028 #ifndef OLENA_TRANSFORMS_FFT_HH
00029 # define OLENA_TRANSFORMS_FFT_HH
00030
00031
00032
00033
00034 # include <oln/config/system.hh>
00035
00036 # ifdef HAVE_FFTW
00037
00038 # include <fftw.h>
00039 # include <rfftw.h>
00040
00041 # include <ntg/basics.hh>
00042 # include <ntg/all.hh>
00043 # include <oln/basics2d.hh>
00044 # include <oln/morpher/piece_morpher.hh>
00045
00046 namespace oln {
00047
00048 namespace internal {
00049
00051 enum fft_dispatch { fft_cplx, fft_real };
00052
00056 template <class T>
00057 struct fft_trait
00058 {
00059 void ensure() { ntg_is_a(T, ntg::real)::ensure(); }
00060
00061 static const fft_dispatch which = fft_real;
00062 typedef fftw_real fftw_input;
00063 };
00064
00070 template <ntg::cplx_representation R, class T>
00071 struct fft_trait<ntg::cplx<R, T> >
00072 {
00073 static const fft_dispatch which = fft_cplx;
00074 typedef fftw_complex fftw_input;
00075 };
00076
00083 template <class T, ntg::cplx_representation R>
00084 class _fft
00085 {
00086 public:
00092 const image2d<ntg::cplx<R, ntg::float_d> > transformed_image() const
00093 {
00094 return trans_im;
00095 }
00096
00102 image2d<ntg::cplx<R, ntg::float_d> >& transformed_image()
00103 {
00104 return trans_im;
00105 }
00106
00117 template <class T1>
00118 image2d<T1> transformed_image_magn(bool ordered = true) const
00119 {
00120 ntg_is_a(T1, ntg::real)::ensure();
00121
00122 image2d<T1> new_im(trans_im.size());
00123
00124 if (ordered)
00125 for (int row = 0; row < new_im.nrows(); ++row)
00126 for (int col = 0; col < new_im.ncols(); ++col)
00127 new_im(row, col) =
00128 trans_im((row + trans_im.nrows() / 2) % trans_im.nrows(),
00129 (col + trans_im.ncols() / 2) % trans_im.ncols()).magn();
00130 else
00131 {
00132 typename image2d<ntg::cplx<R, ntg::float_d> >::iter it(trans_im);
00133
00134 for_all(it)
00135 new_im[it] = trans_im[it].magn();
00136 }
00137
00138 return new_im;
00139 }
00140
00149 image2d<T> transformed_image_magn(bool ordered = true) const
00150 {
00151 return transformed_image_magn<T>(ordered);
00152 }
00153
00165 template <class T1>
00166 image2d<T1> transformed_image_clipped_magn(const ntg::float_d clip,
00167 bool ordered = true) const
00168 {
00169 ntg_is_a(T1, ntg::real)::ensure();
00170
00171 image2d<T1> new_im(trans_im.size());
00172 ntg::range<ntg::float_d, ntg::bounded_u<0, 1>, ntg::saturate> c(clip);
00173 ntg::float_d max = 0;
00174 typename image2d<ntg::cplx<R, ntg::float_d> >::iter_type it(trans_im);
00175
00176 for_all(it)
00177 if (max < trans_im[it].magn())
00178 max = trans_im[it].magn();
00179
00180 if (ordered)
00181 for (int row = 0; row < new_im.nrows(); ++row)
00182 for (int col = 0; col < new_im.ncols(); ++col)
00183 {
00184 if (trans_im((row + trans_im.nrows() / 2) % trans_im.nrows(),
00185 (col + trans_im.ncols() / 2) %
00186 trans_im.ncols()).magn() >= max * c)
00187 new_im(row, col) = ntg_max_val(T);
00188 else
00189 new_im(row, col) =
00190 ntg_max_val(T) *
00191 trans_im((row + trans_im.nrows() / 2) % trans_im.nrows(),
00192 (col + trans_im.ncols() / 2) %
00193 trans_im.ncols()).magn() / (max * c);
00194 }
00195 else
00196 for_all(it)
00197 {
00198 if (trans_im[it].magn() >= max * c)
00199 new_im[it] = ntg_max_val(T);
00200 else
00201 new_im[it] = ntg_max_val(T) *
00202 trans_im[it].magn() / (max * c);
00203 }
00204
00205 return new_im;
00206 }
00207
00217 image2d<T> transformed_image_clipped_magn(const ntg::float_d clip,
00218 bool ordered = true) const
00219 {
00220 return transformed_image_clipped_magn<T>(clip, ordered);
00221 }
00222
00233 template <class T1>
00234 image2d<T1> transformed_image_clipped_magn(bool ordered = true) const
00235 {
00236 return transformed_image_clipped_magn<T1>(1, ordered);
00237 }
00238
00249 image2d<T> transformed_image_clipped_magn(bool ordered = true) const
00250 {
00251 return transformed_image_clipped_magn<T>(1, ordered);
00252 }
00253
00254
00265 template <class T1>
00266 image2d<T1> transformed_image_log_magn(const ntg::range<ntg::float_d,
00267 ntg::bounded_u<0, 1000>,
00268 ntg::saturate> a,
00269 const ntg::range<ntg::float_d,
00270 ntg::bounded_u<0, 1000>,
00271 ntg::saturate> b,
00272 bool ordered = true) const
00273 {
00274 ntg_is_a(T1, ntg::real)::ensure();
00275
00276 image2d<T1> new_im(trans_im.size());
00277 ntg::float_d max = 0;
00278 typename image2d<ntg::cplx<R, ntg::float_d> >::iter_type it(trans_im);
00279
00280 for_all(it)
00281 if (max < trans_im[it].magn())
00282 max = trans_im[it].magn();
00283
00284 if (ordered)
00285 for (int row = 0; row < new_im.nrows(); ++row)
00286 for (int col = 0; col < new_im.ncols(); ++col)
00287 new_im(row, col) =
00288 log(a + b * trans_im((row + trans_im.nrows() / 2) %
00289 trans_im.nrows(),
00290 (col + trans_im.ncols() / 2) %
00291 trans_im.ncols()).magn()) /
00292 log(a + b * max) * ntg_max_val(T);
00293 else
00294 for_all(it)
00295 new_im[it] = log(a + b * trans_im[it].magn()) /
00296 log(a + b * max) * ntg_max_val(T);
00297
00298 return new_im;
00299 }
00300
00301
00312 image2d<T> transformed_image_log_magn(const ntg::range<ntg::float_d,
00313 ntg::bounded_u<0, 1000>,
00314 ntg::saturate> a,
00315 const ntg::range<ntg::float_d,
00316 ntg::bounded_u<0, 1000>,
00317 ntg::saturate> b,
00318 bool ordered = true) const
00319 {
00320 return transformed_image_log_magn<T>(a, b, ordered);
00321 }
00322
00333 template <class T1>
00334 image2d<T1> transformed_image_log_magn(bool ordered = true) const
00335 {
00336 return transformed_image_log_magn<T1>(1, 100, ordered);
00337 }
00338
00347 image2d<T> transformed_image_log_magn(bool ordered = true) const
00348 {
00349 return transformed_image_log_magn<T>(1, 100, ordered);
00350 }
00351
00357 ~_fft()
00358 {
00359 delete [] in;
00360 delete [] out;
00361 fftwnd_destroy_plan(p);
00362 fftwnd_destroy_plan(p_inv);
00363 }
00364
00365 protected:
00366
00367 typename fft_trait<T>::fftw_input *in;
00368 fftw_complex *out;
00369 fftwnd_plan p;
00370 fftwnd_plan p_inv;
00371 image2d<ntg::cplx<R, ntg::float_d> > trans_im;
00372
00373 };
00374
00375 }
00376
00377 namespace transforms {
00378
00386 template <class T,
00387 ntg::cplx_representation R = ntg::polar,
00388 internal::fft_dispatch which = internal::fft_trait<T>::which >
00389 class fft;
00390
00397 template <class T, ntg::cplx_representation R>
00398 class fft<T, R, internal::fft_real> : public internal::_fft<T, R>
00399 {
00400
00401 public:
00402
00410 fft(const image2d<T>& original_im)
00411 {
00412 this->in = new fftw_real[original_im.nrows() * original_im.ncols()];
00413 this->out = new fftw_complex[original_im.nrows() *
00414 (original_im.ncols() / 2 + 1)];
00415
00416 for (int row = 0; row < original_im.nrows(); ++row)
00417 for (int col = 0; col < original_im.ncols(); ++col)
00418 this->in[row * original_im.ncols() + col] = original_im(row, col);
00419
00420 this->p = rfftw2d_create_plan(original_im.nrows(), original_im.ncols(),
00421 FFTW_REAL_TO_COMPLEX, FFTW_ESTIMATE);
00422 this->p_inv = rfftw2d_create_plan(original_im.nrows(), original_im.ncols(),
00423 FFTW_COMPLEX_TO_REAL, FFTW_ESTIMATE);
00424
00425 this->trans_im = image2d<ntg::cplx<R, ntg::float_d> >(original_im.size());
00426 }
00427
00431 image2d<ntg::cplx<R, ntg::float_d> > transform()
00432 {
00433 rfftwnd_one_real_to_complex(this->p, this->in, this->out);
00434
00435 unsigned denom = this->trans_im.nrows() * this->trans_im.ncols();
00436 int i = 0;
00437 for (int row = 0; row < this->trans_im.nrows(); ++row)
00438 for (int col = 0; col <= this->trans_im.ncols() / 2; ++col)
00439 {
00440 this->out[i].re /= denom;
00441 this->out[i].im /= denom;
00442 this->trans_im(row, col) = ntg::cplx<ntg::rect, ntg::float_d>(this->out[i].re,
00443 this->out[i].im);
00444 ++i;
00445 }
00446 for (int row = 0; row < this->trans_im.nrows(); ++row)
00447 for (int col = this->trans_im.ncols() - 1; col > this->trans_im.ncols() / 2; --col)
00448 this->trans_im(row, col) = this->trans_im(this->trans_im.nrows() - row - 1,
00449 this->trans_im.ncols() - col - 1);
00450 return this->trans_im;
00451 }
00452
00458 template <class T1>
00459 image2d<T1> transform_inv()
00460 {
00461 ntg_is_a(T1, ntg::real)::ensure();
00462
00463 for (int row = 0; row < this->trans_im.nrows(); ++row)
00464 for (int col = 0; col <= this->trans_im.ncols() / 2; ++col)
00465 {
00466 this->out[row * (this->trans_im.ncols() / 2 + 1) + col].re =
00467 this->trans_im(row, col).real();
00468 this->out[row * (this->trans_im.ncols() / 2 + 1) + col].im =
00469 this->trans_im(row, col).imag();
00470 }
00471 rfftwnd_one_complex_to_real(this->p_inv, this->out, this->in);
00472
00473 image2d<T1> new_im(this->trans_im.size());
00474 int i = 0;
00475 for (int row = 0; row < this->trans_im.nrows(); ++row)
00476 for (int col = 0; col < this->trans_im.ncols(); ++col)
00477 {
00478 new_im(row, col) = (this->in[i] >= ntg_inf_val(T1) ?
00479 (this->in[i] <= ntg_sup_val(T1) ?
00480 this->in [i] :
00481 ntg_sup_val(T1)) :
00482 ntg_inf_val(T1));
00483 ++i;
00484 }
00485 return new_im;
00486 }
00487
00510 template <class T1>
00511 image2d<T1> shift_transform_inv()
00512 {
00513 image2d<T1> t = transform_inv<T1>();
00514 image2d<T1> st(t.size());
00515
00516
00517
00518 typedef morpher::piece_morpher< image2d<T1> > piece_t;
00519 piece_t t_1(t, dpoint2d(0, 0),
00520 image2d_size((t.size().nrows() - 1) / 2,
00521 (t.size().ncols() - 1) / 2,
00522 t.border()));
00523 piece_t t_1_dest(st, dpoint2d(t.nrows() - t_1.nrows(),
00524 t.ncols() - t_1.ncols()),
00525 image2d_size(t_1.nrows(), t_1.ncols(),
00526 t.border()));
00527 piece_t t_2(t, dpoint2d(0, (t.size().ncols() - 1) / 2),
00528 image2d_size((t.size().nrows() - 1) / 2,
00529 t.size().ncols() - (t.size().ncols() - 1) / 2,
00530 t.border()));
00531 piece_t t_2_dest(st, dpoint2d(t.nrows() - t_2.nrows(), 0),
00532 image2d_size(t_2.nrows(), t_2.ncols(),
00533 t.border()));
00534 piece_t t_3(t, dpoint2d((t.size().nrows() - 1) / 2, 0),
00535 image2d_size(t.size().nrows() - (t.size().nrows() - 1) / 2,
00536 (t.size().ncols() - 1) / 2,
00537 t.border()));
00538 piece_t t_3_dest(st, dpoint2d(0, t.ncols() - t_3.ncols()),
00539 image2d_size(t_3.nrows(), t_3.ncols(),
00540 t.border()));
00541 piece_t t_4(t, dpoint2d((t.size().nrows() - 1) / 2,
00542 (t.size().ncols() - 1) / 2),
00543 image2d_size(t.size().nrows() - (t.size().nrows() - 1) / 2,
00544 t.size().ncols() - (t.size().ncols() - 1) / 2,
00545 t.border()));
00546 piece_t t_4_dest(st, dpoint2d(0, 0),
00547 image2d_size(t_4.nrows(), t_4.ncols(),
00548 t.border()));
00549
00550 oln_iter_type(piece_t) i1(t_1);
00551 for_all(i1)
00552 t_1_dest[i1] = t_1[i1];
00553 oln_iter_type(piece_t) i2(t_2);
00554 for_all(i2)
00555 t_2_dest[i2] = t_2[i2];
00556 oln_iter_type(piece_t) i3(t_3);
00557 for_all(i3)
00558 t_3_dest[i3] = t_3[i3];
00559 oln_iter_type(piece_t) i4(t_4);
00560 for_all(i4)
00561 t_4_dest[i4] = t_4[i4];
00562
00563 return st;
00564 }
00565
00569 image2d<T> transform_inv()
00570 {
00571 return transform_inv<T>();
00572 }
00573
00578 image2d<T> shift_transform_inv()
00579 {
00580 return shift_transform_inv<T>();
00581 }
00582
00583 };
00584
00591 template <class T, ntg::cplx_representation R>
00592 class fft<T, R, internal::fft_cplx> : public internal::_fft<T, R>
00593 {
00594
00595 public:
00596
00604 template <ntg::cplx_representation R1>
00605 fft(const image2d<ntg::cplx<R1, T> >& original_im)
00606 {
00607 this->in = new fftw_complex[original_im.nrows() * original_im.ncols()];
00608 this->out = new fftw_complex[original_im.nrows() * original_im.ncols()];
00609
00610 for (int row = 0; row < original_im.nrows(); ++row)
00611 for (int col = 0; col < original_im.ncols(); ++col)
00612 {
00613 this->in[row * original_im.ncols() + col].re =
00614 original_im(row, col).real();
00615 this->in[row * original_im.ncols() + col].im =
00616 original_im(row, col).imag();
00617 }
00618
00619 this->p = fftw2d_create_plan(original_im.nrows(), original_im.ncols(),
00620 FFTW_FORWARD, FFTW_ESTIMATE);
00621 this->p_inv = fftw2d_create_plan(original_im.nrows(), original_im.ncols(),
00622 FFTW_BACKWARD, FFTW_ESTIMATE);
00623
00624 this->trans_im = image2d<ntg::cplx<ntg::rect, ntg::float_d> >(original_im.size());
00625 }
00626
00630 image2d<ntg::cplx<R, ntg::float_d> > transform()
00631 {
00632 fftwnd_one(this->p, this->in, this->out);
00633
00634 unsigned denom = this->trans_im.nrows() * this->trans_im.ncols();
00635 int i = 0;
00636 for (int row = 0; row < this->trans_im.nrows(); ++row)
00637 for (int col = 0; col < this->trans_im.ncols(); ++col)
00638 {
00639 this->out[i].re /= denom;
00640 this->out[i].im /= denom;
00641 this->trans_im(row, col) = ntg::cplx<ntg::rect, ntg::float_d>(this->out[i].re,
00642 this->out[i].im);
00643 ++i;
00644 }
00645 return this->trans_im;
00646 }
00647
00653 template <class T1>
00654 image2d<ntg::cplx<R, T1> > transform_inv()
00655 {
00656 for (int row = 0; row < this->trans_im.nrows(); ++row)
00657 for (int col = 0; col < this->trans_im.ncols(); ++col)
00658 {
00659 this->out[row * this->trans_im.ncols() + col].re = this->trans_im(row, col).real();
00660 this->out[row * this->trans_im.ncols() + col].im = this->trans_im(row, col).imag();
00661 }
00662 fftwnd_one(this->p_inv, this->out, this->in);
00663
00664 image2d<ntg::cplx<R, T1> > new_im(this->trans_im.size());
00665 int i = 0;
00666 for (int row = 0; row < this->trans_im.nrows(); ++row)
00667 for (int col = 0; col < this->trans_im.ncols(); ++col)
00668 {
00669 new_im(row, col) = this->in[i];
00670 ++i;
00671 }
00672 return new_im;
00673 }
00674
00678 image2d<ntg::cplx<R, T> > transform_inv()
00679 {
00680 return transform_inv<T>();
00681 }
00682
00683 };
00684
00685 }
00686
00687 }
00688
00689 # endif // HAVE_FFTW
00690
00691 #endif // ! OLENA_TRANSFORM_FFT_HH