NFFT Logo 3.1.0 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 3110 2009-03-13 16:32:18Z 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 how to use 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 a number of M arbitrary nodes x_j in [-1,1]. */
00053   const int N = 8, M = 10;
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)), *gamma = malloc((N+2)*sizeof(double));
00064 
00065   /* alpha[0] and beta[0] are not referenced. */
00066   alpha[0] = beta[0] = 0.0;
00067   /* gamma[0] contains P_0(x) which is a constant. */
00068   gamma[0] = 1.0;
00069 
00070   /* Actual three-term recurrence coefficients for Legendre polynomials */
00071   {
00072     int k;
00073     for (k = 0; k <= N; k++)
00074     {
00075       alpha[k+1] = ((double)(2*k+1))/((double)(k+1));
00076       beta[k+1] = 0.0;
00077       gamma[k+1] = -((double)(k))/((double)(k+1));
00078     }
00079   }
00080 
00081   printf(
00082     "Computing a fast polynomial transform (FPT) and a fast nonequispaced \n"
00083     "cosine transform (NFCT) to evaluate\n\n"
00084     "  f_j = a_0 P_0(x) + a_1 P_1(x) + ... + a_N P_N(x), j=1,2,...,M,\n\n"
00085     "with N=%d, M=%d, x_j random nodes in [-1,1], and a_k, k=0,1,...,N random\n"
00086     "Fourier coefficients in [-1,1], where P_k, k=0,1,...,N are the Legendre "
00087     "polynomials.",
00088     N,M
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     /* Plan for NFCT */
00105     nfct_plan p;
00106 
00107     /* Init transform plan. */
00108     nfct_init_1d(&p, N+1, M);
00109 
00110     /* random nodes */
00111     {
00112       int j;
00113       printf("\n1) Random nodes x_j, j=1,2,...,M:\n");
00114       for (j = 0; j < M; j++)
00115       {
00116         p.x[j] = acos(2.0*drand48() - 1.0)/KPI2;
00117         /* for debugging: use acos(2.0*((double)j)/((double)(M-1))-1.0)/KPI2 */
00118         printf("   x_%-2d = %+5.3lE\n",j+1,cos(KPI2*p.x[j]));
00119       }
00120     }
00121 
00122     if(p.nfct_flags & PRE_PSI)
00123       nfct_precompute_psi(&p);
00124 
00125     {
00126       /* Array for Fourier coefficients. */
00127       double _Complex *a = malloc((N+1)*sizeof(double _Complex));
00128       double _Complex *b = malloc((N+1)*sizeof(double _Complex));
00129 
00130       /* random Fourier coefficients */
00131       {
00132         int k;
00133         printf("\n2) Random Fourier coefficients a_k, k=0,1,...,N:\n");
00134         for (k = 0; k <= N; k++)
00135         {
00136           a[k] = 2.0*drand48() - 1.0; /* for debugging: use k+1 */
00137           printf("   a_%-2d = %+5.3lE\n",k,creal(a[k]));
00138         }
00139       }
00140 
00141       /* fast polynomial transform */
00142       fpt_trafo(set,0,a,b,N,0U);
00143 
00144       {
00145         int k;
00146         for (k = 0; k <= N; k++)
00147           p.f_hat[k] = creal(b[k]);
00148       }
00149 
00150       /* nonequispaced discrete cosine transform */
00151       nfct_trafo(&p);
00152 
00153       {
00154         int j;
00155         printf("\n3) Function values f_j, j=1,1,...,M:\n");
00156         for (j = 0; j < M; j++)
00157           printf("   f_%-2d = %+5.3lE\n",j+1,p.f[j]);
00158       }
00159 
00160       /* cleanup */
00161       free(a);
00162       free(b);
00163     }
00164 
00165     /* cleanup */
00166     nfct_finalize(&p);
00167   }
00168 
00169   /* cleanup */
00170   fpt_finalize(set);
00171   free(alpha);
00172   free(beta);
00173   free(gamma);
00174 
00175   return EXIT_SUCCESS;
00176 }

Generated on 19 Mar 2009 by Doxygen 1.5.3