My Project
so3_test.c
Go to the documentation of this file.
1 // S03 package to perform Wigner transform on the rotation group SO(3)
2 // Copyright (C) 2013 Martin Büttner and Jason McEwen
3 // See LICENSE.txt for license details
4 
27 #include <stdio.h>
28 #include <stdlib.h>
29 #include <complex.h>
30 #include <time.h>
31 #include <fftw3.h>
32 
33 #include "so3.h"
34 
35 #define NREPEAT 2
36 #define MIN(a,b) ((a < b) ? (a) : (b))
37 #define MAX(a,b) ((a > b) ? (a) : (b))
38 
39 double get_max_error(complex double *expected, complex double *actual, int n);
40 double ran2_dp(int idum);
41 void so3_test_gen_flmn_complex(complex double *flmn, const so3_parameters_t *parameters, int seed);
42 void so3_test_gen_flmn_real(complex double *flmn, const so3_parameters_t *parameters, int seed);
43 
44 int main(int argc, char **argv)
45 {
46  so3_parameters_t parameters = {};
47  int L, N, L0;
48  complex double *flmn_orig, *flmn_syn_direct, *flmn_syn_ssht;
49  complex double *f_direct, *f_ssht;
50  double *f_real_direct, *f_real_ssht;
51  int seed;
52  clock_t time_start, time_end, time_start_ssht, time_end_ssht, time_start_direct, time_end_direct;
53  int i, sampling_scheme, n_order, storage_mode, n_mode, real, routine, steerable;
54  int flmn_size;
55 
56  const char *n_order_str[SO3_N_ORDER_SIZE];
57  const char *storage_mode_str[SO3_STORAGE_SIZE];
58  const char *n_mode_str[SO3_N_MODE_SIZE];
59  const char *reality_str[2];
60  //const char *routine_str[2];
61  const char *sampling_str[2];
62  const char *steerable_str[2];
63 
64  n_order_str[SO3_N_ORDER_ZERO_FIRST] = "n = 0 first";
65  n_order_str[SO3_N_ORDER_NEGATIVE_FIRST] = "n = -N+1 first";
66 
67  storage_mode_str[SO3_STORAGE_PADDED] = "PADDED storage";
68  storage_mode_str[SO3_STORAGE_COMPACT] = "COMPACT storage";
69 
70  n_mode_str[SO3_N_MODE_ALL] = "ALL n";
71  n_mode_str[SO3_N_MODE_EVEN] = "only EVEN n";
72  n_mode_str[SO3_N_MODE_ODD] = "only ODD n";
73  n_mode_str[SO3_N_MODE_MAXIMUM] = "only |n| = N-1";
74  n_mode_str[SO3_N_MODE_L] = "only |n| = el";
75 
76  reality_str[0] = "COMPLEX";
77  reality_str[1] = "REAL";
78 
79  //routine_str[0] = "using SSHT";
80  //routine_str[1] = "not using SSHT";
81 
82  sampling_str[0] = "MW";
83  sampling_str[1] = "MWSS";
84 
85  steerable_str[0] = "NOT STEERABLE";
86  steerable_str[1] = "STEERABLE";
87 
88  // one element for each combination of sampling scheme, storage mode, n-mode and
89  // real or complex signal.
90  //double errors[2][SO3_SAMPLING_SIZE][SO3_N_MODE_SIZE][SO3_STORAGE_SIZE][SO3_N_MODE_SIZE][2][2];
91  //double durations_forward[2][SO3_SAMPLING_SIZE][SO3_N_MODE_SIZE][SO3_STORAGE_SIZE][SO3_N_MODE_SIZE][2][2];
92  //double durations_inverse[2][SO3_SAMPLING_SIZE][SO3_N_MODE_SIZE][SO3_STORAGE_SIZE][SO3_N_MODE_SIZE][2][2];
93  double errors[2][SO3_SAMPLING_SIZE][SO3_N_MODE_SIZE][SO3_STORAGE_SIZE][SO3_N_MODE_SIZE][2];
94  double durations_forward[2][SO3_SAMPLING_SIZE][SO3_N_MODE_SIZE][SO3_STORAGE_SIZE][SO3_N_MODE_SIZE][2];
95  double durations_inverse[2][SO3_SAMPLING_SIZE][SO3_N_MODE_SIZE][SO3_STORAGE_SIZE][SO3_N_MODE_SIZE][2];
96 
97  //fftw_init_threads();
98  //int nthreads = omp_get_max_threads();
99  //printf("Using %d threads.\n", nthreads);
100  //fftw_plan_with_nthreads(nthreads);
101 
102  // Parse command line arguments
103  L = N = 4;
104  int show_arrays = 2;
105  L0 = 0;
106  seed = 1;
107  if (argc > 1)
108  {
109  L = atoi(argv[1]);
110  if (argc > 2)
111  N = atoi(argv[2]);
112  else
113  N = L;
114  }
115 
116  if (argc > 3)
117  L0 = atoi(argv[3]);
118 
119  if (argc > 4)
120  seed = atoi(argv[4]);
121 
122  parameters.L0 = L0;
123  parameters.L = L;
124  parameters.N = N;
125  parameters.verbosity = 0;
126 
127  // (2*N-1)*L*L is the largest number of flmn ever needed. For more
128  // compact storage modes, only part of the memory will be used.
129  flmn_orig = malloc((2*N-1)*L*L * sizeof *flmn_orig);
130  SO3_ERROR_MEM_ALLOC_CHECK(flmn_orig);
131  flmn_syn_direct = malloc((2*N-1)*L*L * sizeof *flmn_syn_direct);
132  SO3_ERROR_MEM_ALLOC_CHECK(flmn_syn_direct);
133  flmn_syn_ssht = malloc((2*N-1)*L*L * sizeof *flmn_syn_ssht);
134  SO3_ERROR_MEM_ALLOC_CHECK(flmn_syn_ssht);
135 
136  // We only need (2*L) * (L+1) * (2*N-1) samples for MW symmetric sampling.
137  // For the usual MW sampling, only part of the memory will be used.
138  f_direct = malloc((2*L)*(L+1)*(2*N-1) * sizeof *f_direct);
139  SO3_ERROR_MEM_ALLOC_CHECK(f_direct);
140  f_real_direct = malloc((2*L)*(L+1)*(2*N-1) * sizeof *f_real_direct);
141  SO3_ERROR_MEM_ALLOC_CHECK(f_real_direct);
142  f_ssht = malloc((2*L)*(L+1)*(2*N-1) * sizeof *f_ssht);
143  SO3_ERROR_MEM_ALLOC_CHECK(f_ssht);
144  f_real_ssht = malloc((2*L)*(L+1)*(2*N-1) * sizeof *f_real_ssht);
145  SO3_ERROR_MEM_ALLOC_CHECK(f_real_ssht);
146 
147  // Write program name.
148  printf("\n");
149  printf("SO3 test program (C implementation)\n");
150  printf("================================================================\n");
151 
152  // steerable == 0 --> don't use steerable
153  // steerable == 1 --> use steerable
154  //for (n_mode = 0; n_mode < SO3_N_MODE_SIZE; ++ n_mode)
155  for (n_mode = 0; n_mode < 2; ++ n_mode)
156  {
157  parameters.n_mode = n_mode;
158 
159  for (real = 0; real < 1; ++real)
160  {
161  parameters.reality = real;
162 
163  // MWSS not yet supported by direct routines
164  for (sampling_scheme = 0; sampling_scheme < 1; ++sampling_scheme)
165  {
166  parameters.sampling_scheme = sampling_scheme;
167 
168  // For real signals, the n_order does not matter, so skip the second option
169  // in that case.
170  //for (n_order = 0; n_order < SO3_N_ORDER_SIZE - real; ++n_order)
171  for (n_order = 0; n_order < 1; ++n_order)
172  {
173  parameters.n_order = n_order;
174 
175  //for (storage_mode = 0; storage_mode < SO3_STORAGE_SIZE; ++storage_mode)
176  for (storage_mode = 0; storage_mode < 1; ++storage_mode)
177  {
178  parameters.storage = storage_mode;
179 
180  for (steerable = 0; steerable < 2; ++steerable)
181  {
182  parameters.steerable = steerable;
183 
184  durations_inverse[steerable][sampling_scheme][n_order][storage_mode][n_mode][real] = 0.0;
185  durations_forward[steerable][sampling_scheme][n_order][storage_mode][n_mode][real] = 0.0;
186  errors[steerable][sampling_scheme][n_order][storage_mode][n_mode][real] = 0.0;
187 
188  printf("\n");
189  printf("Testing a %s signal with %s with %s sampling using %s with %s. N-Mode: %s. Running %d times: ",
190  reality_str[real],
191  steerable_str[steerable],
192  sampling_str[sampling_scheme],
193  storage_mode_str[storage_mode],
194  n_order_str[n_order],
195  n_mode_str[n_mode],
196  NREPEAT);
197  printf("\n");
198 
199 
200  printf("----------------------------------------------------------------------------------------------------------------------------------------------------------\n");
201  printf("Tot Variation:");
202  printf(" |-| ");
203  printf(" Inv Variation:");
204  printf(" |-| ");
205  printf(" Direct Error: ");
206  printf("|-| ");
207  printf(" SSHT Error:");
208  printf(" |-| ");
209  printf(" Inv T Direct:");
210  printf(" |-| ");
211  printf(" For T Direct:");
212  printf(" |-| ");
213  printf(" Inv T SSHT:");
214  printf(" |-| ");
215  printf(" For T SSHT:\n");
216  printf("----------------------------------------------------------------------------------------------------------------------------------------------------------\n");
217 
218  for (i = 0; i <NREPEAT; ++i)
219  {
220  int j;
221  double duration, tot, tot_inverse, error_direct, error_ssht;
222  double inverse_duration_ssht, forward_duration_ssht;
223  double inverse_duration_direct, forward_duration_direct;
224 
225 
226 
227  // Reset output array
228  for (j = 0; j < (2*N-1)*L*L; ++j)
229  {
230  flmn_syn_direct[j] = 0.0;
231  flmn_syn_ssht[j] = 0.0;
232  }
233 
234  //count = 0.0;
235  tot = 0.0;
236  tot_inverse = 0.0;
237  error_direct = 0.0;
238  error_ssht = 0.0;
239 
240  //parameters.steerable = steerable;
241 
242  if (real) so3_test_gen_flmn_real(flmn_orig, &parameters, seed);
243  else so3_test_gen_flmn_complex(flmn_orig, &parameters, seed);
244 
245  time_start = clock();
246 
247  time_start_direct = clock();
248  if (real) so3_core_inverse_direct_real(f_real_direct, flmn_orig, &parameters);
249  else so3_core_inverse_direct(f_direct, flmn_orig, &parameters);
250  time_end_direct = clock();
251 
252  inverse_duration_direct = (time_end_direct - time_start_direct)/ (double)CLOCKS_PER_SEC;
253 
254  time_start_ssht = clock();
255  if (real) so3_core_inverse_via_ssht_real(f_real_ssht, flmn_orig, &parameters);
256  else so3_core_inverse_via_ssht(f_ssht, flmn_orig, &parameters);
257  time_end_ssht = clock();
258 
259  inverse_duration_ssht = (time_end_ssht - time_start_ssht)/ (double)CLOCKS_PER_SEC;
260 
261  time_end = clock();
262 
263  duration = (time_end - time_start) / (double)CLOCKS_PER_SEC;
264  if (!i || duration < durations_inverse[steerable][sampling_scheme][n_order][storage_mode][n_mode][real])
265  durations_inverse[steerable][sampling_scheme][n_order][storage_mode][n_mode][real] = duration;
266 
267  //parameters.steerable = 0;
268 
269  time_start = clock();
270 
271  time_start_direct = clock();
272  if (real) so3_core_forward_direct_real(flmn_syn_direct, f_real_direct, &parameters);
273  else so3_core_forward_direct(flmn_syn_direct, f_direct, &parameters);
274  time_end_direct = clock();
275 
276  forward_duration_direct = (time_end_direct - time_start_direct)/ (double)CLOCKS_PER_SEC;
277 
278  time_start_ssht = clock();
279  if (real) so3_core_forward_via_ssht_real(flmn_syn_ssht, f_real_ssht, &parameters);
280  else so3_core_forward_via_ssht(flmn_syn_ssht, f_ssht, &parameters);
281  time_end_ssht = clock();
282 
283  forward_duration_ssht = (time_end_ssht - time_start_ssht)/ (double)CLOCKS_PER_SEC;
284 
285  time_end = clock();
286 
287  duration = (time_end - time_start) / (double)CLOCKS_PER_SEC;
288  if (!i || duration < durations_forward[steerable][sampling_scheme][n_order][storage_mode][n_mode][real])
289  durations_forward[steerable][sampling_scheme][n_order][storage_mode][n_mode][real] = duration;
290 
291  flmn_size = so3_sampling_flmn_size(&parameters);
292  errors[steerable][sampling_scheme][n_order][storage_mode][n_mode][real] += get_max_error(flmn_orig, flmn_syn_ssht, flmn_size)/NREPEAT;
293  //printf("%.10e,",*flmn_syn);
294  //printf("%.10e",*flmn_orig);
295  //printf(" ---");
296 
297 
298  //for (j = 0; j < (2*N-1)*L*L; ++j)
299  // tot = 0.0;
300  // tot_inverse = 0.0;
301  // error_direct = 0.0;
302  // error_ssht = 0.0;
303  // tot += cabs(flmn_syn_direct[j] - flmn_syn_ssht[j]);
304  // if (real) tot_inverse += cabs(f_real_direct[j] - f_real_ssht[j]);
305  // else tot_inverse += cabs(f_direct[j] - f_ssht[j]);
306 
307  tot = get_max_error(flmn_syn_direct, flmn_syn_ssht, flmn_size);
308 
309  if (real) tot_inverse = get_max_error(f_real_direct, f_real_ssht, flmn_size);
310  else tot_inverse = get_max_error(f_direct, f_ssht, flmn_size);
311 
312  error_direct = get_max_error(flmn_orig, flmn_syn_direct, flmn_size);
313  error_ssht = get_max_error(flmn_orig, flmn_syn_ssht, flmn_size);
314 
315  printf("%.7e", tot);
316  printf(" |-| ");
317  printf("%.7e", tot_inverse);
318  printf(" |-| ");
319  printf("%.7e", error_direct);
320  printf(" |-| ");
321  printf("%.7e", error_ssht);
322  printf(" |-| ");
323  printf("%.7e", inverse_duration_direct);
324  printf(" |-| ");
325  printf("%.7e", forward_duration_direct);
326  printf(" |-| ");
327  printf("%.7e", inverse_duration_ssht);
328  printf(" |-| ");
329  printf("%.7e\n", forward_duration_ssht);
330  //for (j = 0; j < (2*N-1)*L*L; ++j)
331  // if (real) tot = flmn_syn[j];
332  //else tot = cabs(flmn_syn[j]);
333 
334  // if (tot != 0.0) count = count + tot;
335  //printf("%.30e", count);
336 
337  }
338 
339  if (show_arrays == 1)
340  {
341  int c1, c2, c3;
342  c1 = c2 = c3 = 0;
343 
344  printf("\n");
345  printf("----------------------------------------------------------------------------------------------------------------------------------------------------------\n");
346  printf("\n");
347  printf("~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ SSHT flmn real component ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~\n");
348  printf("----------------------------------------------------------------------------------------------------------------------------------------------------------\n");
349  for (int k = 0; k < (2*N-1)*L*L; ++k)
350  {
351  if (k < 10)
352  {
353  printf("|| %d ||",k);
354  }
355  if (k > 9 && k < 100)
356  {
357  printf("|| %d ||",k);
358  }
359  if (k > 99)
360  {
361  printf("|| %d ||",k);
362  }
363 
364  if (creal(flmn_syn_ssht[k]) >= 0.0)
365  {
366  if (creal(flmn_syn_ssht[k]) == 0.0) printf(" --------- ");
367  else printf(" %.3e ", creal(flmn_syn_ssht[k]));
368  }
369  else
370  {
371  printf("%.3e ", creal(flmn_syn_ssht[k]));
372  }
373  ++c1;
374 
375  if (c1 == 7)
376  {
377  printf("\n");
378  c1 = 0;
379  }
380  }
381 
382  printf("\n");
383  printf("----------------------------------------------------------------------------------------------------------------------------------------------------------\n");
384  printf("\n");
385  printf("~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ DIRECT flmn real component ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~\n");
386  printf("----------------------------------------------------------------------------------------------------------------------------------------------------------\n");
387  for (int h = 0; h < (2*N-1)*L*L; ++h)
388  {
389  if (h < 10)
390  {
391  printf("|| %d ||",h);
392  }
393  if (h > 9 && h < 100)
394  {
395  printf("|| %d ||",h);
396  }
397  if (h > 99)
398  {
399  printf("|| %d ||",h);
400  }
401 
402  if (creal(flmn_syn_direct[h]) >= 0.0)
403  {
404  if (creal(flmn_syn_direct[h]) == 0.0) printf(" --------- ");
405  else printf(" %.3e ", creal(flmn_syn_direct[h]));
406  }
407  else
408  {
409  printf("%.3e ", creal(flmn_syn_direct[h]));
410  }
411 
412  ++c2;
413 
414  if (c2 == 7)
415  {
416  printf("\n");
417  c2 = 0;
418  }
419  }
420 
421  printf("\n");
422  printf("----------------------------------------------------------------------------------------------------------------------------------------------------------\n");
423  printf("\n");
424  printf("~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ TRUE flmn real component~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~\n");
425  printf("----------------------------------------------------------------------------------------------------------------------------------------------------------\n");
426  for (int r = 0; r < (2*N-1)*L*L; ++r)
427  {
428  if (r < 10)
429  {
430  printf("|| %d ||",r);
431  }
432  if (r > 9 && r < 100)
433  {
434  printf("|| %d ||",r);
435  }
436  if (r > 99)
437  {
438  printf("|| %d ||",r);
439  }
440 
441  if (creal(flmn_orig[r]) >= 0.0)
442  {
443  if (creal(flmn_orig[r]) == 0.0) printf(" --------- ");
444  else printf(" %.3e ", creal(flmn_orig[r]));
445  }
446  else
447  {
448  printf("%.3e ", creal(flmn_orig[r]));
449  }
450  ++c3;
451 
452  if (c3 == 7)
453  {
454  printf("\n");
455  c3 = 0;
456  }
457  }
458  }
459 
460  if (show_arrays == 2)
461  {
462  int c1, c2, c3;
463  c1 = c2 = c3 = 0;
464 
465  printf("\n");
466  printf("----------------------------------------------------------------------------------------------------------------------------------------------------------\n");
467  printf("\n");
468  printf("~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ SSHT real part of f ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~\n");
469  printf("----------------------------------------------------------------------------------------------------------------------------------------------------------\n");
470  for (int k = 0; k < (2*N-1)*L*L; ++k)
471  {
472  if (k < 10)
473  {
474  printf("|| %d ||",k);
475  }
476  if (k > 9 && k < 100)
477  {
478  printf("|| %d ||",k);
479  }
480  if (k > 99)
481  {
482  printf("|| %d ||",k);
483  }
484 
485  if (creal(f_ssht[k]) >= 0.0)
486  {
487  if (creal(f_ssht[k]) == 0.0) printf(" --------- ");
488  else printf(" %.3e ", creal(f_ssht[k]));
489  }
490  else
491  {
492  printf("%.3e ", creal(f_ssht[k]));
493  }
494  ++c1;
495 
496  if (c1 == 7)
497  {
498  printf("\n");
499  c1 = 0;
500  }
501  }
502 
503  printf("\n");
504  printf("----------------------------------------------------------------------------------------------------------------------------------------------------------\n");
505  printf("\n");
506  printf("~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ DIRECT real part of f ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~\n");
507  printf("----------------------------------------------------------------------------------------------------------------------------------------------------------\n");
508  for (int h = 0; h < (2*N-1)*L*L; ++h)
509  {
510  if (h < 10)
511  {
512  printf("|| %d ||",h);
513  }
514  if (h > 9 && h < 100)
515  {
516  printf("|| %d ||",h);
517  }
518  if (h > 99)
519  {
520  printf("|| %d ||",h);
521  }
522 
523  if (creal(f_direct[h]) >= 0.0)
524  {
525  if (creal(f_direct[h]) == 0.0) printf(" --------- ");
526  else printf(" %.3e ", creal(f_direct[h]));
527  }
528  else
529  {
530  printf("%.3e ", creal(f_direct[h]));
531  }
532 
533  ++c2;
534 
535  if (c2 == 7)
536  {
537  printf("\n");
538  c2 = 0;
539  }
540  }
541 
542  printf("\n");
543  printf("----------------------------------------------------------------------------------------------------------------------------------------------------------\n");
544  printf("\n");
545  printf("~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ Difference between real part of f's ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~\n");
546  printf("----------------------------------------------------------------------------------------------------------------------------------------------------------\n");
547  for (int r = 0; r < (2*N-1)*L*L; ++r)
548  {
549  if (r < 10)
550  {
551  printf("|| %d ||",r);
552  }
553  if (r > 9 && r < 100)
554  {
555  printf("|| %d ||",r);
556  }
557  if (r > 99)
558  {
559  printf("|| %d ||",r);
560  }
561 
562  if (creal(f_direct[r] - f_ssht[r]) >= 0.0)
563  {
564  if ((creal(f_direct[r]) - creal(f_ssht[r])) <= 0.0000000000001) printf(" --------- ");
565  else printf(" %.3e ", (creal(f_direct[r]) - creal(f_ssht[r])));
566  }
567  else
568  {
569  if ((creal(f_direct[r]) - creal(f_ssht[r])) >= -0.0000000000001) printf(" --------- ");
570  else printf("%.3e ", (creal(f_direct[r]) - creal(f_ssht[r])));
571  }
572  ++c3;
573 
574  if (c3 == 7)
575  {
576  printf("\n");
577  c3 = 0;
578  }
579  }
580  }
581 
582 
583 
584  printf("\n");
585  printf("----------------------------------------------------------------------------------------------------------------------------------------------------------\n");
586  }
587  }
588  }
589  }
590  }
591  }
592 
593 
594  free(flmn_orig);
595  free(flmn_syn_ssht);
596  free(flmn_syn_direct);
597  free(f_direct);
598  free(f_ssht);
599  free(f_real_direct);
600  free(f_real_ssht);
601 
602  // =========================================================================
603  // Summarise results
604 
605  printf("================================================================\n");
606  printf("Summary\n\n");
607  printf("Runs = %40d\n", NREPEAT);
608  printf("L0 = %40d\n", L0);
609  printf("L = %40d\n", L);
610  printf("N = %40d\n\n", N);
611 
612  for (steerable = 1; steerable < 2; ++steerable)
613  {
614 
615  printf("Results for %s...\n", steerable_str[steerable]);
616  // real == 0 --> complex signal
617  // real == 1 --> real signal
618 
619  for (real = 0; real < 2; ++real)
620  {
621  //printf(" ...with %s signals...\n", reality_str[real]);
622 
623  // MWSS not yet supported by direct routines
624  for (sampling_scheme = 0; sampling_scheme < 1; ++sampling_scheme)
625  {
626  //printf(" ...using %s sampling...\n", sampling_str[sampling_scheme]);
627 
628  //for (storage_mode = 0; storage_mode < SO3_STORAGE_SIZE; ++storage_mode)
629  for (storage_mode = 0; storage_mode < 1; ++storage_mode)
630  {
631  //printf(" ...with %s...\n", storage_mode_str[storage_mode]);
632 
633  // For real signals, the n_order does not matter, so skip the second option
634  // in that case.
635  for (n_order = 0; n_order < SO3_N_ORDER_SIZE - real; ++n_order)
636  {
637  //printf(" ...using %s...\n", n_order_str[n_order]);
638 
639  for (n_mode = 0; n_mode < SO3_N_MODE_SIZE; ++ n_mode)
640  {
641  //printf(" ...and %s...\n", n_mode_str[n_mode]);
642  //printf(" Minimum time for forward transform: %fs\n", durations_forward[steerable][sampling_scheme][n_order][storage_mode][n_mode][real]);
643  //printf(" Minimum time for inverse transform: %fs\n", durations_inverse[steerable][sampling_scheme][n_order][storage_mode][n_mode][real]);
644  //printf(" Average max errors for round-trip: %e\n", errors[steerable][sampling_scheme][n_order][storage_mode][n_mode][real]);
645  }
646  }
647  printf("\n");
648  }
649  }
650  }
651  }
652 
653 
654  return 0;
655 }
656 
657 double get_max_error(complex double *expected, complex double *actual, int n)
658 {
659  int i;
660  double error, maxError = 0;
661 
662  for (i = 0; i < n; ++i)
663  {
664  error = cabs(expected[i] - actual[i]);
665  maxError = MAX(error, maxError);
666  }
667 
668  return maxError;
669 }
670 
688  complex double *flmn,
689  const so3_parameters_t *parameters,
690  int seed)
691 {
692  int L0, L, N;
693  int i, el, m, n, n_start, n_stop, n_inc, ind;
694 
695  L0 = parameters->L0;
696  L = parameters->L;
697  N = parameters->N;
698 
699  for (i = 0; i < (2*N-1)*L*L; ++i)
700  flmn[i] = 0.0;
701 
702  switch (parameters->n_mode)
703  {
704  case SO3_N_MODE_ALL:
705  case SO3_N_MODE_L:
706  n_start = -N+1;
707  n_stop = N-1;
708  n_inc = 1;
709  break;
710  case SO3_N_MODE_EVEN:
711  n_start = ((N-1) % 2 == 0) ? -N+1 : -N+2;
712  n_stop = ((N-1) % 2 == 0) ? N-1 : N-2;
713  n_inc = 2;
714  break;
715  case SO3_N_MODE_ODD:
716  n_start = ((N-1) % 2 != 0) ? -N+1 : -N+2;
717  n_stop = ((N-1) % 2 != 0) ? N-1 : N-2;
718  n_inc = 2;
719  break;
720  case SO3_N_MODE_MAXIMUM:
721  n_start = -N+1;
722  n_stop = N-1;
723  if (N > 1)
724  n_inc = 2*N - 2;
725  else
726  n_inc = 1;
727  break;
728  default:
729  SO3_ERROR_GENERIC("Invalid n-mode.");
730  }
731 
732  for (n = n_start; n <= n_stop; n += n_inc)
733  {
734  for (el = MAX(L0, abs(n)); el < L; ++el)
735  {
736  if (parameters->n_mode == SO3_N_MODE_L && el != abs(n))
737  break;
738 
739  for (m = -el; m <= el; ++m)
740  {
741  so3_sampling_elmn2ind(&ind, el, m, n, parameters);
742  flmn[ind] = (2.0*ran2_dp(seed) - 1.0) + I * (2.0*ran2_dp(seed) - 1.0);
743  }
744  }
745  }
746 }
747 
768  complex double *flmn,
769  const so3_parameters_t *parameters,
770  int seed)
771 {
772  int L0, L, N;
773  int i, el, m, n, n_start, n_stop, n_inc, ind;
774  double real, imag;
775 
776  L0 = parameters->L0;
777  L = parameters->L;
778  N = parameters->N;
779 
780  for (i = 0; i < (2*N-1)*L*L; ++i)
781  flmn[i] = 0.0;
782 
783  switch (parameters->n_mode)
784  {
785  case SO3_N_MODE_ALL:
786  case SO3_N_MODE_L:
787  n_start = 0;
788  n_stop = N-1;
789  n_inc = 1;
790  break;
791  case SO3_N_MODE_EVEN:
792  n_start = 0;
793  n_stop = ((N-1) % 2 == 0) ? N-1 : N-2;
794  n_inc = 2;
795  break;
796  case SO3_N_MODE_ODD:
797  n_start = 1;
798  n_stop = ((N-1) % 2 != 0) ? N-1 : N-2;
799  n_inc = 2;
800  break;
801  case SO3_N_MODE_MAXIMUM:
802  n_start = N-1;
803  n_stop = N-1;
804  n_inc = 1;
805  break;
806  default:
807  SO3_ERROR_GENERIC("Invalid n-mode.");
808  }
809 
810  for (n = n_start; n <= n_stop; n += n_inc)
811  {
812  if (n == 0)
813  {
814  // Fill fl00 with random real values
815  for (el = L0; el < L; ++el)
816  {
817  if (parameters->n_mode == SO3_N_MODE_L && el != 0)
818  break;
819 
820  so3_sampling_elmn2ind_real(&ind, el, 0, 0, parameters);
821  flmn[ind] = (2.0*ran2_dp(seed) - 1.0);
822  }
823 
824  // Fill fl+-m0 with conjugated random values
825 
826  for (el = L0; el < L; ++el)
827  {
828  if (parameters->n_mode == SO3_N_MODE_L && el != 0)
829  break;
830 
831  for (m = 1; m <= el; ++m)
832  {
833  real = (2.0*ran2_dp(seed) - 1.0);
834  imag = (2.0*ran2_dp(seed) - 1.0);
835  so3_sampling_elmn2ind_real(&ind, el, m, 0, parameters);
836  flmn[ind] = real + imag * I;
837  so3_sampling_elmn2ind_real(&ind, el, -m, 0, parameters);
838  flmn[ind] = real - imag * I;
839  if (m % 2)
840  flmn[ind] = - real + imag * I;
841  else
842  flmn[ind] = real - imag * I;
843  }
844  }
845  }
846  else
847  {
848  for (el = MAX(L0, n); el < L; ++el)
849  {
850  if (parameters->n_mode == SO3_N_MODE_L && el != n)
851  break;
852 
853  for (m = -el; m <= el; ++m)
854  {
855 
856  so3_sampling_elmn2ind_real(&ind, el, m, n, parameters);
857  flmn[ind] = (2.0*ran2_dp(seed) - 1.0) + I * (2.0*ran2_dp(seed) - 1.0);
858  }
859  }
860  }
861  }
862 }
863 
877 double ran2_dp(int idum)
878 {
879  int IM1=2147483563,IM2=2147483399,IMM1=IM1-1,
880  IA1=40014,IA2=40692,IQ1=53668,IQ2=52774,IR1=12211,IR2=3791,
881  NTAB=32,NDIV=1+IMM1/NTAB;
882 
883  double AM=1./IM1,EPS=1.2e-7,RNMX=1.-EPS;
884  int j,k;
885  static int iv[32],iy,idum2 = 123456789;
886  // N.B. in C static variables are initialised to 0 by default.
887 
888  if (idum <= 0) {
889  idum= (-idum>1 ? -idum : 1); // max(-idum,1);
890  idum2=idum;
891  for(j=NTAB+8;j>=1;j--) {
892  k=idum/IQ1;
893  idum=IA1*(idum-k*IQ1)-k*IR1;
894  if (idum < 0) idum=idum+IM1;
895  if (j < NTAB) iv[j-1]=idum;
896  }
897  iy=iv[0];
898  }
899  k=idum/IQ1;
900  idum=IA1*(idum-k*IQ1)-k*IR1;
901  if (idum < 0) idum=idum+IM1;
902  k=idum2/IQ2;
903  idum2=IA2*(idum2-k*IQ2)-k*IR2;
904  if (idum2 < 0) idum2=idum2+IM2;
905  j=1+iy/NDIV;
906  iy=iv[j-1]-idum2;
907  iv[j-1]=idum;
908  if(iy < 1)iy=iy+IMM1;
909  return (AM*iy < RNMX ? AM*iy : RNMX); // min(AM*iy,RNMX);
910 }
so3_core_inverse_direct
void so3_core_inverse_direct(complex double *f, const complex double *flmn, const so3_parameters_t *parameters)
Definition: so3_core.c:951
so3_test_gen_flmn_complex
void so3_test_gen_flmn_complex(complex double *flmn, const so3_parameters_t *parameters, int seed)
Definition: so3_test.c:687
so3_sampling_elmn2ind_real
void so3_sampling_elmn2ind_real(int *ind, int el, int m, int n, const so3_parameters_t *parameters)
Definition: so3_sampling.c:569
so3_core_forward_direct_real
void so3_core_forward_direct_real(complex double *flmn, const double *f, const so3_parameters_t *parameters)
Definition: so3_core.c:2255
MAX
#define MAX(a, b)
Definition: so3_test.c:37
so3_core_inverse_via_ssht
void so3_core_inverse_via_ssht(complex double *f, const complex double *flmn, const so3_parameters_t *parameters)
Definition: so3_core.c:48
so3_core_inverse_direct_real
void so3_core_inverse_direct_real(double *f, const complex double *flmn, const so3_parameters_t *parameters)
Definition: so3_core.c:1850
so3_core_forward_via_ssht_real
void so3_core_forward_via_ssht_real(complex double *flmn, const double *f, const so3_parameters_t *parameters)
Definition: so3_core.c:700
so3_core_forward_via_ssht
void so3_core_forward_via_ssht(complex double *flmn, const complex double *f, const so3_parameters_t *parameters)
Definition: so3_core.c:253
main
int main(int argc, char **argv)
Definition: so3_test.c:44
so3_sampling_flmn_size
int so3_sampling_flmn_size(const so3_parameters_t *parameters)
Definition: so3_sampling.c:314
so3_test_gen_flmn_real
void so3_test_gen_flmn_real(complex double *flmn, const so3_parameters_t *parameters, int seed)
Definition: so3_test.c:767
so3_core_inverse_via_ssht_real
void so3_core_inverse_via_ssht_real(double *f, const complex double *flmn, const so3_parameters_t *parameters)
Definition: so3_core.c:475
so3_core_forward_direct
void so3_core_forward_direct(complex double *flmn, const complex double *f, const so3_parameters_t *parameters)
Definition: so3_core.c:1354
ran2_dp
double ran2_dp(int idum)
Definition: so3_test.c:877
get_max_error
double get_max_error(complex double *expected, complex double *actual, int n)
Definition: so3_test.c:657
NREPEAT
#define NREPEAT
Definition: so3_test.c:35
so3_sampling_elmn2ind
void so3_sampling_elmn2ind(int *ind, int el, int m, int n, const so3_parameters_t *parameters)
Definition: so3_sampling.c:367