00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022 #include <stdlib.h>
00023 #include <stdio.h>
00024 #include <math.h>
00025
00026
00027 #include <complex.h>
00028
00029
00030 #include "nfft3.h"
00031
00032
00033
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
00043 #define KPI2 6.2831853071795864769252867665590057683943387987502
00044
00045 int main(void)
00046 {
00047
00048
00049
00050
00051
00052
00053 const int N = 8, M = 10;
00054
00055
00056
00057
00058
00059 fpt_set set = fpt_init(1,lrint(ceil(log2((double)N))),0U);
00060
00061
00062 double *alpha = malloc((N+2)*sizeof(double)),
00063 *beta = malloc((N+2)*sizeof(double)), *gamma = malloc((N+2)*sizeof(double));
00064
00065
00066 alpha[0] = beta[0] = 0.0;
00067
00068 gamma[0] = 1.0;
00069
00070
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
00092 srand48(616642);
00093
00094
00095
00096
00097
00098
00099
00100
00101 fpt_precompute(set,0,alpha,beta,gamma,0,1000.0);
00102
00103 {
00104
00105 nfct_plan p;
00106
00107
00108 nfct_init_1d(&p, N+1, M);
00109
00110
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
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
00127 double _Complex *a = malloc((N+1)*sizeof(double _Complex));
00128 double _Complex *b = malloc((N+1)*sizeof(double _Complex));
00129
00130
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;
00137 printf(" a_%-2d = %+5.3lE\n",k,creal(a[k]));
00138 }
00139 }
00140
00141
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
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
00161 free(a);
00162 free(b);
00163 }
00164
00165
00166 nfct_finalize(&p);
00167 }
00168
00169
00170 fpt_finalize(set);
00171 free(alpha);
00172 free(beta);
00173 free(gamma);
00174
00175 return EXIT_SUCCESS;
00176 }