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 <complex.h>
00024
00025 #include "nfft3util.h"
00026 #include "nfft3.h"
00027
00037 void reconstruct(char* filename,int N,int M,int Z,int iteration, int weight, fftw_complex *mem)
00038 {
00039 int j,k,l,z;
00040 double real,imag;
00041 nfft_plan my_plan;
00042 solver_plan_complex my_iplan;
00043 FILE* fin;
00044 int my_N[2],my_n[2];
00045 double tmp, epsilon=0.0000003;
00046
00047
00048 unsigned infft_flags = CGNR | PRECOMPUTE_DAMP;
00049
00050
00051 my_N[0]=N;my_n[0]=ceil(N*1.2);
00052 my_N[1]=N; my_n[1]=ceil(N*1.2);
00053 nfft_init_guru(&my_plan, 2, my_N, M/Z, my_n, 6, PRE_PHI_HUT| PRE_PSI|
00054 MALLOC_X| MALLOC_F_HAT| MALLOC_F|
00055 FFTW_INIT| FFT_OUT_OF_PLACE,
00056 FFTW_MEASURE| FFTW_DESTROY_INPUT);
00057
00058
00059 if(my_plan.nfft_flags & PRE_LIN_PSI)
00060 nfft_precompute_lin_psi(&my_plan);
00061
00062
00063 if (weight)
00064 infft_flags = infft_flags | PRECOMPUTE_WEIGHT;
00065
00066
00067 solver_init_advanced_complex(&my_iplan,(mv_plan_complex*)(&my_plan), infft_flags );
00068
00069
00070 if(my_iplan.flags & PRECOMPUTE_WEIGHT)
00071 {
00072 fin=fopen("weights.dat","r");
00073 for(j=0;j<my_plan.M_total;j++)
00074 {
00075 fscanf(fin,"%le ",&my_iplan.w[j]);
00076 }
00077 fclose(fin);
00078 }
00079
00080
00081 if(my_iplan.flags & PRECOMPUTE_DAMP)
00082 {
00083 for(j=0;j<N;j++){
00084 for(k=0;k<N;k++) {
00085 int j2= j-N/2;
00086 int k2= k-N/2;
00087 double r=sqrt(j2*j2+k2*k2);
00088 if(r>(double) N/2)
00089 my_iplan.w_hat[j*N+k]=0.0;
00090 else
00091 my_iplan.w_hat[j*N+k]=1.0;
00092 }
00093 }
00094 }
00095
00096
00097 fin=fopen(filename,"r");
00098
00099
00100 for(z=0;z<Z;z++) {
00101
00102
00103 for(j=0;j<my_plan.M_total;j++)
00104 {
00105 fscanf(fin,"%le %le %le %le %le ",&my_plan.x[2*j+0],&my_plan.x[2*j+1], &tmp,
00106 &real,&imag);
00107 my_iplan.y[j] = real + _Complex_I*imag;
00108 }
00109
00110
00111 if(z==0 && my_plan.nfft_flags & PRE_PSI)
00112 nfft_precompute_psi(&my_plan);
00113
00114
00115 if(z==0 && my_plan.nfft_flags & PRE_FULL_PSI)
00116 nfft_precompute_full_psi(&my_plan);
00117
00118
00119 for(k=0;k<my_plan.N_total;k++)
00120 my_iplan.f_hat_iter[k]=0.0;
00121
00122
00123 solver_before_loop_complex(&my_iplan);
00124 for(l=0;l<iteration;l++)
00125 {
00126
00127 if(my_iplan.dot_r_iter<epsilon)
00128 break;
00129 fprintf(stderr,"%e, %i of %i\n",sqrt(my_iplan.dot_r_iter),
00130 iteration*z+l+1,iteration*Z);
00131 solver_loop_one_step_complex(&my_iplan);
00132 }
00133 for(k=0;k<my_plan.N_total;k++) {
00134
00135
00136 mem[(Z*N*N/2+z*N*N+ k)%(Z*N*N)] = my_iplan.f_hat_iter[k];
00137 }
00138 }
00139
00140 fclose(fin);
00141
00142
00143 solver_finalize_complex(&my_iplan);
00144
00145
00146 nfft_finalize(&my_plan);
00147 }
00148
00153 void print(int N,int M,int Z, fftw_complex *mem)
00154 {
00155 int i,j;
00156 FILE* fout_real;
00157 FILE* fout_imag;
00158 fout_real=fopen("output_real.dat","w");
00159 fout_imag=fopen("output_imag.dat","w");
00160
00161 for(i=0;i<Z;i++) {
00162 for (j=0;j<N*N;j++) {
00163 fprintf(fout_real,"%le ",creal(mem[(Z*N*N/2+i*N*N+ j)%(Z*N*N)]) /Z);
00164 fprintf(fout_imag,"%le ",cimag(mem[(Z*N*N/2+i*N*N+ j)%(Z*N*N)]) /Z);
00165 }
00166 fprintf(fout_real,"\n");
00167 fprintf(fout_imag,"\n");
00168 }
00169
00170 fclose(fout_real);
00171 fclose(fout_imag);
00172 }
00173
00174 int main(int argc, char **argv)
00175 {
00176 fftw_complex *mem;
00177 fftw_plan plan;
00178 int N,M,Z;
00179
00180 if (argc <= 6) {
00181 printf("usage: ./reconstruct FILENAME N M Z ITER WEIGHTS\n");
00182 return 1;
00183 }
00184
00185 N=atoi(argv[2]);
00186 M=atoi(argv[3]);
00187 Z=atoi(argv[4]);
00188
00189
00190
00191 mem = (fftw_complex*) nfft_malloc(sizeof(fftw_complex) * atoi(argv[2]) * atoi(argv[2]) * atoi(argv[4]));
00192
00193
00194 plan = fftw_plan_many_dft(1, &Z, N*N,
00195 mem, NULL,
00196 N*N, 1,
00197 mem, NULL,
00198 N*N,1 ,
00199 FFTW_BACKWARD, FFTW_MEASURE);
00200
00201
00202 reconstruct(argv[1],N,M,Z,atoi(argv[5]),atoi(argv[6]),mem);
00203
00204
00205 fftw_execute(plan);
00206
00207
00208 print(N,M,Z, mem);
00209
00210
00211 nfft_free(mem);
00212 fftw_destroy_plan(plan);
00213 return 1;
00214 }
00215