NFFT  3.4.1
radon.c
Go to the documentation of this file.
1 /*
2  * Copyright (c) 2002, 2017 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 
39 #include <stdio.h>
40 #include <math.h>
41 #include <stdlib.h>
42 #include <string.h>
43 #include <complex.h>
44 
45 #define NFFT_PRECISION_DOUBLE
46 
47 #include "nfft3mp.h"
48 
50 /*#define KERNEL(r) 1.0 */
51 #define KERNEL(r) (NFFT_K(1.0)-NFFT_M(fabs)((NFFT_R)(r))/((NFFT_R)S/2))
52 
56 static int polar_grid(int T, int S, NFFT_R *x, NFFT_R *w)
57 {
58  int t, r;
59  NFFT_R W = (NFFT_R) T * (((NFFT_R) S / NFFT_K(2.0)) * ((NFFT_R) S / NFFT_K(2.0)) + NFFT_K(1.0) / NFFT_K(4.0));
60 
61  for (t = -T / 2; t < T / 2; t++)
62  {
63  for (r = -S / 2; r < S / 2; r++)
64  {
65  x[2 * ((t + T / 2) * S + (r + S / 2)) + 0] = (NFFT_R) r / (NFFT_R)(S) * NFFT_M(cos)(NFFT_KPI * (NFFT_R)(t) / (NFFT_R)(T));
66  x[2 * ((t + T / 2) * S + (r + S / 2)) + 1] = (NFFT_R) r / (NFFT_R)(S) * NFFT_M(sin)(NFFT_KPI * (NFFT_R)(t) / (NFFT_R)(T));
67  if (r == 0)
68  w[(t + T / 2) * S + (r + S / 2)] = NFFT_K(1.0) / NFFT_K(4.0) / W;
69  else
70  w[(t + T / 2) * S + (r + S / 2)] = NFFT_M(fabs)((NFFT_R) r) / W;
71  }
72  }
73 
74  return 0;
75 }
76 
80 static int linogram_grid(int T, int S, NFFT_R *x, NFFT_R *w)
81 {
82  int t, r;
83  NFFT_R W = (NFFT_R) T * (((NFFT_R) S / NFFT_K(2.0)) * ((NFFT_R) S / NFFT_K(2.0)) + NFFT_K(1.0) / NFFT_K(4.0));
84 
85  for (t = -T / 2; t < T / 2; t++)
86  {
87  for (r = -S / 2; r < S / 2; r++)
88  {
89  if (t < 0)
90  {
91  x[2 * ((t + T / 2) * S + (r + S / 2)) + 0] = (NFFT_R) r / (NFFT_R)(S);
92  x[2 * ((t + T / 2) * S + (r + S / 2)) + 1] = NFFT_K(4.0) * ((NFFT_R)(t) + (NFFT_R)(T) / NFFT_K(4.0)) / (NFFT_R)(T) * (NFFT_R)(r)
93  / (NFFT_R)(S);
94  }
95  else
96  {
97  x[2 * ((t + T / 2) * S + (r + S / 2)) + 0] = -NFFT_K(4.0) * ((NFFT_R)(t) - (NFFT_R)(T) / NFFT_K(4.0)) / (NFFT_R)(T)
98  * (NFFT_R)(r) / (NFFT_R)(S);
99  x[2 * ((t + T / 2) * S + (r + S / 2)) + 1] = (NFFT_R) r / (NFFT_R)(S);
100  }
101  if (r == 0)
102  w[(t + T / 2) * S + (r + S / 2)] = NFFT_K(1.0) / NFFT_K(4.0) / W;
103  else
104  w[(t + T / 2) * S + (r + S / 2)] = NFFT_M(fabs)((NFFT_R) r) / W;
105  }
106  }
107 
108  return 0;
109 }
110 
114 static int Radon_trafo(int (*gridfcn)(), int T, int S, NFFT_R *f, int NN, NFFT_R *Rf)
115 {
116  int j, k;
117  NFFT(plan) my_nfft_plan;
119  NFFT_C *fft;
120  FFTW(plan) my_fftw_plan;
122  int t, r;
123  NFFT_R *x, *w;
125  int N[2], n[2];
126  int M = T * S;
127 
128  N[0] = NN;
129  n[0] = 2 * N[0];
130  N[1] = NN;
131  n[1] = 2 * N[1];
132 
133  fft = (NFFT_C *) NFFT(malloc)((size_t)(S) * sizeof(NFFT_C));
134  my_fftw_plan = FFTW(plan_dft_1d)(S, fft, fft, FFTW_BACKWARD, FFTW_MEASURE);
135 
136  x = (NFFT_R *) NFFT(malloc)((size_t)(2 * T * S) * (sizeof(NFFT_R)));
137  if (x == NULL)
138  return EXIT_FAILURE;
139 
140  w = (NFFT_R *) NFFT(malloc)((size_t)(T * S) * (sizeof(NFFT_R)));
141  if (w == NULL)
142  return EXIT_FAILURE;
143 
145  NFFT(init_guru)(&my_nfft_plan, 2, N, M, n, 4,
148  FFTW_MEASURE | FFTW_DESTROY_INPUT);
149 
151  gridfcn(T, S, 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.flags & PRE_LIN_PSI)
160  NFFT(precompute_lin_psi)(&my_nfft_plan);
161 
162  if (my_nfft_plan.flags & PRE_PSI)
163  NFFT(precompute_psi)(&my_nfft_plan);
164 
165  if (my_nfft_plan.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 * NFFT_K(0.0);
171 
173  NFFT(trafo)(&my_nfft_plan);
174 
176  for (t = 0; t < T; t++)
177  {
178  fft[0] = NFFT_K(0.0);
179  for (r = -S / 2 + 1; r < S / 2; r++)
180  fft[r + S / 2] = KERNEL(r) * my_nfft_plan.f[t * S + (r + S / 2)];
181 
182  NFFT(fftshift_complex_int)(fft, 1, &S);
183  FFTW(execute)(my_fftw_plan);
184  NFFT(fftshift_complex_int)(fft, 1, &S);
185 
186  for (r = 0; r < S; r++)
187  Rf[t * S + r] = NFFT_M(creal)(fft[r]) / (NFFT_R)(S);
188 
189  /* for(r=0; r<R/2; r++)
190  Rf[t*R+(r+R/2)] = creal(cexp(-I*NFFT_KPI*r)*fft[r]);
191  for(r=0; r<R/2; r++)
192  Rf[t*R+r] = creal(cexp(-I*NFFT_KPI*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, S;
211  FILE *fp;
212  int N;
213  NFFT_R *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(EXIT_FAILURE);
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  S = atoi(argv[4]);
234  /*printf("N=%d, %s grid with T=%d, R=%d. \n",N,argv[1],T,R);*/
235 
236  f = (NFFT_R *) NFFT(malloc)((size_t)(N * N) * (sizeof(NFFT_R)));
237  Rf = (NFFT_R *) NFFT(malloc)((size_t)(T * S) * (sizeof(NFFT_R)));
238 
240  fp = fopen("input_data.bin", "rb");
241  if (fp == NULL)
242  return EXIT_FAILURE;
243  fread(f, sizeof(NFFT_R), (size_t)(N * N), fp);
244  fclose(fp);
245 
247  Radon_trafo(gridfcn, T, S, f, N, Rf);
248 
250  fp = fopen("sinogram_data.bin", "wb+");
251  if (fp == NULL)
252  return EXIT_FAILURE;
253  fwrite(Rf, sizeof(NFFT_R), (size_t)(T * S), fp);
254  fclose(fp);
255 
257  NFFT(free)(f);
258  NFFT(free)(Rf);
259 
260  return EXIT_SUCCESS;
261 }
static int Radon_trafo(int(*gridfcn)(), int T, int S, NFFT_R *f, int NN, NFFT_R *Rf)
computes the NFFT-based discrete Radon transform of f on the grid given by gridfcn() with T angles an...
Definition: radon.c:114
#define MALLOC_X
Definition: nfft3.h:199
#define MALLOC_F_HAT
Definition: nfft3.h:200
static int linogram_grid(int T, int S, NFFT_R *x, NFFT_R *w)
generates the points x with weights w for the linogram grid with T slopes and R offsets ...
Definition: radon.c:80
int main(int argc, char **argv)
simple test program for the discrete Radon transform
Definition: radon.c:207
static void fft(int N, int M, int Z, fftw_complex *mem)
fft makes an 1D-ftt for every knot through all layers
#define FFTW_INIT
Definition: nfft3.h:203
#define MALLOC_F
Definition: nfft3.h:201
#define FFT_OUT_OF_PLACE
Definition: nfft3.h:202
#define PRE_LIN_PSI
Definition: nfft3.h:195
#define PRE_PSI
Definition: nfft3.h:197
#define PRE_FULL_PSI
Definition: nfft3.h:198
#define PRE_PHI_HUT
Definition: nfft3.h:193
static int polar_grid(int T, int S, NFFT_R *x, NFFT_R *w)
generates the points x with weights w for the polar grid with T angles and R offsets ...
Definition: radon.c:56
#define KERNEL(r)
define weights of kernel function for discrete Radon transform
Definition: radon.c:51