NFFT Logo 3.2.3
reconstruct_data_2d.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_2d.c 3896 2012-10-10 12:19:26Z tovo $ */
20 #include "config.h"
21 
22 #include <math.h>
23 #include <stdlib.h>
24 #ifdef HAVE_COMPLEX_H
25 #include <complex.h>
26 #endif
27 
28 #include "nfft3util.h"
29 #include "nfft3.h"
30 #include "infft.h"
31 
41 static void reconstruct(char* filename,int N,int M,int iteration, int weight)
42 {
43  int j,k,l; /* some variables */
44  ticks t0, t1;
45  double real,imag,t; /* to read the real and imag part of a complex number */
46  nfft_plan my_plan; /* plan for the two dimensional nfft */
47  solver_plan_complex my_iplan; /* plan for the two dimensional infft */
48  FILE* fin; /* input file */
49  FILE* fout_real; /* output file */
50  FILE* fout_imag; /* output file */
51  int my_N[2],my_n[2]; /* to init the nfft */
52  double epsilon=0.0000003; /* epsilon is a the break criterium for
53  the iteration */
54  unsigned infft_flags = CGNR | PRECOMPUTE_DAMP; /* flags for the infft*/
55  int m = 6;
56  double alpha = 2.0;
57  /* initialise my_plan */
58  my_N[0]=N; my_n[0]=ceil(N*alpha);
59  my_N[1]=N; my_n[1]=ceil(N*alpha);
60  nfft_init_guru(&my_plan, 2, my_N, M, my_n, m, PRE_PHI_HUT| PRE_PSI|
61  MALLOC_X| MALLOC_F_HAT| MALLOC_F|
62  FFTW_INIT| FFT_OUT_OF_PLACE,
63  FFTW_MEASURE| FFTW_DESTROY_INPUT);
64 
65  /* precompute lin psi if set */
66  if(my_plan.nfft_flags & PRE_LIN_PSI)
67  nfft_precompute_lin_psi(&my_plan);
68 
69  /* set the flags for the infft*/
70  if (weight)
71  infft_flags = infft_flags | PRECOMPUTE_WEIGHT;
72 
73  /* initialise my_iplan, advanced */
74  solver_init_advanced_complex(&my_iplan,(nfft_mv_plan_complex*)&my_plan, infft_flags );
75 
76  /* get the weights */
77  if(my_iplan.flags & PRECOMPUTE_WEIGHT)
78  {
79  fin=fopen("weights.dat","r");
80  for(j=0;j<my_plan.M_total;j++)
81  {
82  fscanf(fin,"%le ",&my_iplan.w[j]);
83  }
84  fclose(fin);
85  }
86 
87  /* get the damping factors */
88  if(my_iplan.flags & PRECOMPUTE_DAMP)
89  {
90  for(j=0;j<N;j++){
91  for(k=0;k<N;k++) {
92  int j2= j-N/2;
93  int k2= k-N/2;
94  double r=sqrt(j2*j2+k2*k2);
95  if(r>(double) N/2)
96  my_iplan.w_hat[j*N+k]=0.0;
97  else
98  my_iplan.w_hat[j*N+k]=1.0;
99  }
100  }
101  }
102 
103  /* open the input file */
104  fin=fopen(filename,"r");
105 
106  /* read x,y,freal and fimag from the knots */
107  for(j=0;j<my_plan.M_total;j++)
108  {
109  fscanf(fin,"%le %le %le %le ",&my_plan.x[2*j+0],&my_plan.x[2*j+1],
110  &real,&imag);
111  my_iplan.y[j] = real + _Complex_I*imag;
112  }
113 
114  fclose(fin);
115 
116  /* precompute psi */
117  if(my_plan.nfft_flags & PRE_PSI)
118  nfft_precompute_psi(&my_plan);
119 
120  /* precompute full psi */
121  if(my_plan.nfft_flags & PRE_FULL_PSI)
122  nfft_precompute_full_psi(&my_plan);
123 
124  /* init some guess */
125  for(k=0;k<my_plan.N_total;k++)
126  my_iplan.f_hat_iter[k]=0.0;
127 
128  t0 = getticks();
129 
130  /* inverse trafo */
131  solver_before_loop_complex(&my_iplan);
132  for(l=0;l<iteration;l++)
133  {
134  /* break if dot_r_iter is smaller than epsilon*/
135  if(my_iplan.dot_r_iter<epsilon)
136  break;
137  fprintf(stderr,"%e, %i of %i\n",sqrt(my_iplan.dot_r_iter),
138  l+1,iteration);
139  solver_loop_one_step_complex(&my_iplan);
140  }
141 
142 
143  t1 = getticks();
144  t=nfft_elapsed_seconds(t1,t0);
145 
146  fout_real=fopen("output_real.dat","w");
147  fout_imag=fopen("output_imag.dat","w");
148 
149  for(k=0;k<my_plan.N_total;k++) {
150  fprintf(fout_real,"%le ", creal(my_iplan.f_hat_iter[k]));
151  fprintf(fout_imag,"%le ", cimag(my_iplan.f_hat_iter[k]));
152  }
153 
154  fclose(fout_real);
155  fclose(fout_imag);
156 
157  /* finalize the infft */
158  solver_finalize_complex(&my_iplan);
159 
160  /* finalize the nfft */
161  nfft_finalize(&my_plan);
162 }
163 
164 int main(int argc, char **argv)
165 {
166  if (argc <= 5) {
167  printf("usage: ./reconstruct_data_2d FILENAME N M ITER WEIGHTS\n");
168  return 1;
169  }
170 
171  reconstruct(argv[1],atoi(argv[2]),atoi(argv[3]),atoi(argv[4]),atoi(argv[5]));
172 
173  return 1;
174 }
175 /* \} */

Generated on Tue Apr 30 2013 by Doxygen 1.8.1