NFFT Logo 3.1.2 API Reference

reconstruct_data_2d.c

00001 /*
00002  * Copyright (c) 2002, 2009 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 3198 2009-05-27 14:16:50Z keiner $ */
00020 
00021 #include <math.h>
00022 #include <stdlib.h>
00023 #include <complex.h>
00024 
00025 #include "nfft3util.h"
00026 #include "nfft3.h"
00027 
00037 void reconstruct(char* filename,int N,int M,int iteration, int weight)
00038 {
00039   int j,k,l;                    /* some variables  */
00040   double real,imag,t;           /* to read the real and imag part of a complex number */
00041   nfft_plan my_plan;            /* plan for the two dimensional nfft  */
00042   solver_plan_complex my_iplan; /* plan for the two dimensional infft */
00043   FILE* fin;                    /* input file                         */
00044   FILE* fout_real;              /* output file                        */
00045   FILE* fout_imag;              /* output file                        */
00046   int my_N[2],my_n[2];          /* to init the nfft */
00047   double epsilon=0.0000003;     /* epsilon is a the break criterium for
00048                                    the iteration */
00049   unsigned infft_flags = CGNR | PRECOMPUTE_DAMP;  /* flags for the infft*/
00050   int m = 6;
00051   double alpha = 2.0;
00052   /* initialise my_plan */
00053   my_N[0]=N; my_n[0]=ceil(N*alpha);
00054   my_N[1]=N; my_n[1]=ceil(N*alpha);
00055   nfft_init_guru(&my_plan, 2, my_N, M, my_n, m, PRE_PHI_HUT| PRE_PSI|
00056                          MALLOC_X| MALLOC_F_HAT| MALLOC_F|
00057                          FFTW_INIT| FFT_OUT_OF_PLACE,
00058                          FFTW_MEASURE| FFTW_DESTROY_INPUT);
00059 
00060   /* precompute lin psi if set */
00061   if(my_plan.nfft_flags & PRE_LIN_PSI)
00062     nfft_precompute_lin_psi(&my_plan);
00063 
00064   /* set the flags for the infft*/
00065   if (weight)
00066     infft_flags = infft_flags | PRECOMPUTE_WEIGHT;
00067 
00068   /* initialise my_iplan, advanced */
00069   solver_init_advanced_complex(&my_iplan,(mv_plan_complex*)&my_plan, infft_flags );
00070 
00071   /* get the weights */
00072   if(my_iplan.flags & PRECOMPUTE_WEIGHT)
00073   {
00074     fin=fopen("weights.dat","r");
00075     for(j=0;j<my_plan.M_total;j++)
00076     {
00077         fscanf(fin,"%le ",&my_iplan.w[j]);
00078     }
00079     fclose(fin);
00080   }
00081 
00082   /* get the damping factors */
00083   if(my_iplan.flags & PRECOMPUTE_DAMP)
00084   {
00085     for(j=0;j<N;j++){
00086       for(k=0;k<N;k++) {
00087         int j2= j-N/2;
00088         int k2= k-N/2;
00089         double r=sqrt(j2*j2+k2*k2);
00090         if(r>(double) N/2)
00091           my_iplan.w_hat[j*N+k]=0.0;
00092         else
00093           my_iplan.w_hat[j*N+k]=1.0;
00094       }
00095     }
00096   }
00097 
00098   /* open the input file */
00099   fin=fopen(filename,"r");
00100 
00101   /* read x,y,freal and fimag from the knots */
00102   for(j=0;j<my_plan.M_total;j++)
00103   {
00104     fscanf(fin,"%le %le %le %le ",&my_plan.x[2*j+0],&my_plan.x[2*j+1],
00105     &real,&imag);
00106     my_iplan.y[j] = real + _Complex_I*imag;
00107   }
00108 
00109   fclose(fin);
00110 
00111   /* precompute psi */
00112   if(my_plan.nfft_flags & PRE_PSI)
00113     nfft_precompute_psi(&my_plan);
00114 
00115   /* precompute full psi */
00116   if(my_plan.nfft_flags & PRE_FULL_PSI)
00117       nfft_precompute_full_psi(&my_plan);
00118 
00119   /* init some guess */
00120   for(k=0;k<my_plan.N_total;k++)
00121     my_iplan.f_hat_iter[k]=0.0;
00122 
00123   t=nfft_second();
00124 
00125   /* inverse trafo */
00126   solver_before_loop_complex(&my_iplan);
00127   for(l=0;l<iteration;l++)
00128   {
00129     /* break if dot_r_iter is smaller than epsilon*/
00130     if(my_iplan.dot_r_iter<epsilon)
00131       break;
00132     fprintf(stderr,"%e,  %i of %i\n",sqrt(my_iplan.dot_r_iter),
00133     l+1,iteration);
00134     solver_loop_one_step_complex(&my_iplan);
00135   }
00136 
00137 
00138   t=nfft_second()-t;
00139 #ifdef HAVE_MALLINFO
00140   fprintf(stderr,"time: %e seconds mem: %i \n",t,nfft_total_used_memory());
00141 #else
00142   fprintf(stderr,"time: %e seconds mem: mallinfo not available\n",t);
00143 #endif
00144 
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 16 Sep 2009 by Doxygen 1.5.3