NFFT  3.3.1
simple_test.c
00001 /*
00002  * Copyright (c) 2002, 2016 Jens Keiner, Stefan Kunis, Daniel Potts
00003  *
00004  * This program is free software; you can redistribute it and/or modify it under
00005  * the terms of the GNU General Public License as published by the Free Software
00006  * Foundation; either version 2 of the License, or (at your option) any later
00007  * version.
00008  *
00009  * This program is distributed in the hope that it will be useful, but WITHOUT
00010  * ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
00011  * FOR A PARTICULAR PURPOSE.  See the GNU General Public License for more
00012  * details.
00013  *
00014  * You should have received a copy of the GNU General Public License along with
00015  * this program; if not, write to the Free Software Foundation, Inc., 51
00016  * Franklin Street, Fifth Floor, Boston, MA 02110-1301, USA.
00017  */
00018 
00019 /* standard headers */
00020 #include <stdlib.h>
00021 #include <stdio.h>
00022 #include <math.h>
00023 
00024 /* It is important to include complex.h before nfft3.h and fftw3.h. */
00025 #include <complex.h>
00026 
00027 #include <fftw3.h>
00028 
00029 /* NFFT3 header */
00030 #include "nfft3.h"
00031 
00032 int main(void)
00033 {
00034   /* This example shows the use of the fast polynomial transform to evaluate a
00035    * finite expansion in Legendre polynomials,
00036    *
00037    *   f(x) = a_0 P_0(x) + a_1 P_1(x) + ... + a_N P_N(x)                     (1)
00038    *
00039    * at the Chebyshev nodes x_j = cos(j*pi/N), j=0,1,...,N. */
00040   const int N = 8;
00041 
00042   /* An fpt_set is a data structure that contains precomputed data for a number
00043    * of different polynomial transforms. Here, we need only one transform. the
00044    * second parameter (t) is the exponent of the maximum transform size desired
00045    * (2^t), i.e., t = 3 means that N in (1) can be at most N = 8. */
00046   fpt_set set = fpt_init(1,lrint(ceil(log2((double)N))),0U);
00047 
00048   /* Three-term recurrence coefficients for Legendre polynomials */
00049   double *alpha = malloc((N+2)*sizeof(double)),
00050     *beta = malloc((N+2)*sizeof(double)),
00051     *gamma = malloc((N+2)*sizeof(double));
00052 
00053   /* alpha[0] and beta[0] are not referenced. */
00054   alpha[0] = beta[0] = 0.0;
00055   /* gamma[0] contains the value of P_0(x) (which is a constant). */
00056   gamma[0] = 1.0;
00057 
00058   /* Actual three-term recurrence coefficients for Legendre polynomials */
00059   {
00060     int k;
00061     for (k = 0; k <= N; k++)
00062     {
00063       alpha[k+1] = ((double)(2*k+1))/((double)(k+1));
00064       beta[k+1] = 0.0;
00065       gamma[k+1] = -((double)(k))/((double)(k+1));
00066     }
00067   }
00068 
00069   printf(
00070     "Computing a fast polynomial transform (FPT) and a fast discrete cosine \n"
00071     "transform (DCT) to evaluate\n\n"
00072     "  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"
00073     "with N=%d, x_j = cos(j*pi/N), j=0,1,...N, the Chebyshev nodes, a_k,\n"
00074     "k=0,1,...,N, random Fourier coefficients in [-1,1]x[-1,1]*I, and P_k,\n"
00075     "k=0,1,...,N, the Legendre polynomials.",N
00076   );
00077 
00078   /* Random seed, makes things reproducible. */
00079   nfft_srand48(314);
00080 
00081   /* The function fpt_repcompute actually does the precomputation for a single
00082    * transform. It needs arrays alpha, beta, and gamma, containing the three-
00083    * term recurrence coefficients, here of the Legendre polynomials. The format
00084    * is explained above. The sixth parameter (k_start) is where the index in the
00085    * linear combination (1) starts, here k_start=0. The seventh parameter
00086    * (kappa) is the threshold which has an influence on the accuracy of the fast
00087    * polynomial transform. Usually, kappa = 1000 is a good choice. */
00088   fpt_precompute(set,0,alpha,beta,gamma,0,1000.0);
00089 
00090 
00091   {
00092     /* Arrays for Fourier coefficients and function values. */
00093     double _Complex *a = malloc((N+1)*sizeof(double _Complex));
00094     double _Complex *b = malloc((N+1)*sizeof(double _Complex));
00095     double *f = malloc((N+1)*sizeof(double _Complex));
00096 
00097     /* Plan for discrete cosine transform */
00098     const int NP1 = N + 1;
00099     fftw_r2r_kind kind = FFTW_REDFT00;
00100     fftw_plan p = fftw_plan_many_r2r(1, &NP1, 1, (double*)b, NULL, 2, 1,
00101       (double*)f, NULL, 1, 1, &kind, 0U);
00102 
00103     /* random Fourier coefficients */
00104     {
00105       int k;
00106       printf("\n2) Random Fourier coefficients a_k, k=0,1,...,N:\n");
00107       for (k = 0; k <= N; k++)
00108       {
00109         a[k] = 2.0*nfft_drand48() - 1.0; /* for debugging: use k+1 */
00110         printf("   a_%-2d = %+5.3lE\n",k,creal(a[k]));
00111       }
00112     }
00113 
00114     /* fast polynomial transform */
00115     fpt_trafo(set,0,a,b,N,0U);
00116 
00117     /* Renormalize coefficients b_j, j=1,2,...,N-1 owing to how FFTW defines a
00118      * DCT-I; see
00119      * http://www.fftw.org/fftw3_doc/1d-Real_002deven-DFTs-_0028DCTs_0029.html
00120      * for details */
00121     {
00122       int j;
00123       for (j = 1; j < N; j++)
00124         b[j] *= 0.5;
00125     }
00126 
00127     /* discrete cosine transform */
00128     fftw_execute(p);
00129 
00130     {
00131       int j;
00132       printf("\n3) Function values f_j, j=1,1,...,M:\n");
00133       for (j = 0; j <= N; j++)
00134         printf("   f_%-2d = %+5.3lE\n",j,f[j]);
00135     }
00136 
00137     /* cleanup */
00138     free(a);
00139     free(b);
00140     free(f);
00141 
00142     /* cleanup */
00143     fftw_destroy_plan(p);
00144   }
00145 
00146   /* cleanup */
00147   fpt_finalize(set);
00148   free(alpha);
00149   free(beta);
00150   free(gamma);
00151 
00152   return EXIT_SUCCESS;
00153 }