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