NFFT Logo 3.2.3
fpt/simple_test.c
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: simple_test.c 3896 2012-10-10 12:19:26Z tovo $ */
20 
21 #include "config.h"
22 
23 /* standard headers */
24 #include <stdlib.h>
25 #include <stdio.h>
26 #include <math.h>
27 
28 /* It is important to include complex.h before nfft3.h. */
29 #ifdef HAVE_COMPLEX_H
30 #include <complex.h>
31 #endif
32 
33 #include <fftw3.h>
34 
35 /* NFFT3 header */
36 #include "nfft3.h"
37 #include "nfft3util.h"
38 
39 /* Two times Pi */
40 #define KPI2 6.2831853071795864769252867665590057683943387987502
41 
42 int main(void)
43 {
44  /* This example shows the use of the fast polynomial transform to evaluate a
45  * finite expansion in Legendre polynomials,
46  *
47  * f(x) = a_0 P_0(x) + a_1 P_1(x) + ... + a_N P_N(x) (1)
48  *
49  * at the Chebyshev nodes x_j = cos(j*pi/N), j=0,1,...,N. */
50  const int N = 8;
51 
52  /* An fpt_set is a data structure that contains precomputed data for a number
53  * of different polynomial transforms. Here, we need only one transform. the
54  * second parameter (t) is the exponent of the maximum transform size desired
55  * (2^t), i.e., t = 3 means that N in (1) can be at most N = 8. */
56  fpt_set set = fpt_init(1,lrint(ceil(log2((double)N))),0U);
57 
58  /* Three-term recurrence coefficients for Legendre polynomials */
59  double *alpha = malloc((N+2)*sizeof(double)),
60  *beta = malloc((N+2)*sizeof(double)),
61  *gamma = malloc((N+2)*sizeof(double));
62 
63  /* alpha[0] and beta[0] are not referenced. */
64  alpha[0] = beta[0] = 0.0;
65  /* gamma[0] contains the value of P_0(x) (which is a constant). */
66  gamma[0] = 1.0;
67 
68  /* Actual three-term recurrence coefficients for Legendre polynomials */
69  {
70  int k;
71  for (k = 0; k <= N; k++)
72  {
73  alpha[k+1] = ((double)(2*k+1))/((double)(k+1));
74  beta[k+1] = 0.0;
75  gamma[k+1] = -((double)(k))/((double)(k+1));
76  }
77  }
78 
79  printf(
80  "Computing a fast polynomial transform (FPT) and a fast discrete cosine \n"
81  "transform (DCT) to evaluate\n\n"
82  " f_j = a_0 P_0(x_j) + a_1 P_1(x_j) + ... + a_N P_N(x_j), j=0,1,...,N,\n\n"
83  "with N=%d, x_j = cos(j*pi/N), j=0,1,...N, the Chebyshev nodes, a_k,\n"
84  "k=0,1,...,N, random Fourier coefficients in [-1,1]x[-1,1]*I, and P_k,\n"
85  "k=0,1,...,N, the Legendre polynomials.",N
86  );
87 
88  /* Random seed, makes things reproducible. */
89  nfft_srand48(314);
90 
91  /* The function fpt_repcompute actually does the precomputation for a single
92  * transform. It needs arrays alpha, beta, and gamma, containing the three-
93  * term recurrence coefficients, here of the Legendre polynomials. The format
94  * is explained above. The sixth parameter (k_start) is where the index in the
95  * linear combination (1) starts, here k_start=0. The seventh parameter
96  * (kappa) is the threshold which has an influence on the accuracy of the fast
97  * polynomial transform. Usually, kappa = 1000 is a good choice. */
98  fpt_precompute(set,0,alpha,beta,gamma,0,1000.0);
99 
100 
101  {
102  /* Arrays for Fourier coefficients and function values. */
103  double _Complex *a = malloc((N+1)*sizeof(double _Complex));
104  double _Complex *b = malloc((N+1)*sizeof(double _Complex));
105  double *f = malloc((N+1)*sizeof(double _Complex));
106 
107  /* Plan for discrete cosine transform */
108  const int NP1 = N + 1;
109  fftw_r2r_kind kind = FFTW_REDFT00;
110  fftw_plan p = fftw_plan_many_r2r(1, &NP1, 1, (double*)b, NULL, 2, 1,
111  (double*)f, NULL, 1, 1, &kind, 0U);
112 
113  /* random Fourier coefficients */
114  {
115  int k;
116  printf("\n2) Random Fourier coefficients a_k, k=0,1,...,N:\n");
117  for (k = 0; k <= N; k++)
118  {
119  a[k] = 2.0*nfft_drand48() - 1.0; /* for debugging: use k+1 */
120  printf(" a_%-2d = %+5.3lE\n",k,creal(a[k]));
121  }
122  }
123 
124  /* fast polynomial transform */
125  fpt_trafo(set,0,a,b,N,0U);
126 
127  /* Renormalize coefficients b_j, j=1,2,...,N-1 owing to how FFTW defines a
128  * DCT-I; see
129  * http://www.fftw.org/fftw3_doc/1d-Real_002deven-DFTs-_0028DCTs_0029.html
130  * for details */
131  {
132  int j;
133  for (j = 1; j < N; j++)
134  b[j] *= 0.5;
135  }
136 
137  /* discrete cosine transform */
138  fftw_execute(p);
139 
140  {
141  int j;
142  printf("\n3) Function values f_j, j=1,1,...,M:\n");
143  for (j = 0; j <= N; j++)
144  printf(" f_%-2d = %+5.3lE\n",j,f[j]);
145  }
146 
147  /* cleanup */
148  free(a);
149  free(b);
150  free(f);
151 
152  /* cleanup */
153  fftw_destroy_plan(p);
154  }
155 
156  /* cleanup */
157  fpt_finalize(set);
158  free(alpha);
159  free(beta);
160  free(gamma);
161 
162  return EXIT_SUCCESS;
163 }

Generated on Tue Apr 30 2013 by Doxygen 1.8.1