NFFT Logo 3.1.2 API Reference

simple_test.c

00001 /*
00002  * Copyright (c) 2002, 2009 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 /* $Id: simple_test.c 3153 2009-04-07 11:05:56Z keiner $ */
00020 
00021 /* standard headers */
00022 #include <stdlib.h>
00023 #include <stdio.h>
00024 #include <math.h>
00025 
00026 /* It is important to include complex.h before nfft3.h. */
00027 #include <complex.h>
00028 
00029 /* NFFT3 header */
00030 #include "nfft3.h"
00031 
00032 /* We declare drand48() and srand48() ourselves, if they are is not declared in
00033  * math.h (e.g. on SuSE 9.3), grrr. */
00034 #include "config.h"
00035 #if HAVE_DECL_DRAND48 == 0
00036   extern double drand48(void);
00037 #endif
00038 #if HAVE_DECL_SRAND48 == 0
00039   extern void srand48(long int);
00040 #endif
00041 
00042 /* Two times Pi */
00043 #define KPI2 6.2831853071795864769252867665590057683943387987502
00044 
00045 int main(void)
00046 {
00047   /* This example shows the use of the fast polynomial transform to evaluate a
00048    * finite expansion in Legendre polynomials,
00049    *
00050    *   f(x) = a_0 P_0(x) + a_1 P_1(x) + ... + a_N P_N(x)                     (1)
00051    *
00052    * at the Chebyshev nodes x_j = cos(j*pi/N), j=0,1,...,N. */
00053   const int N = 8;
00054 
00055   /* An fpt_set is a data structure that contains precomputed data for a number
00056    * of different polynomial transforms. Here, we need only one transform. the
00057    * second parameter (t) is the exponent of the maximum transform size desired
00058    * (2^t), i.e., t = 3 means that N in (1) can be at most N = 8. */
00059   fpt_set set = fpt_init(1,lrint(ceil(log2((double)N))),0U);
00060 
00061   /* Three-term recurrence coefficients for Legendre polynomials */
00062   double *alpha = malloc((N+2)*sizeof(double)),
00063     *beta = malloc((N+2)*sizeof(double)),
00064     *gamma = malloc((N+2)*sizeof(double));
00065 
00066   /* alpha[0] and beta[0] are not referenced. */
00067   alpha[0] = beta[0] = 0.0;
00068   /* gamma[0] contains the value of P_0(x) (which is a constant). */
00069   gamma[0] = 1.0;
00070 
00071   /* Actual three-term recurrence coefficients for Legendre polynomials */
00072   {
00073     int k;
00074     for (k = 0; k <= N; k++)
00075     {
00076       alpha[k+1] = ((double)(2*k+1))/((double)(k+1));
00077       beta[k+1] = 0.0;
00078       gamma[k+1] = -((double)(k))/((double)(k+1));
00079     }
00080   }
00081 
00082   printf(
00083     "Computing a fast polynomial transform (FPT) and a fast discrete cosine \n"
00084     "transform (DCT) to evaluate\n\n"
00085     "  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"
00086     "with N=%d, x_j = cos(j*pi/N), j=0,1,...N, the Chebyshev nodes, a_k,\n"
00087     "k=0,1,...,N, random Fourier coefficients in [-1,1]x[-1,1]*I, and P_k,\n"
00088     "k=0,1,...,N, the Legendre polynomials.",N
00089   );
00090 
00091   /* Random seed, makes things reproducible. */
00092   srand48(616642);
00093 
00094   /* The function fpt_repcompute actually does the precomputation for a single
00095    * transform. It needs arrays alpha, beta, and gamma, containing the three-
00096    * term recurrence coefficients, here of the Legendre polynomials. The format
00097    * is explained above. The sixth parameter (k_start) is where the index in the
00098    * linear combination (1) starts, here k_start=0. The seventh parameter
00099    * (kappa) is the threshold which has an influence on the accuracy of the fast
00100    * polynomial transform. Usually, kappa = 1000 is a good choice. */
00101   fpt_precompute(set,0,alpha,beta,gamma,0,1000.0);
00102 
00103 
00104   {
00105     /* Arrays for Fourier coefficients and function values. */
00106     double _Complex *a = malloc((N+1)*sizeof(double _Complex));
00107     double _Complex *b = malloc((N+1)*sizeof(double _Complex));
00108     double *f = malloc((N+1)*sizeof(double _Complex));
00109 
00110     /* Plan for discrete cosine transform */
00111     const int NP1 = N + 1;
00112     fftw_r2r_kind kind = FFTW_REDFT00;
00113     fftw_plan p = fftw_plan_many_r2r(1, &NP1, 1, (double*)b, NULL, 2, 1,
00114       (double*)f, NULL, 1, 1, &kind, 0U);
00115 
00116     /* random Fourier coefficients */
00117     {
00118       int k;
00119       printf("\n2) Random Fourier coefficients a_k, k=0,1,...,N:\n");
00120       for (k = 0; k <= N; k++)
00121       {
00122         a[k] = 2.0*drand48() - 1.0; /* for debugging: use k+1 */
00123         printf("   a_%-2d = %+5.3lE\n",k,creal(a[k]));
00124       }
00125     }
00126 
00127     /* fast polynomial transform */
00128     fpt_trafo(set,0,a,b,N,0U);
00129 
00130     /* Renormalize coefficients b_j, j=1,2,...,N-1 owing to how FFTW defines a
00131      * DCT-I; see
00132      * http://www.fftw.org/fftw3_doc/1d-Real_002deven-DFTs-_0028DCTs_0029.html
00133      * for details */
00134     {
00135       int j;
00136       for (j = 1; j < N; j++)
00137         b[j] *= 0.5;
00138     }
00139 
00140     /* discrete cosine transform */
00141     fftw_execute(p);
00142 
00143     {
00144       int j;
00145       printf("\n3) Function values f_j, j=1,1,...,M:\n");
00146       for (j = 0; j <= N; j++)
00147         printf("   f_%-2d = %+5.3lE\n",j,f[j]);
00148     }
00149 
00150     /* cleanup */
00151     free(a);
00152     free(b);
00153     free(f);
00154 
00155     /* cleanup */
00156     fftw_destroy_plan(p);
00157   }
00158 
00159   /* cleanup */
00160   fpt_finalize(set);
00161   free(alpha);
00162   free(beta);
00163   free(gamma);
00164 
00165   return EXIT_SUCCESS;
00166 }

Generated on 16 Sep 2009 by Doxygen 1.5.3