NFFT Logo 3.2.3
radon.c
Go to the documentation of this file.
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: radon.c 3896 2012-10-10 12:19:26Z tovo $ */
20 
40 #include "config.h"
41 
42 #include <stdio.h>
43 #include <math.h>
44 #include <stdlib.h>
45 #include <string.h>
46 #ifdef HAVE_COMPLEX_H
47 #include <complex.h>
48 #endif
49 
50 #include "nfft3util.h"
51 #include "nfft3.h"
52 
54 /*#define KERNEL(r) 1.0 */
55 #define KERNEL(r) (1.0-fabs((double)(r))/((double)R/2))
56 
60 static int polar_grid(int T, int R, double *x, double *w)
61 {
62  int t, r;
63  double W=(double)T*(((double)R/2.0)*((double)R/2.0)+1.0/4.0);
64 
65  for(t=-T/2; t<T/2; t++)
66  {
67  for(r=-R/2; r<R/2; r++)
68  {
69  x[2*((t+T/2)*R+(r+R/2))+0] = (double)r/R*cos(PI*t/T);
70  x[2*((t+T/2)*R+(r+R/2))+1] = (double)r/R*sin(PI*t/T);
71  if (r==0)
72  w[(t+T/2)*R+(r+R/2)] = 1.0/4.0/W;
73  else
74  w[(t+T/2)*R+(r+R/2)] = fabs((double)r)/W;
75  }
76  }
77 
78  return 0;
79 }
80 
84 static int linogram_grid(int T, int R, double *x, double *w)
85 {
86  int t, r;
87  double W=(double)T*(((double)R/2.0)*((double)R/2.0)+1.0/4.0);
88 
89  for(t=-T/2; t<T/2; t++)
90  {
91  for(r=-R/2; r<R/2; r++)
92  {
93  if(t<0)
94  {
95  x[2*((t+T/2)*R+(r+R/2))+0] = (double)r/R;
96  x[2*((t+T/2)*R+(r+R/2))+1] = (double)4*(t+T/4)/T*r/R;
97  }
98  else
99  {
100  x[2*((t+T/2)*R+(r+R/2))+0] = -(double)4*(t-T/4)/T*r/R;
101  x[2*((t+T/2)*R+(r+R/2))+1] = (double)r/R;
102  }
103  if (r==0)
104  w[(t+T/2)*R+(r+R/2)] = 1.0/4.0/W;
105  else
106  w[(t+T/2)*R+(r+R/2)] = fabs((double)r)/W;
107  }
108  }
109 
110  return 0;
111 }
112 
116 int Radon_trafo(int (*gridfcn)(), int T, int R, double *f, int NN, double *Rf)
117 {
118  int j,k;
119  nfft_plan my_nfft_plan;
121  fftw_complex *fft;
122  fftw_plan my_fftw_plan;
124  int t,r;
125  double *x, *w;
127  int N[2],n[2];
128  int M=T*R;
129 
130  N[0]=NN; n[0]=2*N[0];
131  N[1]=NN; n[1]=2*N[1];
132 
133  fft = (fftw_complex *)nfft_malloc(R*sizeof(fftw_complex));
134  my_fftw_plan = fftw_plan_dft_1d(R,fft,fft,FFTW_BACKWARD,FFTW_MEASURE);
135 
136  x = (double *)nfft_malloc(2*T*R*(sizeof(double)));
137  if (x==NULL)
138  return -1;
139 
140  w = (double *)nfft_malloc(T*R*(sizeof(double)));
141  if (w==NULL)
142  return -1;
143 
145  nfft_init_guru(&my_nfft_plan, 2, N, M, n, 4,
146  PRE_PHI_HUT| PRE_PSI| MALLOC_X | MALLOC_F_HAT| MALLOC_F| FFTW_INIT | FFT_OUT_OF_PLACE,
147  FFTW_MEASURE| FFTW_DESTROY_INPUT);
148 
149 
151  gridfcn(T,R,x,w);
152  for(j=0;j<my_nfft_plan.M_total;j++)
153  {
154  my_nfft_plan.x[2*j+0] = x[2*j+0];
155  my_nfft_plan.x[2*j+1] = x[2*j+1];
156  }
157 
159  if(my_nfft_plan.nfft_flags & PRE_LIN_PSI)
160  nfft_precompute_lin_psi(&my_nfft_plan);
161 
162  if(my_nfft_plan.nfft_flags & PRE_PSI)
163  nfft_precompute_psi(&my_nfft_plan);
164 
165  if(my_nfft_plan.nfft_flags & PRE_FULL_PSI)
166  nfft_precompute_full_psi(&my_nfft_plan);
167 
169  for(k=0;k<my_nfft_plan.N_total;k++)
170  my_nfft_plan.f_hat[k] = f[k] + _Complex_I*0.0;
171 
173  nfft_trafo(&my_nfft_plan);
174 
176  for(t=0; t<T; t++)
177  {
178  fft[0]=0.0;
179  for(r=-R/2+1; r<R/2; r++)
180  fft[r+R/2] = KERNEL(r)*my_nfft_plan.f[t*R+(r+R/2)];
181 
182  nfft_fftshift_complex(fft, 1, &R);
183  fftw_execute(my_fftw_plan);
184  nfft_fftshift_complex(fft, 1, &R);
185 
186  for(r=0; r<R; r++)
187  Rf[t*R+r] = creal(fft[r])/R;
188 
189 /* for(r=0; r<R/2; r++)
190  Rf[t*R+(r+R/2)] = creal(cexp(-I*PI*r)*fft[r]);
191  for(r=0; r<R/2; r++)
192  Rf[t*R+r] = creal(cexp(-I*PI*r)*fft[r+R/2]);
193  */
194  }
195 
197  fftw_destroy_plan(my_fftw_plan);
198  nfft_free(fft);
199  nfft_finalize(&my_nfft_plan);
200  nfft_free(x);
201  nfft_free(w);
202  return 0;
203 }
204 
207 int main(int argc,char **argv)
208 {
209  int (*gridfcn)();
210  int T, R;
211  FILE *fp;
212  int N;
213  double *f, *Rf;
214 
215  if( argc!=5 )
216  {
217  printf("radon gridfcn N T R\n");
218  printf("\n");
219  printf("gridfcn \"polar\" or \"linogram\" \n");
220  printf("N image size NxN \n");
221  printf("T number of slopes \n");
222  printf("R number of offsets \n");
223  exit(-1);
224  }
225 
226  if (strcmp(argv[1],"polar") == 0)
227  gridfcn = polar_grid;
228  else
229  gridfcn = linogram_grid;
230 
231  N = atoi(argv[2]);
232  T = atoi(argv[3]);
233  R = atoi(argv[4]);
234  /*printf("N=%d, %s grid with T=%d, R=%d. \n",N,argv[1],T,R);*/
235 
236  f = (double *)nfft_malloc(N*N*(sizeof(double)));
237  Rf = (double *)nfft_malloc(T*R*(sizeof(double)));
238 
240  fp=fopen("input_data.bin","rb");
241  if (fp==NULL)
242  return(-1);
243  fread(f,sizeof(double),N*N,fp);
244  fclose(fp);
245 
247  Radon_trafo(gridfcn,T,R,f,N,Rf);
248 
250  fp=fopen("sinogram_data.bin","wb+");
251  if (fp==NULL)
252  return(-1);
253  fwrite(Rf,sizeof(double),T*R,fp);
254  fclose(fp);
255 
257  nfft_free(f);
258  nfft_free(Rf);
259 
260  return EXIT_SUCCESS;
261 }

Generated on Tue Apr 30 2013 by Doxygen 1.8.1