![]() |
NFFT
3.3.2
|
00001 /* 00002 * Copyright (c) 2002, 2016 Jens Keiner, Stefan Kunis, Daniel Potts 00003 * 00004 * This program is free software; you can redistribute it and/or modify it under 00005 * the terms of the GNU General Public License as published by the Free Software 00006 * Foundation; either version 2 of the License, or (at your option) any later 00007 * version. 00008 * 00009 * This program is distributed in the hope that it will be useful, but WITHOUT 00010 * ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS 00011 * FOR A PARTICULAR PURPOSE. See the GNU General Public License for more 00012 * details. 00013 * 00014 * You should have received a copy of the GNU General Public License along with 00015 * this program; if not, write to the Free Software Foundation, Inc., 51 00016 * Franklin Street, Fifth Floor, Boston, MA 02110-1301, USA. 00017 */ 00018 00029 #include <math.h> 00030 #include <stdlib.h> 00031 #include <complex.h> 00032 00033 #define NFFT_PRECISION_DOUBLE 00034 00035 #include "nfft3mp.h" 00036 00043 NFFT_R GLOBAL_elapsed_time; 00044 00058 static int mpolar_grid(int T, int S, NFFT_R *x, NFFT_R *w) 00059 { 00060 int t, r; 00061 NFFT_R W; 00062 int R2 = 2 * (int)(NFFT_M(lrint)(NFFT_M(ceil)(NFFT_M(sqrt)(NFFT_K(2.0)) * (NFFT_R)(S) / NFFT_K(2.0)))); 00063 NFFT_R xx, yy; 00064 int M = 0; 00065 00066 for (t = -T / 2; t < T / 2; t++) 00067 { 00068 for (r = -R2 / 2; r < R2 / 2; r++) 00069 { 00070 xx = (NFFT_R) (r) / (NFFT_R)(S) * NFFT_M(cos)(NFFT_KPI * (NFFT_R)(t) / (NFFT_R)(T)); 00071 yy = (NFFT_R) (r) / (NFFT_R)(S) * NFFT_M(sin)(NFFT_KPI * (NFFT_R)(t) / (NFFT_R)(T)); 00072 00073 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)) 00074 & ((-NFFT_K(0.5) - NFFT_K(1.0) / (NFFT_R) S) <= yy) 00075 & (yy <= (NFFT_K(0.5) + NFFT_K(1.0) / (NFFT_R) S))) 00076 { 00077 x[2 * M + 0] = xx; 00078 x[2 * M + 1] = yy; 00079 00080 if (r == 0) 00081 w[M] = NFFT_K(1.0) / NFFT_K(4.0); 00082 else 00083 w[M] = NFFT_M(fabs)((NFFT_R) r); 00084 00085 M++; 00086 } 00087 } 00088 } 00089 00091 W = NFFT_K(0.0); 00092 for (t = 0; t < M; t++) 00093 W += w[t]; 00094 00095 for (t = 0; t < M; t++) 00096 w[t] /= W; 00097 00098 return M; 00099 } 00100 00102 static int mpolar_dft(NFFT_C *f_hat, int NN, NFFT_C *f, int T, int S, int m) 00103 { 00104 double t0, t1; 00105 int j, k; 00106 NFFT(plan) my_nfft_plan; 00108 NFFT_R *x, *w; 00110 int N[2], n[2]; 00111 int M; 00113 N[0] = NN; 00114 n[0] = 2 * N[0]; 00115 N[1] = NN; 00116 n[1] = 2 * N[1]; 00118 x = (NFFT_R *) NFFT(malloc)((size_t)(5 * (T / 2) * S) * (sizeof(NFFT_R))); 00119 if (x == NULL) 00120 return EXIT_FAILURE; 00121 00122 w = (NFFT_R *) NFFT(malloc)((size_t)(5 * (T * S) / 4) * (sizeof(NFFT_R))); 00123 if (w == NULL) 00124 return EXIT_FAILURE; 00125 00127 M = mpolar_grid(T, S, x, w); 00128 NFFT(init_guru)(&my_nfft_plan, 2, N, M, n, m, 00129 PRE_PHI_HUT | PRE_PSI | MALLOC_X | MALLOC_F_HAT | MALLOC_F | FFTW_INIT 00130 | FFT_OUT_OF_PLACE, 00131 FFTW_MEASURE | FFTW_DESTROY_INPUT); 00132 00134 for (j = 0; j < my_nfft_plan.M_total; j++) 00135 { 00136 my_nfft_plan.x[2 * j + 0] = x[2 * j + 0]; 00137 my_nfft_plan.x[2 * j + 1] = x[2 * j + 1]; 00138 } 00139 00141 for (k = 0; k < my_nfft_plan.N_total; k++) 00142 my_nfft_plan.f_hat[k] = f_hat[k]; 00143 00144 t0 = NFFT(clock_gettime_seconds)(); 00145 00147 NFFT(trafo_direct)(&my_nfft_plan); 00148 00149 t1 = NFFT(clock_gettime_seconds)(); 00150 GLOBAL_elapsed_time = (t1 - t0); 00151 00153 for (j = 0; j < my_nfft_plan.M_total; j++) 00154 f[j] = my_nfft_plan.f[j]; 00155 00157 NFFT(finalize)(&my_nfft_plan); 00158 NFFT(free)(x); 00159 NFFT(free)(w); 00160 00161 return EXIT_SUCCESS; 00162 } 00163 00165 static int mpolar_fft(NFFT_C *f_hat, int NN, NFFT_C *f, int T, int S, int m) 00166 { 00167 double t0, t1; 00168 int j, k; 00169 NFFT(plan) my_nfft_plan; 00171 NFFT_R *x, *w; 00173 int N[2], n[2]; 00174 int M; 00176 N[0] = NN; 00177 n[0] = 2 * N[0]; 00178 N[1] = NN; 00179 n[1] = 2 * N[1]; 00181 x = (NFFT_R *) NFFT(malloc)((size_t)(5 * T * S / 2) * (sizeof(NFFT_R))); 00182 if (x == NULL) 00183 return EXIT_FAILURE; 00184 00185 w = (NFFT_R *) NFFT(malloc)((size_t)(5 * T * S / 4) * (sizeof(NFFT_R))); 00186 if (w == NULL) 00187 return EXIT_FAILURE; 00188 00190 M = mpolar_grid(T, S, x, w); 00191 NFFT(init_guru)(&my_nfft_plan, 2, N, M, n, m, 00192 PRE_PHI_HUT | PRE_PSI | MALLOC_X | MALLOC_F_HAT | MALLOC_F | FFTW_INIT 00193 | FFT_OUT_OF_PLACE, 00194 FFTW_MEASURE | FFTW_DESTROY_INPUT); 00195 00198 for (j = 0; j < my_nfft_plan.M_total; j++) 00199 { 00200 my_nfft_plan.x[2 * j + 0] = x[2 * j + 0]; 00201 my_nfft_plan.x[2 * j + 1] = x[2 * j + 1]; 00202 } 00203 00205 if (my_nfft_plan.flags & PRE_LIN_PSI) 00206 NFFT(precompute_lin_psi)(&my_nfft_plan); 00207 00208 if (my_nfft_plan.flags & PRE_PSI) 00209 NFFT(precompute_psi)(&my_nfft_plan); 00210 00211 if (my_nfft_plan.flags & PRE_FULL_PSI) 00212 NFFT(precompute_full_psi)(&my_nfft_plan); 00213 00215 for (k = 0; k < my_nfft_plan.N_total; k++) 00216 my_nfft_plan.f_hat[k] = f_hat[k]; 00217 00218 t0 = NFFT(clock_gettime_seconds)(); 00219 00221 NFFT(trafo)(&my_nfft_plan); 00222 00223 t1 = NFFT(clock_gettime_seconds)(); 00224 GLOBAL_elapsed_time = (t1 - t0); 00225 00227 for (j = 0; j < my_nfft_plan.M_total; j++) 00228 f[j] = my_nfft_plan.f[j]; 00229 00231 NFFT(finalize)(&my_nfft_plan); 00232 NFFT(free)(x); 00233 NFFT(free)(w); 00234 00235 return EXIT_SUCCESS; 00236 } 00237 00239 static int inverse_mpolar_fft(NFFT_C *f, int T, int S, NFFT_C *f_hat, int NN, int max_i, 00240 int m) 00241 { 00242 double t0, t1; 00243 int j, k; 00244 NFFT(plan) my_nfft_plan; 00245 SOLVER(plan_complex) my_infft_plan; 00247 NFFT_R *x, *w; 00248 int l; 00250 int N[2], n[2]; 00251 int M; 00253 N[0] = NN; 00254 n[0] = 2 * N[0]; 00255 N[1] = NN; 00256 n[1] = 2 * N[1]; 00258 x = (NFFT_R *) NFFT(malloc)((size_t)(5 * T * S / 2) * (sizeof(NFFT_R))); 00259 if (x == NULL) 00260 return EXIT_FAILURE; 00261 00262 w = (NFFT_R *) NFFT(malloc)((size_t)(5 * T * S / 4) * (sizeof(NFFT_R))); 00263 if (w == NULL) 00264 return EXIT_FAILURE; 00265 00267 M = mpolar_grid(T, S, x, w); 00268 NFFT(init_guru)(&my_nfft_plan, 2, N, M, n, m, 00269 PRE_PHI_HUT | PRE_PSI | MALLOC_X | MALLOC_F_HAT | MALLOC_F | FFTW_INIT 00270 | FFT_OUT_OF_PLACE, 00271 FFTW_MEASURE | FFTW_DESTROY_INPUT); 00272 00274 SOLVER(init_advanced_complex)(&my_infft_plan, 00275 (NFFT(mv_plan_complex)*) (&my_nfft_plan), CGNR | PRECOMPUTE_WEIGHT); 00276 00278 for (j = 0; j < my_nfft_plan.M_total; j++) 00279 { 00280 my_nfft_plan.x[2 * j + 0] = x[2 * j + 0]; 00281 my_nfft_plan.x[2 * j + 1] = x[2 * j + 1]; 00282 my_infft_plan.y[j] = f[j]; 00283 my_infft_plan.w[j] = w[j]; 00284 } 00285 00287 if (my_nfft_plan.flags & PRE_LIN_PSI) 00288 NFFT(precompute_lin_psi)(&my_nfft_plan); 00289 00290 if (my_nfft_plan.flags & PRE_PSI) 00291 NFFT(precompute_psi)(&my_nfft_plan); 00292 00293 if (my_nfft_plan.flags & PRE_FULL_PSI) 00294 NFFT(precompute_full_psi)(&my_nfft_plan); 00295 00297 if (my_infft_plan.flags & PRECOMPUTE_DAMP) 00298 for (j = 0; j < my_nfft_plan.N[0]; j++) 00299 for (k = 0; k < my_nfft_plan.N[1]; k++) 00300 { 00301 my_infft_plan.w_hat[j * my_nfft_plan.N[1] + k] = ( 00302 NFFT_M(sqrt)( 00303 NFFT_M(pow)((NFFT_R)(j - my_nfft_plan.N[0] / 2), NFFT_K(2.0)) 00304 + NFFT_M(pow)((NFFT_R)(k - my_nfft_plan.N[1] / 2), NFFT_K(2.0))) 00305 > (NFFT_R)(my_nfft_plan.N[0] / 2) ? NFFT_K(0.0) : NFFT_K(1.0)); 00306 } 00307 00309 for (k = 0; k < my_nfft_plan.N_total; k++) 00310 my_infft_plan.f_hat_iter[k] = NFFT_K(0.0) + _Complex_I * NFFT_K(0.0); 00311 00312 t0 = NFFT(clock_gettime_seconds)(); 00313 00315 SOLVER(before_loop_complex)(&my_infft_plan); 00316 00317 if (max_i < 1) 00318 { 00319 l = 1; 00320 for (k = 0; k < my_nfft_plan.N_total; k++) 00321 my_infft_plan.f_hat_iter[k] = my_infft_plan.p_hat_iter[k]; 00322 } 00323 else 00324 { 00325 for (l = 1; l <= max_i; l++) 00326 { 00327 SOLVER(loop_one_step_complex)(&my_infft_plan); 00328 } 00329 } 00330 00331 t1 = NFFT(clock_gettime_seconds)(); 00332 GLOBAL_elapsed_time = (t1 - t0); 00333 00335 for (k = 0; k < my_nfft_plan.N_total; k++) 00336 f_hat[k] = my_infft_plan.f_hat_iter[k]; 00337 00339 SOLVER(finalize_complex)(&my_infft_plan); 00340 NFFT(finalize)(&my_nfft_plan); 00341 NFFT(free)(x); 00342 NFFT(free)(w); 00343 00344 return EXIT_SUCCESS; 00345 } 00346 00348 static int comparison_fft(FILE *fp, int N, int T, int S) 00349 { 00350 double t0, t1; 00351 FFTW(plan) my_fftw_plan; 00352 NFFT_C *f_hat, *f; 00353 int m, k; 00354 NFFT_R t_fft, t_dft_mpolar; 00355 00356 f_hat = (NFFT_C *) NFFT(malloc)(sizeof(NFFT_C) * (size_t)(N * N)); 00357 f = (NFFT_C *) NFFT(malloc)(sizeof(NFFT_C) * (size_t)((T * S / 4) * 5)); 00358 00359 my_fftw_plan = FFTW(plan_dft_2d)(N, N, f_hat, f, FFTW_BACKWARD, FFTW_MEASURE); 00360 00361 for (k = 0; k < N * N; k++) 00362 f_hat[k] = NFFT(drand48)() + _Complex_I * NFFT(drand48)(); 00363 00364 t0 = NFFT(clock_gettime_seconds)(); 00365 for (m = 0; m < 65536 / N; m++) 00366 { 00367 FFTW(execute)(my_fftw_plan); 00368 /* touch */ 00369 f_hat[2] = NFFT_K(2.0) * f_hat[0]; 00370 } 00371 t1 = NFFT(clock_gettime_seconds)(); 00372 GLOBAL_elapsed_time = (t1 - t0); 00373 t_fft = (NFFT_R)(N) * GLOBAL_elapsed_time / NFFT_K(65536.0); 00374 00375 if (N < 256) 00376 { 00377 mpolar_dft(f_hat, N, f, T, S, 1); 00378 t_dft_mpolar = GLOBAL_elapsed_time; 00379 } 00380 00381 for (m = 3; m <= 9; m += 3) 00382 { 00383 if ((m == 3) && (N < 256)) 00384 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); 00385 else if (m == 3) 00386 fprintf(fp, "%d\t&\t&\t%1.1" NFFT__FES__ "&\t &\t%d\t", N, t_fft, m); 00387 else 00388 fprintf(fp, " \t&\t&\t &\t &\t%d\t", m); 00389 00390 printf("N=%d\tt_fft=%1.1" NFFT__FES__ "\tt_dft_mpolar=%1.1" NFFT__FES__ "\tm=%d\t", N, t_fft, 00391 t_dft_mpolar, m); 00392 00393 mpolar_fft(f_hat, N, f, T, S, m); 00394 fprintf(fp, "%1.1" NFFT__FES__ "&\t", GLOBAL_elapsed_time); 00395 printf("t_mpolar=%1.1" NFFT__FES__ "\t", GLOBAL_elapsed_time); 00396 inverse_mpolar_fft(f, T, S, f_hat, N, 2 * m, m); 00397 if (m == 9) 00398 fprintf(fp, "%1.1" NFFT__FES__ "\\\\\\hline\n", GLOBAL_elapsed_time); 00399 else 00400 fprintf(fp, "%1.1" NFFT__FES__ "\\\\\n", GLOBAL_elapsed_time); 00401 printf("t_impolar=%1.1" NFFT__FES__ "\n", GLOBAL_elapsed_time); 00402 } 00403 00404 fflush(fp); 00405 00406 NFFT(free)(f); 00407 NFFT(free)(f_hat); 00408 00409 return EXIT_SUCCESS; 00410 } 00411 00413 int main(int argc, char **argv) 00414 { 00415 int N; 00416 int T, S; 00417 int M; 00418 NFFT_R *x, *w; 00419 NFFT_C *f_hat, *f, *f_direct, *f_tilde; 00420 int k; 00421 int max_i; 00422 int m; 00423 NFFT_R temp1, temp2, E_max = NFFT_K(0.0); 00424 FILE *fp1, *fp2; 00425 char filename[30]; 00426 int logN; 00427 00428 if (argc != 4) 00429 { 00430 printf("mpolar_fft_test N T R \n"); 00431 printf("\n"); 00432 printf("N mpolar FFT of size NxN \n"); 00433 printf("T number of slopes \n"); 00434 printf("R number of offsets \n"); 00435 00437 printf("\nHence, comparison FFTW, mpolar FFT and inverse mpolar FFT\n"); 00438 fp1 = fopen("mpolar_comparison_fft.dat", "w"); 00439 if (fp1 == NULL) 00440 return (-1); 00441 for (logN = 4; logN <= 8; logN++) 00442 comparison_fft(fp1, (int)(1U << logN), 3 * (int)(1U << logN), 00443 3 * (int)(1U << (logN - 1))); 00444 fclose(fp1); 00445 00446 exit(EXIT_FAILURE); 00447 } 00448 00449 N = atoi(argv[1]); 00450 T = atoi(argv[2]); 00451 S = atoi(argv[3]); 00452 printf("N=%d, modified polar grid with T=%d, R=%d => ", N, T, S); 00453 00454 x = (NFFT_R *) NFFT(malloc)((size_t)(5 * T * S / 2) * (sizeof(NFFT_R))); 00455 w = (NFFT_R *) NFFT(malloc)((size_t)(5 * T * S / 4) * (sizeof(NFFT_R))); 00456 00457 f_hat = (NFFT_C *) NFFT(malloc)(sizeof(NFFT_C) * (size_t)(N * N)); 00458 f = (NFFT_C *) NFFT(malloc)(sizeof(NFFT_C) * (size_t)(1.25 * T * S)); /* 4/pi*log(1+sqrt(2)) = 1.122... < 1.25 */ 00459 f_direct = (NFFT_C *) NFFT(malloc)(sizeof(NFFT_C) * (size_t)(1.25 * T * S)); 00460 f_tilde = (NFFT_C *) NFFT(malloc)(sizeof(NFFT_C) * (size_t)(N * N)); 00461 00463 M = mpolar_grid(T, S, x, w); 00464 printf("M=%d.\n", M); 00465 00467 fp1 = fopen("input_data_r.dat", "r"); 00468 fp2 = fopen("input_data_i.dat", "r"); 00469 if ((fp1 == NULL) || (fp2 == NULL)) 00470 return (-1); 00471 for (k = 0; k < N * N; k++) 00472 { 00473 fscanf(fp1, NFFT__FR__ " ", &temp1); 00474 fscanf(fp2, NFFT__FR__ " ", &temp2); 00475 f_hat[k] = temp1 + _Complex_I * temp2; 00476 } 00477 fclose(fp1); 00478 fclose(fp2); 00479 00481 mpolar_dft(f_hat, N, f_direct, T, S, 1); 00482 // mpolar_fft(f_hat,N,f_direct,T,R,12); 00483 00485 printf("\nTest of the mpolar FFT: \n"); 00486 fp1 = fopen("mpolar_fft_error.dat", "w+"); 00487 for (m = 1; m <= 12; m++) 00488 { 00490 mpolar_fft(f_hat, N, f, T, S, m); 00491 00493 E_max = NFFT(error_l_infty_complex)(f_direct, f, M); 00494 printf("m=%2d: E_max = %" NFFT__FES__ "\n", m, E_max); 00495 fprintf(fp1, "%" NFFT__FES__ "\n", E_max); 00496 } 00497 fclose(fp1); 00498 00500 for (m = 3; m <= 9; m += 3) 00501 { 00502 printf("\nTest of the inverse mpolar FFT for m=%d: \n", m); 00503 sprintf(filename, "mpolar_ifft_error%d.dat", m); 00504 fp1 = fopen(filename, "w+"); 00505 for (max_i = 0; max_i <= 20; max_i += 2) 00506 { 00508 inverse_mpolar_fft(f_direct, T, S, f_tilde, N, max_i, m); 00509 00511 E_max = NFFT(error_l_infty_complex)(f_hat, f_tilde, N * N); 00512 printf("%3d iterations: E_max = %" NFFT__FES__ "\n", max_i, E_max); 00513 fprintf(fp1, "%" NFFT__FES__ "\n", E_max); 00514 } 00515 fclose(fp1); 00516 } 00517 00519 NFFT(free)(x); 00520 NFFT(free)(w); 00521 NFFT(free)(f_hat); 00522 NFFT(free)(f); 00523 NFFT(free)(f_direct); 00524 NFFT(free)(f_tilde); 00525 00526 return EXIT_SUCCESS; 00527 } 00528 /* \} */