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

reconstruct_data_3d.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 Z,int iteration, int weight)
00016 {
00017   int j,k,z,l;                  /* 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   FILE* fout_real;              /* output file (real part) */
00023   FILE* fout_imag;              /* output file (imag part) */
00024   int my_N[3],my_n[3];          /* to init the nfft */
00025   double epsilon=0.0000003;     /* tmp to read the obsolent z from 700.acs
00026                                    epsilon is a the break criterion for
00027                                    the iteration */
00028   unsigned infft_flags = CGNR | PRECOMPUTE_DAMP;  /* flags for the infft */
00029                                    
00030   /* initialise my_plan, specific. 
00031      we don't precompute psi */
00032   my_N[0]=Z; my_n[0]=ceil(Z*1.2);
00033   my_N[1]=N; my_n[1]=ceil(N*1.2);
00034   my_N[2]=N; my_n[2]=ceil(N*1.2);
00035   nfft_init_guru(&my_plan, 3, my_N, M, my_n, 6,
00036                       PRE_PHI_HUT| PRE_PSI |MALLOC_X| MALLOC_F_HAT|
00037                       MALLOC_F| FFTW_INIT| FFT_OUT_OF_PLACE,
00038                       FFTW_MEASURE| FFTW_DESTROY_INPUT);
00039  
00040   /* precompute lin psi */
00041   if(my_plan.nfft_flags & PRE_LIN_PSI)
00042     nfft_precompute_lin_psi(&my_plan);
00043                     
00044   if (weight)
00045     infft_flags = infft_flags | PRECOMPUTE_WEIGHT;
00046   
00047   /* initialise my_iplan, advanced */
00048   infft_init_advanced(&my_iplan,&my_plan, infft_flags );
00049  
00050   /* get the weights */
00051   if(my_iplan.flags & PRECOMPUTE_WEIGHT)
00052   {
00053     fin=fopen("weights.dat","r");
00054     for(j=0;j<M;j++)
00055     {
00056       fscanf(fin,"%le ",&my_iplan.w[j]);
00057     }
00058     fclose(fin);
00059   }
00060   
00061   /* get the damping factors */
00062   if(my_iplan.flags & PRECOMPUTE_DAMP)
00063   {
00064     for(j=0;j<N;j++){
00065       for(k=0;k<N;k++) {
00066         for(z=0;z<N;z++) {
00067         int j2= j-N/2;
00068         int k2= k-N/2;
00069         int z2= z-N/2;
00070         double r=sqrt(j2*j2+k2*k2+z2*z2);
00071         if(r>(double) N/2) 
00072           my_iplan.w_hat[z*N*N+j*N+k]=0.0;
00073         else
00074           my_iplan.w_hat[z*N*N+j*N+k]=1.0;
00075         }
00076       }   
00077     }
00078   }
00079 
00080   /* open the input file */
00081   fin=fopen(filename,"r");
00082   
00083   /* open the output files */
00084   fout_real=fopen("output_real.dat","w");
00085   fout_imag=fopen("output_imag.dat","w");
00086 
00087   /* read x,y,freal and fimag from the knots */
00088   for(j=0;j<M;j++)
00089   {
00090     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],
00091     &real,&imag);
00092     my_iplan.y[j] = real + I*imag;
00093   }
00094 
00095   /* precompute psi */
00096   if(my_plan.nfft_flags & PRE_PSI)
00097     nfft_precompute_psi(&my_plan);      
00098   
00099   /* precompute full psi */
00100   if(my_plan.nfft_flags & PRE_FULL_PSI)
00101     nfft_precompute_full_psi(&my_plan);
00102 
00103   /* init some guess */
00104   for(k=0;k<my_plan.N_total;k++)
00105     my_iplan.f_hat_iter[k]=0.0;
00106   
00107   /* inverse trafo */
00108   infft_before_loop(&my_iplan);
00109   for(l=0;l<iteration;l++)
00110   {
00111     /* break if dot_r_iter is smaller than epsilon*/
00112     if(my_iplan.dot_r_iter<epsilon)
00113       break;
00114     fprintf(stderr,"%e,  %i of %i\n",sqrt(my_iplan.dot_r_iter),
00115     l+1,iteration);
00116     infft_loop_one_step(&my_iplan);
00117   }
00118   
00119   for(l=0;l<Z;l++)
00120   {
00121     for(k=0;k<N*N;k++)
00122     {
00123       /* write every Layer in the files */
00124       fprintf(fout_real,"%le ",creal(my_iplan.f_hat_iter[ k+N*N*l ]));
00125       fprintf(fout_imag,"%le ",cimag(my_iplan.f_hat_iter[ k+N*N*l ]));
00126     }  
00127     fprintf(fout_real,"\n");
00128     fprintf(fout_imag,"\n");
00129   }
00130   
00131   fclose(fout_real);
00132   fclose(fout_imag);
00133 
00134   infft_finalize(&my_iplan);
00135   nfft_finalize(&my_plan);
00136 }
00137 
00138 int main(int argc, char **argv)
00139 {
00140   if (argc <= 6) { 
00141     printf("usage: ./reconstruct3D FILENAME N M Z ITER WEIGHTS\n");
00142     return 1;
00143   }
00144 
00145   reconstruct(argv[1],atoi(argv[2]),atoi(argv[3]),atoi(argv[4]),atoi(argv[5]),atoi(argv[6]));
00146   return 1;
00147 }
00148 /* \} */

Generated on 1 Nov 2006 by Doxygen 1.4.4