File: linogram_fft_test.c.in

package info (click to toggle)
nfft 3.5.3-5
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid, trixie
  • size: 8,972 kB
  • sloc: ansic: 44,442; lisp: 1,697; makefile: 1,322; sh: 744; sed: 5
file content (501 lines) | stat: -rw-r--r-- 15,336 bytes parent folder | download
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
/*
 * Copyright (c) 2002, 2017 Jens Keiner, Stefan Kunis, Daniel Potts
 *
 * This program is free software; you can redistribute it and/or modify it under
 * the terms of the GNU General Public License as published by the Free Software
 * Foundation; either version 2 of the License, or (at your option) any later
 * version.
 *
 * This program is distributed in the hope that it will be useful, but WITHOUT
 * ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
 * FOR A PARTICULAR PURPOSE.  See the GNU General Public License for more
 * details.
 *
 * You should have received a copy of the GNU General Public License along with
 * this program; if not, write to the Free Software Foundation, Inc., 51
 * Franklin Street, Fifth Floor, Boston, MA 02110-1301, USA.
 */

/**
 * \file polarFFT/linogram_fft_test.c
 * \brief NFFT-based pseudo-polar FFT and inverse.
 *
 * Computes the NFFT-based pseudo-polar FFT and its inverse.
 * \author Markus Fenn
 * \date 2006
 */

#include <math.h>
#include <stdlib.h>
#include <complex.h>

#define @NFFT_PRECISION_MACRO@

#include "nfft3mp.h"

/**
 * \defgroup applications_polarFFT_linogramm linogram_fft_test
 * \ingroup applications_polarFFT
 * \{
 */

NFFT_R GLOBAL_elapsed_time;

/** Generates the points x with weights w
 *  for the linogram grid with T slopes and R offsets.
 */
static int linogram_grid(int T, int rr, NFFT_R *x, NFFT_R *w)
{
  int t, r;
  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));

  for (t = -T / 2; t < T / 2; t++)
  {
    for (r = -rr / 2; r < rr / 2; r++)
    {
      if (t < 0)
      {
        x[2 * ((t + T / 2) * rr + (r + rr / 2)) + 0] = (NFFT_R) (r) / (NFFT_R)(rr);
        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)
            * (NFFT_R)(r) / (NFFT_R)(rr);
      }
      else
      {
        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)
            * (NFFT_R)(r) / (NFFT_R)(rr);
        x[2 * ((t + T / 2) * rr + (r + rr / 2)) + 1] = (NFFT_R) r / (NFFT_R)(rr);
      }
      if (r == 0)
        w[(t + T / 2) * rr + (r + rr / 2)] = NFFT_K(1.0) / NFFT_K(4.0) / W;
      else
        w[(t + T / 2) * rr + (r + rr / 2)] = NFFT_M(fabs)((NFFT_R) r) / W;
    }
  }

  return T * rr; /** return the number of knots        */
}

/** discrete pseudo-polar FFT */
static int linogram_dft(NFFT_C *f_hat, int NN, NFFT_C *f, int T, int rr, int m)
{
  double t0, t1;
  int j, k; /**< index for nodes and freqencies   */
  NFFT(plan) my_nfft_plan; /**< plan for the nfft-2D             */

  NFFT_R *x, *w; /**< knots and associated weights     */

  int N[2], n[2];
  int M = T * rr; /**< number of knots                  */

  N[0] = NN;
  n[0] = 2 * N[0]; /**< oversampling factor sigma=2      */
  N[1] = NN;
  n[1] = 2 * N[1]; /**< oversampling factor sigma=2      */

  x = (NFFT_R *) NFFT(malloc)((size_t)(2 * T * rr) * (sizeof(NFFT_R)));
  if (x == NULL)
    return EXIT_FAILURE;

  w = (NFFT_R *) NFFT(malloc)((size_t)(T * rr) * (sizeof(NFFT_R)));
  if (w == NULL)
    return EXIT_FAILURE;

  /** init two dimensional NFFT plan */
  NFFT(init_guru)(&my_nfft_plan, 2, N, M, n, m,
      PRE_PHI_HUT | PRE_PSI | MALLOC_X | MALLOC_F_HAT | MALLOC_F | FFTW_INIT,
      FFTW_MEASURE);

  /** init nodes from linogram grid*/
  linogram_grid(T, rr, x, w);
  for (j = 0; j < my_nfft_plan.M_total; j++)
  {
    my_nfft_plan.x[2 * j + 0] = x[2 * j + 0];
    my_nfft_plan.x[2 * j + 1] = x[2 * j + 1];
  }

  /** init Fourier coefficients from given image */
  for (k = 0; k < my_nfft_plan.N_total; k++)
    my_nfft_plan.f_hat[k] = f_hat[k];

  /** NFFT-2D */
  t0 = NFFT(clock_gettime_seconds)();
  NFFT(trafo_direct)(&my_nfft_plan);
  t1 = NFFT(clock_gettime_seconds)();
  GLOBAL_elapsed_time = (t1 - t0);

  /** copy result */
  for (j = 0; j < my_nfft_plan.M_total; j++)
    f[j] = my_nfft_plan.f[j];

  /** finalise the plans and free the variables */
  NFFT(finalize)(&my_nfft_plan);
  NFFT(free)(x);
  NFFT(free)(w);

  return EXIT_SUCCESS;
}

