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

reconstruct_data_2d.c

00001 #include <math.h>
00002 #include <stdlib.h>
00003 #include "util.h"
00004 #include "nfft3.h"
00005 
00015 void reconstruct(char* filename,int N,int M,int iteration, int weight)
00016 {
00017   int j,k,l;                    /* some variables  */
00018   double real,imag,t;           /* to read the real and imag part of a complex number */
00019   nfft_plan my_plan;            /* plan for the two dimensional nfft  */
00020   infft_plan my_iplan;          /* plan for the two dimensional infft */
00021   FILE* fin;                    /* input file                         */
00022   FILE* fout_real;              /* output file                        */
00023   FILE* fout_imag;              /* output file                        */
00024   int my_N[2],my_n[2];          /* to init the nfft */
00025   double epsilon=0.0000003;     /* epsilon is a the break criterium for
00026                                    the iteration */
00027   unsigned infft_flags = CGNR | PRECOMPUTE_DAMP;  /* flags for the infft*/
00028   int m = 6;
00029   double alpha = 2.0;
00030   /* initialise my_plan */
00031   my_N[0]=N; my_n[0]=ceil(N*alpha);
00032   my_N[1]=N; my_n[1]=ceil(N*alpha);
00033   nfft_init_guru(&my_plan, 2, my_N, M, my_n, m, PRE_PHI_HUT| PRE_PSI|
00034                          MALLOC_X| MALLOC_F_HAT| MALLOC_F| 
00035                          FFTW_INIT| FFT_OUT_OF_PLACE,
00036                          FFTW_MEASURE| FFTW_DESTROY_INPUT);
00037   
00038   /* precompute lin psi if set */
00039   if(my_plan.nfft_flags & PRE_LIN_PSI)
00040     nfft_precompute_lin_psi(&my_plan);
00041   
00042   /* set the flags for the infft*/
00043   if (weight)
00044     infft_flags = infft_flags | PRECOMPUTE_WEIGHT;
00045 
00046   /* initialise my_iplan, advanced */
00047   infft_init_advanced(&my_iplan,&my_plan, infft_flags );
00048 
00049   /* get the weights */
00050   if(my_iplan.flags & PRECOMPUTE_WEIGHT)
00051   {
00052     fin=fopen("weights.dat","r");
00053     for(j=0;j<my_plan.M_total;j++) 
00054     {
00055         fscanf(fin,"%le ",&my_iplan.w[j]);
00056     }
00057     fclose(fin);
00058   }
00059   
00060   /* get the damping factors */
00061   if(my_iplan.flags & PRECOMPUTE_DAMP)
00062   {
00063     for(j=0;j<N;j++){
00064       for(k=0;k<N;k++) {
00065         int j2= j-N/2;
00066         int k2= k-N/2;
00067         double r=sqrt(j2*j2+k2*k2);
00068         if(r>(double) N/2) 
00069           my_iplan.w_hat[j*N+k]=0.0;
00070         else
00071           my_iplan.w_hat[j*N+k]=1.0;
00072       }   
00073     }
00074   }
00075   
00076   /* open the input file */
00077   fin=fopen(filename,"r"); 
00078 
00079   /* read x,y,freal and fimag from the knots */
00080   for(j=0;j<my_plan.M_total;j++)
00081   {
00082     fscanf(fin,"%le %le %le %le ",&my_plan.x[2*j+0],&my_plan.x[2*j+1],
00083     &real,&imag);
00084     my_iplan.y[j] = real + I*imag;
00085   }
00086 
00087   fclose(fin);
00088   
00089   /* precompute psi */
00090   if(my_plan.nfft_flags & PRE_PSI)
00091     nfft_precompute_psi(&my_plan);
00092 
00093   /* precompute full psi */
00094   if(my_plan.nfft_flags & PRE_FULL_PSI)
00095       nfft_precompute_full_psi(&my_plan);
00096 
00097   /* init some guess */
00098   for(k=0;k<my_plan.N_total;k++)
00099     my_iplan.f_hat_iter[k]=0.0;
00100 
00101   t=nfft_second();
00102     
00103   /* inverse trafo */
00104   infft_before_loop(&my_iplan);
00105   for(l=0;l<iteration;l++)
00106   {
00107     /* break if dot_r_iter is smaller than epsilon*/
00108     if(my_iplan.dot_r_iter<epsilon)
00109       break;
00110     fprintf(stderr,"%e,  %i of %i\n",sqrt(my_iplan.dot_r_iter),
00111     l+1,iteration);
00112     infft_loop_one_step(&my_iplan);
00113   }
00114 
00115   
00116   t=nfft_second()-t;
00117 #ifdef HAVE_MALLINFO
00118   fprintf(stderr,"time: %e seconds mem: %i \n",t,nfft_total_used_memory());
00119 #else
00120   fprintf(stderr,"time: %e seconds mem: mallinfo not available\n",t);
00121 #endif
00122 
00123   
00124   fout_real=fopen("output_real.dat","w");
00125   fout_imag=fopen("output_imag.dat","w");
00126 
00127   for(k=0;k<my_plan.N_total;k++) {
00128     fprintf(fout_real,"%le ", creal(my_iplan.f_hat_iter[k]));
00129     fprintf(fout_imag,"%le ", cimag(my_iplan.f_hat_iter[k]));
00130   }
00131 
00132   fclose(fout_real);
00133   fclose(fout_imag);
00134 
00135   /* finalize the infft */
00136   infft_finalize(&my_iplan);
00137   
00138   /* finalize the nfft */
00139   nfft_finalize(&my_plan);
00140 }
00141 
00142 int main(int argc, char **argv)
00143 {
00144   if (argc <= 5) {
00145     printf("usage: ./reconstruct_data_2d FILENAME N M ITER WEIGHTS\n");
00146     return 1;
00147   }
00148   
00149   reconstruct(argv[1],atoi(argv[2]),atoi(argv[3]),atoi(argv[4]),atoi(argv[5]));
00150    
00151   return 1;
00152 }
00153 /* \} */

Generated on 1 Nov 2006 by Doxygen 1.4.4