NFFT Logo 3.2.3
reconstruct_data_inh_2d1d.c
1 /*
2  * Copyright (c) 2002, 2012 Jens Keiner, Stefan Kunis, Daniel Potts
3  *
4  * This program is free software; you can redistribute it and/or modify it under
5  * the terms of the GNU General Public License as published by the Free Software
6  * Foundation; either version 2 of the License, or (at your option) any later
7  * version.
8  *
9  * This program is distributed in the hope that it will be useful, but WITHOUT
10  * ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
11  * FOR A PARTICULAR PURPOSE. See the GNU General Public License for more
12  * details.
13  *
14  * You should have received a copy of the GNU General Public License along with
15  * this program; if not, write to the Free Software Foundation, Inc., 51
16  * Franklin Street, Fifth Floor, Boston, MA 02110-1301, USA.
17  */
18 
19 /* $Id: reconstruct_data_inh_2d1d.c 3896 2012-10-10 12:19:26Z tovo $ */
20 #include "config.h"
21 
22 #include <stdlib.h>
23 #include <math.h>
24 #include <limits.h>
25 #ifdef HAVE_COMPLEX_H
26 #include <complex.h>
27 #endif
28 
29 #include "nfft3.h"
30 #include "nfft3util.h"
31 #include "infft.h"
32 
39 static void reconstruct(char* filename,int N,int M,int iteration , int weight)
40 {
41  int j,k,l;
42  double time,min_time,max_time,min_inh,max_inh;
43  ticks t0, t1;
44  double t,real,imag;
45  double w,epsilon=0.0000003; /* epsilon is a the break criterium for
46  the iteration */;
47  mri_inh_2d1d_plan my_plan;
48  solver_plan_complex my_iplan;
49  FILE* fp,*fw,*fout_real,*fout_imag,*finh,*ftime;
50  int my_N[3],my_n[3];
51  int flags = PRE_PHI_HUT| PRE_PSI |MALLOC_X| MALLOC_F_HAT|
52  MALLOC_F| FFTW_INIT| FFT_OUT_OF_PLACE;
53  unsigned infft_flags = CGNR | PRECOMPUTE_DAMP;
54 
55  double Ts;
56  double W,T;
57  int N3;
58  int m=2;
59  double sigma = 1.25;
60 
61  ftime=fopen("readout_time.dat","r");
62  finh=fopen("inh.dat","r");
63 
64  min_time=INT_MAX; max_time=INT_MIN;
65  for(j=0;j<M;j++)
66  {
67  fscanf(ftime,"%le ",&time);
68  if(time<min_time)
69  min_time = time;
70  if(time>max_time)
71  max_time = time;
72  }
73 
74  fclose(ftime);
75 
76  Ts=(min_time+max_time)/2.0;
77 
78 
79  min_inh=INT_MAX; max_inh=INT_MIN;
80  for(j=0;j<N*N;j++)
81  {
82  fscanf(finh,"%le ",&w);
83  if(w<min_inh)
84  min_inh = w;
85  if(w>max_inh)
86  max_inh = w;
87  }
88  fclose(finh);
89 
90  N3=ceil((NFFT_MAX(fabs(min_inh),fabs(max_inh))*(max_time-min_time)/2.0+(m)/(2*sigma))*4*sigma);
91  /* N3 has to be even */
92  if(N3%2!=0)
93  N3++;
94 
95  T=((max_time-min_time)/2.0)/(0.5-((double) (m))/N3);
96  W=N3/T;
97 
98  my_N[0]=N; my_n[0]=ceil(N*sigma);
99  my_N[1]=N; my_n[1]=ceil(N*sigma);
100  my_N[2]=N3; my_n[2]=N3;
101 
102  /* initialise nfft */
103  mri_inh_2d1d_init_guru(&my_plan, my_N, M, my_n, m, sigma, flags,
104  FFTW_MEASURE| FFTW_DESTROY_INPUT);
105 
106 
107  /* precompute lin psi if set */
108  if(my_plan.plan.nfft_flags & PRE_LIN_PSI)
109  nfft_precompute_lin_psi(&my_plan.plan);
110 
111  if (weight)
112  infft_flags = infft_flags | PRECOMPUTE_WEIGHT;
113 
114  /* initialise my_iplan, advanced */
115  solver_init_advanced_complex(&my_iplan,(nfft_mv_plan_complex*)(&my_plan), infft_flags );
116 
117  /* get the weights */
118  if(my_iplan.flags & PRECOMPUTE_WEIGHT)
119  {
120  fw=fopen("weights.dat","r");
121  for(j=0;j<my_plan.M_total;j++)
122  {
123  fscanf(fw,"%le ",&my_iplan.w[j]);
124  }
125  fclose(fw);
126  }
127 
128  /* get the damping factors */
129  if(my_iplan.flags & PRECOMPUTE_DAMP)
130  {
131  for(j=0;j<N;j++){
132  for(k=0;k<N;k++) {
133  int j2= j-N/2;
134  int k2= k-N/2;
135  double r=sqrt(j2*j2+k2*k2);
136  if(r>(double) N/2)
137  my_iplan.w_hat[j*N+k]=0.0;
138  else
139  my_iplan.w_hat[j*N+k]=1.0;
140  }
141  }
142  }
143 
144  fp=fopen(filename,"r");
145  ftime=fopen("readout_time.dat","r");
146 
147  for(j=0;j<my_plan.M_total;j++)
148  {
149  fscanf(fp,"%le %le %le %le",&my_plan.plan.x[2*j+0],&my_plan.plan.x[2*j+1],&real,&imag);
150  my_iplan.y[j]=real+ _Complex_I*imag;
151  fscanf(ftime,"%le ",&my_plan.t[j]);
152 
153  my_plan.t[j] = (my_plan.t[j]-Ts)/T;
154  }
155  fclose(fp);
156  fclose(ftime);
157 
158 
159  finh=fopen("inh.dat","r");
160  for(j=0;j<N*N;j++)
161  {
162  fscanf(finh,"%le ",&my_plan.w[j]);
163  my_plan.w[j]/=W;
164  }
165  fclose(finh);
166 
167 
168  if(my_plan.plan.nfft_flags & PRE_PSI) {
169  nfft_precompute_psi(&my_plan.plan);
170  }
171  if(my_plan.plan.nfft_flags & PRE_FULL_PSI) {
172  nfft_precompute_full_psi(&my_plan.plan);
173  }
174 
175  /* init some guess */
176  for(j=0;j<my_plan.N_total;j++)
177  {
178  my_iplan.f_hat_iter[j]=0.0;
179  }
180 
181  t0 = getticks();
182 
183  /* inverse trafo */
184  solver_before_loop_complex(&my_iplan);
185  for(l=0;l<iteration;l++)
186  {
187  /* break if dot_r_iter is smaller than epsilon*/
188  if(my_iplan.dot_r_iter<epsilon)
189  break;
190  fprintf(stderr,"%e, %i of %i\n",sqrt(my_iplan.dot_r_iter),
191  l+1,iteration);
192  solver_loop_one_step_complex(&my_iplan);
193  }
194 
195  t1 = getticks();
196  t = nfft_elapsed_seconds(t1,t0);
197 
198  fout_real=fopen("output_real.dat","w");
199  fout_imag=fopen("output_imag.dat","w");
200 
201  for (j=0;j<N*N;j++) {
202  /* Verschiebung wieder herausrechnen */
203  my_iplan.f_hat_iter[j]*=cexp(-2.0*_Complex_I*PI*Ts*my_plan.w[j]*W);
204 
205  fprintf(fout_real,"%le ",creal(my_iplan.f_hat_iter[j]));
206  fprintf(fout_imag,"%le ",cimag(my_iplan.f_hat_iter[j]));
207  }
208 
209  fclose(fout_real);
210  fclose(fout_imag);
211  solver_finalize_complex(&my_iplan);
212  mri_inh_2d1d_finalize(&my_plan);
213 }
214 
215 
216 int main(int argc, char **argv)
217 {
218  if (argc <= 5) {
219 
220  printf("usage: ./reconstruct_data_inh_2d1d FILENAME N M ITER WEIGHTS\n");
221  return 1;
222  }
223 
224  reconstruct(argv[1],atoi(argv[2]),atoi(argv[3]),atoi(argv[4]),atoi(argv[5]));
225 
226  return 1;
227 }
228 /* \} */

Generated on Tue Apr 30 2013 by Doxygen 1.8.1