18 #ifndef itkFFTWCommon_h
19 #define itkFFTWCommon_h
21 #if defined(ITK_USE_FFTWF) || defined(ITK_USE_FFTWD)
22 # if defined(ITK_USE_CUFFTW)
28 # if !defined(FFTW_WISDOM_ONLY)
31 # error "FFTW 3.3.2 or later is required so that FFTW_WISDOM_ONLY is defined."
53 template <
typename TPixel>
63 #if defined(ITK_USE_FFTWF)
75 # ifdef ITK_USE_CUFFTW
87 bool canDestroyInput =
false)
89 return Plan_dft_c2r(1, &n, in, out, flags, threads, canDestroyInput);
99 bool canDestroyInput =
false)
101 auto * sizes =
new int[2];
104 PlanType plan = Plan_dft_c2r(2, sizes, in, out, flags, threads, canDestroyInput);
117 bool canDestroyInput =
false)
119 auto * sizes =
new int[3];
123 PlanType plan = Plan_dft_c2r(3, sizes, in, out, flags, threads, canDestroyInput);
135 bool canDestroyInput =
false)
137 # ifndef ITK_USE_CUFFTW
139 fftwf_plan_with_nthreads(threads);
145 unsigned roflags = flags;
146 if (!(flags & FFTW_ESTIMATE))
148 roflags = flags | FFTW_WISDOM_ONLY;
150 PlanType plan = fftwf_plan_dft_c2r(rank, n, in, out, roflags);
157 plan = fftwf_plan_dft_c2r(rank, n, in, out, flags);
163 for (
int i = 0; i < rank; i++)
168 fftwf_plan_dft_c2r(rank, n, din, out, flags);
171 plan = fftwf_plan_dft_c2r(rank, n, in, out, roflags);
173 # ifndef ITK_USE_CUFFTW
177 itkAssertOrThrowMacro(plan !=
nullptr,
"PLAN_CREATION_FAILED ");
188 bool canDestroyInput =
false)
190 return Plan_dft_r2c(1, &n, in, out, flags, threads, canDestroyInput);
200 bool canDestroyInput =
false)
202 auto * sizes =
new int[2];
205 PlanType plan = Plan_dft_r2c(2, sizes, in, out, flags, threads, canDestroyInput);
218 bool canDestroyInput =
false)
220 auto * sizes =
new int[3];
224 PlanType plan = Plan_dft_r2c(3, sizes, in, out, flags, threads, canDestroyInput);
236 bool canDestroyInput =
false)
239 # ifndef ITK_USE_CUFFTW
241 fftwf_plan_with_nthreads(threads);
247 unsigned roflags = flags;
248 if (!(flags & FFTW_ESTIMATE))
250 roflags = flags | FFTW_WISDOM_ONLY;
252 PlanType plan = fftwf_plan_dft_r2c(rank, n, in, out, roflags);
259 plan = fftwf_plan_dft_r2c(rank, n, in, out, flags);
265 for (
int i = 0; i < rank; i++)
270 fftwf_plan_dft_r2c(rank, n, din, out, flags);
273 plan = fftwf_plan_dft_r2c(rank, n, in, out, roflags);
275 # ifndef ITK_USE_CUFFTW
279 itkAssertOrThrowMacro(plan !=
nullptr,
"PLAN_CREATION_FAILED ");
290 bool canDestroyInput =
false)
292 return Plan_dft(1, &n, in, out, sign, flags, threads, canDestroyInput);
303 bool canDestroyInput =
false)
305 auto * sizes =
new int[2];
308 PlanType plan = Plan_dft(2, sizes, in, out, sign, flags, threads, canDestroyInput);
322 bool canDestroyInput =
false)
324 auto * sizes =
new int[3];
328 PlanType plan = Plan_dft(3, sizes, in, out, sign, flags, threads, canDestroyInput);
341 bool canDestroyInput =
false)
343 # ifndef ITK_USE_CUFFTW
345 fftwf_plan_with_nthreads(threads);
351 unsigned roflags = flags;
352 if (!(flags & FFTW_ESTIMATE))
354 roflags = flags | FFTW_WISDOM_ONLY;
356 PlanType plan = fftwf_plan_dft(rank, n, in, out, sign, roflags);
363 plan = fftwf_plan_dft(rank, n, in, out, sign, flags);
369 for (
int i = 0; i < rank; i++)
374 fftwf_plan_dft(rank, n, din, out, sign, flags);
377 plan = fftwf_plan_dft(rank, n, in, out, sign, roflags);
379 # ifndef ITK_USE_CUFFTW
383 itkAssertOrThrowMacro(plan !=
nullptr,
"PLAN_CREATION_FAILED ");
396 # ifndef ITK_USE_CUFFTW
399 fftwf_destroy_plan(p);
403 #endif // ITK_USE_FFTWF
406 #if defined(ITK_USE_FFTWD)
417 # ifdef ITK_USE_CUFFTW
429 bool canDestroyInput =
false)
431 return Plan_dft_c2r(1, &n, in, out, flags, threads, canDestroyInput);
441 bool canDestroyInput =
false)
443 auto * sizes =
new int[2];
446 PlanType plan = Plan_dft_c2r(2, sizes, in, out, flags, threads, canDestroyInput);
459 bool canDestroyInput =
false)
461 auto * sizes =
new int[3];
465 PlanType plan = Plan_dft_c2r(3, sizes, in, out, flags, threads, canDestroyInput);
477 bool canDestroyInput =
false)
479 # ifndef ITK_USE_CUFFTW
481 fftw_plan_with_nthreads(threads);
487 unsigned roflags = flags;
488 if (!(flags & FFTW_ESTIMATE))
490 roflags = flags | FFTW_WISDOM_ONLY;
492 PlanType plan = fftw_plan_dft_c2r(rank, n, in, out, roflags);
499 plan = fftw_plan_dft_c2r(rank, n, in, out, flags);
505 for (
int i = 0; i < rank; i++)
510 fftw_plan_dft_c2r(rank, n, din, out, flags);
513 plan = fftw_plan_dft_c2r(rank, n, in, out, roflags);
515 # ifndef ITK_USE_CUFFTW
519 itkAssertOrThrowMacro(plan !=
nullptr,
"PLAN_CREATION_FAILED ");
530 bool canDestroyInput =
false)
532 return Plan_dft_r2c(1, &n, in, out, flags, threads, canDestroyInput);
542 bool canDestroyInput =
false)
544 auto * sizes =
new int[2];
547 PlanType plan = Plan_dft_r2c(2, sizes, in, out, flags, threads, canDestroyInput);
560 bool canDestroyInput =
false)
562 auto * sizes =
new int[3];
566 PlanType plan = Plan_dft_r2c(3, sizes, in, out, flags, threads, canDestroyInput);
578 bool canDestroyInput =
false)
580 # ifndef ITK_USE_CUFFTW
582 fftw_plan_with_nthreads(threads);
588 unsigned roflags = flags;
589 if (!(flags & FFTW_ESTIMATE))
591 roflags = flags | FFTW_WISDOM_ONLY;
593 PlanType plan = fftw_plan_dft_r2c(rank, n, in, out, roflags);
600 plan = fftw_plan_dft_r2c(rank, n, in, out, flags);
606 for (
int i = 0; i < rank; i++)
611 fftw_plan_dft_r2c(rank, n, din, out, flags);
614 plan = fftw_plan_dft_r2c(rank, n, in, out, roflags);
616 # ifndef ITK_USE_CUFFTW
620 itkAssertOrThrowMacro(plan !=
nullptr,
"PLAN_CREATION_FAILED ");
631 bool canDestroyInput =
false)
633 return Plan_dft(1, &n, in, out, sign, flags, threads, canDestroyInput);
644 bool canDestroyInput =
false)
646 auto * sizes =
new int[2];
649 PlanType plan = Plan_dft(2, sizes, in, out, sign, flags, threads, canDestroyInput);
663 bool canDestroyInput =
false)
665 auto * sizes =
new int[3];
669 PlanType plan = Plan_dft(3, sizes, in, out, sign, flags, threads, canDestroyInput);
682 bool canDestroyInput =
false)
684 # ifndef ITK_USE_CUFFTW
686 fftw_plan_with_nthreads(threads);
692 unsigned roflags = flags;
693 if (!(flags & FFTW_ESTIMATE))
695 roflags = flags | FFTW_WISDOM_ONLY;
697 PlanType plan = fftw_plan_dft(rank, n, in, out, sign, roflags);
704 plan = fftw_plan_dft(rank, n, in, out, sign, flags);
710 for (
int i = 0; i < rank; i++)
715 fftw_plan_dft(rank, n, din, out, sign, flags);
718 plan = fftw_plan_dft(rank, n, in, out, sign, roflags);
720 # ifndef ITK_USE_CUFFTW
724 itkAssertOrThrowMacro(plan !=
nullptr,
"PLAN_CREATION_FAILED ");
737 # ifndef ITK_USE_CUFFTW
740 fftw_destroy_plan(p);