NFFT Logo 3.2.3
reconstruct_data_3d.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_3d.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 
40 static void reconstruct(char* filename,int N,int M,int Z,int iteration, int weight)
41 {
42  int j,k,z,l; /* 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  FILE* fout_real; /* output file (real part) */
48  FILE* fout_imag; /* output file (imag part) */
49  int my_N[3],my_n[3]; /* to init the nfft */
50  double epsilon=0.0000003; /* tmp to read the obsolent z from 700.acs
51  epsilon is a the break criterion for
52  the iteration */
53  unsigned infft_flags = CGNR | PRECOMPUTE_DAMP; /* flags for the infft */
54 
55  /* initialise my_plan, specific.
56  we don't precompute psi */
57  my_N[0]=Z; my_n[0]=ceil(Z*1.2);
58  my_N[1]=N; my_n[1]=ceil(N*1.2);
59  my_N[2]=N; my_n[2]=ceil(N*1.2);
60  nfft_init_guru(&my_plan, 3, my_N, M, my_n, 6,
61  PRE_PHI_HUT| PRE_PSI |MALLOC_X| MALLOC_F_HAT|
62  MALLOC_F| FFTW_INIT| FFT_OUT_OF_PLACE,
63  FFTW_MEASURE| FFTW_DESTROY_INPUT);
64 
65  /* precompute lin psi */
66  if(my_plan.nfft_flags & PRE_LIN_PSI)
67  nfft_precompute_lin_psi(&my_plan);
68 
69  if (weight)
70  infft_flags = infft_flags | PRECOMPUTE_WEIGHT;
71 
72  /* initialise my_iplan, advanced */
73  solver_init_advanced_complex(&my_iplan,(nfft_mv_plan_complex*)(&my_plan), infft_flags );
74 
75  /* get the weights */
76  if(my_iplan.flags & PRECOMPUTE_WEIGHT)
77  {
78  fin=fopen("weights.dat","r");
79  for(j=0;j<M;j++)
80  {
81  fscanf(fin,"%le ",&my_iplan.w[j]);
82  }
83  fclose(fin);
84  }
85 
86  /* get the damping factors */
87  if(my_iplan.flags & PRECOMPUTE_DAMP)
88  {
89  for(j=0;j<N;j++){
90  for(k=0;k<N;k++) {
91  for(z=0;z<N;z++) {
92  int j2= j-N/2;
93  int k2= k-N/2;
94  int z2= z-N/2;
95  double r=sqrt(j2*j2+k2*k2+z2*z2);
96  if(r>(double) N/2)
97  my_iplan.w_hat[z*N*N+j*N+k]=0.0;
98  else
99  my_iplan.w_hat[z*N*N+j*N+k]=1.0;
100  }
101  }
102  }
103  }
104 
105  /* open the input file */
106  fin=fopen(filename,"r");
107 
108  /* open the output files */
109  fout_real=fopen("output_real.dat","w");
110  fout_imag=fopen("output_imag.dat","w");
111 
112  /* read x,y,freal and fimag from the knots */
113  for(j=0;j<M;j++)
114  {
115  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],
116  &real,&imag);
117  my_iplan.y[j] = real + _Complex_I*imag;
118  }
119 
120  /* precompute psi */
121  if(my_plan.nfft_flags & PRE_PSI)
122  nfft_precompute_psi(&my_plan);
123 
124  /* precompute full psi */
125  if(my_plan.nfft_flags & PRE_FULL_PSI)
126  nfft_precompute_full_psi(&my_plan);
127 
128  /* init some guess */
129  for(k=0;k<my_plan.N_total;k++)
130  my_iplan.f_hat_iter[k]=0.0;
131 
132  /* inverse trafo */
133  solver_before_loop_complex(&my_iplan);
134  for(l=0;l<iteration;l++)
135  {
136  /* break if dot_r_iter is smaller than epsilon*/
137  if(my_iplan.dot_r_iter<epsilon)
138  break;
139  fprintf(stderr,"%e, %i of %i\n",sqrt(my_iplan.dot_r_iter),
140  l+1,iteration);
141  solver_loop_one_step_complex(&my_iplan);
142  }
143 
144  for(l=0;l<Z;l++)
145  {
146  for(k=0;k<N*N;k++)
147  {
148  /* write every Layer in the files */
149  fprintf(fout_real,"%le ",creal(my_iplan.f_hat_iter[ k+N*N*l ]));
150  fprintf(fout_imag,"%le ",cimag(my_iplan.f_hat_iter[ k+N*N*l ]));
151  }
152  fprintf(fout_real,"\n");
153  fprintf(fout_imag,"\n");
154  }
155 
156  fclose(fout_real);
157  fclose(fout_imag);
158 
159  solver_finalize_complex(&my_iplan);
160  nfft_finalize(&my_plan);
161 }
162 
163 int main(int argc, char **argv)
164 {
165  if (argc <= 6) {
166  printf("usage: ./reconstruct3D FILENAME N M Z ITER WEIGHTS\n");
167  return 1;
168  }
169 
170  reconstruct(argv[1],atoi(argv[2]),atoi(argv[3]),atoi(argv[4]),atoi(argv[5]),atoi(argv[6]));
171  return 1;
172 }
173 /* \} */

Generated on Tue Apr 30 2013 by Doxygen 1.8.1