/** NFFT-based pseudo-polar FFT */
static int linogram_fft(NFFT_C *f_hat, int NN, NFFT_C *f, int T, int rr, int m)
{
  double t0, t1;
  int j, k; /**< index for nodes and freqencies   */
  NFFT(plan) my_nfft_plan; /**< plan for the nfft-2D             */

  NFFT_R *x, *w; /**< knots and associated weights     */

  int N[2], n[2];
  int M = T * rr; /**< number of knots                  */

  N[0] = NN;
  n[0] = 2 * N[0]; /**< oversampling factor sigma=2      */
  N[1] = NN;
  n[1] = 2 * N[1]; /**< oversampling factor sigma=2      */

  x = (NFFT_R *) NFFT(malloc)((size_t)(2 * T * rr) * (sizeof(NFFT_R)));
  if (x == NULL)
    return EXIT_FAILURE;

  w = (NFFT_R *) NFFT(malloc)((size_t)(T * rr) * (sizeof(NFFT_R)));
  if (w == NULL)
    return EXIT_FAILURE;

  /** init two dimensional NFFT plan */
  NFFT(init_guru)(&my_nfft_plan, 2, N, M, n, m,
      PRE_PHI_HUT | PRE_PSI | MALLOC_X | MALLOC_F_HAT | MALLOC_F | FFTW_INIT
          | FFT_OUT_OF_PLACE,
      FFTW_MEASURE | FFTW_DESTROY_INPUT);

  /** init nodes from linogram grid*/
  linogram_grid(T, rr, x, w);
  for (j = 0; j < my_nfft_plan.M_total; j++)
  {
    my_nfft_plan.x[2 * j + 0] = x[2 * j + 0];
    my_nfft_plan.x[2 * j + 1] = x[2 * j + 1];
  }

  /** precompute psi, the entries of the matrix B */
  if (my_nfft_plan.flags & PRE_LIN_PSI)
    NFFT(precompute_lin_psi)(&my_nfft_plan);

  if (my_nfft_plan.flags & PRE_PSI)
    NFFT(precompute_psi)(&my_nfft_plan);

  if (my_nfft_plan.flags & PRE_FULL_PSI)
    NFFT(precompute_full_psi)(&my_nfft_plan);

  /** init Fourier coefficients from given image */
  for (k = 0; k < my_nfft_plan.N_total; k++)
    my_nfft_plan.f_hat[k] = f_hat[k];

  /** NFFT-2D */
  t0 = NFFT(clock_gettime_seconds)();
  NFFT(trafo)(&my_nfft_plan);
  t1 = NFFT(clock_gettime_seconds)();
  GLOBAL_elapsed_time = (t1 - t0);

  /** copy result */
  for (j = 0; j < my_nfft_plan.M_total; j++)
    f[j] = my_nfft_plan.f[j];

  /** finalise the plans and free the variables */
  NFFT(finalize)(&my_nfft_plan);
  NFFT(free)(x);
  NFFT(free)(w);

  return EXIT_SUCCESS;
}

