/** * Initialisation of a transform plan, guru. * * \arg ths The pointer to a fpt plan * \arg N The number of source nodes * \arg M The number of target nodes * \arg sigma The parameter of the Gaussian * \arg n The polynomial expansion degree * \arg p the periodisation length, at least 1 * \arg m The spatial cut-off of the nfft * \arg flags FGT flags to use * * \author Stefan Kunis */ void fgt_init_guru(fgt_plan *ths, int N, int M, double _Complex sigma, int n, double p, int m, unsigned flags) { int j,n_fftw; fftw_plan fplan; ths->M = M; ths->N = N; ths->sigma = sigma; ths->flags = flags; ths->x = (double*)nfft_malloc(ths->N*sizeof(double)); ths->y = (double*)nfft_malloc(ths->M*sizeof(double)); ths->alpha = (double _Complex*)nfft_malloc(ths->N*sizeof(double _Complex)); ths->f = (double _Complex*)nfft_malloc(ths->M*sizeof(double _Complex)); ths->n = n; ths->p = p; ths->b = (double _Complex*)nfft_malloc(ths->n*sizeof(double _Complex)); ths->nplan1 = (nfft_plan*) nfft_malloc(sizeof(nfft_plan)); ths->nplan2 = (nfft_plan*) nfft_malloc(sizeof(nfft_plan)); n_fftw=X(next_power_of_2)(2*ths->n); nfft_init_guru(ths->nplan1, 1, &(ths->n), ths->N, &n_fftw, m, PRE_PHI_HUT| PRE_PSI| MALLOC_X| MALLOC_F_HAT| FFTW_INIT, FFTW_MEASURE); nfft_init_guru(ths->nplan2, 1, &(ths->n), ths->M, &n_fftw, m, PRE_PHI_HUT| PRE_PSI| MALLOC_X| FFTW_INIT, FFTW_MEASURE); ths->nplan1->f = ths->alpha; ths->nplan2->f_hat = ths->nplan1->f_hat; ths->nplan2->f = ths->f; if(ths->flags & FGT_APPROX_B) { fplan = fftw_plan_dft_1d(ths->n, ths->b, ths->b, FFTW_FORWARD, FFTW_MEASURE); for(j=0; j<ths->n; j++) ths->b[j] = cexp(-ths->p*ths->p*ths->sigma*(j-ths->n/2)*(j-ths->n/2)/ ((double)ths->n*ths->n)) / ths->n; nfft_fftshift_complex(ths->b, 1, &ths->n); fftw_execute(fplan); nfft_fftshift_complex(ths->b, 1, &ths->n); fftw_destroy_plan(fplan); } else { for(j=0; j<ths->n; j++) ths->b[j] = 1.0/ths->p * csqrt(PI/ths->sigma)* cexp(-PI*PI*(j-ths->n/2)*(j-ths->n/2)/ (ths->p*ths->p*ths->sigma)); } }
/** precomputation for fastsum */ void fastsum_precompute(fastsum_plan *ths) { int j,k,t; int n_total; ticks t0, t1; ths->MEASURE_TIME_t[0] = 0.0; ths->MEASURE_TIME_t[1] = 0.0; ths->MEASURE_TIME_t[2] = 0.0; ths->MEASURE_TIME_t[3] = 0.0; #ifdef MEASURE_TIME t0 = getticks(); #endif if (ths->flags & NEARFIELD_BOXES) { BuildBox(ths); } else { /** sort source knots */ BuildTree(ths->d,0,ths->x,ths->alpha,ths->N_total); } #ifdef MEASURE_TIME t1 = getticks(); ths->MEASURE_TIME_t[3] += nfft_elapsed_seconds(t1,t0); #endif #ifdef MEASURE_TIME t0 = getticks(); #endif /** precompute spline values for near field*/ if (!(ths->flags & EXACT_NEARFIELD)) { if (ths->d==1) #pragma omp parallel for default(shared) private(k) for (k=-ths->Ad/2-2; k <= ths->Ad/2+2; k++) ths->Add[k+ths->Ad/2+2] = regkern1(ths->k, ths->eps_I*(double)k/ths->Ad*2, ths->p, ths->kernel_param, ths->eps_I, ths->eps_B); else #pragma omp parallel for default(shared) private(k) for (k=0; k <= ths->Ad+2; k++) ths->Add[k] = regkern3(ths->k, ths->eps_I*(double)k/ths->Ad, ths->p, ths->kernel_param, ths->eps_I, ths->eps_B); } #ifdef MEASURE_TIME t1 = getticks(); ths->MEASURE_TIME_t[0] += nfft_elapsed_seconds(t1,t0); #endif #ifdef MEASURE_TIME t0 = getticks(); #endif /** init NFFT plan for transposed transform in first step*/ for (k=0; k<ths->mv1.M_total; k++) for (t=0; t<ths->mv1.d; t++) ths->mv1.x[ths->mv1.d*k+t] = - ths->x[ths->mv1.d*k+t]; /* note the factor -1 for transposed transform instead of adjoint*/ /** precompute psi, the entries of the matrix B */ if(ths->mv1.nfft_flags & PRE_LIN_PSI) nfft_precompute_lin_psi(&(ths->mv1)); if(ths->mv1.nfft_flags & PRE_PSI) nfft_precompute_psi(&(ths->mv1)); if(ths->mv1.nfft_flags & PRE_FULL_PSI) nfft_precompute_full_psi(&(ths->mv1)); #ifdef MEASURE_TIME t1 = getticks(); ths->MEASURE_TIME_t[1] += nfft_elapsed_seconds(t1,t0); #endif /** init Fourier coefficients */ for(k=0; k<ths->mv1.M_total;k++) ths->mv1.f[k] = ths->alpha[k]; #ifdef MEASURE_TIME t0 = getticks(); #endif /** init NFFT plan for transform in third step*/ for (j=0; j<ths->mv2.M_total; j++) for (t=0; t<ths->mv2.d; t++) ths->mv2.x[ths->mv2.d*j+t] = - ths->y[ths->mv2.d*j+t]; /* note the factor -1 for conjugated transform instead of standard*/ /** precompute psi, the entries of the matrix B */ if(ths->mv2.nfft_flags & PRE_LIN_PSI) nfft_precompute_lin_psi(&(ths->mv2)); if(ths->mv2.nfft_flags & PRE_PSI) nfft_precompute_psi(&(ths->mv2)); if(ths->mv2.nfft_flags & PRE_FULL_PSI) nfft_precompute_full_psi(&(ths->mv2)); #ifdef MEASURE_TIME t1 = getticks(); ths->MEASURE_TIME_t[2] += nfft_elapsed_seconds(t1,t0); #endif #ifdef MEASURE_TIME t0 = getticks(); #endif /** precompute Fourier coefficients of regularised kernel*/ n_total = 1; for (t=0; t<ths->d; t++) n_total *= ths->n; #pragma omp parallel for default(shared) private(j,k,t) for (j=0; j<n_total; j++) { if (ths->d==1) ths->b[j] = regkern1(ths->k, (double)j / (ths->n) - 0.5, ths->p, ths->kernel_param, ths->eps_I, ths->eps_B)/n_total; else { k=j; ths->b[j]=0.0; for (t=0; t<ths->d; t++) { ths->b[j] += ((double)(k % (ths->n)) / (ths->n) - 0.5) * ((double)(k % (ths->n)) / (ths->n) - 0.5); k = k / (ths->n); } ths->b[j] = regkern3(ths->k, sqrt(ths->b[j]), ths->p, ths->kernel_param, ths->eps_I, ths->eps_B)/n_total; } } nfft_fftshift_complex(ths->b, ths->mv1.d, ths->mv1.N); fftw_execute(ths->fft_plan); nfft_fftshift_complex(ths->b, ths->mv1.d, ths->mv1.N); #ifdef MEASURE_TIME t1 = getticks(); ths->MEASURE_TIME_t[0] += nfft_elapsed_seconds(t1,t0); #endif }
/** 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; }