32 #define NFFT_PRECISION_DOUBLE
42 NFFT_R GLOBAL_elapsed_time;
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));
52 for (t = -T / 2; t < T / 2; t++)
54 for (r = -rr / 2; r < rr / 2; r++)
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);
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);
69 w[(t + T / 2) * rr + (r + rr / 2)] = NFFT_K(1.0) / NFFT_K(4.0) / W;
71 w[(t + T / 2) * rr + (r + rr / 2)] = NFFT_M(fabs)((NFFT_R) r) / W;
79 static int linogram_dft(NFFT_C *f_hat,
int NN, NFFT_C *f,
int T,
int rr,
int m)
83 NFFT(plan) my_nfft_plan;
95 x = (NFFT_R *) NFFT(malloc)((size_t)(2 * T * rr) * (
sizeof(NFFT_R)));
99 w = (NFFT_R *) NFFT(malloc)((size_t)(T * rr) * (
sizeof(NFFT_R)));
104 NFFT(init_guru)(&my_nfft_plan, 2, N, M, n, m,
110 for (j = 0; j < my_nfft_plan.M_total; j++)
112 my_nfft_plan.x[2 * j + 0] = x[2 * j + 0];
113 my_nfft_plan.x[2 * j + 1] = x[2 * j + 1];
117 for (k = 0; k < my_nfft_plan.N_total; k++)
118 my_nfft_plan.f_hat[k] = f_hat[k];
121 t0 = NFFT(clock_gettime_seconds)();
122 NFFT(trafo_direct)(&my_nfft_plan);
123 t1 = NFFT(clock_gettime_seconds)();
124 GLOBAL_elapsed_time = (t1 - t0);
127 for (j = 0; j < my_nfft_plan.M_total; j++)
128 f[j] = my_nfft_plan.f[j];
131 NFFT(finalize)(&my_nfft_plan);
139 static int linogram_fft(NFFT_C *f_hat,
int NN, NFFT_C *f,
int T,
int rr,
int m)
143 NFFT(plan) my_nfft_plan;
155 x = (NFFT_R *) NFFT(malloc)((size_t)(2 * T * rr) * (
sizeof(NFFT_R)));
159 w = (NFFT_R *) NFFT(malloc)((size_t)(T * rr) * (
sizeof(NFFT_R)));
164 NFFT(init_guru)(&my_nfft_plan, 2, N, M, n, m,
167 FFTW_MEASURE | FFTW_DESTROY_INPUT);
171 for (j = 0; j < my_nfft_plan.M_total; j++)
173 my_nfft_plan.x[2 * j + 0] = x[2 * j + 0];
174 my_nfft_plan.x[2 * j + 1] = x[2 * j + 1];
179 NFFT(precompute_lin_psi)(&my_nfft_plan);
181 if (my_nfft_plan.flags &
PRE_PSI)
182 NFFT(precompute_psi)(&my_nfft_plan);
185 NFFT(precompute_full_psi)(&my_nfft_plan);
188 for (k = 0; k < my_nfft_plan.N_total; k++)
189 my_nfft_plan.f_hat[k] = f_hat[k];
192 t0 = NFFT(clock_gettime_seconds)();
193 NFFT(trafo)(&my_nfft_plan);
194 t1 = NFFT(clock_gettime_seconds)();
195 GLOBAL_elapsed_time = (t1 - t0);
198 for (j = 0; j < my_nfft_plan.M_total; j++)
199 f[j] = my_nfft_plan.f[j];
202 NFFT(finalize)(&my_nfft_plan);
215 NFFT(plan) my_nfft_plan;
216 SOLVER(plan_complex) my_infft_plan;
229 x = (NFFT_R *) NFFT(malloc)((size_t)(2 * T * rr) * (
sizeof(NFFT_R)));
233 w = (NFFT_R *) NFFT(malloc)((size_t)(T * rr) * (
sizeof(NFFT_R)));
238 NFFT(init_guru)(&my_nfft_plan, 2, N, M, n, m,
241 FFTW_MEASURE | FFTW_DESTROY_INPUT);
244 SOLVER(init_advanced_complex)(&my_infft_plan,
249 for (j = 0; j < my_nfft_plan.M_total; j++)
251 my_nfft_plan.x[2 * j + 0] = x[2 * j + 0];
252 my_nfft_plan.x[2 * j + 1] = x[2 * j + 1];
253 my_infft_plan.y[j] = f[j];
254 my_infft_plan.w[j] = w[j];
259 NFFT(precompute_lin_psi)(&my_nfft_plan);
261 if (my_nfft_plan.flags &
PRE_PSI)
262 NFFT(precompute_psi)(&my_nfft_plan);
265 NFFT(precompute_full_psi)(&my_nfft_plan);
269 for (j = 0; j < my_nfft_plan.N[0]; j++)
270 for (k = 0; k < my_nfft_plan.N[1]; k++)
272 my_infft_plan.w_hat[j * my_nfft_plan.N[1] + k] = (
274 NFFT_M(pow)((NFFT_R)(j - my_nfft_plan.N[0] / 2), NFFT_K(2.0))
275 + NFFT_M(pow)((NFFT_R)(k - my_nfft_plan.N[1] / 2), NFFT_K(2.0)))
276 > (NFFT_R)(my_nfft_plan.N[0] / 2) ? 0 : 1);
280 for (k = 0; k < my_nfft_plan.N_total; k++)
281 my_infft_plan.f_hat_iter[k] = NFFT_K(0.0) + _Complex_I * NFFT_K(0.0);
283 t0 = NFFT(clock_gettime_seconds)();
285 SOLVER(before_loop_complex)(&my_infft_plan);
290 for (k = 0; k < my_nfft_plan.N_total; k++)
291 my_infft_plan.f_hat_iter[k] = my_infft_plan.p_hat_iter[k];
295 for (l = 1; l <=
max_i; l++)
297 SOLVER(loop_one_step_complex)(&my_infft_plan);
301 t1 = NFFT(clock_gettime_seconds)();
302 GLOBAL_elapsed_time = (t1 - t0);
305 for (k = 0; k < my_nfft_plan.N_total; k++)
306 f_hat[k] = my_infft_plan.f_hat_iter[k];
309 SOLVER(finalize_complex)(&my_infft_plan);
310 NFFT(finalize)(&my_nfft_plan);
321 FFTW(plan) my_fftw_plan;
324 NFFT_R t_fft, t_dft_linogram;
326 f_hat = (NFFT_C *) NFFT(malloc)(
sizeof(NFFT_C) * (
size_t)(N * N));
327 f = (NFFT_C *) NFFT(malloc)(
sizeof(NFFT_C) * (
size_t)((T * rr / 4) * 5));
329 my_fftw_plan = FFTW(plan_dft_2d)(N, N, f_hat, f, FFTW_BACKWARD, FFTW_MEASURE);
331 for (k = 0; k < N * N; k++)
332 f_hat[k] = NFFT(drand48)() + _Complex_I * NFFT(drand48)();
334 t0 = NFFT(clock_gettime_seconds)();
335 for (m = 0; m < 65536 / N; m++)
337 FFTW(execute)(my_fftw_plan);
339 f_hat[2] = NFFT_K(2.0) * f_hat[0];
341 t1 = NFFT(clock_gettime_seconds)();
342 GLOBAL_elapsed_time = (t1 - t0);
343 t_fft = (NFFT_R)(N) * GLOBAL_elapsed_time / NFFT_K(65536.0);
348 t_dft_linogram = GLOBAL_elapsed_time;
351 for (m = 3; m <= 9; m += 3)
353 if ((m == 3) && (N < 256))
354 fprintf(fp,
"%d\t&\t&\t%1.1" NFFT__FES__
"&\t%1.1" NFFT__FES__
"&\t%d\t", N, t_fft, t_dft_linogram,
357 fprintf(fp,
"%d\t&\t&\t%1.1" NFFT__FES__
"&\t &\t%d\t", N, t_fft, m);
359 fprintf(fp,
" \t&\t&\t &\t &\t%d\t", m);
361 printf(
"N=%d\tt_fft=%1.1" NFFT__FES__
"\tt_dft_linogram=%1.1" NFFT__FES__
"\tm=%d\t", N, t_fft,
365 fprintf(fp,
"%1.1" NFFT__FES__
"&\t", GLOBAL_elapsed_time);
366 printf(
"t_linogram=%1.1" NFFT__FES__
"\t", GLOBAL_elapsed_time);
369 fprintf(fp,
"%1.1" NFFT__FES__
"\\\\\\hline\n", GLOBAL_elapsed_time);
371 fprintf(fp,
"%1.1" NFFT__FES__
"\\\\\n", GLOBAL_elapsed_time);
372 printf(
"t_ilinogram=%1.1" NFFT__FES__
"\n", GLOBAL_elapsed_time);
384 int main(
int argc,
char **argv)
390 NFFT_C *f_hat, *f, *f_direct, *f_tilde;
394 NFFT_R temp1, temp2, E_max = NFFT_K(0.0);
401 printf(
"linogram_fft_test N T R \n");
403 printf(
"N linogram FFT of size NxN \n");
404 printf(
"T number of slopes \n");
405 printf(
"R number of offsets \n");
408 printf(
"\nHence, comparison FFTW, linogram FFT and inverse linogram FFT\n");
409 fp1 = fopen(
"linogram_comparison_fft.dat",
"w");
412 for (logN = 4; logN <= 8; logN++)
414 3 * (
int)(1U << (logN - 1)));
423 printf(
"N=%d, linogram grid with T=%d, R=%d => ", N, T, rr);
425 x = (NFFT_R *) NFFT(malloc)((size_t)(5 * T * rr / 2) * (
sizeof(NFFT_R)));
426 w = (NFFT_R *) NFFT(malloc)((size_t)(5 * T * rr / 4) * (
sizeof(NFFT_R)));
428 f_hat = (NFFT_C *) NFFT(malloc)(
sizeof(NFFT_C) * (
size_t)(N * N));
429 f = (NFFT_C *) NFFT(malloc)(
sizeof(NFFT_C) * (
size_t)(5 * T * rr / 4));
430 f_direct = (NFFT_C *) NFFT(malloc)(
sizeof(NFFT_C) * (
size_t)(5 * T * rr / 4));
431 f_tilde = (NFFT_C *) NFFT(malloc)(
sizeof(NFFT_C) * (
size_t)(N * N));
435 printf(
"M=%d.\n", M);
438 fp1 = fopen(
"input_data_r.dat",
"r");
439 fp2 = fopen(
"input_data_i.dat",
"r");
440 if ((fp1 == NULL) || (fp2 == NULL))
442 for (k = 0; k < N * N; k++)
444 fscanf(fp1, NFFT__FR__
" ", &temp1);
445 fscanf(fp2, NFFT__FR__
" ", &temp2);
446 f_hat[k] = temp1 + _Complex_I * temp2;
456 printf(
"\nTest of the linogram FFT: \n");
457 fp1 = fopen(
"linogram_fft_error.dat",
"w+");
458 for (m = 1; m <= 12; m++)
464 E_max = NFFT(error_l_infty_complex)(f_direct, f, M);
467 printf(
"m=%2d: E_max = %" NFFT__FES__
"\n", m, E_max);
468 fprintf(fp1,
"%" NFFT__FES__
"\n", E_max);
473 for (m = 3; m <= 9; m += 3)
475 printf(
"\nTest of the inverse linogram FFT for m=%d: \n", m);
476 sprintf(filename,
"linogram_ifft_error%d.dat", m);
477 fp1 = fopen(filename,
"w+");
484 E_max = NFFT(error_l_infty_complex)(f_hat, f_tilde, N * N);
485 printf(
"%3d iterations: E_max = %" NFFT__FES__
"\n",
max_i, E_max);
486 fprintf(fp1,
"%" NFFT__FES__
"\n", E_max);
496 NFFT(free)(f_direct);
static int max_i(int a, int b)
max
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.
static int linogram_fft(NFFT_C *f_hat, int NN, NFFT_C *f, int T, int rr, int m)
NFFT-based pseudo-polar FFT.
int main(int argc, char **argv)
test program for various parameters
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.
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_dft(NFFT_C *f_hat, int NN, NFFT_C *f, int T, int rr, int m)
discrete pseudo-polar FFT
#define PRECOMPUTE_WEIGHT