/** NFFT-based inverse pseudo-polar FFT */
static int inverse_linogram_fft(NFFT_C *f, int T, int rr, NFFT_C *f_hat, int NN,
    int max_i, int m)
{
  double t0, t1;
  int j, k; /**< index for nodes and freqencies   */
  NFFT(plan) my_nfft_plan; /**< plan for the nfft-2D             */
  SOLVER(plan_complex) my_infft_plan; /**< plan for the inverse nfft        */

  NFFT_R *x, *w; /**< knots and associated weights     */
  int l; /**< index for iterations             */

  int N[2], n[2];
  int M = T * rr; /**< number of knots                  */

  N[0] = NN;
  n[0] = 2 * N[0]; /**< oversampling factor sigma=2      */
  N[1] = NN;
  n[1] = 2 * N[1]; /**< oversampling factor sigma=2      */

  x = (NFFT_R *) NFFT(malloc)((size_t)(2 * T * rr) * (sizeof(NFFT_R)));
  if (x == NULL)
    return EXIT_FAILURE;

  w = (NFFT_R *) NFFT(malloc)((size_t)(T * rr) * (sizeof(NFFT_R)));
  if (w == NULL)
    return EXIT_FAILURE;

  /** init two dimensional NFFT plan */
  NFFT(init_guru)(&my_nfft_plan, 2, N, M, n, m,
      PRE_PHI_HUT | PRE_PSI | MALLOC_X | MALLOC_F_HAT | MALLOC_F | FFTW_INIT
          | FFT_OUT_OF_PLACE,
      FFTW_MEASURE | FFTW_DESTROY_INPUT);

  /** init two dimensional infft plan */
  SOLVER(init_advanced_complex)(&my_infft_plan,
      (NFFT(mv_plan_complex)*) (&my_nfft_plan), CGNR | PRECOMPUTE_WEIGHT);

  /** init nodes, given samples and weights */
  linogram_grid(T, rr, x, w);
  for (j = 0; j < my_nfft_plan.M_total; j++)
  {
    my_nfft_plan.x[2 * j + 0] = x[2 * j + 0];
    my_nfft_plan.x[2 * j + 1] = x[2 * j + 1];
    my_infft_plan.y[j] = f[j];
    my_infft_plan.w[j] = w[j];
  }

  /** precompute psi, the entries of the matrix B */
  if (my_nfft_plan.flags & PRE_LIN_PSI)
    NFFT(precompute_lin_psi)(&my_nfft_plan);

  if (my_nfft_plan.flags & PRE_PSI)
    NFFT(precompute_psi)(&my_nfft_plan);

  if (my_nfft_plan.flags & PRE_FULL_PSI)
    NFFT(precompute_full_psi)(&my_nfft_plan);

  /** initialise damping factors */
  if (my_infft_plan.flags & PRECOMPUTE_DAMP)
    for (j = 0; j < my_nfft_plan.N[0]; j++)
      for (k = 0; k < my_nfft_plan.N[1]; k++)
      {
        my_infft_plan.w_hat[j * my_nfft_plan.N[1] + k] = (
            NFFT_M(sqrt)(
                NFFT_M(pow)((NFFT_R)(j - my_nfft_plan.N[0] / 2), NFFT_K(2.0))
                    + NFFT_M(pow)((NFFT_R)(k - my_nfft_plan.N[1] / 2), NFFT_K(2.0)))
                > (NFFT_R)(my_nfft_plan.N[0] / 2) ? 0 : 1);
      }

  /** initialise some guess f_hat_0 */
  for (k = 0; k < my_nfft_plan.N_total; k++)
    my_infft_plan.f_hat_iter[k] = NFFT_K(0.0) + _Complex_I * NFFT_K(0.0);

  t0 = NFFT(clock_gettime_seconds)();
  /** solve the system */
  SOLVER(before_loop_complex)(&my_infft_plan);

  if (max_i < 1)
  {
    l = 1;
    for (k = 0; k < my_nfft_plan.N_total; k++)
      my_infft_plan.f_hat_iter[k] = my_infft_plan.p_hat_iter[k];
  }
  else
  {
    for (l = 1; l <= max_i; l++)
    {
      SOLVER(loop_one_step_complex)(&my_infft_plan);
    }
  }

  t1 = NFFT(clock_gettime_seconds)();
  GLOBAL_elapsed_time = (t1 - t0);

  /** copy result */
  for (k = 0; k < my_nfft_plan.N_total; k++)
    f_hat[k] = my_infft_plan.f_hat_iter[k];

  /** finalise the plans and free the variables */
  SOLVER(finalize_complex)(&my_infft_plan);
  NFFT(finalize)(&my_nfft_plan);
  NFFT(free)(x);
  NFFT(free)(w);

  return EXIT_SUCCESS;
}

