NFFT  3.4.1
linogram_fft_test.c
Go to the documentation of this file.
1 /*
2  * Copyright (c) 2002, 2017 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 
28 #include <math.h>
29 #include <stdlib.h>
30 #include <complex.h>
31 
32 #define NFFT_PRECISION_DOUBLE
33 
34 #include "nfft3mp.h"
35 
42 NFFT_R GLOBAL_elapsed_time;
43 
47 static int linogram_grid(int T, int rr, NFFT_R *x, NFFT_R *w)
48 {
49  int t, r;
50  NFFT_R W = (NFFT_R) T * (((NFFT_R) rr / NFFT_K(2.0)) * ((NFFT_R) rr / NFFT_K(2.0)) + NFFT_K(1.0) / NFFT_K(4.0));
51 
52  for (t = -T / 2; t < T / 2; t++)
53  {
54  for (r = -rr / 2; r < rr / 2; r++)
55  {
56  if (t < 0)
57  {
58  x[2 * ((t + T / 2) * rr + (r + rr / 2)) + 0] = (NFFT_R) (r) / (NFFT_R)(rr);
59  x[2 * ((t + T / 2) * rr + (r + rr / 2)) + 1] = NFFT_K(4.0) * ((NFFT_R)(t) + (NFFT_R)(T) / NFFT_K(4.0)) / (NFFT_R)(T)
60  * (NFFT_R)(r) / (NFFT_R)(rr);
61  }
62  else
63  {
64  x[2 * ((t + T / 2) * rr + (r + rr / 2)) + 0] = -NFFT_K(4.0) * ((NFFT_R)(t) - (NFFT_R)(T) / NFFT_K(4.0)) / (NFFT_R)(T)
65  * (NFFT_R)(r) / (NFFT_R)(rr);
66  x[2 * ((t + T / 2) * rr + (r + rr / 2)) + 1] = (NFFT_R) r / (NFFT_R)(rr);
67  }
68  if (r == 0)
69  w[(t + T / 2) * rr + (r + rr / 2)] = NFFT_K(1.0) / NFFT_K(4.0) / W;
70  else
71  w[(t + T / 2) * rr + (r + rr / 2)] = NFFT_M(fabs)((NFFT_R) r) / W;
72  }
73  }
74 
75  return T * rr;
76 }
77 
79 static int linogram_dft(NFFT_C *f_hat, int NN, NFFT_C *f, int T, int rr, int m)
80 {
81  double t0, t1;
82  int j, k;
83  NFFT(plan) my_nfft_plan;
85  NFFT_R *x, *w;
87  int N[2], n[2];
88  int M = T * rr;
90  N[0] = NN;
91  n[0] = 2 * N[0];
92  N[1] = NN;
93  n[1] = 2 * N[1];
95  x = (NFFT_R *) NFFT(malloc)((size_t)(2 * T * rr) * (sizeof(NFFT_R)));
96  if (x == NULL)
97  return EXIT_FAILURE;
98 
99  w = (NFFT_R *) NFFT(malloc)((size_t)(T * rr) * (sizeof(NFFT_R)));
100  if (w == NULL)
101  return EXIT_FAILURE;
102 
104  NFFT(init_guru)(&my_nfft_plan, 2, N, M, n, m,
107  FFTW_MEASURE | FFTW_DESTROY_INPUT);
108 
110  linogram_grid(T, rr, x, w);
111  for (j = 0; j < my_nfft_plan.M_total; j++)
112  {
113  my_nfft_plan.x[2 * j + 0] = x[2 * j + 0];
114  my_nfft_plan.x[2 * j + 1] = x[2 * j + 1];
115  }
116 
118  for (k = 0; k < my_nfft_plan.N_total; k++)
119  my_nfft_plan.f_hat[k] = f_hat[k];
120 
122  t0 = NFFT(clock_gettime_seconds)();
123  NFFT(trafo_direct)(&my_nfft_plan);
124  t1 = NFFT(clock_gettime_seconds)();
125  GLOBAL_elapsed_time = (t1 - t0);
126 
128  for (j = 0; j < my_nfft_plan.M_total; j++)
129  f[j] = my_nfft_plan.f[j];
130 
132  NFFT(finalize)(&my_nfft_plan);
133  NFFT(free)(x);
134  NFFT(free)(w);
135 
136  return EXIT_SUCCESS;
137 }
138 
140 static int linogram_fft(NFFT_C *f_hat, int NN, NFFT_C *f, int T, int rr, int m)
141 {
142  double t0, t1;
143  int j, k;
144  NFFT(plan) my_nfft_plan;
146  NFFT_R *x, *w;
148  int N[2], n[2];
149  int M = T * rr;
151  N[0] = NN;
152  n[0] = 2 * N[0];
153  N[1] = NN;
154  n[1] = 2 * N[1];
156  x = (NFFT_R *) NFFT(malloc)((size_t)(2 * T * rr) * (sizeof(NFFT_R)));
157  if (x == NULL)
158  return EXIT_FAILURE;
159 
160  w = (NFFT_R *) NFFT(malloc)((size_t)(T * rr) * (sizeof(NFFT_R)));
161  if (w == NULL)
162  return EXIT_FAILURE;
163 
165  NFFT(init_guru)(&my_nfft_plan, 2, N, M, n, m,
168  FFTW_MEASURE | FFTW_DESTROY_INPUT);
169 
171  linogram_grid(T, rr, x, w);
172  for (j = 0; j < my_nfft_plan.M_total; j++)
173  {
174  my_nfft_plan.x[2 * j + 0] = x[2 * j + 0];
175  my_nfft_plan.x[2 * j + 1] = x[2 * j + 1];
176  }
177 
179  if (my_nfft_plan.flags & PRE_LIN_PSI)
180  NFFT(precompute_lin_psi)(&my_nfft_plan);
181 
182  if (my_nfft_plan.flags & PRE_PSI)
183  NFFT(precompute_psi)(&my_nfft_plan);
184 
185  if (my_nfft_plan.flags & PRE_FULL_PSI)
186  NFFT(precompute_full_psi)(&my_nfft_plan);
187 
189  for (k = 0; k < my_nfft_plan.N_total; k++)
190  my_nfft_plan.f_hat[k] = f_hat[k];
191 
193  t0 = NFFT(clock_gettime_seconds)();
194  NFFT(trafo)(&my_nfft_plan);
195  t1 = NFFT(clock_gettime_seconds)();
196  GLOBAL_elapsed_time = (t1 - t0);
197 
199  for (j = 0; j < my_nfft_plan.M_total; j++)
200  f[j] = my_nfft_plan.f[j];
201 
203  NFFT(finalize)(&my_nfft_plan);
204  NFFT(free)(x);
205  NFFT(free)(w);
206 
207  return EXIT_SUCCESS;
208 }
209 
211 static int inverse_linogram_fft(NFFT_C *f, int T, int rr, NFFT_C *f_hat, int NN,
212  int max_i, int m)
213 {
214  double t0, t1;
215  int j, k;
216  NFFT(plan) my_nfft_plan;
217  SOLVER(plan_complex) my_infft_plan;
219  NFFT_R *x, *w;
220  int l;
222  int N[2], n[2];
223  int M = T * rr;
225  N[0] = NN;
226  n[0] = 2 * N[0];
227  N[1] = NN;
228  n[1] = 2 * N[1];
230  x = (NFFT_R *) NFFT(malloc)((size_t)(2 * T * rr) * (sizeof(NFFT_R)));
231  if (x == NULL)
232  return EXIT_FAILURE;
233 
234  w = (NFFT_R *) NFFT(malloc)((size_t)(T * rr) * (sizeof(NFFT_R)));
235  if (w == NULL)
236  return EXIT_FAILURE;
237 
239  NFFT(init_guru)(&my_nfft_plan, 2, N, M, n, m,
242  FFTW_MEASURE | FFTW_DESTROY_INPUT);
243 
245  SOLVER(init_advanced_complex)(&my_infft_plan,
246  (NFFT(mv_plan_complex)*) (&my_nfft_plan), CGNR | PRECOMPUTE_WEIGHT);
247 
249  linogram_grid(T, rr, x, w);
250  for (j = 0; j < my_nfft_plan.M_total; j++)
251  {
252  my_nfft_plan.x[2 * j + 0] = x[2 * j + 0];
253  my_nfft_plan.x[2 * j + 1] = x[2 * j + 1];
254  my_infft_plan.y[j] = f[j];
255  my_infft_plan.w[j] = w[j];
256  }
257 
259  if (my_nfft_plan.flags & PRE_LIN_PSI)
260  NFFT(precompute_lin_psi)(&my_nfft_plan);
261 
262  if (my_nfft_plan.flags & PRE_PSI)
263  NFFT(precompute_psi)(&my_nfft_plan);
264 
265  if (my_nfft_plan.flags & PRE_FULL_PSI)
266  NFFT(precompute_full_psi)(&my_nfft_plan);
267 
269  if (my_infft_plan.flags & PRECOMPUTE_DAMP)
270  for (j = 0; j < my_nfft_plan.N[0]; j++)
271  for (k = 0; k < my_nfft_plan.N[1]; k++)
272  {
273  my_infft_plan.w_hat[j * my_nfft_plan.N[1] + k] = (
274  NFFT_M(sqrt)(
275  NFFT_M(pow)((NFFT_R)(j - my_nfft_plan.N[0] / 2), NFFT_K(2.0))
276  + NFFT_M(pow)((NFFT_R)(k - my_nfft_plan.N[1] / 2), NFFT_K(2.0)))
277  > (NFFT_R)(my_nfft_plan.N[0] / 2) ? 0 : 1);
278  }
279 
281  for (k = 0; k < my_nfft_plan.N_total; k++)
282  my_infft_plan.f_hat_iter[k] = NFFT_K(0.0) + _Complex_I * NFFT_K(0.0);
283 
284  t0 = NFFT(clock_gettime_seconds)();
286  SOLVER(before_loop_complex)(&my_infft_plan);
287 
288  if (max_i < 1)
289  {
290  l = 1;
291  for (k = 0; k < my_nfft_plan.N_total; k++)
292  my_infft_plan.f_hat_iter[k] = my_infft_plan.p_hat_iter[k];
293  }
294  else
295  {
296  for (l = 1; l <= max_i; l++)
297  {
298  SOLVER(loop_one_step_complex)(&my_infft_plan);
299  }
300  }
301 
302  t1 = NFFT(clock_gettime_seconds)();
303  GLOBAL_elapsed_time = (t1 - t0);
304 
306  for (k = 0; k < my_nfft_plan.N_total; k++)
307  f_hat[k] = my_infft_plan.f_hat_iter[k];
308 
310  SOLVER(finalize_complex)(&my_infft_plan);
311  NFFT(finalize)(&my_nfft_plan);
312  NFFT(free)(x);
313  NFFT(free)(w);
314 
315  return EXIT_SUCCESS;
316 }
317 
319 static int comparison_fft(FILE *fp, int N, int T, int rr)
320 {
321  double t0, t1;
322  FFTW(plan) my_fftw_plan;
323  NFFT_C *f_hat, *f;
324  int m, k;
325  NFFT_R t_fft, t_dft_linogram;
326 
327  f_hat = (NFFT_C *) NFFT(malloc)(sizeof(NFFT_C) * (size_t)(N * N));
328  f = (NFFT_C *) NFFT(malloc)(sizeof(NFFT_C) * (size_t)((T * rr / 4) * 5));
329 
330  my_fftw_plan = FFTW(plan_dft_2d)(N, N, f_hat, f, FFTW_BACKWARD, FFTW_MEASURE);
331 
332  for (k = 0; k < N * N; k++)
333  f_hat[k] = NFFT(drand48)() + _Complex_I * NFFT(drand48)();
334 
335  t0 = NFFT(clock_gettime_seconds)();
336  for (m = 0; m < 65536 / N; m++)
337  {
338  FFTW(execute)(my_fftw_plan);
339  /* touch */
340  f_hat[2] = NFFT_K(2.0) * f_hat[0];
341  }
342  t1 = NFFT(clock_gettime_seconds)();
343  GLOBAL_elapsed_time = (t1 - t0);
344  t_fft = (NFFT_R)(N) * GLOBAL_elapsed_time / NFFT_K(65536.0);
345 
346  if (N < 256)
347  {
348  linogram_dft(f_hat, N, f, T, rr, m);
349  t_dft_linogram = GLOBAL_elapsed_time;
350  }
351 
352  for (m = 3; m <= 9; m += 3)
353  {
354  if ((m == 3) && (N < 256))
355  fprintf(fp, "%d\t&\t&\t%1.1" NFFT__FES__ "&\t%1.1" NFFT__FES__ "&\t%d\t", N, t_fft, t_dft_linogram,
356  m);
357  else if (m == 3)
358  fprintf(fp, "%d\t&\t&\t%1.1" NFFT__FES__ "&\t &\t%d\t", N, t_fft, m);
359  else
360  fprintf(fp, " \t&\t&\t &\t &\t%d\t", m);
361 
362  printf("N=%d\tt_fft=%1.1" NFFT__FES__ "\tt_dft_linogram=%1.1" NFFT__FES__ "\tm=%d\t", N, t_fft,
363  t_dft_linogram, m);
364 
365  linogram_fft(f_hat, N, f, T, rr, m);
366  fprintf(fp, "%1.1" NFFT__FES__ "&\t", GLOBAL_elapsed_time);
367  printf("t_linogram=%1.1" NFFT__FES__ "\t", GLOBAL_elapsed_time);
368  inverse_linogram_fft(f, T, rr, f_hat, N, m + 3, m);
369  if (m == 9)
370  fprintf(fp, "%1.1" NFFT__FES__ "\\\\\\hline\n", GLOBAL_elapsed_time);
371  else
372  fprintf(fp, "%1.1" NFFT__FES__ "\\\\\n", GLOBAL_elapsed_time);
373  printf("t_ilinogram=%1.1" NFFT__FES__ "\n", GLOBAL_elapsed_time);
374  }
375 
376  fflush(fp);
377 
378  NFFT(free)(f);
379  NFFT(free)(f_hat);
380 
381  return EXIT_SUCCESS;
382 }
383 
385 int main(int argc, char **argv)
386 {
387  int N;
388  int T, rr;
389  int M;
390  NFFT_R *x, *w;
391  NFFT_C *f_hat, *f, *f_direct, *f_tilde;
392  int k;
393  int max_i;
394  int m;
395  NFFT_R temp1, temp2, E_max = NFFT_K(0.0);
396  FILE *fp1, *fp2;
397  char filename[30];
398  int logN;
399 
400  if (argc != 4)
401  {
402  printf("linogram_fft_test N T R \n");
403  printf("\n");
404  printf("N linogram FFT of size NxN \n");
405  printf("T number of slopes \n");
406  printf("R number of offsets \n");
407 
409  printf("\nHence, comparison FFTW, linogram FFT and inverse linogram FFT\n");
410  fp1 = fopen("linogram_comparison_fft.dat", "w");
411  if (fp1 == NULL)
412  return (-1);
413  for (logN = 4; logN <= 8; logN++)
414  comparison_fft(fp1, (int)(1U << logN), 3 * (int)(1U << logN),
415  3 * (int)(1U << (logN - 1)));
416  fclose(fp1);
417 
418  return EXIT_FAILURE;
419  }
420 
421  N = atoi(argv[1]);
422  T = atoi(argv[2]);
423  rr = atoi(argv[3]);
424  printf("N=%d, linogram grid with T=%d, R=%d => ", N, T, rr);
425 
426  x = (NFFT_R *) NFFT(malloc)((size_t)(5 * T * rr / 2) * (sizeof(NFFT_R)));
427  w = (NFFT_R *) NFFT(malloc)((size_t)(5 * T * rr / 4) * (sizeof(NFFT_R)));
428 
429  f_hat = (NFFT_C *) NFFT(malloc)(sizeof(NFFT_C) * (size_t)(N * N));
430  f = (NFFT_C *) NFFT(malloc)(sizeof(NFFT_C) * (size_t)(5 * T * rr / 4)); /* 4/pi*log(1+sqrt(2)) = 1.122... < 1.25 */
431  f_direct = (NFFT_C *) NFFT(malloc)(sizeof(NFFT_C) * (size_t)(5 * T * rr / 4));
432  f_tilde = (NFFT_C *) NFFT(malloc)(sizeof(NFFT_C) * (size_t)(N * N));
433 
435  M = linogram_grid(T, rr, x, w);
436  printf("M=%d.\n", M);
437 
439  fp1 = fopen("input_data_r.dat", "r");
440  fp2 = fopen("input_data_i.dat", "r");
441  if ((fp1 == NULL) || (fp2 == NULL))
442  return EXIT_FAILURE;
443  for (k = 0; k < N * N; k++)
444  {
445  fscanf(fp1, NFFT__FR__ " ", &temp1);
446  fscanf(fp2, NFFT__FR__ " ", &temp2);
447  f_hat[k] = temp1 + _Complex_I * temp2;
448  }
449  fclose(fp1);
450  fclose(fp2);
451 
453  linogram_dft(f_hat, N, f_direct, T, rr, 1);
454  // linogram_fft(f_hat,N,f_direct,T,R,12);
455 
457  printf("\nTest of the linogram FFT: \n");
458  fp1 = fopen("linogram_fft_error.dat", "w+");
459  for (m = 1; m <= 12; m++)
460  {
462  linogram_fft(f_hat, N, f, T, rr, m);
463 
465  E_max = NFFT(error_l_infty_complex)(f_direct, f, M);
466  //E_max=NFFT(error_l_2_complex)(f_direct,f,M);
467 
468  printf("m=%2d: E_max = %" NFFT__FES__ "\n", m, E_max);
469  fprintf(fp1, "%" NFFT__FES__ "\n", E_max);
470  }
471  fclose(fp1);
472 
474  for (m = 3; m <= 9; m += 3)
475  {
476  printf("\nTest of the inverse linogram FFT for m=%d: \n", m);
477  sprintf(filename, "linogram_ifft_error%d.dat", m);
478  fp1 = fopen(filename, "w+");
479  for (max_i = 0; max_i <= 20; max_i += 2)
480  {
482  inverse_linogram_fft(f_direct, T, rr, f_tilde, N, max_i, m);
483 
485  E_max = NFFT(error_l_infty_complex)(f_hat, f_tilde, N * N);
486  printf("%3d iterations: E_max = %" NFFT__FES__ "\n", max_i, E_max);
487  fprintf(fp1, "%" NFFT__FES__ "\n", E_max);
488  }
489  fclose(fp1);
490  }
491 
493  NFFT(free)(x);
494  NFFT(free)(w);
495  NFFT(free)(f_hat);
496  NFFT(free)(f);
497  NFFT(free)(f_direct);
498  NFFT(free)(f_tilde);
499 
500  return EXIT_SUCCESS;
501 }
502 /* \} */
static int max_i(int a, int b)
max
Definition: fastsum.c:46
#define PRECOMPUTE_DAMP
Definition: nfft3.h:792
#define MALLOC_X
Definition: nfft3.h:199
#define MALLOC_F_HAT
Definition: nfft3.h:200
#define PRECOMPUTE_WEIGHT
Definition: nfft3.h:791
static int inverse_linogram_fft(NFFT_C *f, int T, int rr, NFFT_C *f_hat, int NN, int max_i, int m)
NFFT-based inverse pseudo-polar FFT.
#define FFTW_INIT
Definition: nfft3.h:203
#define MALLOC_F
Definition: nfft3.h:201
#define FFT_OUT_OF_PLACE
Definition: nfft3.h:202
#define PRE_LIN_PSI
Definition: nfft3.h:195
#define PRE_PSI
Definition: nfft3.h:197
#define CGNR
Definition: nfft3.h:788
int main(int argc, char **argv)
test program for various parameters
static int linogram_dft(NFFT_C *f_hat, int NN, NFFT_C *f, int T, int rr, int m)
discrete pseudo-polar FFT
#define PRE_FULL_PSI
Definition: nfft3.h:198
#define PRE_PHI_HUT
Definition: nfft3.h:193
static int comparison_fft(FILE *fp, int N, int T, int rr)
Comparison of the FFTW, linogram FFT, and inverse linogram FFT.
static int linogram_fft(NFFT_C *f_hat, int NN, NFFT_C *f, int T, int rr, int m)
NFFT-based pseudo-polar FFT.
static int linogram_grid(int T, int rr, NFFT_R *x, NFFT_R *w)
Generates the points x with weights w for the linogram grid with T slopes and R offsets.