/** computes the inverse discrete Radon transform of Rf * on the grid given by gridfcn() with T angles and R offsets * by a NFFT-based CG-type algorithm */ int Inverse_Radon_trafo(int (*gridfcn)(), int T, int S, double *Rf, int NN, double *f, int max_i) { 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 */ fftw_complex *fft; /**< variable for the fftw-1Ds */ fftw_plan my_fftw_plan; /**< plan for the fftw-1Ds */ int t,r; /**< index for directions and offsets */ double *x, *w; /**< knots and associated weights */ int l; /**< index for iterations */ int N[2],n[2]; int M=T*S; N[0]=NN; n[0]=2*N[0]; N[1]=NN; n[1]=2*N[1]; fft = (fftw_complex *)nfft_malloc(S*sizeof(fftw_complex)); my_fftw_plan = fftw_plan_dft_1d(S,fft,fft,FFTW_FORWARD,FFTW_MEASURE); x = (double *)nfft_malloc(2*T*S*(sizeof(double))); if (x==NULL) return -1; w = (double *)nfft_malloc(T*S*(sizeof(double))); if (w==NULL) return -1; /** init two dimensional NFFT plan */ nfft_init_guru(&my_nfft_plan, 2, N, M, n, 4, 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 and weights of grid*/ gridfcn(T,S,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]; if (j%S) my_infft_plan.w[j] = w[j]; else my_infft_plan.w[j] = 0.0; } /** precompute psi, the entries of the matrix B */ if(my_nfft_plan.nfft_flags & PRE_LIN_PSI) nfft_precompute_lin_psi(&my_nfft_plan); if(my_nfft_plan.nfft_flags & PRE_PSI) nfft_precompute_psi(&my_nfft_plan); if(my_nfft_plan.nfft_flags & PRE_FULL_PSI) nfft_precompute_full_psi(&my_nfft_plan); /** compute 1D-ffts and init given samples and weights */ for(t=0; t<T; t++) { /* for(r=0; r<R/2; r++) fft[r] = cexp(I*KPI*r)*Rf[t*R+(r+R/2)]; for(r=0; r<R/2; r++) fft[r+R/2] = cexp(I*KPI*r)*Rf[t*R+r]; */ for(r=0; r<S; r++) fft[r] = Rf[t*S+r] + _Complex_I*0.0; nfft_fftshift_complex(fft, 1, &S); fftw_execute(my_fftw_plan); nfft_fftshift_complex(fft, 1, &S); my_infft_plan.y[t*S] = 0.0; for(r=-S/2+1; r<S/2; r++) my_infft_plan.y[t*S+(r+S/2)] = fft[r+S/2]/KERNEL(r); } /** initialise some guess f_hat_0 */ for(k=0;k<my_nfft_plan.N_total;k++) my_infft_plan.f_hat_iter[k] = 0.0 + _Complex_I*0.0; /** 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); /*if (sqrt(my_infft_plan.dot_r_iter)<=1e-12) break;*/ } } /*printf("after %d iteration(s): weighted 2-norm of original residual vector = %g\n",l-1,sqrt(my_infft_plan.dot_r_iter));*/ /** copy result */ for(k=0;k<my_nfft_plan.N_total;k++) f[k] = creal(my_infft_plan.f_hat_iter[k]); /** finalise the plans and free the variables */ fftw_destroy_plan(my_fftw_plan); nfft_free(fft); solver_finalize_complex(&my_infft_plan); nfft_finalize(&my_nfft_plan); nfft_free(x); nfft_free(w); return 0; }
/** computes the inverse discrete Radon transform of Rf * on the grid given by gridfcn() with T angles and R offsets * by a NFFT-based CG-type algorithm */ static int inverse_radon_trafo(int (*gridfcn)(), int T, int S, NFFT_R *Rf, int NN, NFFT_R *f, int max_i) { 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_C *fft; /**< variable for the fftw-1Ds */ FFTW(plan) my_fftw_plan; /**< plan for the fftw-1Ds */ int t, r; /**< index for directions and offsets */ NFFT_R *x, *w; /**< knots and associated weights */ int l; /**< index for iterations */ int N[2], n[2]; int M = T * S; N[0] = NN; n[0] = 2 * N[0]; N[1] = NN; n[1] = 2 * N[1]; fft = (NFFT_C *) NFFT(malloc)((size_t)(S) * sizeof(NFFT_C)); my_fftw_plan = FFTW(plan_dft_1d)(S, fft, fft, FFTW_FORWARD, FFTW_MEASURE); x = (NFFT_R *) NFFT(malloc)((size_t)(2 * T * S) * (sizeof(NFFT_R))); if (x == NULL) return EXIT_FAILURE; w = (NFFT_R *) NFFT(malloc)((size_t)(T * S) * (sizeof(NFFT_R))); if (w == NULL) return EXIT_FAILURE; /** init two dimensional NFFT plan */ NFFT(init_guru)(&my_nfft_plan, 2, N, M, n, 4, 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 and weights of grid*/ gridfcn(T, S, 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]; if (j % S) my_infft_plan.w[j] = w[j]; else my_infft_plan.w[j] = NFFT_K(0.0); } /** 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); /** compute 1D-ffts and init given samples and weights */ for (t = 0; t < T; t++) { /* for(r=0; r<R/2; r++) fft[r] = cexp(I*NFFT_KPI*r)*Rf[t*R+(r+R/2)]; for(r=0; r<R/2; r++) fft[r+R/2] = cexp(I*NFFT_KPI*r)*Rf[t*R+r]; */ for (r = 0; r < S; r++) fft[r] = Rf[t * S + r] + _Complex_I * NFFT_K(0.0); NFFT(fftshift_complex_int)(fft, 1, &S); FFTW(execute)(my_fftw_plan); NFFT(fftshift_complex_int)(fft, 1, &S); my_infft_plan.y[t * S] = NFFT_K(0.0); for (r = -S / 2 + 1; r < S / 2; r++) my_infft_plan.y[t * S + (r + S / 2)] = fft[r + S / 2] / KERNEL(r); } /** 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); /** 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); /*if (sqrt(my_infft_plan.dot_r_iter)<=1e-12) break;*/ } } /*printf("after %d iteration(s): weighted 2-norm of original residual vector = %g\n",l-1,sqrt(my_infft_plan.dot_r_iter));*/ /** copy result */ for (k = 0; k < my_nfft_plan.N_total; k++) f[k] = NFFT_M(creal)(my_infft_plan.f_hat_iter[k]); /** finalise the plans and free the variables */ FFTW(destroy_plan)(my_fftw_plan); NFFT(free)(fft); SOLVER(finalize_complex)(&my_infft_plan); NFFT(finalize)(&my_nfft_plan); NFFT(free)(x); NFFT(free)(w); return 0; }