/** Comparison of the FFTW, linogram FFT, and inverse linogram FFT */
static int comparison_fft(FILE *fp, int N, int T, int rr)
{
  double t0, t1;
  FFTW(plan) my_fftw_plan;
  NFFT_C *f_hat, *f;
  int m, k;
  NFFT_R t_fft, t_dft_linogram;

  f_hat = (NFFT_C *) NFFT(malloc)(sizeof(NFFT_C) * (size_t)(N * N));
  f = (NFFT_C *) NFFT(malloc)(sizeof(NFFT_C) * (size_t)((T * rr / 4) * 5));

  my_fftw_plan = FFTW(plan_dft_2d)(N, N, f_hat, f, FFTW_BACKWARD, FFTW_MEASURE);

  for (k = 0; k < N * N; k++)
    f_hat[k] = NFFT(drand48)() + _Complex_I * NFFT(drand48)();

  t0 = NFFT(clock_gettime_seconds)();
  for (m = 0; m < 65536 / N; m++)
  {
    FFTW(execute)(my_fftw_plan);
    /* touch */
    f_hat[2] = NFFT_K(2.0) * f_hat[0];
  }
  t1 = NFFT(clock_gettime_seconds)();
  GLOBAL_elapsed_time = (t1 - t0);
  t_fft = (NFFT_R)(N) * GLOBAL_elapsed_time / NFFT_K(65536.0);

  if (N < 256)
  {
    linogram_dft(f_hat, N, f, T, rr, m);
    t_dft_linogram = GLOBAL_elapsed_time;
  }

  for (m = 3; m <= 9; m += 3)
  {
    if ((m == 3) && (N < 256))
      fprintf(fp, "%d\t&\t&\t%1.1" NFFT__FES__ "&\t%1.1" NFFT__FES__ "&\t%d\t", N, t_fft, t_dft_linogram,
          m);
    else if (m == 3)
      fprintf(fp, "%d\t&\t&\t%1.1" NFFT__FES__ "&\t       &\t%d\t", N, t_fft, m);
    else
      fprintf(fp, "  \t&\t&\t       &\t       &\t%d\t", m);

    printf("N=%d\tt_fft=%1.1" NFFT__FES__ "\tt_dft_linogram=%1.1" NFFT__FES__ "\tm=%d\t", N, t_fft,
        t_dft_linogram, m);

    linogram_fft(f_hat, N, f, T, rr, m);
    fprintf(fp, "%1.1" NFFT__FES__ "&\t", GLOBAL_elapsed_time);
    printf("t_linogram=%1.1" NFFT__FES__ "\t", GLOBAL_elapsed_time);
    inverse_linogram_fft(f, T, rr, f_hat, N, m + 3, m);
    if (m == 9)
      fprintf(fp, "%1.1" NFFT__FES__ "\\\\\\hline\n", GLOBAL_elapsed_time);
    else
      fprintf(fp, "%1.1" NFFT__FES__ "\\\\\n", GLOBAL_elapsed_time);
    printf("t_ilinogram=%1.1" NFFT__FES__ "\n", GLOBAL_elapsed_time);
  }

  fflush(fp);

  NFFT(free)(f);
  NFFT(free)(f_hat);

  return EXIT_SUCCESS;
}

