NFFT  3.3.0
mpolar_fft_test.c
Go to the documentation of this file.
1 /*
2  * Copyright (c) 2002, 2015 Jens Keiner, Stefan Kunis, Daniel Potts
3  *
4  * This program is free software; you can redistribute it and/or modify it under
5  * the terms of the GNU General Public License as published by the Free Software
6  * Foundation; either version 2 of the License, or (at your option) any later
7  * version.
8  *
9  * This program is distributed in the hope that it will be useful, but WITHOUT
10  * ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
11  * FOR A PARTICULAR PURPOSE. See the GNU General Public License for more
12  * details.
13  *
14  * You should have received a copy of the GNU General Public License along with
15  * this program; if not, write to the Free Software Foundation, Inc., 51
16  * Franklin Street, Fifth Floor, Boston, MA 02110-1301, USA.
17  */
18 
19 /* $Id$ */
20 
31 #include <math.h>
32 #include <stdlib.h>
33 #include <complex.h>
34 
35 #define NFFT_PRECISION_DOUBLE
36 
37 #include "nfft3mp.h"
38 
45 NFFT_R GLOBAL_elapsed_time;
46 
60 static int mpolar_grid(int T, int S, NFFT_R *x, NFFT_R *w)
61 {
62  int t, r;
63  NFFT_R W;
64  int R2 = 2 * (int)(NFFT_M(lrint)(NFFT_M(ceil)(NFFT_M(sqrt)(NFFT_K(2.0)) * (NFFT_R)(S) / NFFT_K(2.0))));
65  NFFT_R xx, yy;
66  int M = 0;
67 
68  for (t = -T / 2; t < T / 2; t++)
69  {
70  for (r = -R2 / 2; r < R2 / 2; r++)
71  {
72  xx = (NFFT_R) (r) / (NFFT_R)(S) * NFFT_M(cos)(NFFT_KPI * (NFFT_R)(t) / (NFFT_R)(T));
73  yy = (NFFT_R) (r) / (NFFT_R)(S) * NFFT_M(sin)(NFFT_KPI * (NFFT_R)(t) / (NFFT_R)(T));
74 
75  if (((-NFFT_K(0.5) - NFFT_K(1.0) / (NFFT_R) S) <= xx) & (xx <= (NFFT_K(0.5) + NFFT_K(1.0) / (NFFT_R) S))
76  & ((-NFFT_K(0.5) - NFFT_K(1.0) / (NFFT_R) S) <= yy)
77  & (yy <= (NFFT_K(0.5) + NFFT_K(1.0) / (NFFT_R) S)))
78  {
79  x[2 * M + 0] = xx;
80  x[2 * M + 1] = yy;
81 
82  if (r == 0)
83  w[M] = NFFT_K(1.0) / NFFT_K(4.0);
84  else
85  w[M] = NFFT_M(fabs)((NFFT_R) r);
86 
87  M++;
88  }
89  }
90  }
91 
93  W = NFFT_K(0.0);
94  for (t = 0; t < M; t++)
95  W += w[t];
96 
97  for (t = 0; t < M; t++)
98  w[t] /= W;
99 
100  return M;
101 }
102 
104 static int mpolar_dft(NFFT_C *f_hat, int NN, NFFT_C *f, int T, int S, int m)
105 {
106  double t0, t1;
107  int j, k;
108  NFFT(plan) my_nfft_plan;
110  NFFT_R *x, *w;
112  int N[2], n[2];
113  int M;
115  N[0] = NN;
116  n[0] = 2 * N[0];
117  N[1] = NN;
118  n[1] = 2 * N[1];
120  x = (NFFT_R *) NFFT(malloc)((size_t)(5 * (T / 2) * S) * (sizeof(NFFT_R)));
121  if (x == NULL)
122  return EXIT_FAILURE;
123 
124  w = (NFFT_R *) NFFT(malloc)((size_t)(5 * (T * S) / 4) * (sizeof(NFFT_R)));
125  if (w == NULL)
126  return EXIT_FAILURE;
127 
129  M = mpolar_grid(T, S, x, w);
130  NFFT(init_guru)(&my_nfft_plan, 2, N, M, n, m,
131  PRE_PHI_HUT | PRE_PSI | MALLOC_X | MALLOC_F_HAT | MALLOC_F | FFTW_INIT
132  | FFT_OUT_OF_PLACE,
133  FFTW_MEASURE | FFTW_DESTROY_INPUT);
134 
136  for (j = 0; j < my_nfft_plan.M_total; j++)
137  {
138  my_nfft_plan.x[2 * j + 0] = x[2 * j + 0];
139  my_nfft_plan.x[2 * j + 1] = x[2 * j + 1];
140  }
141 
143  for (k = 0; k < my_nfft_plan.N_total; k++)
144  my_nfft_plan.f_hat[k] = f_hat[k];
145 
146  t0 = NFFT(clock_gettime_seconds)();
147 
149  NFFT(trafo_direct)(&my_nfft_plan);
150 
151  t1 = NFFT(clock_gettime_seconds)();
152  GLOBAL_elapsed_time = (t1 - t0);
153 
155  for (j = 0; j < my_nfft_plan.M_total; j++)
156  f[j] = my_nfft_plan.f[j];
157 
159  NFFT(finalize)(&my_nfft_plan);
160  NFFT(free)(x);
161  NFFT(free)(w);
162 
163  return EXIT_SUCCESS;
164 }
165 
167 static int mpolar_fft(NFFT_C *f_hat, int NN, NFFT_C *f, int T, int S, int m)
168 {
169  double t0, t1;
170  int j, k;
171  NFFT(plan) my_nfft_plan;
173  NFFT_R *x, *w;
175  int N[2], n[2];
176  int M;
178  N[0] = NN;
179  n[0] = 2 * N[0];
180  N[1] = NN;
181  n[1] = 2 * N[1];
183  x = (NFFT_R *) NFFT(malloc)((size_t)(5 * T * S / 2) * (sizeof(NFFT_R)));
184  if (x == NULL)
185  return EXIT_FAILURE;
186 
187  w = (NFFT_R *) NFFT(malloc)((size_t)(5 * T * S / 4) * (sizeof(NFFT_R)));
188  if (w == NULL)
189  return EXIT_FAILURE;
190 
192  M = mpolar_grid(T, S, x, w);
193  NFFT(init_guru)(&my_nfft_plan, 2, N, M, n, m,
194  PRE_PHI_HUT | PRE_PSI | MALLOC_X | MALLOC_F_HAT | MALLOC_F | FFTW_INIT
195  | FFT_OUT_OF_PLACE,
196  FFTW_MEASURE | FFTW_DESTROY_INPUT);
197 
200  for (j = 0; j < my_nfft_plan.M_total; j++)
201  {
202  my_nfft_plan.x[2 * j + 0] = x[2 * j + 0];
203  my_nfft_plan.x[2 * j + 1] = x[2 * j + 1];
204  }
205 
207  if (my_nfft_plan.flags & PRE_LIN_PSI)
208  NFFT(precompute_lin_psi)(&my_nfft_plan);
209 
210  if (my_nfft_plan.flags & PRE_PSI)
211  NFFT(precompute_psi)(&my_nfft_plan);
212 
213  if (my_nfft_plan.flags & PRE_FULL_PSI)
214  NFFT(precompute_full_psi)(&my_nfft_plan);
215 
217  for (k = 0; k < my_nfft_plan.N_total; k++)
218  my_nfft_plan.f_hat[k] = f_hat[k];
219 
220  t0 = NFFT(clock_gettime_seconds)();
221 
223  NFFT(trafo)(&my_nfft_plan);
224 
225  t1 = NFFT(clock_gettime_seconds)();
226  GLOBAL_elapsed_time = (t1 - t0);
227 
229  for (j = 0; j < my_nfft_plan.M_total; j++)
230  f[j] = my_nfft_plan.f[j];
231 
233  NFFT(finalize)(&my_nfft_plan);
234  NFFT(free)(x);
235  NFFT(free)(w);
236 
237  return EXIT_SUCCESS;
238 }
239 
241 static int inverse_mpolar_fft(NFFT_C *f, int T, int S, NFFT_C *f_hat, int NN, int max_i,
242  int m)
243 {
244  double t0, t1;
245  int j, k;
246  NFFT(plan) my_nfft_plan;
247  SOLVER(plan_complex) my_infft_plan;
249  NFFT_R *x, *w;
250  int l;
252  int N[2], n[2];
253  int M;
255  N[0] = NN;
256  n[0] = 2 * N[0];
257  N[1] = NN;
258  n[1] = 2 * N[1];
260  x = (NFFT_R *) NFFT(malloc)((size_t)(5 * T * S / 2) * (sizeof(NFFT_R)));
261  if (x == NULL)
262  return EXIT_FAILURE;
263 
264  w = (NFFT_R *) NFFT(malloc)((size_t)(5 * T * S / 4) * (sizeof(NFFT_R)));
265  if (w == NULL)
266  return EXIT_FAILURE;
267 
269  M = mpolar_grid(T, S, x, w);
270  NFFT(init_guru)(&my_nfft_plan, 2, N, M, n, m,
271  PRE_PHI_HUT | PRE_PSI | MALLOC_X | MALLOC_F_HAT | MALLOC_F | FFTW_INIT
272  | FFT_OUT_OF_PLACE,
273  FFTW_MEASURE | FFTW_DESTROY_INPUT);
274 
276  SOLVER(init_advanced_complex)(&my_infft_plan,
277  (NFFT(mv_plan_complex)*) (&my_nfft_plan), CGNR | PRECOMPUTE_WEIGHT);
278 
280  for (j = 0; j < my_nfft_plan.M_total; j++)
281  {
282  my_nfft_plan.x[2 * j + 0] = x[2 * j + 0];
283  my_nfft_plan.x[2 * j + 1] = x[2 * j + 1];
284  my_infft_plan.y[j] = f[j];
285  my_infft_plan.w[j] = w[j];
286  }
287 
289  if (my_nfft_plan.flags & PRE_LIN_PSI)
290  NFFT(precompute_lin_psi)(&my_nfft_plan);
291 
292  if (my_nfft_plan.flags & PRE_PSI)
293  NFFT(precompute_psi)(&my_nfft_plan);
294 
295  if (my_nfft_plan.flags & PRE_FULL_PSI)
296  NFFT(precompute_full_psi)(&my_nfft_plan);
297 
299  if (my_infft_plan.flags & PRECOMPUTE_DAMP)
300  for (j = 0; j < my_nfft_plan.N[0]; j++)
301  for (k = 0; k < my_nfft_plan.N[1]; k++)
302  {
303  my_infft_plan.w_hat[j * my_nfft_plan.N[1] + k] = (
304  NFFT_M(sqrt)(
305  NFFT_M(pow)((NFFT_R)(j - my_nfft_plan.N[0] / 2), NFFT_K(2.0))
306  + NFFT_M(pow)((NFFT_R)(k - my_nfft_plan.N[1] / 2), NFFT_K(2.0)))
307  > (NFFT_R)(my_nfft_plan.N[0] / 2) ? NFFT_K(0.0) : NFFT_K(1.0));
308  }
309 
311  for (k = 0; k < my_nfft_plan.N_total; k++)
312  my_infft_plan.f_hat_iter[k] = NFFT_K(0.0) + _Complex_I * NFFT_K(0.0);
313 
314  t0 = NFFT(clock_gettime_seconds)();
315 
317  SOLVER(before_loop_complex)(&my_infft_plan);
318 
319  if (max_i < 1)
320  {
321  l = 1;
322  for (k = 0; k < my_nfft_plan.N_total; k++)
323  my_infft_plan.f_hat_iter[k] = my_infft_plan.p_hat_iter[k];
324  }
325  else
326  {
327  for (l = 1; l <= max_i; l++)
328  {
329  SOLVER(loop_one_step_complex)(&my_infft_plan);
330  }
331  }
332 
333  t1 = NFFT(clock_gettime_seconds)();
334  GLOBAL_elapsed_time = (t1 - t0);
335 
337  for (k = 0; k < my_nfft_plan.N_total; k++)
338  f_hat[k] = my_infft_plan.f_hat_iter[k];
339 
341  SOLVER(finalize_complex)(&my_infft_plan);
342  NFFT(finalize)(&my_nfft_plan);
343  NFFT(free)(x);
344  NFFT(free)(w);
345 
346  return EXIT_SUCCESS;
347 }
348 
350 static int comparison_fft(FILE *fp, int N, int T, int S)
351 {
352  double t0, t1;
353  FFTW(plan) my_fftw_plan;
354  NFFT_C *f_hat, *f;
355  int m, k;
356  NFFT_R t_fft, t_dft_mpolar;
357 
358  f_hat = (NFFT_C *) NFFT(malloc)(sizeof(NFFT_C) * (size_t)(N * N));
359  f = (NFFT_C *) NFFT(malloc)(sizeof(NFFT_C) * (size_t)((T * S / 4) * 5));
360 
361  my_fftw_plan = FFTW(plan_dft_2d)(N, N, f_hat, f, FFTW_BACKWARD, FFTW_MEASURE);
362 
363  for (k = 0; k < N * N; k++)
364  f_hat[k] = NFFT(drand48)() + _Complex_I * NFFT(drand48)();
365 
366  t0 = NFFT(clock_gettime_seconds)();
367  for (m = 0; m < 65536 / N; m++)
368  {
369  FFTW(execute)(my_fftw_plan);
370  /* touch */
371  f_hat[2] = NFFT_K(2.0) * f_hat[0];
372  }
373  t1 = NFFT(clock_gettime_seconds)();
374  GLOBAL_elapsed_time = (t1 - t0);
375  t_fft = (NFFT_R)(N) * GLOBAL_elapsed_time / NFFT_K(65536.0);
376 
377  if (N < 256)
378  {
379  mpolar_dft(f_hat, N, f, T, S, 1);
380  t_dft_mpolar = GLOBAL_elapsed_time;
381  }
382 
383  for (m = 3; m <= 9; m += 3)
384  {
385  if ((m == 3) && (N < 256))
386  fprintf(fp, "%d\t&\t&\t%1.1" NFFT__FES__ "&\t%1.1" NFFT__FES__ "&\t%d\t", N, t_fft, t_dft_mpolar, m);
387  else if (m == 3)
388  fprintf(fp, "%d\t&\t&\t%1.1" NFFT__FES__ "&\t &\t%d\t", N, t_fft, m);
389  else
390  fprintf(fp, " \t&\t&\t &\t &\t%d\t", m);
391 
392  printf("N=%d\tt_fft=%1.1" NFFT__FES__ "\tt_dft_mpolar=%1.1" NFFT__FES__ "\tm=%d\t", N, t_fft,
393  t_dft_mpolar, m);
394 
395  mpolar_fft(f_hat, N, f, T, S, m);
396  fprintf(fp, "%1.1" NFFT__FES__ "&\t", GLOBAL_elapsed_time);
397  printf("t_mpolar=%1.1" NFFT__FES__ "\t", GLOBAL_elapsed_time);
398  inverse_mpolar_fft(f, T, S, f_hat, N, 2 * m, m);
399  if (m == 9)
400  fprintf(fp, "%1.1" NFFT__FES__ "\\\\\\hline\n", GLOBAL_elapsed_time);
401  else
402  fprintf(fp, "%1.1" NFFT__FES__ "\\\\\n", GLOBAL_elapsed_time);
403  printf("t_impolar=%1.1" NFFT__FES__ "\n", GLOBAL_elapsed_time);
404  }
405 
406  fflush(fp);
407 
408  NFFT(free)(f);
409  NFFT(free)(f_hat);
410 
411  return EXIT_SUCCESS;
412 }
413 
415 int main(int argc, char **argv)
416 {
417  int N;
418  int T, S;
419  int M;
420  NFFT_R *x, *w;
421  NFFT_C *f_hat, *f, *f_direct, *f_tilde;
422  int k;
423  int max_i;
424  int m;
425  NFFT_R temp1, temp2, E_max = NFFT_K(0.0);
426  FILE *fp1, *fp2;
427  char filename[30];
428  int logN;
429 
430  if (argc != 4)
431  {
432  printf("mpolar_fft_test N T R \n");
433  printf("\n");
434  printf("N mpolar FFT of size NxN \n");
435  printf("T number of slopes \n");
436  printf("R number of offsets \n");
437 
439  printf("\nHence, comparison FFTW, mpolar FFT and inverse mpolar FFT\n");
440  fp1 = fopen("mpolar_comparison_fft.dat", "w");
441  if (fp1 == NULL)
442  return (-1);
443  for (logN = 4; logN <= 8; logN++)
444  comparison_fft(fp1, (int)(1U << logN), 3 * (int)(1U << logN),
445  3 * (int)(1U << (logN - 1)));
446  fclose(fp1);
447 
448  exit(EXIT_FAILURE);
449  }
450 
451  N = atoi(argv[1]);
452  T = atoi(argv[2]);
453  S = atoi(argv[3]);
454  printf("N=%d, modified polar grid with T=%d, R=%d => ", N, T, S);
455 
456  x = (NFFT_R *) NFFT(malloc)((size_t)(5 * T * S / 2) * (sizeof(NFFT_R)));
457  w = (NFFT_R *) NFFT(malloc)((size_t)(5 * T * S / 4) * (sizeof(NFFT_R)));
458 
459  f_hat = (NFFT_C *) NFFT(malloc)(sizeof(NFFT_C) * (size_t)(N * N));
460  f = (NFFT_C *) NFFT(malloc)(sizeof(NFFT_C) * (size_t)(1.25 * T * S)); /* 4/pi*log(1+sqrt(2)) = 1.122... < 1.25 */
461  f_direct = (NFFT_C *) NFFT(malloc)(sizeof(NFFT_C) * (size_t)(1.25 * T * S));
462  f_tilde = (NFFT_C *) NFFT(malloc)(sizeof(NFFT_C) * (size_t)(N * N));
463 
465  M = mpolar_grid(T, S, x, w);
466  printf("M=%d.\n", M);
467 
469  fp1 = fopen("input_data_r.dat", "r");
470  fp2 = fopen("input_data_i.dat", "r");
471  if ((fp1 == NULL) || (fp2 == NULL))
472  return (-1);
473  for (k = 0; k < N * N; k++)
474  {
475  fscanf(fp1, NFFT__FR__ " ", &temp1);
476  fscanf(fp2, NFFT__FR__ " ", &temp2);
477  f_hat[k] = temp1 + _Complex_I * temp2;
478  }
479  fclose(fp1);
480  fclose(fp2);
481 
483  mpolar_dft(f_hat, N, f_direct, T, S, 1);
484  // mpolar_fft(f_hat,N,f_direct,T,R,12);
485 
487  printf("\nTest of the mpolar FFT: \n");
488  fp1 = fopen("mpolar_fft_error.dat", "w+");
489  for (m = 1; m <= 12; m++)
490  {
492  mpolar_fft(f_hat, N, f, T, S, m);
493 
495  E_max = NFFT(error_l_infty_complex)(f_direct, f, M);
496  printf("m=%2d: E_max = %" NFFT__FES__ "\n", m, E_max);
497  fprintf(fp1, "%" NFFT__FES__ "\n", E_max);
498  }
499  fclose(fp1);
500 
502  for (m = 3; m <= 9; m += 3)
503  {
504  printf("\nTest of the inverse mpolar FFT for m=%d: \n", m);
505  sprintf(filename, "mpolar_ifft_error%d.dat", m);
506  fp1 = fopen(filename, "w+");
507  for (max_i = 0; max_i <= 20; max_i += 2)
508  {
510  inverse_mpolar_fft(f_direct, T, S, f_tilde, N, max_i, m);
511 
513  E_max = NFFT(error_l_infty_complex)(f_hat, f_tilde, N * N);
514  printf("%3d iterations: E_max = %" NFFT__FES__ "\n", max_i, E_max);
515  fprintf(fp1, "%" NFFT__FES__ "\n", E_max);
516  }
517  fclose(fp1);
518  }
519 
521  NFFT(free)(x);
522  NFFT(free)(w);
523  NFFT(free)(f_hat);
524  NFFT(free)(f);
525  NFFT(free)(f_direct);
526  NFFT(free)(f_tilde);
527 
528  return EXIT_SUCCESS;
529 }
530 /* \} */
static int max_i(int a, int b)
max
Definition: fastsum.c:48
static int inverse_mpolar_fft(NFFT_C *f, int T, int S, NFFT_C *f_hat, int NN, int max_i, int m)
inverse NFFT-based mpolar FFT
static int mpolar_fft(NFFT_C *f_hat, int NN, NFFT_C *f, int T, int S, int m)
NFFT-based mpolar FFT.
static int mpolar_dft(NFFT_C *f_hat, int NN, NFFT_C *f, int T, int S, int m)
discrete mpolar FFT
static int mpolar_grid(int T, int S, NFFT_R *x, NFFT_R *w)
Generates the points with weights for the modified polar grid with angles and offsets...
static int comparison_fft(FILE *fp, int N, int T, int S)
Comparison of the FFTW, mpolar FFT, and inverse mpolar FFT.