00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021 #include <stdlib.h>
00022 #include <math.h>
00023 #include <limits.h>
00024 #include <complex.h>
00025
00026 #include "nfft3util.h"
00027 #include "nfft3.h"
00028
00038 void reconstruct(char* filename,int N,int M,int iteration, int weight)
00039 {
00040 int j,k,l;
00041 nnfft_plan my_plan;
00042 solver_plan_complex my_iplan;
00043 FILE* fin;
00044 FILE* finh;
00045 FILE* ftime;
00046 FILE* fout_real;
00047 FILE* fout_imag;
00048 int my_N[3],my_n[3];
00049 double t,epsilon=0.0000003;
00050
00051 unsigned infft_flags = CGNR | PRECOMPUTE_DAMP;
00052 double time,min_time,max_time,min_inh,max_inh;
00053 double real,imag;
00054 double *w;
00055
00056 double Ts;
00057 double W;
00058 int N3;
00059 int m=2;
00060 double sigma = 1.25;
00061
00062 w = (double*)nfft_malloc(N*N*sizeof(double));
00063
00064 ftime=fopen("readout_time.dat","r");
00065 finh=fopen("inh.dat","r");
00066
00067 min_time=INT_MAX; max_time=INT_MIN;
00068 for(j=0;j<M;j++)
00069 {
00070 fscanf(ftime,"%le ",&time);
00071 if(time<min_time)
00072 min_time = time;
00073 if(time>max_time)
00074 max_time = time;
00075 }
00076
00077 fclose(ftime);
00078
00079 Ts=(min_time+max_time)/2.0;
00080
00081 min_inh=INT_MAX; max_inh=INT_MIN;
00082 for(j=0;j<N*N;j++)
00083 {
00084 fscanf(finh,"%le ",&w[j]);
00085 if(w[j]<min_inh)
00086 min_inh = w[j];
00087 if(w[j]>max_inh)
00088 max_inh = w[j];
00089 }
00090 fclose(finh);
00091
00092 N3=ceil((NFFT_MAX(fabs(min_inh),fabs(max_inh))*(max_time-min_time)/2.0)*4);
00093
00094
00095 W=NFFT_MAX(fabs(min_inh),fabs(max_inh))*2.0;
00096
00097 fprintf(stderr,"3: %i %e %e %e %e %e %e\n",N3,W,min_inh,max_inh,min_time,max_time,Ts);
00098
00099
00100 my_N[0]=N;my_n[0]=ceil(N*sigma);
00101 my_N[1]=N; my_n[1]=ceil(N*sigma);
00102 my_N[2]=N3; my_n[2]=ceil(N3*sigma);
00103 nnfft_init_guru(&my_plan, 3, N*N, M, my_N,my_n,m,
00104 PRE_PSI| PRE_PHI_HUT| MALLOC_X| MALLOC_V| MALLOC_F_HAT| MALLOC_F );
00105
00106
00107 if(my_plan.nnfft_flags & PRE_LIN_PSI)
00108 nnfft_precompute_lin_psi(&my_plan);
00109
00110
00111 if (weight)
00112 infft_flags = infft_flags | PRECOMPUTE_WEIGHT;
00113
00114
00115 solver_init_advanced_complex(&my_iplan,(mv_plan_complex*)(&my_plan), infft_flags );
00116
00117
00118 if(my_iplan.flags & PRECOMPUTE_WEIGHT)
00119 {
00120 fin=fopen("weights.dat","r");
00121 for(j=0;j<my_plan.M_total;j++)
00122 {
00123 fscanf(fin,"%le ",&my_iplan.w[j]);
00124 }
00125 fclose(fin);
00126 }
00127
00128
00129 if(my_iplan.flags & PRECOMPUTE_DAMP)
00130 {
00131 for(j=0;j<N;j++){
00132 for(k=0;k<N;k++) {
00133 int j2= j-N/2;
00134 int k2= k-N/2;
00135 double r=sqrt(j2*j2+k2*k2);
00136 if(r>(double) N/2)
00137 my_iplan.w_hat[j*N+k]=0.0;
00138 else
00139 my_iplan.w_hat[j*N+k]=1.0;
00140 }
00141 }
00142 }
00143
00144
00145 fin=fopen(filename,"r");
00146 ftime=fopen("readout_time.dat","r");
00147
00148 for(j=0;j<my_plan.M_total;j++)
00149 {
00150 fscanf(fin,"%le %le %le %le ",&my_plan.x[3*j+0],&my_plan.x[3*j+1],&real,&imag);
00151 my_iplan.y[j]=real+ _Complex_I*imag;
00152 fscanf(ftime,"%le ",&my_plan.x[3*j+2]);
00153
00154 my_plan.x[3*j+2] = (my_plan.x[3*j+2]-Ts)*W/N3;
00155 }
00156
00157 for(j=0;j<N;j++)
00158 {
00159 for(l=0;l<N;l++)
00160 {
00161 my_plan.v[3*(N*j+l)+0]=(((double) j) -(((double) N)/2.0))/((double) N);
00162 my_plan.v[3*(N*j+l)+1]=(((double) l) -(((double) N)/2.0))/((double) N);
00163 my_plan.v[3*(N*j+l)+2] = w[N*j+l]/W ;
00164 }
00165 }
00166
00167
00168 if(my_plan.nnfft_flags & PRE_PSI) {
00169 nnfft_precompute_psi(&my_plan);
00170 if(my_plan.nnfft_flags & PRE_FULL_PSI)
00171 nnfft_precompute_full_psi(&my_plan);
00172 }
00173
00174 if(my_plan.nnfft_flags & PRE_PHI_HUT)
00175 nnfft_precompute_phi_hut(&my_plan);
00176
00177
00178 for(k=0;k<my_plan.N_total;k++)
00179 {
00180 my_iplan.f_hat_iter[k]=0.0;
00181 }
00182
00183 t=nfft_second();
00184
00185
00186 solver_before_loop_complex(&my_iplan);
00187 for(l=0;l<iteration;l++)
00188 {
00189
00190 if(my_iplan.dot_r_iter<epsilon)
00191 break;
00192 fprintf(stderr,"%e, %i of %i\n",sqrt(my_iplan.dot_r_iter),
00193 l+1,iteration);
00194 solver_loop_one_step_complex(&my_iplan);
00195 }
00196
00197 t=nfft_second()-t;
00198 #ifdef HAVE_TOTAL_USED_MEMORY
00199 fprintf(stderr,"time: %e seconds mem: %i \n",t,nfft_total_used_memory());
00200 #else
00201 fprintf(stderr,"time: %e seconds mem: mallinfo not available\n",t);
00202 #endif
00203
00204 fout_real=fopen("output_real.dat","w");
00205 fout_imag=fopen("output_imag.dat","w");
00206
00207 for(k=0;k<my_plan.N_total;k++) {
00208
00209 my_iplan.f_hat_iter[k]*=cexp(2.0*_Complex_I*PI*Ts*w[k]);
00210
00211 fprintf(fout_real,"%le ", creal(my_iplan.f_hat_iter[k]));
00212 fprintf(fout_imag,"%le ", cimag(my_iplan.f_hat_iter[k]));
00213 }
00214
00215
00216 fclose(fout_real);
00217 fclose(fout_imag);
00218
00219
00220
00221 solver_finalize_complex(&my_iplan);
00222
00223
00224 nnfft_finalize(&my_plan);
00225
00226 nfft_free(w);
00227 }
00228
00229 int main(int argc, char **argv)
00230 {
00231 if (argc <= 5) {
00232 printf("usage: ./reconstruct_data_inh_nnfft FILENAME N M ITER WEIGHTS\n");
00233 return 1;
00234 }
00235
00236 reconstruct(argv[1],atoi(argv[2]),atoi(argv[3]),atoi(argv[4]),atoi(argv[5]));
00237
00238 return 1;
00239 }
00240