NFFT Logo 3.2.2
reconstruct_data_2d.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_2d.c 3896 2012-10-10 12:19:26Z tovo $ */
00020 #include "config.h"
00021 
00022 #include <math.h>
00023 #include <stdlib.h>
00024 #ifdef HAVE_COMPLEX_H
00025 #include <complex.h>
00026 #endif
00027 
00028 #include "nfft3util.h"
00029 #include "nfft3.h"
00030 #include "infft.h"
00031 
00041 static void reconstruct(char* filename,int N,int M,int iteration, int weight)
00042 {
00043   int j,k,l;                    /* some variables  */
00044   ticks t0, t1;
00045   double real,imag,t;           /* to read the real and imag part of a complex number */
00046   nfft_plan my_plan;            /* plan for the two dimensional nfft  */
00047   solver_plan_complex my_iplan; /* plan for the two dimensional infft */
00048   FILE* fin;                    /* input file                         */
00049   FILE* fout_real;              /* output file                        */
00050   FILE* fout_imag;              /* output file                        */
00051   int my_N[2],my_n[2];          /* to init the nfft */
00052   double epsilon=0.0000003;     /* epsilon is a the break criterium for
00053                                    the iteration */
00054   unsigned infft_flags = CGNR | PRECOMPUTE_DAMP;  /* flags for the infft*/
00055   int m = 6;
00056   double alpha = 2.0;
00057   /* initialise my_plan */
00058   my_N[0]=N; my_n[0]=ceil(N*alpha);
00059   my_N[1]=N; my_n[1]=ceil(N*alpha);
00060   nfft_init_guru(&my_plan, 2, my_N, M, my_n, m, PRE_PHI_HUT| PRE_PSI|
00061                          MALLOC_X| MALLOC_F_HAT| MALLOC_F|
00062                          FFTW_INIT| FFT_OUT_OF_PLACE,
00063                          FFTW_MEASURE| FFTW_DESTROY_INPUT);
00064 
00065   /* precompute lin psi if set */
00066   if(my_plan.nfft_flags & PRE_LIN_PSI)
00067     nfft_precompute_lin_psi(&my_plan);
00068 
00069   /* set the flags for the infft*/
00070   if (weight)
00071     infft_flags = infft_flags | PRECOMPUTE_WEIGHT;
00072 
00073   /* initialise my_iplan, advanced */
00074   solver_init_advanced_complex(&my_iplan,(nfft_mv_plan_complex*)&my_plan, infft_flags );
00075 
00076   /* get the weights */
00077   if(my_iplan.flags & PRECOMPUTE_WEIGHT)
00078   {
00079     fin=fopen("weights.dat","r");
00080     for(j=0;j<my_plan.M_total;j++)
00081     {
00082         fscanf(fin,"%le ",&my_iplan.w[j]);
00083     }
00084     fclose(fin);
00085   }
00086 
00087   /* get the damping factors */
00088   if(my_iplan.flags & PRECOMPUTE_DAMP)
00089   {
00090     for(j=0;j<N;j++){
00091       for(k=0;k<N;k++) {
00092         int j2= j-N/2;
00093         int k2= k-N/2;
00094         double r=sqrt(j2*j2+k2*k2);
00095         if(r>(double) N/2)
00096           my_iplan.w_hat[j*N+k]=0.0;
00097         else
00098           my_iplan.w_hat[j*N+k]=1.0;
00099       }
00100     }
00101   }
00102 
00103   /* open the input file */
00104   fin=fopen(filename,"r");
00105 
00106   /* read x,y,freal and fimag from the knots */
00107   for(j=0;j<my_plan.M_total;j++)
00108   {
00109     fscanf(fin,"%le %le %le %le ",&my_plan.x[2*j+0],&my_plan.x[2*j+1],
00110     &real,&imag);
00111     my_iplan.y[j] = real + _Complex_I*imag;
00112   }
00113 
00114   fclose(fin);
00115 
00116   /* precompute psi */
00117   if(my_plan.nfft_flags & PRE_PSI)
00118     nfft_precompute_psi(&my_plan);
00119 
00120   /* precompute full psi */
00121   if(my_plan.nfft_flags & PRE_FULL_PSI)
00122       nfft_precompute_full_psi(&my_plan);
00123 
00124   /* init some guess */
00125   for(k=0;k<my_plan.N_total;k++)
00126     my_iplan.f_hat_iter[k]=0.0;
00127 
00128   t0 = getticks();
00129 
00130   /* inverse trafo */
00131   solver_before_loop_complex(&my_iplan);
00132   for(l=0;l<iteration;l++)
00133   {
00134     /* break if dot_r_iter is smaller than epsilon*/
00135     if(my_iplan.dot_r_iter<epsilon)
00136       break;
00137     fprintf(stderr,"%e,  %i of %i\n",sqrt(my_iplan.dot_r_iter),
00138     l+1,iteration);
00139     solver_loop_one_step_complex(&my_iplan);
00140   }
00141 
00142 
00143   t1 = getticks();
00144   t=nfft_elapsed_seconds(t1,t0);
00145 
00146   fout_real=fopen("output_real.dat","w");
00147   fout_imag=fopen("output_imag.dat","w");
00148 
00149   for(k=0;k<my_plan.N_total;k++) {
00150     fprintf(fout_real,"%le ", creal(my_iplan.f_hat_iter[k]));
00151     fprintf(fout_imag,"%le ", cimag(my_iplan.f_hat_iter[k]));
00152   }
00153 
00154   fclose(fout_real);
00155   fclose(fout_imag);
00156 
00157   /* finalize the infft */
00158   solver_finalize_complex(&my_iplan);
00159 
00160   /* finalize the nfft */
00161   nfft_finalize(&my_plan);
00162 }
00163 
00164 int main(int argc, char **argv)
00165 {
00166   if (argc <= 5) {
00167     printf("usage: ./reconstruct_data_2d FILENAME N M ITER WEIGHTS\n");
00168     return 1;
00169   }
00170 
00171   reconstruct(argv[1],atoi(argv[2]),atoi(argv[3]),atoi(argv[4]),atoi(argv[5]));
00172 
00173   return 1;
00174 }
00175 /* \} */

Generated on Fri Oct 12 2012 by Doxygen 1.8.0-20120409