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

reconstruct_data_2d1d.c

00001 #include <stdlib.h>
00002 #include <math.h>
00003 #include "util.h"
00004 #include "nfft3.h"
00005 
00015 void reconstruct(char* filename,int N,int M,int Z,int iteration, int weight, fftw_complex *mem)
00016 {
00017   int j,k,l,z;                  /* some variables  */
00018   double real,imag;             /* 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   int my_N[2],my_n[2];          /* to init the nfft */
00023   double tmp, epsilon=0.0000003;/* tmp to read the obsolent z from the input file
00024                                    epsilon is the break criterium for
00025                                    the iteration */
00026   unsigned infft_flags = CGNR | PRECOMPUTE_DAMP; /* flags for the infft */
00027                                    
00028   /* initialise my_plan */
00029   my_N[0]=N;my_n[0]=ceil(N*1.2);
00030   my_N[1]=N; my_n[1]=ceil(N*1.2);
00031   nfft_init_guru(&my_plan, 2, my_N, M/Z, my_n, 6, PRE_PHI_HUT| PRE_PSI|
00032                          MALLOC_X| MALLOC_F_HAT| MALLOC_F| 
00033                         FFTW_INIT| FFT_OUT_OF_PLACE,
00034                         FFTW_MEASURE| FFTW_DESTROY_INPUT);
00035   
00036   /* precompute lin psi if set */
00037   if(my_plan.nfft_flags & PRE_LIN_PSI)
00038     nfft_precompute_lin_psi(&my_plan);
00039   
00040   /* set the flags for the infft*/
00041   if (weight)
00042     infft_flags = infft_flags | PRECOMPUTE_WEIGHT;
00043   
00044   /* initialise my_iplan, advanced */
00045   infft_init_advanced(&my_iplan,&my_plan, infft_flags );
00046 
00047   /* get the weights */
00048   if(my_iplan.flags & PRECOMPUTE_WEIGHT)
00049   {
00050     fin=fopen("weights.dat","r");
00051     for(j=0;j<my_plan.M_total;j++) 
00052     {
00053         fscanf(fin,"%le ",&my_iplan.w[j]);
00054     }
00055     fclose(fin);
00056   }
00057   
00058   /* get the damping factors */
00059   if(my_iplan.flags & PRECOMPUTE_DAMP)
00060   {
00061     for(j=0;j<N;j++){
00062       for(k=0;k<N;k++) {
00063         int j2= j-N/2;
00064         int k2= k-N/2;
00065         double r=sqrt(j2*j2+k2*k2);
00066         if(r>(double) N/2) 
00067           my_iplan.w_hat[j*N+k]=0.0;
00068         else
00069           my_iplan.w_hat[j*N+k]=1.0;
00070       }   
00071     }
00072   }
00073   
00074   /* open the input file */
00075   fin=fopen(filename,"r"); 
00076 
00077   /* For every Layer*/
00078   for(z=0;z<Z;z++) {
00079 
00080     /* read x,y,freal and fimag from the knots */
00081     for(j=0;j<my_plan.M_total;j++)
00082     {
00083       fscanf(fin,"%le %le %le %le %le ",&my_plan.x[2*j+0],&my_plan.x[2*j+1], &tmp,
00084       &real,&imag);
00085       my_iplan.y[j] = real + I*imag;
00086     }
00087     
00088     /* precompute psi if set just one time because the knots equal each plane */
00089     if(z==0 && my_plan.nfft_flags & PRE_PSI) 
00090       nfft_precompute_psi(&my_plan);
00091       
00092     /* precompute full psi if set just one time because the knots equal each plane */
00093     if(z==0 && my_plan.nfft_flags & PRE_FULL_PSI)
00094       nfft_precompute_full_psi(&my_plan);
00095 
00096     /* init some guess */
00097     for(k=0;k<my_plan.N_total;k++)
00098       my_iplan.f_hat_iter[k]=0.0;
00099  
00100     /* inverse trafo */
00101     infft_before_loop(&my_iplan);
00102     for(l=0;l<iteration;l++)
00103     {
00104       /* break if dot_r_iter is smaller than epsilon*/
00105       if(my_iplan.dot_r_iter<epsilon)
00106       break;
00107       fprintf(stderr,"%e,  %i of %i\n",sqrt(my_iplan.dot_r_iter),
00108       iteration*z+l+1,iteration*Z);
00109       infft_loop_one_step(&my_iplan);
00110     }
00111     for(k=0;k<my_plan.N_total;k++) {
00112       /* write every slice in the memory.
00113       here we make an fftshift direct */
00114       mem[(Z*N*N/2+z*N*N+ k)%(Z*N*N)] = my_iplan.f_hat_iter[k];
00115     }
00116   }
00117 
00118   fclose(fin);
00119   
00120   /* finalize the infft */
00121   infft_finalize(&my_iplan);
00122   
00123   /* finalize the nfft */
00124   nfft_finalize(&my_plan);
00125 }
00126 
00131 void print(int N,int M,int Z, fftw_complex *mem)
00132 {
00133   int i,j;
00134   FILE* fout_real;
00135   FILE* fout_imag;
00136   fout_real=fopen("output_real.dat","w");
00137   fout_imag=fopen("output_imag.dat","w");
00138 
00139   for(i=0;i<Z;i++) {
00140     for (j=0;j<N*N;j++) {
00141       fprintf(fout_real,"%le ",creal(mem[(Z*N*N/2+i*N*N+ j)%(Z*N*N)]) /Z);
00142       fprintf(fout_imag,"%le ",cimag(mem[(Z*N*N/2+i*N*N+ j)%(Z*N*N)]) /Z);
00143     }
00144     fprintf(fout_real,"\n");
00145     fprintf(fout_imag,"\n");
00146   }
00147 
00148   fclose(fout_real);
00149   fclose(fout_imag);
00150 }
00151 
00152 int main(int argc, char **argv)
00153 {
00154   fftw_complex *mem;
00155   fftw_plan plan;
00156   int N,M,Z;  
00157 
00158   if (argc <= 6) {
00159     printf("usage: ./reconstruct FILENAME N M Z ITER WEIGHTS\n");
00160     return 1;
00161   }
00162 
00163   N=atoi(argv[2]);
00164   M=atoi(argv[3]);
00165   Z=atoi(argv[4]);
00166 
00167   /* Allocate memory to hold every layer in memory after the
00168   2D-infft */
00169   mem = (fftw_complex*) fftw_malloc(sizeof(fftw_complex) * atoi(argv[2]) * atoi(argv[2]) * atoi(argv[4]));
00170 
00171   /* Create plan for the 1d-ifft */
00172   plan = fftw_plan_many_dft(1, &Z, N*N,
00173                                   mem, NULL,
00174                                   N*N, 1,
00175                                   mem, NULL,
00176                                   N*N,1 ,
00177                                   FFTW_BACKWARD, FFTW_MEASURE);
00178   
00179   /* execute the 2d-infft's */
00180   reconstruct(argv[1],N,M,Z,atoi(argv[5]),atoi(argv[6]),mem);
00181 
00182   /* execute the 1d-fft's */
00183   fftw_execute(plan);
00184 
00185   /* write the memory back in files */
00186   print(N,M,Z, mem);
00187 
00188   /* free memory */
00189   fftw_free(mem);
00190   fftw_destroy_plan(plan);
00191   return 1;
00192 }
00193 /* \} */

Generated on 1 Nov 2006 by Doxygen 1.4.4