62 NFFT(init_guru)((NFFT(plan)*) ths, 1, &N, M, &n, m,
64 FFTW_ESTIMATE | FFTW_PRESERVE_INPUT);
66 ths->idx0 = (INT*) NFFT(malloc)((size_t)(M) *
sizeof(INT));
67 ths->deltax0 = (R*) NFFT(malloc)((size_t)(M) *
sizeof(R));
81 NFFT(plan)* cths = (NFFT(plan)*) ths;
83 for (j = 0; j < cths->M_total; j++)
85 ths->idx0[j] = (LRINT(ROUND((cths->x[j] + K(0.5)) * (R)(cths->n[0])))
86 + cths->n[0] / 2) % cths->n[0];
87 ths->deltax0[j] = cths->x[j]
88 - (ROUND((cths->x[j] + K(0.5)) * (R)(cths->n[0])) / (R)(cths->n[0]) - K(0.5));
101 NFFT(free)(ths->deltax0);
102 NFFT(free)(ths->idx0);
104 NFFT(finalize)((NFFT(plan)*) ths);
124 NFFT(plan) *cths = (NFFT(plan)*) ths;
126 for (j = 0, f = cths->f; j < cths->M_total; j++)
129 for (k = 0; k < cths->n_total; k++)
130 cths->g1[k] = K(0.0);
132 for (k = -cths->N_total / 2, g1 = cths->g1 + cths->n_total
133 - cths->N_total / 2, f_hat = cths->f_hat; k < 0; k++)
134 (*g1++) = CPOW(-K2PI * II * (R)(k), (R)(cths->m)) * (*f_hat++);
136 cths->g1[0] = cths->f_hat[cths->N_total / 2];
138 for (k = 1, g1 = cths->g1 + 1, f_hat = cths->f_hat + cths->N_total / 2 + 1;
139 k < cths->N_total / 2; k++)
140 (*g1++) = CPOW(-K2PI * II * (R)(k), (R)(cths->m)) * (*f_hat++);
142 for (l = cths->m - 1; l >= 0; l--)
144 for (k = -cths->N_total / 2, g1 = cths->g1 + cths->n_total
145 - cths->N_total / 2; k < 0; k++)
146 (*g1++) /= (-K2PI * II * (R)(k));
148 for (k = 1, g1 = cths->g1 + 1; k < cths->N_total / 2; k++)
149 (*g1++) /= (-K2PI * II * (R)(k));
151 FFTW(execute)(cths->my_fftw_plan1);
153 ll = (l == 0 ? 1 : l);
155 for (j = 0, f = cths->f, deltax = ths->deltax0, idx = ths->idx0;
156 j < cths->M_total; j++, f++)
157 (*f) = ((*f) * (*deltax++) + cths->g2[*idx++]) / (R)(ll);
175 int m_taylor,
unsigned test_accuracy)
178 R t_ndft, t_nfft, t_taylor, t;
185 printf(
"%d\t%d\t", N, M);
189 NFFT(init_guru)(&np, 1, &N, M, &n, m,
191 FFTW_ESTIMATE | FFTW_DESTROY_INPUT);
195 np.f_hat = tp.p.f_hat;
200 swapndft = (C*) NFFT(malloc)((size_t)(M) *
sizeof(C));
203 NFFT(vrand_shifted_unit_double)(np.x, np.M_total);
210 NFFT(precompute_one_psi)(&np);
213 NFFT(vrand_unit_complex)(np.f_hat, np.N_total);
218 CSWAP(np.f, swapndft);
222 while (t_ndft < K(0.01))
226 NFFT(trafo_direct)(&np);
228 t = NFFT(elapsed_seconds)(t1, t0);
233 CSWAP(np.f, swapndft);
234 printf(
"%.2" __FES__
"\t", t_ndft);
242 while (t_nfft < K(0.01))
248 t = NFFT(elapsed_seconds)(t1, t0);
253 printf(
"%.2" __FES__
"\t%d\t%.2" __FES__
"\t", ((R)(n)) / ((R)(N)), m, t_nfft);
256 printf(
"%.2" __FES__
"\t", NFFT(error_l_infty_complex)(swapndft, np.f, np.M_total));
263 while (t_taylor < K(0.01))
269 t = NFFT(elapsed_seconds)(t1, t0);
274 printf(
"%.2" __FES__
"\t%d\t%.2" __FES__
"\t", ((R)(n_taylor)) / ((R)(N)), m_taylor, t_taylor);
277 printf(
"%.2" __FES__
"\n", NFFT(error_l_infty_complex)(swapndft, np.f, np.M_total));
285 NFFT(free)(swapndft);
291 int main(
int argc,
char **argv)
298 "taylor_nfft type first last trials sigma_nfft sigma_taylor.\n");
302 fprintf(stderr,
"Testing the Nfft & a Taylor expansion based version.\n\n");
303 fprintf(stderr,
"Columns: N, M, t_ndft, sigma_nfft, m_nfft, t_nfft, e_nfft");
304 fprintf(stderr,
", sigma_taylor, m_taylor, t_taylor, e_taylor\n");
307 if (atoi(argv[1]) == 0)
309 fprintf(stderr,
"Fixed target accuracy, timings.\n\n");
310 int arg2 = atoi(argv[2]);
311 int arg3 = atoi(argv[3]);
312 int arg4 = atoi(argv[4]);
313 for (l = arg2; l <= arg3; l++)
315 int N = (int)(1U << l);
316 int M = (int)(1U << l);
317 int arg5 = (int)(atof(argv[5]) * N);
318 int arg6 = (int)(atof(argv[6]) * N);
319 for (trial = 0; trial < arg4; trial++)
327 if (atoi(argv[1]) == 1)
329 int arg2 = atoi(argv[2]);
330 int arg3 = atoi(argv[3]);
331 int arg4 = atoi(argv[4]);
332 int N = atoi(argv[7]);
333 int arg5 = (int) (atof(argv[5]) * N);
334 int arg6 = (int) (atof(argv[6]) * N);
335 fprintf(stderr,
"Fixed N=M=%d, error vs. m.\n\n", N);
336 for (m = arg2; m <= arg3; m++)
338 for (trial = 0; trial < arg4; trial++)
static void taylor_finalize(taylor_plan *ths)
Destroys a transform plan.
static void taylor_precompute(taylor_plan *ths)
Precomputation of weights and indices in Taylor expansion.
static void taylor_trafo(taylor_plan *ths)
Executes a Taylor-NFFT, see equation (1.1) in [Guide], computes fast and approximate by means of a Ta...
static void taylor_init(taylor_plan *ths, int N, int M, int n, int m)
Initialisation of a transform plan.
static void taylor_time_accuracy(int N, int M, int n, int m, int n_taylor, int m_taylor, unsigned test_accuracy)
Compares NDFT, NFFT, and Taylor-NFFT.
Header file for the nfft3 library.
#define CSWAP(x, y)
Swap two vectors.