/** test program for various parameters */
int main(int argc, char **argv)
{
  int N; /**< linogram FFT size NxN              */
  int T, rr; /**< number of directions/offsets     */
  int M; /**< number of knots of linogram grid   */
  NFFT_R *x, *w; /**< knots and associated weights     */
  NFFT_C *f_hat, *f, *f_direct, *f_tilde;
  int k;
  int max_i; /**< number of iterations             */
  int m;
  NFFT_R temp1, temp2, E_max = NFFT_K(0.0);
  FILE *fp1, *fp2;
  char filename[30];
  int logN;

  if (argc != 4)
  {
    printf("linogram_fft_test N T R \n");
    printf("\n");
    printf("N          linogram FFT of size NxN    \n");
    printf("T          number of slopes          \n");
    printf("R          number of offsets         \n");

    /** Hence, comparison of the FFTW, linogram FFT, and inverse linogram FFT */
    printf("\nHence, comparison FFTW, linogram FFT and inverse linogram FFT\n");
    fp1 = fopen("linogram_comparison_fft.dat", "w");
    if (fp1 == NULL)
      return (-1);
    for (logN = 4; logN <= 8; logN++)
      comparison_fft(fp1, (int)(1U << logN), 3 * (int)(1U << logN),
          3 * (int)(1U << (logN - 1)));
    fclose(fp1);

    return EXIT_FAILURE;
  }

  N = atoi(argv[1]);
  T = atoi(argv[2]);
  rr = atoi(argv[3]);
  printf("N=%d, linogram grid with T=%d, R=%d => ", N, T, rr);

  x = (NFFT_R *) NFFT(malloc)((size_t)(5 * T * rr / 2) * (sizeof(NFFT_R)));
  w = (NFFT_R *) NFFT(malloc)((size_t)(5 * T * rr / 4) * (sizeof(NFFT_R)));

  f_hat = (NFFT_C *) NFFT(malloc)(sizeof(NFFT_C) * (size_t)(N * N));
  f = (NFFT_C *) NFFT(malloc)(sizeof(NFFT_C) * (size_t)(5 * T * rr / 4)); /* 4/pi*log(1+sqrt(2)) = 1.122... < 1.25 */
  f_direct = (NFFT_C *) NFFT(malloc)(sizeof(NFFT_C) * (size_t)(5 * T * rr / 4));
  f_tilde = (NFFT_C *) NFFT(malloc)(sizeof(NFFT_C) * (size_t)(N * N));

  /** generate knots of linogram grid */
  M = linogram_grid(T, rr, x, w);
  printf("M=%d.\n", M);

  /** load data */
  fp1 = fopen("input_data_r.dat", "r");
  fp2 = fopen("input_data_i.dat", "r");
  if ((fp1 == NULL) || (fp2 == NULL))
    return EXIT_FAILURE;
  for (k = 0; k < N * N; k++)
  {
    fscanf(fp1, NFFT__FR__ " ", &temp1);
    fscanf(fp2, NFFT__FR__ " ", &temp2);
    f_hat[k] = temp1 + _Complex_I * temp2;
  }
  fclose(fp1);
  fclose(fp2);

  /** direct linogram FFT */
  linogram_dft(f_hat, N, f_direct, T, rr, 1);
  //  linogram_fft(f_hat,N,f_direct,T,R,12);

  /** Test of the linogram FFT with different m */
  printf("\nTest of the linogram FFT: \n");
  fp1 = fopen("linogram_fft_error.dat", "w+");
  for (m = 1; m <= 12; m++)
  {
    /** fast linogram FFT */
    linogram_fft(f_hat, N, f, T, rr, m);

    /** error of fast linogram FFT */
    E_max = NFFT(error_l_infty_complex)(f_direct, f, M);
    //E_max=NFFT(error_l_2_complex)(f_direct,f,M);

    printf("m=%2d: E_max = %" NFFT__FES__ "\n", m, E_max);
    fprintf(fp1, "%" NFFT__FES__ "\n", E_max);
  }
  fclose(fp1);

  /** Test of the inverse linogram FFT for different m in dependece of the iteration number*/
  for (m = 3; m <= 9; m += 3)
  {
    printf("\nTest of the inverse linogram FFT for m=%d: \n", m);
    sprintf(filename, "linogram_ifft_error%d.dat", m);
    fp1 = fopen(filename, "w+");
    for (max_i = 0; max_i <= 20; max_i += 2)
    {
      /** inverse linogram FFT */
      inverse_linogram_fft(f_direct, T, rr, f_tilde, N, max_i, m);

      /** compute maximum error */
      E_max = NFFT(error_l_infty_complex)(f_hat, f_tilde, N * N);
      printf("%3d iterations: E_max = %" NFFT__FES__ "\n", max_i, E_max);
      fprintf(fp1, "%" NFFT__FES__ "\n", E_max);
    }
    fclose(fp1);
  }

  /** free the variables */
  NFFT(free)(x);
  NFFT(free)(w);
  NFFT(free)(f_hat);
  NFFT(free)(f);
  NFFT(free)(f_direct);
  NFFT(free)(f_tilde);

  return EXIT_SUCCESS;
}
/* \} */