NFFT Logo 3.2.3
reconstruct_data_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_2d1d.c 3896 2012-10-10 12:19:26Z tovo $ */
20 #include "config.h"
21 
22 #include <stdlib.h>
23 #include <math.h>
24 #ifdef HAVE_COMPLEX_H
25 #include <complex.h>
26 #endif
27 
28 #include "nfft3util.h"
29 #include "nfft3.h"
30 
40 static void reconstruct(char* filename,int N,int M,int Z,int iteration, int weight, fftw_complex *mem)
41 {
42  int j,k,l,z; /* some variables */
43  double real,imag; /* to read the real and imag part of a complex number */
44  nfft_plan my_plan; /* plan for the two dimensional nfft */
45  solver_plan_complex my_iplan; /* plan for the two dimensional infft */
46  FILE* fin; /* input file */
47  int my_N[2],my_n[2]; /* to init the nfft */
48  double tmp, epsilon=0.0000003;/* tmp to read the obsolent z from the input file
49  epsilon is the break criterium for
50  the iteration */
51  unsigned infft_flags = CGNR | PRECOMPUTE_DAMP; /* flags for the infft */
52 
53  /* initialise my_plan */
54  my_N[0]=N;my_n[0]=ceil(N*1.2);
55  my_N[1]=N; my_n[1]=ceil(N*1.2);
56  nfft_init_guru(&my_plan, 2, my_N, M/Z, my_n, 6, PRE_PHI_HUT| PRE_PSI|
57  MALLOC_X| MALLOC_F_HAT| MALLOC_F|
58  FFTW_INIT| FFT_OUT_OF_PLACE,
59  FFTW_MEASURE| FFTW_DESTROY_INPUT);
60 
61  /* precompute lin psi if set */
62  if(my_plan.nfft_flags & PRE_LIN_PSI)
63  nfft_precompute_lin_psi(&my_plan);
64 
65  /* set the flags for the infft*/
66  if (weight)
67  infft_flags = infft_flags | PRECOMPUTE_WEIGHT;
68 
69  /* initialise my_iplan, advanced */
70  solver_init_advanced_complex(&my_iplan,(nfft_mv_plan_complex*)(&my_plan), infft_flags );
71 
72  /* get the weights */
73  if(my_iplan.flags & PRECOMPUTE_WEIGHT)
74  {
75  fin=fopen("weights.dat","r");
76  for(j=0;j<my_plan.M_total;j++)
77  {
78  fscanf(fin,"%le ",&my_iplan.w[j]);
79  }
80  fclose(fin);
81  }
82 
83  /* get the damping factors */
84  if(my_iplan.flags & PRECOMPUTE_DAMP)
85  {
86  for(j=0;j<N;j++){
87  for(k=0;k<N;k++) {
88  int j2= j-N/2;
89  int k2= k-N/2;
90  double r=sqrt(j2*j2+k2*k2);
91  if(r>(double) N/2)
92  my_iplan.w_hat[j*N+k]=0.0;
93  else
94  my_iplan.w_hat[j*N+k]=1.0;
95  }
96  }
97  }
98 
99  /* open the input file */
100  fin=fopen(filename,"r");
101 
102  /* For every Layer*/
103  for(z=0;z<Z;z++) {
104 
105  /* read x,y,freal and fimag from the knots */
106  for(j=0;j<my_plan.M_total;j++)
107  {
108  fscanf(fin,"%le %le %le %le %le ",&my_plan.x[2*j+0],&my_plan.x[2*j+1], &tmp,
109  &real,&imag);
110  my_iplan.y[j] = real + _Complex_I*imag;
111  }
112 
113  /* precompute psi if set just one time because the knots equal each plane */
114  if(z==0 && my_plan.nfft_flags & PRE_PSI)
115  nfft_precompute_psi(&my_plan);
116 
117  /* precompute full psi if set just one time because the knots equal each plane */
118  if(z==0 && my_plan.nfft_flags & PRE_FULL_PSI)
119  nfft_precompute_full_psi(&my_plan);
120 
121  /* init some guess */
122  for(k=0;k<my_plan.N_total;k++)
123  my_iplan.f_hat_iter[k]=0.0;
124 
125  /* inverse trafo */
126  solver_before_loop_complex(&my_iplan);
127  for(l=0;l<iteration;l++)
128  {
129  /* break if dot_r_iter is smaller than epsilon*/
130  if(my_iplan.dot_r_iter<epsilon)
131  break;
132  fprintf(stderr,"%e, %i of %i\n",sqrt(my_iplan.dot_r_iter),
133  iteration*z+l+1,iteration*Z);
134  solver_loop_one_step_complex(&my_iplan);
135  }
136  for(k=0;k<my_plan.N_total;k++) {
137  /* write every slice in the memory.
138  here we make an fftshift direct */
139  mem[(Z*N*N/2+z*N*N+ k)%(Z*N*N)] = my_iplan.f_hat_iter[k];
140  }
141  }
142 
143  fclose(fin);
144 
145  /* finalize the infft */
146  solver_finalize_complex(&my_iplan);
147 
148  /* finalize the nfft */
149  nfft_finalize(&my_plan);
150 }
151 
156 static void print(int N,int M,int Z, fftw_complex *mem)
157 {
158  int i,j;
159  FILE* fout_real;
160  FILE* fout_imag;
161  fout_real=fopen("output_real.dat","w");
162  fout_imag=fopen("output_imag.dat","w");
163 
164  for(i=0;i<Z;i++) {
165  for (j=0;j<N*N;j++) {
166  fprintf(fout_real,"%le ",creal(mem[(Z*N*N/2+i*N*N+ j)%(Z*N*N)]) /Z);
167  fprintf(fout_imag,"%le ",cimag(mem[(Z*N*N/2+i*N*N+ j)%(Z*N*N)]) /Z);
168  }
169  fprintf(fout_real,"\n");
170  fprintf(fout_imag,"\n");
171  }
172 
173  fclose(fout_real);
174  fclose(fout_imag);
175 }
176 
177 int main(int argc, char **argv)
178 {
179  fftw_complex *mem;
180  fftw_plan plan;
181  int N,M,Z;
182 
183  if (argc <= 6) {
184  printf("usage: ./reconstruct FILENAME N M Z ITER WEIGHTS\n");
185  return 1;
186  }
187 
188  N=atoi(argv[2]);
189  M=atoi(argv[3]);
190  Z=atoi(argv[4]);
191 
192  /* Allocate memory to hold every layer in memory after the
193  2D-infft */
194  mem = (fftw_complex*) nfft_malloc(sizeof(fftw_complex) * atoi(argv[2]) * atoi(argv[2]) * atoi(argv[4]));
195 
196  /* Create plan for the 1d-ifft */
197  plan = fftw_plan_many_dft(1, &Z, N*N,
198  mem, NULL,
199  N*N, 1,
200  mem, NULL,
201  N*N,1 ,
202  FFTW_BACKWARD, FFTW_MEASURE);
203 
204  /* execute the 2d-infft's */
205  reconstruct(argv[1],N,M,Z,atoi(argv[5]),atoi(argv[6]),mem);
206 
207  /* execute the 1d-fft's */
208  fftw_execute(plan);
209 
210  /* write the memory back in files */
211  print(N,M,Z, mem);
212 
213  /* free memory */
214  nfft_free(mem);
215  fftw_destroy_plan(plan);
216  return 1;
217 }
218 /* \} */

Generated on Tue Apr 30 2013 by Doxygen 1.8.1