S2kit  1.1
Toolkit for working with functions defined on the sphere
test_s2_semi_fly.c
Go to the documentation of this file.
1 
46 #include <math.h>
47 #include <stdio.h>
48 #include <stdlib.h>
49 #include <time.h>
50 
51 #include <fftw3.h>
52 
53 #include "s2kit/FST_semi_fly.h"
54 #include "s2kit/weights.h"
55 #include "s2kit/util.h"
56 
57 #include "util/csecond.h"
58 
59 int main(int argc, char** argv) {
60  if (argc < 3) {
61  fprintf(stdout, "Usage: test_s2_semi_fly bw loops [error_file]\n");
62  exit(0);
63  }
64 
65  int bw = atoi(argv[1]);
66  int loops = atoi(argv[2]);
67 
68  int size = 2 * bw;
69 
70  time_t seed;
71  time(&seed);
72  srand48(seed); // init random generator
73 
74  double* weights = (double*)malloc(sizeof(double) * 4 * bw);
75  double* rdata = (double*)malloc(sizeof(double) * (size * size));
76  double* idata = (double*)malloc(sizeof(double) * (size * size));
77  double* workspace = (double*)malloc(sizeof(double) * ((10 * (bw * bw)) + (24 * bw)));
78 
79  // Make DCT plans. Note that I will be using the GURU interface to execute these plans within the routines
80  fftw_plan DCT_plan = fftw_plan_r2r_1d(size, weights, rdata, FFTW_REDFT10, FFTW_ESTIMATE);
81  fftw_plan inv_DCT_plan = fftw_plan_r2r_1d(size, weights, rdata, FFTW_REDFT01, FFTW_ESTIMATE);
82 
83  // fftw "preamble"
84  // Note that FFT plan places the output in a transposed array
85  int rank = 1;
86  fftw_iodim dims[rank];
87  dims[0].n = size;
88  dims[0].is = 1;
89  dims[0].os = size;
90 
91  int howmany_rank = 1;
92  fftw_iodim howmany_dims[howmany_rank];
93  howmany_dims[0].n = size;
94  howmany_dims[0].is = size;
95  howmany_dims[0].os = 1;
96 
97  fftw_plan FFT_plan = fftw_plan_guru_split_dft(rank, dims, howmany_rank, howmany_dims, rdata, idata, workspace,
98  workspace + (4 * bw * bw), FFTW_ESTIMATE);
99 
100  // Note that FFT plans assumes that I'm working with a transposed array, e.g. the inputs for a length 2*bw transform
101  // are placed every 2*bw apart, the output will be consecutive entries in the array
102 
103  rank = 1;
104  dims[0].n = size;
105  dims[0].is = size;
106  dims[0].os = 1;
107 
108  howmany_rank = 1;
109  howmany_dims[0].n = size;
110  howmany_dims[0].is = 1;
111  howmany_dims[0].os = size;
112 
113  fftw_plan inv_FFT_plan = fftw_plan_guru_split_dft(rank, dims, howmany_rank, howmany_dims, rdata, idata, workspace,
114  workspace + (4 * bw * bw), FFTW_ESTIMATE);
115 
116  GenerateWeightsForDLT(bw, weights);
117 
118  double* rcoeffs = (double*)malloc(sizeof(double) * (bw * bw));
119  double* icoeffs = (double*)malloc(sizeof(double) * (bw * bw));
120 
121  double* rresult = (double*)malloc(sizeof(double) * (bw * bw));
122  double* iresult = (double*)malloc(sizeof(double) * (bw * bw));
123 
124  double* relerror = (double*)malloc(sizeof(double) * loops);
125  double* curmax = (double*)malloc(sizeof(double) * loops);
126 
127  // TODO more tests with dataformat 0 and cutoff
128  int cutoff = (bw / 2) - 10; // = bw // seminaive all orders
129  DataFormat data_format = COMPLEX;
130 
131  double fwd_time = 0.0;
132  double inv_time = 0.0;
133  double total_error = 0.0;
134  double total_relerror = 0.0;
135 
136  fprintf(stdout, "about to enter loop\n\n");
137  for (int i = 0; i < loops; ++i) {
138 
139  // generate spherical harmonic coefficients of a real-valued function
140  for (int m = 0; m < bw; ++m)
141  for (int l = m; l < bw; ++l) {
142  double x = 2.0 * (drand48() - 0.5);
143  double y = 2.0 * (drand48() - 0.5);
144 
145  int index = IndexOfHarmonicCoeff(m, l, bw);
146  rcoeffs[index] = x;
147  icoeffs[index] = y;
148 
149  index = IndexOfHarmonicCoeff(-m, l, bw);
150  rcoeffs[index] = pow(-1.0, m) * x;
151  icoeffs[index] = pow(-1.0, m + 1) * y;
152  }
153 
154  // have to zero out the m=0 coefficients, since those are real
155  for (int m = 0; m < bw; ++m)
156  icoeffs[m] = 0.0;
157 
158  double time_start = csecond();
159  // inverse spherical transform
160  InvFSTSemiFly(rcoeffs, icoeffs, rdata, idata, bw, workspace, data_format, cutoff, &inv_DCT_plan, &inv_FFT_plan);
161 
162  double duration = csecond() - time_start;
163  inv_time += duration;
164  fprintf(stdout, "inv time \t = %.4e\n", duration);
165 
166  time_start = csecond();
167  // forward spherical transform
168  FSTSemiFly(rdata, idata, rresult, iresult, bw, workspace, data_format, cutoff, &DCT_plan, &FFT_plan, weights);
169 
170  duration = csecond() - time_start;
171  fwd_time += duration;
172  fprintf(stdout, "forward time \t = %.4e\n", duration);
173 
174  // compute the error
175  relerror[i] = 0.0;
176  curmax[i] = 0.0;
177 
178  for (int j = 0; j < (bw * bw); ++j) {
179  double rtmp = rresult[j] - rcoeffs[j];
180  double itmp = iresult[j] - icoeffs[j];
181  double origmag = sqrt((rcoeffs[j] * rcoeffs[j]) + (icoeffs[j] * icoeffs[j]));
182  double tmpmag = sqrt((rtmp * rtmp) + (itmp * itmp));
183  relerror[i] = fmax(relerror[i], tmpmag / (origmag + pow(10.0, -50.0)));
184  curmax[i] = fmax(curmax[i], tmpmag);
185  }
186 
187  fprintf(stdout, "r-o error\t = %.12f\n", curmax[i]);
188  fprintf(stdout, "(r-o)/o error\t = %.12f\n\n", relerror[i]);
189 
190  total_error += curmax[i];
191  total_relerror += relerror[i];
192  }
193 
194  double avg_error = total_error / loops;
195  double avg_relerror = total_relerror / loops;
196  double stddev_error = 0.0;
197  double stddev_relerror = 0.0;
198 
199  for (int i = 0; i < loops; ++i) {
200  stddev_error += pow(avg_error - curmax[i], 2.0);
201  stddev_relerror += pow(avg_relerror - relerror[i], 2.0);
202  }
203 
204  if (loops != 1) {
205  stddev_error = sqrt(stddev_error / (loops - 1));
206  stddev_relerror = sqrt(stddev_relerror / (loops - 1));
207  }
208 
209  fprintf(stdout, "Program: test_s2_semi_fly\n");
210  fprintf(stdout, "Bandwidth = %d\n", bw);
211 
212  const char* timing_object = "cpu";
213 #ifdef WALLCLOCK
214  timing_object = "wall";
215 #endif
216 
217  double total_time = inv_time + fwd_time;
218 
219  fprintf(stdout, "Total elapsed %s time:\t\t\t %.4e seconds.\n", timing_object, total_time);
220  fprintf(stdout, "Average %s time forward per iteration:\t %.4e seconds.\n", timing_object, fwd_time / loops);
221  fprintf(stdout, "Average %s time inverse per iteration:\t %.4e seconds.\n", timing_object, inv_time / loops);
222 
223  fprintf(stdout, "Average r-o error:\t\t %.4e\t", total_error / loops);
224  fprintf(stdout, "std dev: %.4e\n", stddev_error);
225  fprintf(stdout, "Average (r-o)/o error:\t\t %.4e\t", total_relerror / loops);
226  fprintf(stdout, "std dev: %.4e\n\n", stddev_relerror);
227 
228  if (argc == 4) {
229  FILE* errorsfp = fopen(argv[3], "w");
230 
231  for (int m = 0; m < bw; ++m) {
232  for (int l = m; l < bw; ++l) {
233  int index = IndexOfHarmonicCoeff(m, l, bw);
234  fprintf(errorsfp, "index = %d\t m = %d\tl = %d\t%.10f %.10f\n", index, m, l,
235  fabs(rcoeffs[index] - rresult[index]), fabs(icoeffs[index] - iresult[index]));
236 
237  index = IndexOfHarmonicCoeff(-m, l, bw);
238  fprintf(errorsfp, "index = %d\t m = %d\tl = %d\t%.10f %.10f\n", index, -m, l,
239  fabs(rcoeffs[index] - rresult[index]), fabs(icoeffs[index] - iresult[index]));
240  }
241  }
242 
243  fclose(errorsfp);
244  }
245 
246  fftw_destroy_plan(inv_FFT_plan);
247  fftw_destroy_plan(FFT_plan);
248  fftw_destroy_plan(inv_DCT_plan);
249  fftw_destroy_plan(DCT_plan);
250 
251  free(weights);
252  free(curmax);
253  free(relerror);
254  free(workspace);
255  free(iresult);
256  free(rresult);
257  free(idata);
258  free(rdata);
259  free(icoeffs);
260  free(rcoeffs);
261 
262  return 0;
263 }
FST_semi_fly.h
weights.h
main
int main(int argc, char **argv)
Definition: test_s2_semi_fly.c:59
IndexOfHarmonicCoeff
int IndexOfHarmonicCoeff(const int, const int, const int)
Gives the position of the coefficient f-hat(m,l) in the one-row array.
Definition: util.c:42
COMPLEX
Definition: util.h:11
FSTSemiFly
void FSTSemiFly(double *, double *, double *, double *, const int, double *, DataFormat, const int, fftw_plan *, fftw_plan *, double *)
Computes spherical harmonic transform using the seminaive and naive algorithms.
Definition: FST_semi_fly.c:67
GenerateWeightsForDLT
void GenerateWeightsForDLT(const int, double *)
Generates weights for both even and odd order Legendre transforms for a given bandwidth.
Definition: weights.c:32
InvFSTSemiFly
void InvFSTSemiFly(double *, double *, double *, double *, const int, double *, DataFormat, const int, fftw_plan *, fftw_plan *)
Computes inverse spherical harmonic transform.
Definition: FST_semi_fly.c:235
DataFormat
DataFormat
Result data format for FST functions.
Definition: util.h:10
csecond.h
csecond
double csecond()
Definition: csecond.c:32
util.h