S2kit  1.1
Toolkit for working with functions defined on the sphere
test_DLT_semi.c
Go to the documentation of this file.
1 
20 #include <math.h>
21 #include <stdio.h>
22 #include <stdlib.h>
23 #include <time.h>
24 
25 #include <fftw3.h>
26 
27 #include "s2kit/cospml.h"
28 #include "s2kit/seminaive.h"
29 #include "s2kit/weights.h"
30 #include "s2kit/chebyshev_nodes.h"
31 
32 #include "util/csecond.h"
33 
34 int main(int argc, char** argv) {
35  if (argc < 4) {
36  fprintf(stdout, "Usage: test_DLT_semi m bw loops\n");
37  exit(0);
38  }
39 
40  int m = atoi(argv[1]);
41  int bw = atoi(argv[2]);
42  int loops = atoi(argv[3]);
43  int n = 2 * bw;
44 
45  double* samples = (double*)malloc(sizeof(double) * n);
46  double* coeffs = (double*)malloc(sizeof(double) * (bw - m));
47  double* new_coeffs = (double*)malloc(sizeof(double) * (bw - m));
48 
49  double* relerror = (double*)malloc(sizeof(double) * loops);
50  double* curmax = (double*)malloc(sizeof(double) * loops);
51 
52  double* eval_args = (double*)malloc(sizeof(double) * n);
53  double* sin_values = (double*)malloc(sizeof(double) * n);
54 
55  // generate sin values (need if order of transform is odd)
56  AcosOfChebyshevNodes(n, eval_args);
57  for (int i = 0; i < n; ++i)
58  sin_values[i] = sin(eval_args[i]);
59 
60  double* workspace = (double*)malloc(sizeof(double) * 9 * bw);
61  double* weights = (double*)malloc(sizeof(double) * 4 * bw);
62 
63  // Note the arrays! Since I'll be using the guru fftw interface,
64  // I don't care what the arrays are - I just want to make sure they're big enough
65  fftw_plan forwardDCT = fftw_plan_r2r_1d(2 * bw, workspace, weights, FFTW_REDFT10, FFTW_ESTIMATE);
66  fftw_plan inverseDCT = fftw_plan_r2r_1d(2 * bw, workspace, weights, FFTW_REDFT01, FFTW_ESTIMATE);
67 
68  GenerateWeightsForDLT(bw, weights);
69 
70  double* cos_pml_table = (double*)malloc(sizeof(double) * TableSize(m, bw));
71  double* transpose_cos_pml_table = (double*)malloc(sizeof(double) * TableSize(m, bw));
72 
73  // precompute cosine series for Pml (Gml) functions, necessary for the forward transform
74  GenerateCosPmlTable(bw, m, cos_pml_table, workspace);
75  // take transpose of the precomputed cosine series, necessary for the inverse transform
76  TransposeCosPmlTable(bw, m, cos_pml_table, transpose_cos_pml_table);
77 
78  double sum_error = 0.0;
79  double sum_relerror = 0.0;
80 
81  double total_time_f = 0.0;
82  double total_time_i = 0.0;
83 
84  long int seed;
85  time(&seed);
86  srand48(seed); // init random generator
87 
88  for (int k = 0; k < loops; ++k) {
89  // generate random coefficients
90  for (int i = 0; i < (bw - m); ++i)
91  coeffs[i] = 2.0 * (drand48() - 0.5);
92 
93  double time_start = csecond();
94  // inverse semi-naive transform
95  InvDLTSemi(coeffs, bw, m, samples, transpose_cos_pml_table, sin_values, workspace, &inverseDCT);
96  total_time_i += (csecond() - time_start);
97 
98  time_start = csecond();
99  // forward semi-naive transform
100  DLTSemi(samples, bw, m, new_coeffs, workspace, cos_pml_table, weights, &forwardDCT);
101  total_time_f += (csecond() - time_start);
102 
103  // tally up the error between the original coefficients and the new ones
104  relerror[k] = 0.0;
105  curmax[k] = 0.0;
106 
107  for (int j = 0; j < bw - m; ++j) {
108  double tmp_error = fabs(coeffs[j] - new_coeffs[j]);
109  double tmp_relerror = tmp_error / (fabs(coeffs[j]) + pow(10.0, -50.0));
110  curmax[k] = fmax(curmax[k], tmp_error);
111  relerror[k] = fmax(relerror[k], tmp_relerror);
112  }
113 
114  sum_error += curmax[k];
115  sum_relerror += relerror[k];
116  }
117 
118  double avg_error = sum_error / loops;
119  double avg_relerror = sum_relerror / loops;
120  double stddev_error = 0.0;
121  double stddev_relerror = 0.0;
122  for (int i = 0; i < loops; ++i) {
123  stddev_error += pow(avg_error - curmax[i], 2.0);
124  stddev_relerror += pow(avg_relerror - relerror[i], 2.0);
125  }
126 
127  if (loops != 1) {
128  stddev_error = sqrt(stddev_error / (loops - 1));
129  stddev_relerror = sqrt(stddev_relerror / (loops - 1));
130  }
131 
132  fprintf(stdout, "bw = %d\tm = %d\n", bw, m);
133  fprintf(stdout, "loops = %d\n", loops);
134 
135  fprintf(stdout, "Average r-o error:\t\t %.4e\t", sum_error / loops);
136  fprintf(stdout, "std dev: %.4e\n", stddev_error);
137  fprintf(stdout, "Average (r-o)/o error:\t\t %.4e\t", sum_relerror / loops);
138  fprintf(stdout, "std dev: %.4e\n\n", stddev_relerror);
139 
140  fprintf(stdout, "average forward time = %.4e\n", total_time_f / loops);
141  fprintf(stdout, "average inverse time = %.4e\n", total_time_i / loops);
142 
143  fftw_destroy_plan(inverseDCT);
144  fftw_destroy_plan(forwardDCT);
145 
146  free(weights);
147  free(workspace);
148  free(curmax);
149  free(relerror);
150  free(sin_values);
151  free(eval_args);
152  free(transpose_cos_pml_table);
153  free(cos_pml_table);
154  free(new_coeffs);
155  free(coeffs);
156  free(samples);
157 
158  return 0;
159 }
TransposeCosPmlTable
void TransposeCosPmlTable(const int, const int, double *, double *)
Definition: cospml.c:301
weights.h
GenerateCosPmlTable
void GenerateCosPmlTable(const int, const int, double *, double *)
Definition: cospml.c:161
chebyshev_nodes.h
InvDLTSemi
void InvDLTSemi(double *, const int, const int, double *, double *, double *, double *, fftw_plan *)
Computes the inverse Legendre transform using the transposed seminaive algorithm.
Definition: seminaive.c:56
GenerateWeightsForDLT
void GenerateWeightsForDLT(const int, double *)
Generates weights for both even and odd order Legendre transforms for a given bandwidth.
Definition: weights.c:32
main
int main(int argc, char **argv)
Definition: test_DLT_semi.c:34
AcosOfChebyshevNodes
void AcosOfChebyshevNodes(const int, double *)
Generates an array of the angular arguments of n Chebyshev nodes.
Definition: chebyshev_nodes.c:16
TableSize
int TableSize(const int, const int)
Computes the number of non-zero entries of coeffs in a table for seminaive transform.
Definition: cospml.c:39
cospml.h
DLTSemi
void DLTSemi(double *, const int, const int, double *, double *, double *, double *, fftw_plan *)
Computes the Legendre transform of data using the seminaive algorithm.
Definition: seminaive.c:153
csecond.h
csecond
double csecond()
Definition: csecond.c:32
seminaive.h