NFFT Logo 3.2.2
reconstruct_data_inh_nnfft.c
00001 /*
00002  * Copyright (c) 2002, 2012 Jens Keiner, Stefan Kunis, Daniel Potts
00003  *
00004  * This program is free software; you can redistribute it and/or modify it under
00005  * the terms of the GNU General Public License as published by the Free Software
00006  * Foundation; either version 2 of the License, or (at your option) any later
00007  * version.
00008  *
00009  * This program is distributed in the hope that it will be useful, but WITHOUT
00010  * ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
00011  * FOR A PARTICULAR PURPOSE.  See the GNU General Public License for more
00012  * details.
00013  *
00014  * You should have received a copy of the GNU General Public License along with
00015  * this program; if not, write to the Free Software Foundation, Inc., 51
00016  * Franklin Street, Fifth Floor, Boston, MA 02110-1301, USA.
00017  */
00018 
00019 /* $Id: reconstruct_data_inh_nnfft.c 3896 2012-10-10 12:19:26Z tovo $ */
00020 #include "config.h"
00021 
00022 #include <stdlib.h>
00023 #include <math.h>
00024 #include <limits.h>
00025 #ifdef HAVE_COMPLEX_H
00026 #include <complex.h>
00027 #endif
00028 
00029 #include "nfft3util.h"
00030 #include "nfft3.h"
00031 #include "infft.h"
00032 
00042 static void reconstruct(char* filename,int N,int M,int iteration, int weight)
00043 {
00044   int j,k,l;                    /* some variables  */
00045   nnfft_plan my_plan;            /* plan for the two dimensional nfft  */
00046   solver_plan_complex my_iplan;          /* plan for the two dimensional infft */
00047   FILE* fin;                    /* input file                         */
00048   FILE* finh;
00049   FILE* ftime;
00050   FILE* fout_real;              /* output file                        */
00051   FILE* fout_imag;              /* output file                        */
00052   int my_N[3],my_n[3];          /* to init the nfft */
00053   ticks t0, t1;
00054   double t,epsilon=0.0000003;     /* epsilon is a the break criterium for
00055                                    the iteration */
00056   unsigned infft_flags = CGNR | PRECOMPUTE_DAMP; /* flags for the infft*/
00057   double time,min_time,max_time,min_inh,max_inh;
00058   double real,imag;
00059   double *w;
00060 
00061   double Ts;
00062   double W;
00063   int N3;
00064   int m=2;
00065   double sigma = 1.25;
00066 
00067   w = (double*)nfft_malloc(N*N*sizeof(double));
00068 
00069   ftime=fopen("readout_time.dat","r");
00070   finh=fopen("inh.dat","r");
00071 
00072   min_time=INT_MAX; max_time=INT_MIN;
00073   for(j=0;j<M;j++)
00074   {
00075     fscanf(ftime,"%le ",&time);
00076     if(time<min_time)
00077       min_time = time;
00078     if(time>max_time)
00079       max_time = time;
00080   }
00081 
00082   fclose(ftime);
00083 
00084   Ts=(min_time+max_time)/2.0;
00085 
00086   min_inh=INT_MAX; max_inh=INT_MIN;
00087   for(j=0;j<N*N;j++)
00088   {
00089     fscanf(finh,"%le ",&w[j]);
00090     if(w[j]<min_inh)
00091       min_inh = w[j];
00092     if(w[j]>max_inh)
00093       max_inh = w[j];
00094   }
00095   fclose(finh);
00096 
00097   N3=ceil((NFFT_MAX(fabs(min_inh),fabs(max_inh))*(max_time-min_time)/2.0)*4);
00098 
00099 
00100   W=NFFT_MAX(fabs(min_inh),fabs(max_inh))*2.0;
00101 
00102   fprintf(stderr,"3:  %i %e %e %e %e %e %e\n",N3,W,min_inh,max_inh,min_time,max_time,Ts);
00103 
00104   /* initialise my_plan */
00105   my_N[0]=N;my_n[0]=ceil(N*sigma);
00106   my_N[1]=N; my_n[1]=ceil(N*sigma);
00107   my_N[2]=N3; my_n[2]=ceil(N3*sigma);
00108   nnfft_init_guru(&my_plan, 3, N*N, M, my_N,my_n,m,
00109         PRE_PSI| PRE_PHI_HUT| MALLOC_X| MALLOC_V| MALLOC_F_HAT| MALLOC_F );
00110 
00111   /* precompute lin psi if set */
00112   if(my_plan.nnfft_flags & PRE_LIN_PSI)
00113     nnfft_precompute_lin_psi(&my_plan);
00114 
00115   /* set the flags for the infft*/
00116   if (weight)
00117     infft_flags = infft_flags | PRECOMPUTE_WEIGHT;
00118 
00119   /* initialise my_iplan, advanced */
00120   solver_init_advanced_complex(&my_iplan,(nfft_mv_plan_complex*)(&my_plan), infft_flags );
00121 
00122   /* get the weights */
00123   if(my_iplan.flags & PRECOMPUTE_WEIGHT)
00124   {
00125     fin=fopen("weights.dat","r");
00126     for(j=0;j<my_plan.M_total;j++)
00127     {
00128         fscanf(fin,"%le ",&my_iplan.w[j]);
00129     }
00130     fclose(fin);
00131   }
00132 
00133   /* get the damping factors */
00134   if(my_iplan.flags & PRECOMPUTE_DAMP)
00135   {
00136     for(j=0;j<N;j++){
00137       for(k=0;k<N;k++) {
00138         int j2= j-N/2;
00139         int k2= k-N/2;
00140         double r=sqrt(j2*j2+k2*k2);
00141         if(r>(double) N/2)
00142           my_iplan.w_hat[j*N+k]=0.0;
00143         else
00144           my_iplan.w_hat[j*N+k]=1.0;
00145       }
00146     }
00147   }
00148 
00149   /* open the input file */
00150   fin=fopen(filename,"r");
00151   ftime=fopen("readout_time.dat","r");
00152 
00153   for(j=0;j<my_plan.M_total;j++)
00154   {
00155     fscanf(fin,"%le %le %le %le ",&my_plan.x[3*j+0],&my_plan.x[3*j+1],&real,&imag);
00156     my_iplan.y[j]=real+ _Complex_I*imag;
00157     fscanf(ftime,"%le ",&my_plan.x[3*j+2]);
00158 
00159     my_plan.x[3*j+2] = (my_plan.x[3*j+2]-Ts)*W/N3;
00160   }
00161 
00162   for(j=0;j<N;j++)
00163     {
00164     for(l=0;l<N;l++)
00165       {
00166         my_plan.v[3*(N*j+l)+0]=(((double) j) -(((double) N)/2.0))/((double) N);
00167         my_plan.v[3*(N*j+l)+1]=(((double) l) -(((double) N)/2.0))/((double) N);
00168         my_plan.v[3*(N*j+l)+2] = w[N*j+l]/W ;
00169       }
00170     }
00171 
00172   /* precompute psi */
00173   if(my_plan.nnfft_flags & PRE_PSI) {
00174     nnfft_precompute_psi(&my_plan);
00175     if(my_plan.nnfft_flags & PRE_FULL_PSI)
00176       nnfft_precompute_full_psi(&my_plan);
00177   }
00178 
00179   if(my_plan.nnfft_flags & PRE_PHI_HUT)
00180     nnfft_precompute_phi_hut(&my_plan);
00181 
00182   /* init some guess */
00183   for(k=0;k<my_plan.N_total;k++)
00184   {
00185     my_iplan.f_hat_iter[k]=0.0;
00186   }
00187 
00188   t0 = getticks();
00189 
00190   /* inverse trafo */
00191   solver_before_loop_complex(&my_iplan);
00192   for(l=0;l<iteration;l++)
00193   {
00194     /* break if dot_r_iter is smaller than epsilon*/
00195     if(my_iplan.dot_r_iter<epsilon)
00196     break;
00197     fprintf(stderr,"%e,  %i of %i\n",sqrt(my_iplan.dot_r_iter),
00198     l+1,iteration);
00199     solver_loop_one_step_complex(&my_iplan);
00200   }
00201 
00202   t1 = getticks();
00203   t = nfft_elapsed_seconds(t1,t0);
00204 
00205   fout_real=fopen("output_real.dat","w");
00206   fout_imag=fopen("output_imag.dat","w");
00207 
00208   for(k=0;k<my_plan.N_total;k++) {
00209 
00210     my_iplan.f_hat_iter[k]*=cexp(2.0*_Complex_I*PI*Ts*w[k]);
00211 
00212     fprintf(fout_real,"%le ", creal(my_iplan.f_hat_iter[k]));
00213     fprintf(fout_imag,"%le ", cimag(my_iplan.f_hat_iter[k]));
00214   }
00215 
00216 
00217   fclose(fout_real);
00218   fclose(fout_imag);
00219 
00220 
00221   /* finalize the infft */
00222   solver_finalize_complex(&my_iplan);
00223 
00224   /* finalize the nfft */
00225   nnfft_finalize(&my_plan);
00226 
00227   nfft_free(w);
00228 }
00229 
00230 int main(int argc, char **argv)
00231 {
00232   if (argc <= 5) {
00233     printf("usage: ./reconstruct_data_inh_nnfft FILENAME N M ITER WEIGHTS\n");
00234     return 1;
00235   }
00236 
00237   reconstruct(argv[1],atoi(argv[2]),atoi(argv[3]),atoi(argv[4]),atoi(argv[5]));
00238 
00239   return 1;
00240 }
00241 /* \} */

Generated on Fri Oct 12 2012 by Doxygen 1.8.0-20120409