NFFT Logo 3.2.2
reconstruct_data_2d1d.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_2d1d.c 3896 2012-10-10 12:19:26Z tovo $ */
00020 #include "config.h"
00021 
00022 #include <stdlib.h>
00023 #include <math.h>
00024 #ifdef HAVE_COMPLEX_H
00025 #include <complex.h>
00026 #endif
00027 
00028 #include "nfft3util.h"
00029 #include "nfft3.h"
00030 
00040 static void reconstruct(char* filename,int N,int M,int Z,int iteration, int weight, fftw_complex *mem)
00041 {
00042   int j,k,l,z;                  /* some variables  */
00043   double real,imag;             /* to read the real and imag part of a complex number */
00044   nfft_plan my_plan;            /* plan for the two dimensional nfft  */
00045   solver_plan_complex my_iplan; /* plan for the two dimensional infft */
00046   FILE* fin;                    /* input file */
00047   int my_N[2],my_n[2];          /* to init the nfft */
00048   double tmp, epsilon=0.0000003;/* tmp to read the obsolent z from the input file
00049                                    epsilon is the break criterium for
00050                                    the iteration */
00051   unsigned infft_flags = CGNR | PRECOMPUTE_DAMP; /* flags for the infft */
00052 
00053   /* initialise my_plan */
00054   my_N[0]=N;my_n[0]=ceil(N*1.2);
00055   my_N[1]=N; my_n[1]=ceil(N*1.2);
00056   nfft_init_guru(&my_plan, 2, my_N, M/Z, my_n, 6, PRE_PHI_HUT| PRE_PSI|
00057                          MALLOC_X| MALLOC_F_HAT| MALLOC_F|
00058                         FFTW_INIT| FFT_OUT_OF_PLACE,
00059                         FFTW_MEASURE| FFTW_DESTROY_INPUT);
00060 
00061   /* precompute lin psi if set */
00062   if(my_plan.nfft_flags & PRE_LIN_PSI)
00063     nfft_precompute_lin_psi(&my_plan);
00064 
00065   /* set the flags for the infft*/
00066   if (weight)
00067     infft_flags = infft_flags | PRECOMPUTE_WEIGHT;
00068 
00069   /* initialise my_iplan, advanced */
00070   solver_init_advanced_complex(&my_iplan,(nfft_mv_plan_complex*)(&my_plan), infft_flags );
00071 
00072   /* get the weights */
00073   if(my_iplan.flags & PRECOMPUTE_WEIGHT)
00074   {
00075     fin=fopen("weights.dat","r");
00076     for(j=0;j<my_plan.M_total;j++)
00077     {
00078         fscanf(fin,"%le ",&my_iplan.w[j]);
00079     }
00080     fclose(fin);
00081   }
00082 
00083   /* get the damping factors */
00084   if(my_iplan.flags & PRECOMPUTE_DAMP)
00085   {
00086     for(j=0;j<N;j++){
00087       for(k=0;k<N;k++) {
00088         int j2= j-N/2;
00089         int k2= k-N/2;
00090         double r=sqrt(j2*j2+k2*k2);
00091         if(r>(double) N/2)
00092           my_iplan.w_hat[j*N+k]=0.0;
00093         else
00094           my_iplan.w_hat[j*N+k]=1.0;
00095       }
00096     }
00097   }
00098 
00099   /* open the input file */
00100   fin=fopen(filename,"r");
00101 
00102   /* For every Layer*/
00103   for(z=0;z<Z;z++) {
00104 
00105     /* read x,y,freal and fimag from the knots */
00106     for(j=0;j<my_plan.M_total;j++)
00107     {
00108       fscanf(fin,"%le %le %le %le %le ",&my_plan.x[2*j+0],&my_plan.x[2*j+1], &tmp,
00109       &real,&imag);
00110       my_iplan.y[j] = real + _Complex_I*imag;
00111     }
00112 
00113     /* precompute psi if set just one time because the knots equal each plane */
00114     if(z==0 && my_plan.nfft_flags & PRE_PSI)
00115       nfft_precompute_psi(&my_plan);
00116 
00117     /* precompute full psi if set just one time because the knots equal each plane */
00118     if(z==0 && my_plan.nfft_flags & PRE_FULL_PSI)
00119       nfft_precompute_full_psi(&my_plan);
00120 
00121     /* init some guess */
00122     for(k=0;k<my_plan.N_total;k++)
00123       my_iplan.f_hat_iter[k]=0.0;
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       iteration*z+l+1,iteration*Z);
00134       solver_loop_one_step_complex(&my_iplan);
00135     }
00136     for(k=0;k<my_plan.N_total;k++) {
00137       /* write every slice in the memory.
00138       here we make an fftshift direct */
00139       mem[(Z*N*N/2+z*N*N+ k)%(Z*N*N)] = my_iplan.f_hat_iter[k];
00140     }
00141   }
00142 
00143   fclose(fin);
00144 
00145   /* finalize the infft */
00146   solver_finalize_complex(&my_iplan);
00147 
00148   /* finalize the nfft */
00149   nfft_finalize(&my_plan);
00150 }
00151 
00156 static void print(int N,int M,int Z, fftw_complex *mem)
00157 {
00158   int i,j;
00159   FILE* fout_real;
00160   FILE* fout_imag;
00161   fout_real=fopen("output_real.dat","w");
00162   fout_imag=fopen("output_imag.dat","w");
00163 
00164   for(i=0;i<Z;i++) {
00165     for (j=0;j<N*N;j++) {
00166       fprintf(fout_real,"%le ",creal(mem[(Z*N*N/2+i*N*N+ j)%(Z*N*N)]) /Z);
00167       fprintf(fout_imag,"%le ",cimag(mem[(Z*N*N/2+i*N*N+ j)%(Z*N*N)]) /Z);
00168     }
00169     fprintf(fout_real,"\n");
00170     fprintf(fout_imag,"\n");
00171   }
00172 
00173   fclose(fout_real);
00174   fclose(fout_imag);
00175 }
00176 
00177 int main(int argc, char **argv)
00178 {
00179   fftw_complex *mem;
00180   fftw_plan plan;
00181   int N,M,Z;
00182 
00183   if (argc <= 6) {
00184     printf("usage: ./reconstruct FILENAME N M Z ITER WEIGHTS\n");
00185     return 1;
00186   }
00187 
00188   N=atoi(argv[2]);
00189   M=atoi(argv[3]);
00190   Z=atoi(argv[4]);
00191 
00192   /* Allocate memory to hold every layer in memory after the
00193   2D-infft */
00194   mem = (fftw_complex*) nfft_malloc(sizeof(fftw_complex) * atoi(argv[2]) * atoi(argv[2]) * atoi(argv[4]));
00195 
00196   /* Create plan for the 1d-ifft */
00197   plan = fftw_plan_many_dft(1, &Z, N*N,
00198                                   mem, NULL,
00199                                   N*N, 1,
00200                                   mem, NULL,
00201                                   N*N,1 ,
00202                                   FFTW_BACKWARD, FFTW_MEASURE);
00203 
00204   /* execute the 2d-infft's */
00205   reconstruct(argv[1],N,M,Z,atoi(argv[5]),atoi(argv[6]),mem);
00206 
00207   /* execute the 1d-fft's */
00208   fftw_execute(plan);
00209 
00210   /* write the memory back in files */
00211   print(N,M,Z, mem);
00212 
00213   /* free memory */
00214   nfft_free(mem);
00215   fftw_destroy_plan(plan);
00216   return 1;
00217 }
00218 /* \} */

Generated on Fri Oct 12 2012 by Doxygen 1.8.0-20120409