![]() |
NFFT
3.3.2
|
00001 /* 00002 * Copyright (c) 2002, 2016 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 #include "config.h" 00019 00020 #include <math.h> 00021 #include <stdlib.h> 00022 #ifdef HAVE_COMPLEX_H 00023 #include <complex.h> 00024 #endif 00025 00026 #include "nfft3.h" 00027 00037 static void reconstruct(char* filename,int N,int M,int Z,int iteration, int weight) 00038 { 00039 int j,k,z,l; /* some variables */ 00040 double real,imag; /* 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 (real part) */ 00045 FILE* fout_imag; /* output file (imag part) */ 00046 int my_N[3],my_n[3]; /* to init the nfft */ 00047 double epsilon=0.0000003; /* tmp to read the obsolent z from 700.acs 00048 epsilon is a the break criterion for 00049 the iteration */ 00050 unsigned infft_flags = CGNR | PRECOMPUTE_DAMP; /* flags for the infft */ 00051 00052 /* initialise my_plan, specific. 00053 we don't precompute psi */ 00054 my_N[0]=Z; my_n[0]=ceil(Z*1.2); 00055 my_N[1]=N; my_n[1]=ceil(N*1.2); 00056 my_N[2]=N; my_n[2]=ceil(N*1.2); 00057 nfft_init_guru(&my_plan, 3, my_N, M, my_n, 6, 00058 PRE_PHI_HUT| PRE_PSI |MALLOC_X| MALLOC_F_HAT| 00059 MALLOC_F| FFTW_INIT| FFT_OUT_OF_PLACE, 00060 FFTW_MEASURE| FFTW_DESTROY_INPUT); 00061 00062 /* precompute lin psi */ 00063 if(my_plan.flags & PRE_LIN_PSI) 00064 nfft_precompute_lin_psi(&my_plan); 00065 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<M;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 for(z=0;z<N;z++) { 00089 int j2= j-N/2; 00090 int k2= k-N/2; 00091 int z2= z-N/2; 00092 double r=sqrt(j2*j2+k2*k2+z2*z2); 00093 if(r>(double) N/2) 00094 my_iplan.w_hat[z*N*N+j*N+k]=0.0; 00095 else 00096 my_iplan.w_hat[z*N*N+j*N+k]=1.0; 00097 } 00098 } 00099 } 00100 } 00101 00102 /* open the input file */ 00103 fin=fopen(filename,"r"); 00104 00105 /* open the output files */ 00106 fout_real=fopen("output_real.dat","w"); 00107 fout_imag=fopen("output_imag.dat","w"); 00108 00109 /* read x,y,freal and fimag from the knots */ 00110 for(j=0;j<M;j++) 00111 { 00112 fscanf(fin,"%le %le %le %le %le ",&my_plan.x[3*j+1],&my_plan.x[3*j+2], &my_plan.x[3*j+0], 00113 &real,&imag); 00114 my_iplan.y[j] = real + _Complex_I*imag; 00115 } 00116 00117 /* precompute psi */ 00118 if(my_plan.flags & PRE_PSI) 00119 nfft_precompute_psi(&my_plan); 00120 00121 /* precompute full psi */ 00122 if(my_plan.flags & PRE_FULL_PSI) 00123 nfft_precompute_full_psi(&my_plan); 00124 00125 /* init some guess */ 00126 for(k=0;k<my_plan.N_total;k++) 00127 my_iplan.f_hat_iter[k]=0.0; 00128 00129 /* inverse trafo */ 00130 solver_before_loop_complex(&my_iplan); 00131 for(l=0;l<iteration;l++) 00132 { 00133 /* break if dot_r_iter is smaller than epsilon*/ 00134 if(my_iplan.dot_r_iter<epsilon) 00135 break; 00136 fprintf(stderr,"%e, %i of %i\n",sqrt(my_iplan.dot_r_iter), 00137 l+1,iteration); 00138 solver_loop_one_step_complex(&my_iplan); 00139 } 00140 00141 for(l=0;l<Z;l++) 00142 { 00143 for(k=0;k<N*N;k++) 00144 { 00145 /* write every Layer in the files */ 00146 fprintf(fout_real,"%le ",creal(my_iplan.f_hat_iter[ k+N*N*l ])); 00147 fprintf(fout_imag,"%le ",cimag(my_iplan.f_hat_iter[ k+N*N*l ])); 00148 } 00149 fprintf(fout_real,"\n"); 00150 fprintf(fout_imag,"\n"); 00151 } 00152 00153 fclose(fout_real); 00154 fclose(fout_imag); 00155 00156 solver_finalize_complex(&my_iplan); 00157 nfft_finalize(&my_plan); 00158 } 00159 00160 int main(int argc, char **argv) 00161 { 00162 if (argc <= 6) { 00163 printf("usage: ./reconstruct3D FILENAME N M Z ITER WEIGHTS\n"); 00164 return 1; 00165 } 00166 00167 reconstruct(argv[1],atoi(argv[2]),atoi(argv[3]),atoi(argv[4]),atoi(argv[5]),atoi(argv[6])); 00168 return 1; 00169 } 00170 /* \} */