/** Reconstruction routine */ static void glacier(int N,int M) { int j,k,k0,k1,l,my_N[2],my_n[2]; double tmp_y; nfft_plan p; solver_plan_complex ip; FILE* fp; /* initialise p */ my_N[0]=N; my_n[0]=X(next_power_of_2)(N); my_N[1]=N; my_n[1]=X(next_power_of_2)(N); nfft_init_guru(&p, 2, my_N, M, my_n, 6, PRE_PHI_HUT| PRE_FULL_PSI| MALLOC_X| MALLOC_F_HAT| MALLOC_F| FFTW_INIT| FFT_OUT_OF_PLACE, FFTW_MEASURE| FFTW_DESTROY_INPUT); /* initialise ip, specific */ solver_init_advanced_complex(&ip,(nfft_mv_plan_complex*)(&p), CGNE| PRECOMPUTE_DAMP); fprintf(stderr,"Using the generic solver!"); /* init nodes */ fp=fopen("input_data.dat","r"); for(j=0;j<p.M_total;j++) { fscanf(fp,"%le %le %le",&p.x[2*j+0],&p.x[2*j+1],&tmp_y); ip.y[j]=tmp_y; } fclose(fp); /* precompute psi */ if(p.nfft_flags & PRE_ONE_PSI) nfft_precompute_one_psi(&p); /* initialise damping factors */ if(ip.flags & PRECOMPUTE_DAMP) for(k0=0;k0<p.N[0];k0++) for(k1=0;k1<p.N[1];k1++) ip.w_hat[k0*p.N[1]+k1]= my_weight(((double)(k0-p.N[0]/2))/p.N[0],0.5,3,0.001)* my_weight(((double)(k1-p.N[1]/2))/p.N[1],0.5,3,0.001); /* init some guess */ for(k=0;k<p.N_total;k++) ip.f_hat_iter[k]=0; /* inverse trafo */ solver_before_loop_complex(&ip); for(l=0;l<40;l++) { fprintf(stderr,"Residual ||r||=%e,\n",sqrt(ip.dot_r_iter)); solver_loop_one_step_complex(&ip); } for(k=0;k<p.N_total;k++) printf("%le %le\n",creal(ip.f_hat_iter[k]),cimag(ip.f_hat_iter[k])); solver_finalize_complex(&ip); nfft_finalize(&p); }
void ipdf_finalize(ipdf_plan *ths){ /* finalise the plans and free the variables */ solver_finalize_complex(&ths->iplan); nfsft_finalize(&ths->plan); }
/** * The main program. * * \param argc The number of arguments * \param argv An array containing the arguments as C-strings * * \return Exit code * * \author Jens Keiner */ int main (int argc, char **argv) { int T; int N; int M; int M2; int t; /* Index variable for testcases */ nfsft_plan plan; /* NFSFT plan */ nfsft_plan plan2; /* NFSFT plan */ solver_plan_complex iplan; /* NFSFT plan */ int j; /* */ int k; /* */ int m; /* */ int use_nfsft; /* */ int use_nfft; /* */ int use_fpt; /* */ int cutoff; /**< The current NFFT cut-off parameter */ double threshold; /**< The current NFSFT threshold parameter */ double re; double im; double a; double *scratch; double xs; double *ys; double *temp; double _Complex *temp2; int qlength; double *qweights; fftw_plan fplan; fpt_set set; int npt; int npt_exp; double *alpha, *beta, *gamma; /* Read the number of testcases. */ fscanf(stdin,"testcases=%d\n",&T); fprintf(stderr,"%d\n",T); /* Process each testcase. */ for (t = 0; t < T; t++) { /* Check if the fast transform shall be used. */ fscanf(stdin,"nfsft=%d\n",&use_nfsft); fprintf(stderr,"%d\n",use_nfsft); if (use_nfsft != NO) { /* Check if the NFFT shall be used. */ fscanf(stdin,"nfft=%d\n",&use_nfft); fprintf(stderr,"%d\n",use_nfsft); if (use_nfft != NO) { /* Read the cut-off parameter. */ fscanf(stdin,"cutoff=%d\n",&cutoff); fprintf(stderr,"%d\n",cutoff); } else { /* TODO remove this */ /* Initialize unused variable with dummy value. */ cutoff = 1; } /* Check if the fast polynomial transform shall be used. */ fscanf(stdin,"fpt=%d\n",&use_fpt); fprintf(stderr,"%d\n",use_fpt); if (use_fpt != NO) { /* Read the NFSFT threshold parameter. */ fscanf(stdin,"threshold=%lf\n",&threshold); fprintf(stderr,"%lf\n",threshold); } else { /* TODO remove this */ /* Initialize unused variable with dummy value. */ threshold = 1000.0; } } else { /* TODO remove this */ /* Set dummy values. */ use_nfft = NO; use_fpt = NO; cutoff = 3; threshold = 1000.0; } /* Read the bandwidth. */ fscanf(stdin,"bandwidth=%d\n",&N); fprintf(stderr,"%d\n",N); /* Do precomputation. */ nfsft_precompute(N,threshold, ((use_nfsft==NO)?(NFSFT_NO_FAST_ALGORITHM):(0U/*NFSFT_NO_DIRECT_ALGORITHM*/)), 0U); /* Read the number of nodes. */ fscanf(stdin,"nodes=%d\n",&M); fprintf(stderr,"%d\n",M); /* */ if ((N+1)*(N+1) > M) { X(next_power_of_2_exp)(N, &npt, &npt_exp); fprintf(stderr, "npt = %d, npt_exp = %d\n", npt, npt_exp); fprintf(stderr,"Optimal interpolation!\n"); scratch = (double*) nfft_malloc(4*sizeof(double)); ys = (double*) nfft_malloc((N+1)*sizeof(double)); temp = (double*) nfft_malloc((2*N+1)*sizeof(double)); temp2 = (double _Complex*) nfft_malloc((N+1)*sizeof(double _Complex)); a = 0.0; for (j = 0; j <= N; j++) { xs = 2.0 + (2.0*j)/(N+1); ys[j] = (2.0-((j == 0)?(1.0):(0.0)))*4.0*nfft_bspline(4,xs,scratch); //fprintf(stdout,"%3d: g(%le) = %le\n",j,xs,ys[j]); a += ys[j]; } //fprintf(stdout,"a = %le\n",a); for (j = 0; j <= N; j++) { ys[j] *= 1.0/a; } qlength = 2*N+1; qweights = (double*) nfft_malloc(qlength*sizeof(double)); fplan = fftw_plan_r2r_1d(N+1, qweights, qweights, FFTW_REDFT00, 0U); for (j = 0; j < N+1; j++) { qweights[j] = -2.0/(4*j*j-1); } fftw_execute(fplan); qweights[0] *= 0.5; for (j = 0; j < N+1; j++) { qweights[j] *= 1.0/(2.0*N+1.0); qweights[2*N+1-1-j] = qweights[j]; } fplan = fftw_plan_r2r_1d(2*N+1, temp, temp, FFTW_REDFT00, 0U); for (j = 0; j <= N; j++) { temp[j] = ((j==0 || j == 2*N)?(1.0):(0.5))*ys[j]; } for (j = N+1; j < 2*N+1; j++) { temp[j] = 0.0; } fftw_execute(fplan); for (j = 0; j < 2*N+1; j++) { temp[j] *= qweights[j]; } fftw_execute(fplan); for (j = 0; j < 2*N+1; j++) { temp[j] *= ((j==0 || j == 2*N)?(1.0):(0.5)); if (j <= N) { temp2[j] = temp[j]; } } set = fpt_init(1, npt_exp, 0U); alpha = (double*) nfft_malloc((N+2)*sizeof(double)); beta = (double*) nfft_malloc((N+2)*sizeof(double)); gamma = (double*) nfft_malloc((N+2)*sizeof(double)); alpha_al_row(alpha, N, 0); beta_al_row(beta, N, 0); gamma_al_row(gamma, N, 0); fpt_precompute(set, 0, alpha, beta, gamma, 0, 1000.0); fpt_transposed(set,0, temp2, temp2, N, 0U); fpt_finalize(set); nfft_free(alpha); nfft_free(beta); nfft_free(gamma); fftw_destroy_plan(fplan); nfft_free(scratch); nfft_free(qweights); nfft_free(ys); nfft_free(temp); } /* Init transform plans. */ nfsft_init_guru(&plan, N, M, ((use_nfft!=0)?(0U):(NFSFT_USE_NDFT)) | ((use_fpt!=0)?(0U):(NFSFT_USE_DPT)) | NFSFT_MALLOC_F | NFSFT_MALLOC_X | NFSFT_MALLOC_F_HAT | NFSFT_NORMALIZED | NFSFT_ZERO_F_HAT, PRE_PHI_HUT | PRE_PSI | FFTW_INIT | FFT_OUT_OF_PLACE, cutoff); if ((N+1)*(N+1) > M) { solver_init_advanced_complex(&iplan, (nfft_mv_plan_complex*)(&plan), CGNE | PRECOMPUTE_DAMP); } else { solver_init_advanced_complex(&iplan, (nfft_mv_plan_complex*)(&plan), CGNR | PRECOMPUTE_WEIGHT | PRECOMPUTE_DAMP); } /* Read the nodes and function values. */ for (j = 0; j < M; j++) { fscanf(stdin,"%le %le %le %le\n",&plan.x[2*j+1],&plan.x[2*j],&re,&im); plan.x[2*j+1] = plan.x[2*j+1]/(2.0*PI); plan.x[2*j] = plan.x[2*j]/(2.0*PI); if (plan.x[2*j] >= 0.5) { plan.x[2*j] = plan.x[2*j] - 1; } iplan.y[j] = re + _Complex_I * im; fprintf(stderr,"%le %le %le %le\n",plan.x[2*j+1],plan.x[2*j], creal(iplan.y[j]),cimag(iplan.y[j])); } /* Read the number of nodes. */ fscanf(stdin,"nodes_eval=%d\n",&M2); fprintf(stderr,"%d\n",M2); /* Init transform plans. */ nfsft_init_guru(&plan2, N, M2, ((use_nfft!=0)?(0U):(NFSFT_USE_NDFT)) | ((use_fpt!=0)?(0U):(NFSFT_USE_DPT)) | NFSFT_MALLOC_F | NFSFT_MALLOC_X | NFSFT_MALLOC_F_HAT | NFSFT_NORMALIZED | NFSFT_ZERO_F_HAT, PRE_PHI_HUT | PRE_PSI | FFTW_INIT | FFT_OUT_OF_PLACE, cutoff); /* Read the nodes and function values. */ for (j = 0; j < M2; j++) { fscanf(stdin,"%le %le\n",&plan2.x[2*j+1],&plan2.x[2*j]); plan2.x[2*j+1] = plan2.x[2*j+1]/(2.0*PI); plan2.x[2*j] = plan2.x[2*j]/(2.0*PI); if (plan2.x[2*j] >= 0.5) { plan2.x[2*j] = plan2.x[2*j] - 1; } fprintf(stderr,"%le %le\n",plan2.x[2*j+1],plan2.x[2*j]); } nfsft_precompute_x(&plan); nfsft_precompute_x(&plan2); /* Frequency weights. */ if ((N+1)*(N+1) > M) { /* Compute Voronoi weights. */ //nfft_voronoi_weights_S2(iplan.w, plan.x, M); /* Print out Voronoi weights. */ /*a = 0.0; for (j = 0; j < plan.M_total; j++) { fprintf(stderr,"%le\n",iplan.w[j]); a += iplan.w[j]; } fprintf(stderr,"sum = %le\n",a);*/ for (j = 0; j < plan.N_total; j++) { iplan.w_hat[j] = 0.0; } for (k = 0; k <= N; k++) { for (j = -k; j <= k; j++) { iplan.w_hat[NFSFT_INDEX(k,j,&plan)] = 1.0/(pow(k+1.0,2.0)); /*temp2[j]*/; } } } else { for (j = 0; j < plan.N_total; j++) { iplan.w_hat[j] = 0.0; } for (k = 0; k <= N; k++) { for (j = -k; j <= k; j++) { iplan.w_hat[NFSFT_INDEX(k,j,&plan)] = 1/(pow(k+1.0,2.5)); } } /* Compute Voronoi weights. */ nfft_voronoi_weights_S2(iplan.w, plan.x, M); /* Print out Voronoi weights. */ a = 0.0; for (j = 0; j < plan.M_total; j++) { fprintf(stderr,"%le\n",iplan.w[j]); a += iplan.w[j]; } fprintf(stderr,"sum = %le\n",a); } fprintf(stderr, "N_total = %d\n", plan.N_total); fprintf(stderr, "M_total = %d\n", plan.M_total); /* init some guess */ for (k = 0; k < plan.N_total; k++) { iplan.f_hat_iter[k] = 0.0; } /* inverse trafo */ solver_before_loop_complex(&iplan); /*for (k = 0; k < plan.M_total; k++) { printf("%le %le\n",creal(iplan.r_iter[k]),cimag(iplan.r_iter[k])); }*/ for (m = 0; m < 29; m++) { fprintf(stderr,"Residual ||r||=%e,\n",sqrt(iplan.dot_r_iter)); solver_loop_one_step_complex(&iplan); } /*NFFT_SWAP_complex(iplan.f_hat_iter, plan.f_hat); nfsft_trafo(&plan); NFFT_SWAP_complex(iplan.f_hat_iter, plan.f_hat); a = 0.0; b = 0.0; for (k = 0; k < plan.M_total; k++) { printf("%le %le %le\n",cabs(iplan.y[k]),cabs(plan.f[k]), cabs(iplan.y[k]-plan.f[k])); a += cabs(iplan.y[k]-plan.f[k])*cabs(iplan.y[k]-plan.f[k]); b += cabs(iplan.y[k])*cabs(iplan.y[k]); } fprintf(stderr,"relative error in 2-norm: %le\n",a/b);*/ NFFT_SWAP_complex(iplan.f_hat_iter, plan2.f_hat); nfsft_trafo(&plan2); NFFT_SWAP_complex(iplan.f_hat_iter, plan2.f_hat); for (k = 0; k < plan2.M_total; k++) { fprintf(stdout,"%le\n",cabs(plan2.f[k])); } solver_finalize_complex(&iplan); nfsft_finalize(&plan); nfsft_finalize(&plan2); /* Delete precomputed data. */ nfsft_forget(); if ((N+1)*(N+1) > M) { nfft_free(temp2); } } /* Process each testcase. */ /* Return exit code for successful run. */ return EXIT_SUCCESS; }
/** Reconstruction routine with cross validation */ static void glacier_cv(int N,int M,int M_cv,unsigned solver_flags) { int j,k,k0,k1,l,my_N[2],my_n[2]; double tmp_y,r; nfft_plan p,cp; solver_plan_complex ip; double _Complex* cp_y; FILE* fp; int M_re=M-M_cv; /* initialise p for reconstruction */ my_N[0]=N; my_n[0]=X(next_power_of_2)(N); my_N[1]=N; my_n[1]=X(next_power_of_2)(N); nfft_init_guru(&p, 2, my_N, M_re, my_n, 6, PRE_PHI_HUT| PRE_FULL_PSI| MALLOC_X| MALLOC_F_HAT| MALLOC_F| FFTW_INIT| FFT_OUT_OF_PLACE, FFTW_MEASURE| FFTW_DESTROY_INPUT); /* initialise ip, specific */ solver_init_advanced_complex(&ip,(nfft_mv_plan_complex*)(&p), solver_flags); /* initialise cp for validation */ cp_y = (double _Complex*) nfft_malloc(M*sizeof(double _Complex)); nfft_init_guru(&cp, 2, my_N, M, my_n, 6, PRE_PHI_HUT| PRE_FULL_PSI| MALLOC_X| MALLOC_F| FFTW_INIT| FFT_OUT_OF_PLACE, FFTW_MEASURE| FFTW_DESTROY_INPUT); cp.f_hat=ip.f_hat_iter; /* set up data in cp and cp_y */ fp=fopen("input_data.dat","r"); for(j=0;j<cp.M_total;j++) { fscanf(fp,"%le %le %le",&cp.x[2*j+0],&cp.x[2*j+1],&tmp_y); cp_y[j]=tmp_y; } fclose(fp); /* copy part of the data to p and ip */ for(j=0;j<p.M_total;j++) { p.x[2*j+0]=cp.x[2*j+0]; p.x[2*j+1]=cp.x[2*j+1]; ip.y[j]=tmp_y; } /* precompute psi */ if(p.nfft_flags & PRE_ONE_PSI) nfft_precompute_one_psi(&p); /* precompute psi */ if(cp.nfft_flags & PRE_ONE_PSI) nfft_precompute_one_psi(&cp); /* initialise damping factors */ if(ip.flags & PRECOMPUTE_DAMP) for(k0=0;k0<p.N[0];k0++) for(k1=0;k1<p.N[1];k1++) ip.w_hat[k0*p.N[1]+k1]= my_weight(((double)(k0-p.N[0]/2))/p.N[0],0.5,3,0.001)* my_weight(((double)(k1-p.N[1]/2))/p.N[1],0.5,3,0.001); /* init some guess */ for(k=0;k<p.N_total;k++) ip.f_hat_iter[k]=0; /* inverse trafo */ solver_before_loop_complex(&ip); // fprintf(stderr,"iteration starts,\t"); for(l=0;l<40;l++) solver_loop_one_step_complex(&ip); //fprintf(stderr,"r=%1.2e, ",sqrt(ip.dot_r_iter)/M_re); NFFT_SWAP_complex(p.f_hat,ip.f_hat_iter); nfft_trafo(&p); NFFT_SWAP_complex(p.f_hat,ip.f_hat_iter); nfft_upd_axpy_complex(p.f,-1,ip.y,M_re); r=sqrt(nfft_dot_complex(p.f,M_re)/nfft_dot_complex(cp_y,M)); fprintf(stderr,"r=%1.2e, ",r); printf("$%1.1e$ & ",r); nfft_trafo(&cp); nfft_upd_axpy_complex(&cp.f[M_re],-1,&cp_y[M_re],M_cv); r=sqrt(nfft_dot_complex(&cp.f[M_re],M_cv)/nfft_dot_complex(cp_y,M)); fprintf(stderr,"r_1=%1.2e\t",r); printf("$%1.1e$ & ",r); nfft_finalize(&cp); solver_finalize_complex(&ip); nfft_finalize(&p); }
static void reconstruct(char* filename,int N,int M,int iteration , int weight) { int j,k,l; double t0, t1; double time,min_time,max_time,min_inh,max_inh; double t,real,imag; double w,epsilon=0.0000003; /* epsilon is a the break criterium for the iteration */; mri_inh_3d_plan my_plan; solver_plan_complex my_iplan; FILE* fp,*fw,*fout_real,*fout_imag,*finh,*ftime; int my_N[3],my_n[3]; int flags = PRE_PHI_HUT| PRE_PSI |MALLOC_X| MALLOC_F_HAT| MALLOC_F| FFTW_INIT| FFT_OUT_OF_PLACE; unsigned infft_flags = CGNR | PRECOMPUTE_DAMP; double Ts; double W; int N3; int m=2; double sigma = 1.25; ftime=fopen("readout_time.dat","r"); finh=fopen("inh.dat","r"); min_time=INT_MAX; max_time=INT_MIN; for(j=0;j<M;j++) { fscanf(ftime,"%le ",&time); if(time<min_time) min_time = time; if(time>max_time) max_time = time; } fclose(ftime); Ts=(min_time+max_time)/2.0; min_inh=INT_MAX; max_inh=INT_MIN; for(j=0;j<N*N;j++) { fscanf(finh,"%le ",&w); if(w<min_inh) min_inh = w; if(w>max_inh) max_inh = w; } fclose(finh); N3=ceil((MAX(fabs(min_inh),fabs(max_inh))*(max_time-min_time)/2.0+m/(2*sigma))*4*sigma); /* N3 has to be even */ if(N3%2!=0) N3++; W= MAX(fabs(min_inh),fabs(max_inh))/(0.5-((double) m)/N3); my_N[0]=N;my_n[0]=ceil(N*sigma); my_N[1]=N; my_n[1]=ceil(N*sigma); my_N[2]=N3; my_n[2]=ceil(N3*sigma); /* initialise nfft */ mri_inh_3d_init_guru(&my_plan, my_N, M, my_n, m, sigma, flags, FFTW_MEASURE| FFTW_DESTROY_INPUT); if (weight) infft_flags = infft_flags | PRECOMPUTE_WEIGHT; /* initialise my_iplan, advanced */ solver_init_advanced_complex(&my_iplan,(nfft_mv_plan_complex*)(&my_plan), infft_flags ); /* get the weights */ if(my_iplan.flags & PRECOMPUTE_WEIGHT) { fw=fopen("weights.dat","r"); for(j=0;j<my_plan.M_total;j++) { fscanf(fw,"%le ",&my_iplan.w[j]); } fclose(fw); } /* get the damping factors */ if(my_iplan.flags & PRECOMPUTE_DAMP) { for(j=0;j<N;j++){ for(k=0;k<N;k++) { int j2= j-N/2; int k2= k-N/2; double r=sqrt(j2*j2+k2*k2); if(r>(double) N/2) my_iplan.w_hat[j*N+k]=0.0; else my_iplan.w_hat[j*N+k]=1.0; } } } fp=fopen(filename,"r"); ftime=fopen("readout_time.dat","r"); for(j=0;j<my_plan.M_total;j++) { fscanf(fp,"%le %le %le %le",&my_plan.plan.x[3*j+0],&my_plan.plan.x[3*j+1],&real,&imag); my_iplan.y[j]=real+ _Complex_I*imag; fscanf(ftime,"%le ",&my_plan.plan.x[3*j+2]); my_plan.plan.x[3*j+2] = (my_plan.plan.x[3*j+2]-Ts)*W/N3; } fclose(fp); fclose(ftime); finh=fopen("inh.dat","r"); for(j=0;j<N*N;j++) { fscanf(finh,"%le ",&my_plan.w[j]); my_plan.w[j]/=W; } fclose(finh); if(my_plan.plan.flags & PRE_PSI) { nfft_precompute_psi(&my_plan.plan); } if(my_plan.plan.flags & PRE_FULL_PSI) { nfft_precompute_full_psi(&my_plan.plan); } /* init some guess */ for(j=0;j<my_plan.N_total;j++) { my_iplan.f_hat_iter[j]=0.0; } t0 = nfft_clock_gettime_seconds(); /* inverse trafo */ solver_before_loop_complex(&my_iplan); for(l=0;l<iteration;l++) { /* break if dot_r_iter is smaller than epsilon*/ if(my_iplan.dot_r_iter<epsilon) break; fprintf(stderr,"%e, %i of %i\n",sqrt(my_iplan.dot_r_iter), l+1,iteration); solver_loop_one_step_complex(&my_iplan); } t1 = nfft_clock_gettime_seconds(); t = t1-t0; fout_real=fopen("output_real.dat","w"); fout_imag=fopen("output_imag.dat","w"); for (j=0;j<N*N;j++) { /* Verschiebung wieder herausrechnen */ my_iplan.f_hat_iter[j]*=cexp(-2.0*_Complex_I*M_PI*Ts*my_plan.w[j]*W); fprintf(fout_real,"%le ",creal(my_iplan.f_hat_iter[j])); fprintf(fout_imag,"%le ",cimag(my_iplan.f_hat_iter[j])); } fclose(fout_real); fclose(fout_imag); solver_finalize_complex(&my_iplan); mri_inh_3d_finalize(&my_plan); }
/** * reconstruct makes an inverse 3d-nfft */ static void reconstruct(char* filename,int N,int M,int Z,int iteration, int weight) { int j,k,z,l; /* some variables */ double real,imag; /* to read the real and imag part of a complex number */ nfft_plan my_plan; /* plan for the two dimensional nfft */ solver_plan_complex my_iplan; /* plan for the two dimensional infft */ FILE* fin; /* input file */ FILE* fout_real; /* output file (real part) */ FILE* fout_imag; /* output file (imag part) */ int my_N[3],my_n[3]; /* to init the nfft */ double epsilon=0.0000003; /* tmp to read the obsolent z from 700.acs epsilon is a the break criterion for the iteration */ unsigned infft_flags = CGNR | PRECOMPUTE_DAMP; /* flags for the infft */ /* initialise my_plan, specific. we don't precompute psi */ my_N[0]=Z; my_n[0]=ceil(Z*1.2); my_N[1]=N; my_n[1]=ceil(N*1.2); my_N[2]=N; my_n[2]=ceil(N*1.2); nfft_init_guru(&my_plan, 3, my_N, M, my_n, 6, PRE_PHI_HUT| PRE_PSI |MALLOC_X| MALLOC_F_HAT| MALLOC_F| FFTW_INIT| FFT_OUT_OF_PLACE, FFTW_MEASURE| FFTW_DESTROY_INPUT); /* precompute lin psi */ if(my_plan.nfft_flags & PRE_LIN_PSI) nfft_precompute_lin_psi(&my_plan); if (weight) infft_flags = infft_flags | PRECOMPUTE_WEIGHT; /* initialise my_iplan, advanced */ solver_init_advanced_complex(&my_iplan,(nfft_mv_plan_complex*)(&my_plan), infft_flags ); /* get the weights */ if(my_iplan.flags & PRECOMPUTE_WEIGHT) { fin=fopen("weights.dat","r"); for(j=0;j<M;j++) { fscanf(fin,"%le ",&my_iplan.w[j]); } fclose(fin); } /* get the damping factors */ if(my_iplan.flags & PRECOMPUTE_DAMP) { for(j=0;j<N;j++){ for(k=0;k<N;k++) { for(z=0;z<N;z++) { int j2= j-N/2; int k2= k-N/2; int z2= z-N/2; double r=sqrt(j2*j2+k2*k2+z2*z2); if(r>(double) N/2) my_iplan.w_hat[z*N*N+j*N+k]=0.0; else my_iplan.w_hat[z*N*N+j*N+k]=1.0; } } } } /* open the input file */ fin=fopen(filename,"r"); /* open the output files */ fout_real=fopen("output_real.dat","w"); fout_imag=fopen("output_imag.dat","w"); /* read x,y,freal and fimag from the knots */ for(j=0;j<M;j++) { fscanf(fin,"%le %le %le %le %le ",&my_plan.x[3*j+1],&my_plan.x[3*j+2], &my_plan.x[3*j+0], &real,&imag); my_iplan.y[j] = real + _Complex_I*imag; } /* precompute psi */ if(my_plan.nfft_flags & PRE_PSI) nfft_precompute_psi(&my_plan); /* precompute full psi */ if(my_plan.nfft_flags & PRE_FULL_PSI) nfft_precompute_full_psi(&my_plan); /* init some guess */ for(k=0;k<my_plan.N_total;k++) my_iplan.f_hat_iter[k]=0.0; /* inverse trafo */ solver_before_loop_complex(&my_iplan); for(l=0;l<iteration;l++) { /* break if dot_r_iter is smaller than epsilon*/ if(my_iplan.dot_r_iter<epsilon) break; fprintf(stderr,"%e, %i of %i\n",sqrt(my_iplan.dot_r_iter), l+1,iteration); solver_loop_one_step_complex(&my_iplan); } for(l=0;l<Z;l++) { for(k=0;k<N*N;k++) { /* write every Layer in the files */ fprintf(fout_real,"%le ",creal(my_iplan.f_hat_iter[ k+N*N*l ])); fprintf(fout_imag,"%le ",cimag(my_iplan.f_hat_iter[ k+N*N*l ])); } fprintf(fout_real,"\n"); fprintf(fout_imag,"\n"); } fclose(fout_real); fclose(fout_imag); solver_finalize_complex(&my_iplan); nfft_finalize(&my_plan); }
/** * reconstruct */ static void reconstruct(char* filename,int N,int M,int iteration, int weight) { int j,k,l; /* some variables */ nnfft_plan my_plan; /* plan for the two dimensional nfft */ solver_plan_complex my_iplan; /* plan for the two dimensional infft */ FILE* fin; /* input file */ FILE* finh; FILE* ftime; FILE* fout_real; /* output file */ FILE* fout_imag; /* output file */ int my_N[3],my_n[3]; /* to init the nfft */ double t0, t1; double t,epsilon=0.0000003; /* epsilon is a the break criterium for the iteration */ unsigned infft_flags = CGNR | PRECOMPUTE_DAMP; /* flags for the infft*/ double time,min_time,max_time,min_inh,max_inh; double real,imag; double *w; double Ts; double W; int N3; int m=2; double sigma = 1.25; w = (double*)nfft_malloc(N*N*sizeof(double)); ftime=fopen("readout_time.dat","r"); finh=fopen("inh.dat","r"); min_time=INT_MAX; max_time=INT_MIN; for(j=0;j<M;j++) { fscanf(ftime,"%le ",&time); if(time<min_time) min_time = time; if(time>max_time) max_time = time; } fclose(ftime); Ts=(min_time+max_time)/2.0; min_inh=INT_MAX; max_inh=INT_MIN; for(j=0;j<N*N;j++) { fscanf(finh,"%le ",&w[j]); if(w[j]<min_inh) min_inh = w[j]; if(w[j]>max_inh) max_inh = w[j]; } fclose(finh); N3=ceil((MAX(fabs(min_inh),fabs(max_inh))*(max_time-min_time)/2.0)*4); W=MAX(fabs(min_inh),fabs(max_inh))*2.0; fprintf(stderr,"3: %i %e %e %e %e %e %e\n",N3,W,min_inh,max_inh,min_time,max_time,Ts); /* initialise my_plan */ my_N[0]=N;my_n[0]=ceil(N*sigma); my_N[1]=N; my_n[1]=ceil(N*sigma); my_N[2]=N3; my_n[2]=ceil(N3*sigma); nnfft_init_guru(&my_plan, 3, N*N, M, my_N,my_n,m, PRE_PSI| PRE_PHI_HUT| MALLOC_X| MALLOC_V| MALLOC_F_HAT| MALLOC_F ); /* precompute lin psi if set */ if(my_plan.nnfft_flags & PRE_LIN_PSI) nnfft_precompute_lin_psi(&my_plan); /* set the flags for the infft*/ if (weight) infft_flags = infft_flags | PRECOMPUTE_WEIGHT; /* initialise my_iplan, advanced */ solver_init_advanced_complex(&my_iplan,(nfft_mv_plan_complex*)(&my_plan), infft_flags ); /* get the weights */ if(my_iplan.flags & PRECOMPUTE_WEIGHT) { fin=fopen("weights.dat","r"); for(j=0;j<my_plan.M_total;j++) { fscanf(fin,"%le ",&my_iplan.w[j]); } fclose(fin); } /* get the damping factors */ if(my_iplan.flags & PRECOMPUTE_DAMP) { for(j=0;j<N;j++){ for(k=0;k<N;k++) { int j2= j-N/2; int k2= k-N/2; double r=sqrt(j2*j2+k2*k2); if(r>(double) N/2) my_iplan.w_hat[j*N+k]=0.0; else my_iplan.w_hat[j*N+k]=1.0; } } } /* open the input file */ fin=fopen(filename,"r"); ftime=fopen("readout_time.dat","r"); for(j=0;j<my_plan.M_total;j++) { fscanf(fin,"%le %le %le %le ",&my_plan.x[3*j+0],&my_plan.x[3*j+1],&real,&imag); my_iplan.y[j]=real+ _Complex_I*imag; fscanf(ftime,"%le ",&my_plan.x[3*j+2]); my_plan.x[3*j+2] = (my_plan.x[3*j+2]-Ts)*W/N3; } for(j=0;j<N;j++) { for(l=0;l<N;l++) { my_plan.v[3*(N*j+l)+0]=(((double) j) -(((double) N)/2.0))/((double) N); my_plan.v[3*(N*j+l)+1]=(((double) l) -(((double) N)/2.0))/((double) N); my_plan.v[3*(N*j+l)+2] = w[N*j+l]/W ; } } /* precompute psi */ if(my_plan.nnfft_flags & PRE_PSI) { nnfft_precompute_psi(&my_plan); if(my_plan.nnfft_flags & PRE_FULL_PSI) nnfft_precompute_full_psi(&my_plan); } if(my_plan.nnfft_flags & PRE_PHI_HUT) nnfft_precompute_phi_hut(&my_plan); /* init some guess */ for(k=0;k<my_plan.N_total;k++) { my_iplan.f_hat_iter[k]=0.0; } t0 = nfft_clock_gettime_seconds(); /* inverse trafo */ solver_before_loop_complex(&my_iplan); for(l=0;l<iteration;l++) { /* break if dot_r_iter is smaller than epsilon*/ if(my_iplan.dot_r_iter<epsilon) break; fprintf(stderr,"%e, %i of %i\n",sqrt(my_iplan.dot_r_iter), l+1,iteration); solver_loop_one_step_complex(&my_iplan); } t1 = nfft_clock_gettime_seconds(); t = t1-t0; fout_real=fopen("output_real.dat","w"); fout_imag=fopen("output_imag.dat","w"); for(k=0;k<my_plan.N_total;k++) { my_iplan.f_hat_iter[k]*=cexp(2.0*_Complex_I*M_PI*Ts*w[k]); fprintf(fout_real,"%le ", creal(my_iplan.f_hat_iter[k])); fprintf(fout_imag,"%le ", cimag(my_iplan.f_hat_iter[k])); } fclose(fout_real); fclose(fout_imag); /* finalize the infft */ solver_finalize_complex(&my_iplan); /* finalize the nfft */ nnfft_finalize(&my_plan); nfft_free(w); }
/** inverse NFFT-based mpolar FFT */ static int inverse_mpolar_fft(fftw_complex *f, int T, int R, fftw_complex *f_hat, int NN, int max_i, int m) { ticks 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 */ double *x, *w; /**< knots and associated weights */ int l; /**< index for iterations */ int N[2],n[2]; int M; /**< 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 = (double *)nfft_malloc(5*T*R/2*(sizeof(double))); if (x==NULL) return -1; w = (double *)nfft_malloc(5*T*R/4*(sizeof(double))); if (w==NULL) return -1; /** init two dimensional NFFT plan */ M=mpolar_grid(T,R,x,w); 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 */ 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.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); /** 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]= (sqrt(pow(j-my_nfft_plan.N[0]/2,2)+pow(k-my_nfft_plan.N[1]/2,2))>(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] = 0.0 + _Complex_I*0.0; t0 = getticks(); /** 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 = getticks(); GLOBAL_elapsed_time = nfft_elapsed_seconds(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; }
/** 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; }
static void simple_test_innfft_1d(void) { int j,k,l,N=8; /**< index for nodes, freqencies, iter*/ nnfft_plan my_plan; /**< plan for the nnfft */ solver_plan_complex my_iplan; /**< plan for the inverse nnfft */ /** initialise an one dimensional plan */ nnfft_init(&my_plan,1 ,8 ,8 ,&N); /** initialise my_iplan */ solver_init_advanced_complex(&my_iplan,(nfft_mv_plan_complex*)(&my_plan),CGNR); /** init pseudo random nodes */ for(j=0;j<my_plan.M_total;j++) my_plan.x[j]=((double)rand())/((double)RAND_MAX)-0.5; /** init pseudo random nodes */ for(k=0;k<my_plan.N_total;k++) my_plan.v[k]=((double)rand())/((double)RAND_MAX)-0.5; /** precompute psi, the entries of the matrix B */ if(my_plan.nnfft_flags & PRE_PSI) nnfft_precompute_psi(&my_plan); if(my_plan.nnfft_flags & PRE_FULL_PSI) nnfft_precompute_full_psi(&my_plan); /** precompute phi_hut, the entries of the matrix D */ if(my_plan.nnfft_flags & PRE_PHI_HUT) nnfft_precompute_phi_hut(&my_plan); /** init pseudo random samples (real) and show them */ for(j=0;j<my_plan.M_total;j++) my_iplan.y[j] = ((double)rand())/((double)RAND_MAX); nfft_vpr_complex(my_iplan.y,my_plan.M_total,"given data, vector given_f"); /** initialise some guess f_hat_0 */ for(k=0;k<my_plan.N_total;k++) my_iplan.f_hat_iter[k] = 0.0; nfft_vpr_complex(my_iplan.f_hat_iter,my_plan.N_total, "approximate solution, vector f_hat_iter"); /** solve the system */ solver_before_loop_complex(&my_iplan); for(l=0;l<8;l++) { printf("iteration l=%d\n",l); solver_loop_one_step_complex(&my_iplan); nfft_vpr_complex(my_iplan.f_hat_iter,my_plan.N_total, "approximate solution, vector f_hat_iter"); CSWAP(my_iplan.f_hat_iter,my_plan.f_hat); nnfft_trafo(&my_plan); nfft_vpr_complex(my_plan.f,my_plan.M_total,"fitting the data, vector f"); CSWAP(my_iplan.f_hat_iter,my_plan.f_hat); } solver_finalize_complex(&my_iplan); nnfft_finalize(&my_plan); }