78 static int polar_grid(
int T,
int R,
double *x,
double *w)
81 double W=(double)T*(((
double)R/2.0)*((double)R/2.0)+1.0/4.0);
83 for(t=-T/2; t<T/2; t++)
85 for(r=-R/2; r<R/2; r++)
87 x[2*((t+T/2)*R+(r+R/2))+0] = (
double)r/R*cos(
PI*t/T);
88 x[2*((t+T/2)*R+(r+R/2))+1] = (
double)r/R*sin(
PI*t/T);
90 w[(t+T/2)*R+(r+R/2)] = 1.0/4.0/W;
92 w[(t+T/2)*R+(r+R/2)] = fabs((
double)r)/W;
100 static int polar_dft(fftw_complex *f_hat,
int NN, fftw_complex *f,
int T,
int R,
int m)
110 N[0]=NN; n[0]=2*N[0];
111 N[1]=NN; n[1]=2*N[1];
122 nfft_init_guru(&my_nfft_plan, 2, N, M, n, m,
123 PRE_PHI_HUT| PRE_PSI| MALLOC_X | MALLOC_F_HAT| MALLOC_F| FFTW_INIT | FFT_OUT_OF_PLACE,
124 FFTW_MEASURE| FFTW_DESTROY_INPUT);
128 for(j=0;j<my_nfft_plan.
M_total;j++)
130 my_nfft_plan.
x[2*j+0] = x[2*j+0];
131 my_nfft_plan.
x[2*j+1] = x[2*j+1];
135 for(k=0;k<my_nfft_plan.
N_total;k++)
136 my_nfft_plan.
f_hat[k] = f_hat[k];
139 nfft_trafo_direct(&my_nfft_plan);
142 for(j=0;j<my_nfft_plan.
M_total;j++)
143 f[j] = my_nfft_plan.
f[j];
146 nfft_finalize(&my_nfft_plan);
154 static int polar_fft(fftw_complex *f_hat,
int NN, fftw_complex *f,
int T,
int R,
int m)
164 N[0]=NN; n[0]=2*N[0];
165 N[1]=NN; n[1]=2*N[1];
176 nfft_init_guru(&my_nfft_plan, 2, N, M, n, m,
177 PRE_PHI_HUT| PRE_PSI| MALLOC_X | MALLOC_F_HAT| MALLOC_F| FFTW_INIT | FFT_OUT_OF_PLACE,
178 FFTW_MEASURE| FFTW_DESTROY_INPUT);
182 for(j=0;j<my_nfft_plan.
M_total;j++)
184 my_nfft_plan.
x[2*j+0] = x[2*j+0];
185 my_nfft_plan.
x[2*j+1] = x[2*j+1];
193 nfft_precompute_psi(&my_nfft_plan);
196 nfft_precompute_full_psi(&my_nfft_plan);
199 for(k=0;k<my_nfft_plan.
N_total;k++)
200 my_nfft_plan.
f_hat[k] = f_hat[k];
206 for(j=0;j<my_nfft_plan.
M_total;j++)
207 f[j] = my_nfft_plan.
f[j];
210 nfft_finalize(&my_nfft_plan);
230 N[0]=NN; n[0]=2*N[0];
231 N[1]=NN; n[1]=2*N[1];
242 nfft_init_guru(&my_nfft_plan, 2, N, M, n, m,
243 PRE_PHI_HUT| PRE_PSI| MALLOC_X | MALLOC_F_HAT| MALLOC_F| FFTW_INIT | FFT_OUT_OF_PLACE,
244 FFTW_MEASURE| FFTW_DESTROY_INPUT);
247 solver_init_advanced_complex(&my_infft_plan,(
nfft_mv_plan_complex*)(&my_nfft_plan), CGNR | PRECOMPUTE_WEIGHT );
251 for(j=0;j<my_nfft_plan.
M_total;j++)
253 my_nfft_plan.
x[2*j+0] = x[2*j+0];
254 my_nfft_plan.
x[2*j+1] = x[2*j+1];
255 my_infft_plan.
y[j] = f[j];
256 my_infft_plan.
w[j] = w[j];
264 nfft_precompute_psi(&my_nfft_plan);
267 nfft_precompute_full_psi(&my_nfft_plan);
270 if(my_infft_plan.
flags & PRECOMPUTE_DAMP)
271 for(j=0;j<my_nfft_plan.
N[0];j++)
272 for(k=0;k<my_nfft_plan.
N[1];k++)
274 my_infft_plan.
w_hat[j*my_nfft_plan.
N[1]+k]=
275 (sqrt(pow((
double)(j-my_nfft_plan.
N[0]/2),2.0)+pow((
double)(k-my_nfft_plan.
N[1]/2),2.0))
276 >((double)(my_nfft_plan.
N[0]/2))?0:1);
280 for(k=0;k<my_nfft_plan.
N_total;k++)
281 my_infft_plan.
f_hat_iter[k] = 0.0 + _Complex_I*0.0;
284 solver_before_loop_complex(&my_infft_plan);
289 for(k=0;k<my_nfft_plan.
N_total;k++)
294 for(l=1;l<=
max_i;l++)
301 for(k=0;k<my_nfft_plan.
N_total;k++)
306 nfft_finalize(&my_nfft_plan);
314 int main(
int argc,
char **argv)
320 fftw_complex *f_hat, *f, *f_direct, *f_tilde;
324 double temp1, temp2, E_max=0.0;
330 printf(
"polar_fft_test N T R \n");
332 printf(
"N polar FFT of size NxN \n");
333 printf(
"T number of slopes \n");
334 printf(
"R number of offsets \n");
341 printf(
"N=%d, polar grid with T=%d, R=%d => ",N,T,R);
343 x = (
double *)
nfft_malloc(2*5*(T/2)*(R/2)*(
sizeof(
double)));
344 w = (
double *)
nfft_malloc(5*(T/2)*(R/2)*(
sizeof(
double)));
346 f_hat = (fftw_complex *)
nfft_malloc(
sizeof(fftw_complex)*N*N);
347 f = (fftw_complex *)
nfft_malloc(
sizeof(fftw_complex)*T*R);
348 f_direct = (fftw_complex *)
nfft_malloc(
sizeof(fftw_complex)*T*R);
349 f_tilde = (fftw_complex *)
nfft_malloc(
sizeof(fftw_complex)*N*N);
355 fp1=fopen(
"input_data_r.dat",
"r");
356 fp2=fopen(
"input_data_i.dat",
"r");
361 fscanf(fp1,
"%le ",&temp1);
362 fscanf(fp2,
"%le ",&temp2);
363 f_hat[k]=temp1+ _Complex_I*temp2;
373 printf(
"\nTest of the polar FFT: \n");
374 fp1=fopen(
"polar_fft_error.dat",
"w+");
375 for (m=1; m<=12; m++)
381 E_max=X(error_l_infty_complex)(f_direct,f,M);
382 printf(
"m=%2d: E_max = %e\n",m,E_max);
383 fprintf(fp1,
"%e\n",E_max);
388 for (m=3; m<=9; m+=3)
390 printf(
"\nTest of the inverse polar FFT for m=%d: \n",m);
391 sprintf(filename,
"polar_ifft_error%d.dat",m);
392 fp1=fopen(filename,
"w+");
393 for (max_i=0; max_i<=100; max_i+=10)
406 E_max=X(error_l_infty_complex)(f_hat,f_tilde,N*N);
407 printf(
"%3d iterations: E_max = %e\n",max_i,E_max);
408 fprintf(fp1,
"%e\n",E_max);
fftw_complex * f_hat
Vector of Fourier coefficients, size is N_total * sizeof( fftw_complex )
int max_i(int a, int b)
max
fftw_complex * p_hat_iter
search direction
double * w
weighting factors
unsigned flags
iteration type
void nfft_precompute_lin_psi(nfft_plan *ths)
create a lookup table, but NOT for each node good idea K=2^xx TODO: estimate K, call from init assume...
fftw_complex * f
Vector of samples, size is M_total * sizeof( fftw_complex )
void nfft_trafo(nfft_plan *ths)
user routines
int M_total
Total number of samples.
Header file for utility functions used by the nfft3 library.
#define PI
Formerly known to be an irrational number.
unsigned nfft_flags
Flags for precomputation, (de)allocation, and FFTW usage, default setting is PRE_PHI_HUT | PRE_PSI | ...
static int polar_dft(fftw_complex *f_hat, int NN, fftw_complex *f, int T, int R, int m)
discrete polar FFT
void * nfft_malloc(size_t n)
static int polar_grid(int T, int R, double *x, double *w)
Generates the points with weights for the polar grid with angles and offsets. ...
fftw_complex * y
right hand side, samples
double * x
Nodes in time/spatial domain, size is doubles.
void solver_loop_one_step_complex(solver_plan_complex *ths)
void solver_loop_one_step
void solver_finalize_complex(solver_plan_complex *ths)
void solver_finalize
int N_total
Total number of Fourier coefficients.
static int polar_fft(fftw_complex *f_hat, int NN, fftw_complex *f, int T, int R, int m)
NFFT-based polar FFT.
static int inverse_polar_fft(fftw_complex *f, int T, int R, fftw_complex *f_hat, int NN, int max_i, int m)
inverse NFFT-based polar FFT
double * w_hat
damping factors
fftw_complex * f_hat_iter
iterative solution