NFFT Logo 3.0 API Reference
Main Page | Modules | Data Structures | Directories | File List | Data Fields | Globals

reconstruct_data_inh_nnfft.c

00001 #include <stdlib.h>
00002 #include <math.h>
00003 #include <limits.h>
00004 #include "util.h"
00005 #include "nfft3.h"
00006 
00016 void reconstruct(char* filename,int N,int M,int iteration, int weight)
00017 { 
00018   int j,k,l;                    /* some variables  */
00019   nnfft_plan my_plan;            /* plan for the two dimensional nfft  */
00020   innfft_plan my_iplan;          /* plan for the two dimensional infft */
00021   FILE* fin;                    /* input file                         */
00022   FILE* finh;
00023   FILE* ftime;
00024   FILE* fout_real;              /* output file                        */
00025   FILE* fout_imag;              /* output file                        */
00026   int my_N[3],my_n[3];          /* to init the nfft */
00027   double t,epsilon=0.0000003;     /* epsilon is a the break criterium for
00028                                    the iteration */
00029   unsigned infft_flags = CGNR | PRECOMPUTE_DAMP; /* flags for the infft*/
00030   double time,min_time,max_time,min_inh,max_inh;
00031   double real,imag;
00032   double *w;
00033 
00034   double Ts;
00035   double W;
00036   int N3;
00037   int m=2;
00038   double sigma = 1.25;
00039   
00040   w = (double*) malloc(N*N*sizeof(double));
00041 
00042   ftime=fopen("readout_time.dat","r");
00043   finh=fopen("inh.dat","r");
00044 
00045   min_time=INT_MAX; max_time=INT_MIN;
00046   for(j=0;j<M;j++)
00047   {
00048     fscanf(ftime,"%le ",&time);
00049     if(time<min_time)
00050       min_time = time;
00051     if(time>max_time)
00052       max_time = time;
00053   }
00054 
00055   fclose(ftime);
00056   
00057   Ts=(min_time+max_time)/2.0;
00058 
00059   min_inh=INT_MAX; max_inh=INT_MIN;
00060   for(j=0;j<N*N;j++)
00061   {
00062     fscanf(finh,"%le ",&w[j]);
00063     if(w[j]<min_inh)
00064       min_inh = w[j];
00065     if(w[j]>max_inh)
00066       max_inh = w[j];
00067   }
00068   fclose(finh);
00069 
00070   N3=ceil((NFFT_MAX(fabs(min_inh),fabs(max_inh))*(max_time-min_time)/2.0)*4);
00071 
00072 
00073   W=NFFT_MAX(fabs(min_inh),fabs(max_inh))*2.0;
00074   
00075   fprintf(stderr,"3:  %i %e %e %e %e %e %e\n",N3,W,min_inh,max_inh,min_time,max_time,Ts);
00076 
00077   /* initialise my_plan */
00078   my_N[0]=N;my_n[0]=ceil(N*sigma);
00079   my_N[1]=N; my_n[1]=ceil(N*sigma);
00080   my_N[2]=N3; my_n[2]=ceil(N3*sigma);
00081   nnfft_init_guru(&my_plan, 3, N*N, M, my_N,my_n,m,
00082         PRE_PSI| PRE_PHI_HUT| MALLOC_X| MALLOC_V| MALLOC_F_HAT| MALLOC_F );
00083         
00084   /* precompute lin psi if set */
00085   if(my_plan.nnfft_flags & PRE_LIN_PSI)
00086     nnfft_precompute_lin_psi(&my_plan);
00087   
00088   /* set the flags for the infft*/
00089   if (weight)
00090     infft_flags = infft_flags | PRECOMPUTE_WEIGHT;
00091 
00092   /* initialise my_iplan, advanced */
00093   innfft_init_advanced(&my_iplan,&my_plan, infft_flags );
00094 
00095   /* get the weights */
00096   if(my_iplan.flags & PRECOMPUTE_WEIGHT)
00097   {
00098     fin=fopen("weights.dat","r");
00099     for(j=0;j<my_plan.M_total;j++)
00100     {
00101         fscanf(fin,"%le ",&my_iplan.w[j]);
00102     }
00103     fclose(fin);
00104   }
00105   
00106   /* get the damping factors */
00107   if(my_iplan.flags & PRECOMPUTE_DAMP)
00108   {
00109     for(j=0;j<N;j++){
00110       for(k=0;k<N;k++) {
00111         int j2= j-N/2;
00112         int k2= k-N/2;
00113         double r=sqrt(j2*j2+k2*k2);
00114         if(r>(double) N/2) 
00115           my_iplan.w_hat[j*N+k]=0.0;
00116         else
00117           my_iplan.w_hat[j*N+k]=1.0;
00118       }   
00119     }
00120   }
00121   
00122   /* open the input file */
00123   fin=fopen(filename,"r"); 
00124   ftime=fopen("readout_time.dat","r");
00125 
00126   for(j=0;j<my_plan.M_total;j++)
00127   {
00128     fscanf(fin,"%le %le %le %le ",&my_plan.x[3*j+0],&my_plan.x[3*j+1],&real,&imag);
00129     my_iplan.y[j]=real+I*imag;
00130     fscanf(ftime,"%le ",&my_plan.x[3*j+2]);
00131 
00132     my_plan.x[3*j+2] = (my_plan.x[3*j+2]-Ts)*W/N3;
00133   }
00134 
00135   for(j=0;j<N;j++)
00136     { 
00137     for(l=0;l<N;l++)
00138       {
00139         my_plan.v[3*(N*j+l)+0]=(((double) j) -(((double) N)/2.0))/((double) N);
00140         my_plan.v[3*(N*j+l)+1]=(((double) l) -(((double) N)/2.0))/((double) N);
00141         my_plan.v[3*(N*j+l)+2] = w[N*j+l]/W ;
00142       }
00143     }
00144  
00145   /* precompute psi */
00146   if(my_plan.nnfft_flags & PRE_PSI) {
00147     nnfft_precompute_psi(&my_plan);
00148     if(my_plan.nnfft_flags & PRE_FULL_PSI)
00149       nnfft_precompute_full_psi(&my_plan);
00150   }
00151   
00152   if(my_plan.nnfft_flags & PRE_PHI_HUT)
00153     nnfft_precompute_phi_hut(&my_plan);
00154     
00155   /* init some guess */
00156   for(k=0;k<my_plan.N_total;k++)
00157   {
00158     my_iplan.f_hat_iter[k]=0.0;
00159   }
00160  
00161   t=nfft_second();
00162   
00163   /* inverse trafo */
00164   innfft_before_loop(&my_iplan);
00165   for(l=0;l<iteration;l++)
00166   {
00167     /* break if dot_r_iter is smaller than epsilon*/
00168     if(my_iplan.dot_r_iter<epsilon)
00169     break;
00170     fprintf(stderr,"%e,  %i of %i\n",sqrt(my_iplan.dot_r_iter),
00171     l+1,iteration);
00172     innfft_loop_one_step(&my_iplan);
00173   }
00174 
00175   t=nfft_second()-t;
00176 #ifdef HAVE_TOTAL_USED_MEMORY
00177 fprintf(stderr,"time: %e seconds mem: %i \n",t,nfft_total_used_memory());
00178 #else
00179 fprintf(stderr,"time: %e seconds mem: mallinfo not available\n",t);
00180 #endif
00181 
00182   fout_real=fopen("output_real.dat","w");
00183   fout_imag=fopen("output_imag.dat","w");
00184 
00185   for(k=0;k<my_plan.N_total;k++) {
00186 
00187     my_iplan.f_hat_iter[k]*=cexp(2.0*I*PI*Ts*w[k]);
00188 
00189     fprintf(fout_real,"%le ", creal(my_iplan.f_hat_iter[k]));
00190     fprintf(fout_imag,"%le ", cimag(my_iplan.f_hat_iter[k]));
00191   }
00192   
00193 
00194   fclose(fout_real);
00195   fclose(fout_imag);
00196 
00197 
00198   /* finalize the infft */
00199   innfft_finalize(&my_iplan);
00200   
00201   /* finalize the nfft */
00202   nnfft_finalize(&my_plan);
00203 
00204   free(w);
00205 }
00206 
00207 int main(int argc, char **argv)
00208 {
00209   if (argc <= 5) {
00210     printf("usage: ./reconstruct_data_inh_nnfft FILENAME N M ITER WEIGHTS\n");
00211     return 1;
00212   }
00213   
00214   reconstruct(argv[1],atoi(argv[2]),atoi(argv[3]),atoi(argv[4]),atoi(argv[5]));
00215    
00216   return 1;
00217 }
00218 /* \} */

Generated on 1 Nov 2006 by Doxygen 1.4.4