My Project
ssht_test.c
Go to the documentation of this file.
1 // SSHT package to perform spin spherical harmonic transforms
2 // Copyright (C) 2011 Jason McEwen
3 // See LICENSE.txt for license details
4 
18 #include <complex.h> // Must be before fftw3.h
19 #include <fftw3.h>
20 #include <omp.h>
21 #include <stdio.h>
22 #include <stdlib.h>
23 #include <time.h>
24 
25 #include "ssht.h"
26 
27 #define NREPEAT 5
28 #define MAX(a, b) ((a) > (b) ? (a) : (b))
29 
30 double ran2_dp(int idum);
31 void ssht_test_gen_flm_complex(complex double* flm, int L, int spin, int seed);
32 void ssht_test_gen_flm_real(complex double* flm, int L, int seed);
33 void ssht_test_gen_lb_flm_complex(complex double* flm, int L_zero, int L, int spin, int seed);
34 void ssht_test_gen_lb_flm_real(complex double* flm, int L0, int L, int seed);
35 
45 int null_test(const complex double *X, int n)
46 {
47  int Y = 0;
48  for(int i = 0; i < n; ++i){ if(cabs(X[i]) != 0.0){ Y = 1; i = n; } }
49  return Y;
50 }
51 
61 int nan_test(const complex double *X, int n)
62 {
63  int Y = 1;
64  for(int i = 0; i < n; ++i){ if(X[i] != X[i]){ Y = 0; i = n; } }
65  return Y;
66 }
67 
68 int main(int argc, char* argv[])
69 {
70 
71  complex double *flm_orig, *flm_syn;
72  complex double *f_mw, *f_mw_ss, *f_gl, *f_dh;
73  double *f_mw_real, *f_mw_ss_real, *f_gl_real, *f_dh_real;
74  complex double *f_mw_lb, *f_mw_lb_ss;
75  double *f_mw_lb_real, *f_mw_lb_ss_real;
76 
77  complex double *f_mw_pole, *f_mw_ss_pole;
78  complex double f_mw_sp, f_mw_ss_np, f_mw_ss_sp;
79  double *f_mw_real_pole, *f_mw_ss_real_pole;
80  double f_mw_real_sp, f_mw_ss_real_sp, f_mw_ss_real_np, phi_sp, phi_np;
81 
82  ssht_dl_method_t dl_method = SSHT_DL_RISBO;
83  int L = 128;
84  int L0 = 32;
85  int spin = 0;
86  int irepeat;
87  int seed = 1;
88  int verbosity = 0;
89  int i;
90  double tmp;
91 
92  clock_t time_start, time_end;
93  double errors_mw[NREPEAT];
94  double errors_mw_lb[NREPEAT];
95  double errors_mw_pole[NREPEAT];
96  double errors_mw_ss[NREPEAT];
97  double errors_mw_lb_ss[NREPEAT];
98  double errors_mw_ss_pole[NREPEAT];
99  double errors_gl[NREPEAT];
100  double errors_dh[NREPEAT];
101  double errors_mw_real[NREPEAT];
102  double errors_mw_lb_real[NREPEAT];
103  double errors_mw_real_pole[NREPEAT];
104  double errors_mw_ss_real[NREPEAT];
105  double errors_mw_lb_ss_real[NREPEAT];
106  double errors_mw_ss_real_pole[NREPEAT];
107  double errors_gl_real[NREPEAT];
108  double errors_dh_real[NREPEAT];
109  double durations_forward_mw[NREPEAT];
110  double durations_inverse_mw[NREPEAT];
111  double durations_forward_mw_lb[NREPEAT];
112  double durations_inverse_mw_lb[NREPEAT];
113  double durations_forward_mw_pole[NREPEAT];
114  double durations_inverse_mw_pole[NREPEAT];
115  double durations_forward_mw_ss[NREPEAT];
116  double durations_inverse_mw_ss[NREPEAT];
117  double durations_forward_mw_lb_ss[NREPEAT];
118  double durations_inverse_mw_lb_ss[NREPEAT];
119  double durations_forward_mw_ss_pole[NREPEAT];
120  double durations_inverse_mw_ss_pole[NREPEAT];
121  double durations_forward_gl[NREPEAT];
122  double durations_inverse_gl[NREPEAT];
123  double durations_forward_dh[NREPEAT];
124  double durations_inverse_dh[NREPEAT];
125  double durations_forward_mw_real[NREPEAT];
126  double durations_inverse_mw_real[NREPEAT];
127  double durations_forward_mw_lb_real[NREPEAT];
128  double durations_inverse_mw_lb_real[NREPEAT];
129  double durations_forward_mw_real_pole[NREPEAT];
130  double durations_inverse_mw_real_pole[NREPEAT];
131  double durations_forward_mw_ss_real[NREPEAT];
132  double durations_inverse_mw_ss_real[NREPEAT];
133  double durations_forward_mw_lb_ss_real[NREPEAT];
134  double durations_inverse_mw_lb_ss_real[NREPEAT];
135  double durations_forward_mw_ss_real_pole[NREPEAT];
136  double durations_inverse_mw_ss_real_pole[NREPEAT];
137  double durations_forward_gl_real[NREPEAT];
138  double durations_inverse_gl_real[NREPEAT];
139  double durations_forward_dh_real[NREPEAT];
140  double durations_inverse_dh_real[NREPEAT];
141 
142  // fftw_init_threads();
143  //int nthreads = omp_get_max_threads();
144  //printf("Using %d threads.\n", nthreads);
145  // fftw_plan_with_nthreads(nthreads);
146 
147  // Parse problem sizes.
148  L = 64;
149  spin = 0;
150  if (argc > 1) {
151  L = atoi(argv[1]);
152  } else {
153  printf("\n");
154  printf("Choosing default L = 128\n");
155  }
156  if (argc > 2) {
157  spin = atoi(argv[2]);
158  } else {
159  printf("\n");
160  printf("Choosing default s = 0\n");
161  }
162  if (argc > 3)
163  L0 = atoi(argv[3]);
164 
165  // Allocate memory.
166  flm_orig = (complex double*)calloc(L * L, sizeof(complex double));
168  flm_syn = (complex double*)calloc(L * L, sizeof(complex double));
170  f_mw = (complex double*)calloc(L * (2 * L - 1), sizeof(complex double));
172  f_mw_lb = (complex double*)calloc(L * (2 * L - 1), sizeof(complex double));
174  f_mw_ss = (complex double*)calloc((L + 1) * (2 * L), sizeof(complex double));
176  f_mw_lb_ss = (complex double*)calloc((L + 1) * (2 * L), sizeof(complex double));
178  f_mw_pole = (complex double*)calloc((L - 1) * (2 * L - 1), sizeof(complex double));
179  SSHT_ERROR_MEM_ALLOC_CHECK(f_mw_pole)
180  f_mw_ss_pole = (complex double*)calloc((L - 1) * (2 * L), sizeof(complex double));
181  SSHT_ERROR_MEM_ALLOC_CHECK(f_mw_ss_pole)
182  f_gl = (complex double*)calloc(L * (2 * L - 1), sizeof(complex double));
184  f_dh = (complex double*)calloc((2 * L) * (2 * L - 1), sizeof(complex double));
186  f_mw_real = (double*)calloc(L * (2 * L - 1), sizeof(double));
187  SSHT_ERROR_MEM_ALLOC_CHECK(f_mw_real)
188  f_mw_lb_real = (double*)calloc(L * (2 * L - 1), sizeof(double));
189  SSHT_ERROR_MEM_ALLOC_CHECK(f_mw_real)
190  f_mw_ss_real = (double*)calloc((L + 1) * (2 * L), sizeof(double));
191  SSHT_ERROR_MEM_ALLOC_CHECK(f_mw_ss_real)
192  f_mw_lb_ss_real = (double*)calloc((L + 1) * (2 * L), sizeof(double));
193  SSHT_ERROR_MEM_ALLOC_CHECK(f_mw_ss_real)
194  f_mw_real_pole = (double*)calloc((L - 1) * (2 * L - 1), sizeof(double));
195  SSHT_ERROR_MEM_ALLOC_CHECK(f_mw_real_pole)
196  f_mw_ss_real_pole = (double*)calloc((L - 1) * (2 * L), sizeof(double));
197  SSHT_ERROR_MEM_ALLOC_CHECK(f_mw_ss_real_pole)
198  f_gl_real = (double*)calloc(L * (2 * L - 1), sizeof(double));
199  SSHT_ERROR_MEM_ALLOC_CHECK(f_gl_real)
200  f_dh_real = (double*)calloc((2 * L) * (2 * L - 1), sizeof(double));
201  SSHT_ERROR_MEM_ALLOC_CHECK(f_dh_real)
202 
203  // Write program name.
204  printf("\n");
205  printf("SSHT test program (C implementation)\n");
206  printf("================================================================\n");
207 
208  // Run algorithm error and timing tests.
209  for (irepeat = 0; irepeat < NREPEAT; irepeat++) {
210 
211  // If spin=0 run tests on algorithms optimised for real spin=0 signal.
212  if (spin == 0) {
213 
214  // =========================================================================
215  // DH real spin=0
216  printf("DH real test no. %d\n", irepeat);
217 
218  ssht_test_gen_flm_real(flm_orig, L, seed);
219  time_start = clock();
220  ssht_core_dh_inverse_sov_real(f_dh_real, flm_orig, L, verbosity);
221  time_end = clock();
222  durations_inverse_dh_real[irepeat] = (time_end - time_start) / (double)CLOCKS_PER_SEC;
223 
224  time_start = clock();
225  ssht_core_dh_forward_sov_real(flm_syn, f_dh_real, L, verbosity);
226  time_end = clock();
227  durations_forward_dh_real[irepeat] = (time_end - time_start) / (double)CLOCKS_PER_SEC;
228 
229  errors_dh_real[irepeat] = 0.0;
230  for (i = 0; i < L * L; i++)
231  errors_dh_real[irepeat] = MAX(cabs(flm_orig[i] - flm_syn[i]), errors_dh_real[irepeat]);
232 
233  printf(" duration_inverse (s) = %40.4f\n",
234  durations_inverse_dh_real[irepeat]);
235  printf(" duration_forward (s) = %40.4f\n",
236  durations_forward_dh_real[irepeat]);
237  printf(" error = %40.5e\n\n",
238  errors_dh_real[irepeat]);
239 
241  printf(" - Null, Nan test result (Input Harmonic Coefficients) : %d, %d\n", null_test(flm_orig, L*L), nan_test(flm_orig, L*L));
242  printf(" - Null, Nan test result (Output Harmonic Coefficients) : %d, %d\n", null_test(flm_syn, L*L), nan_test(flm_syn, L*L));
243 
244  // =========================================================================
245  // GL real spin=0
246  printf("GL real test no. %d\n", irepeat);
247 
248  ssht_test_gen_flm_real(flm_orig, L, seed);
249  time_start = clock();
250  ssht_core_gl_inverse_sov_real(f_gl_real, flm_orig, L, verbosity);
251  time_end = clock();
252  durations_inverse_gl_real[irepeat] = (time_end - time_start) / (double)CLOCKS_PER_SEC;
253 
254  time_start = clock();
255  ssht_core_gl_forward_sov_real(flm_syn, f_gl_real, L, verbosity);
256  time_end = clock();
257  durations_forward_gl_real[irepeat] = (time_end - time_start) / (double)CLOCKS_PER_SEC;
258 
259  errors_gl_real[irepeat] = 0.0;
260  for (i = 0; i < L * L; i++)
261  errors_gl_real[irepeat] = MAX(cabs(flm_orig[i] - flm_syn[i]), errors_gl_real[irepeat]);
262 
263  printf(" duration_inverse (s) = %40.4f\n",
264  durations_inverse_gl_real[irepeat]);
265  printf(" duration_forward (s) = %40.4f\n",
266  durations_forward_gl_real[irepeat]);
267  printf(" error = %40.5e\n\n",
268  errors_gl_real[irepeat]);
269 
271  printf(" - Null, Nan test result (Input Harmonic Coefficients) : %d, %d\n", null_test(flm_orig, L*L), nan_test(flm_orig, L*L));
272  printf(" - Null, Nan test result (Output Harmonic Coefficients) : %d, %d\n", null_test(flm_syn, L*L), nan_test(flm_syn, L*L));
273 
274  // =========================================================================
275  // MW real spin=0
276  printf("MW real test no. %d\n", irepeat);
277 
278  ssht_test_gen_flm_real(flm_orig, L, seed);
279  time_start = clock();
280  ssht_core_mw_inverse_sov_sym_real(f_mw_real, flm_orig, L,
281  dl_method, verbosity);
282  time_end = clock();
283  durations_inverse_mw_real[irepeat] = (time_end - time_start) / (double)CLOCKS_PER_SEC;
284 
285  time_start = clock();
286  ssht_core_mw_forward_sov_conv_sym_real(flm_syn, f_mw_real, L,
287  dl_method, verbosity);
288  time_end = clock();
289  durations_forward_mw_real[irepeat] = (time_end - time_start) / (double)CLOCKS_PER_SEC;
290 
291  errors_mw_real[irepeat] = 0.0;
292  for (i = 0; i < L * L; i++)
293  errors_mw_real[irepeat] = MAX(cabs(flm_orig[i] - flm_syn[i]), errors_mw_real[irepeat]);
294 
295  printf(" duration_inverse (s) = %40.4f\n",
296  durations_inverse_mw_real[irepeat]);
297  printf(" duration_forward (s) = %40.4f\n",
298  durations_forward_mw_real[irepeat]);
299  printf(" error = %40.5e\n\n",
300  errors_mw_real[irepeat]);
301 
303  printf(" - Null, Nan test result (Input Harmonic Coefficients) : %d, %d\n", null_test(flm_orig, L*L), nan_test(flm_orig, L*L));
304  printf(" - Null, Nan test result (Output Harmonic Coefficients) : %d, %d\n", null_test(flm_syn, L*L), nan_test(flm_syn, L*L));
305 
306  // =========================================================================
307  // MW real spin=0 with lower band-limit
308  printf("MW real (lower band-limit) test no. %d\n", irepeat);
309 
310  ssht_test_gen_lb_flm_real(flm_orig, L0, L, seed);
311  time_start = clock();
312  ssht_core_mw_lb_inverse_sov_sym_real(f_mw_lb_real, flm_orig, L0, L,
313  dl_method, verbosity);
314  time_end = clock();
315  durations_inverse_mw_lb_real[irepeat] = (time_end - time_start) / (double)CLOCKS_PER_SEC;
316 
317  time_start = clock();
318  ssht_core_mw_lb_forward_sov_conv_sym_real(flm_syn, f_mw_lb_real, L0, L,
319  dl_method, verbosity);
320  time_end = clock();
321  durations_forward_mw_lb_real[irepeat] = (time_end - time_start) / (double)CLOCKS_PER_SEC;
322 
323  errors_mw_lb_real[irepeat] = 0.0;
324  for (i = 0; i < L * L; i++)
325  errors_mw_lb_real[irepeat] = MAX(cabs(flm_orig[i] - flm_syn[i]), errors_mw_lb_real[irepeat]);
326 
327  printf(" duration_inverse (s) = %40.4f\n",
328  durations_inverse_mw_lb_real[irepeat]);
329  printf(" duration_forward (s) = %40.4f\n",
330  durations_forward_mw_lb_real[irepeat]);
331  printf(" error = %40.5e\n\n",
332  errors_mw_lb_real[irepeat]);
333 
335  printf(" - Null, Nan test result (Input Harmonic Coefficients) : %d, %d\n", null_test(flm_orig, L*L), nan_test(flm_orig, L*L));
336  printf(" - Null, Nan test result (Output Harmonic Coefficients) : %d, %d\n", null_test(flm_syn, L*L), nan_test(flm_syn, L*L));
337 
338  // =========================================================================
339  // MW real spin=0 pole
340  printf("MW real pole test no. %d\n", irepeat);
341 
342  ssht_test_gen_flm_real(flm_orig, L, seed);
343  time_start = clock();
345  &f_mw_real_sp,
346  flm_orig, L,
347  dl_method, verbosity);
348  time_end = clock();
349  durations_inverse_mw_real_pole[irepeat] = (time_end - time_start) / (double)CLOCKS_PER_SEC;
350 
351  time_start = clock();
352  ssht_core_mw_forward_sov_conv_sym_real_pole(flm_syn, f_mw_real_pole,
353  f_mw_real_sp,
354  L,
355  dl_method, verbosity);
356  time_end = clock();
357  durations_forward_mw_real_pole[irepeat] = (time_end - time_start) / (double)CLOCKS_PER_SEC;
358 
359  errors_mw_real_pole[irepeat] = 0.0;
360  for (i = 0; i < L * L; i++)
361  errors_mw_real_pole[irepeat] = MAX(cabs(flm_orig[i] - flm_syn[i]), errors_mw_real_pole[irepeat]);
362 
363  printf(" duration_inverse (s) = %40.4f\n",
364  durations_inverse_mw_real_pole[irepeat]);
365  printf(" duration_forward (s) = %40.4f\n",
366  durations_forward_mw_real_pole[irepeat]);
367  printf(" error = %40.5e\n\n",
368  errors_mw_real_pole[irepeat]);
369 
371  printf(" - Null, Nan test result (Input Harmonic Coefficients) : %d, %d\n", null_test(flm_orig, L*L), nan_test(flm_orig, L*L));
372  printf(" - Null, Nan test result (Output Harmonic Coefficients) : %d, %d\n", null_test(flm_syn, L*L), nan_test(flm_syn, L*L));
373 
374  // =========================================================================
375  // MW SS real spin=0
376  printf("MW SS real test no. %d\n", irepeat);
377 
378  ssht_test_gen_flm_real(flm_orig, L, seed);
379  time_start = clock();
380  ssht_core_mw_inverse_sov_sym_ss_real(f_mw_ss_real, flm_orig, L,
381  dl_method, verbosity);
382  time_end = clock();
383  durations_inverse_mw_ss_real[irepeat] = (time_end - time_start) / (double)CLOCKS_PER_SEC;
384 
385  time_start = clock();
386  ssht_core_mw_forward_sov_conv_sym_ss_real(flm_syn, f_mw_ss_real, L,
387  dl_method, verbosity);
388  time_end = clock();
389  durations_forward_mw_ss_real[irepeat] = (time_end - time_start) / (double)CLOCKS_PER_SEC;
390 
391  errors_mw_ss_real[irepeat] = 0.0;
392  for (i = 0; i < L * L; i++)
393  errors_mw_ss_real[irepeat] = MAX(cabs(flm_orig[i] - flm_syn[i]), errors_mw_ss_real[irepeat]);
394 
395  printf(" duration_inverse (s) = %40.4f\n",
396  durations_inverse_mw_ss_real[irepeat]);
397  printf(" duration_forward (s) = %40.4f\n",
398  durations_forward_mw_ss_real[irepeat]);
399  printf(" error = %40.5e\n\n",
400  errors_mw_ss_real[irepeat]);
401 
403  printf(" - Null, Nan test result (Input Harmonic Coefficients) : %d, %d\n", null_test(flm_orig, L*L), nan_test(flm_orig, L*L));
404  printf(" - Null, Nan test result (Output Harmonic Coefficients) : %d, %d\n", null_test(flm_syn, L*L), nan_test(flm_syn, L*L));
405 
406  // =========================================================================
407  // MW SS real spin=0 with lower band-limit
408  printf("MW SS real (lower band-limit) test no. %d\n", irepeat);
409 
410  ssht_test_gen_lb_flm_real(flm_orig, L0, L, seed);
411  time_start = clock();
412  ssht_core_mw_lb_inverse_sov_sym_ss_real(f_mw_lb_ss_real, flm_orig, L0, L,
413  dl_method, verbosity);
414  time_end = clock();
415  durations_inverse_mw_lb_ss_real[irepeat] = (time_end - time_start) / (double)CLOCKS_PER_SEC;
416 
417  time_start = clock();
418  ssht_core_mw_lb_forward_sov_conv_sym_ss_real(flm_syn, f_mw_lb_ss_real, L0, L,
419  dl_method, verbosity);
420  time_end = clock();
421  durations_forward_mw_lb_ss_real[irepeat] = (time_end - time_start) / (double)CLOCKS_PER_SEC;
422 
423  errors_mw_lb_ss_real[irepeat] = 0.0;
424  for (i = 0; i < L * L; i++)
425  errors_mw_lb_ss_real[irepeat] = MAX(cabs(flm_orig[i] - flm_syn[i]), errors_mw_lb_ss_real[irepeat]);
426 
427  printf(" duration_inverse (s) = %40.4f\n",
428  durations_inverse_mw_lb_ss_real[irepeat]);
429  printf(" duration_forward (s) = %40.4f\n",
430  durations_forward_mw_lb_ss_real[irepeat]);
431  printf(" error = %40.5e\n\n",
432  errors_mw_lb_ss_real[irepeat]);
433 
435  printf(" - Null, Nan test result (Input Harmonic Coefficients) : %d, %d\n", null_test(flm_orig, L*L), nan_test(flm_orig, L*L));
436  printf(" - Null, Nan test result (Output Harmonic Coefficients) : %d, %d\n", null_test(flm_syn, L*L), nan_test(flm_syn, L*L));
437 
438  // =========================================================================
439  // MW SS real spin=0 pole
440  printf("MW SS real pole test no. %d\n", irepeat);
441 
442  ssht_test_gen_flm_real(flm_orig, L, seed);
443  time_start = clock();
445  &f_mw_ss_real_np,
446  &f_mw_ss_real_sp,
447  flm_orig, L,
448  dl_method, verbosity);
449  time_end = clock();
450  durations_inverse_mw_ss_real_pole[irepeat] = (time_end - time_start) / (double)CLOCKS_PER_SEC;
451 
452  time_start = clock();
454  f_mw_ss_real_pole,
455  f_mw_ss_real_np,
456  f_mw_ss_real_sp,
457  L,
458  dl_method, verbosity);
459  time_end = clock();
460  durations_forward_mw_ss_real_pole[irepeat] = (time_end - time_start) / (double)CLOCKS_PER_SEC;
461 
462  errors_mw_ss_real_pole[irepeat] = 0.0;
463  for (i = 0; i < L * L; i++)
464  errors_mw_ss_real_pole[irepeat] = MAX(cabs(flm_orig[i] - flm_syn[i]), errors_mw_ss_real_pole[irepeat]);
465 
466  printf(" duration_inverse (s) = %40.4f\n",
467  durations_inverse_mw_ss_real_pole[irepeat]);
468  printf(" duration_forward (s) = %40.4f\n",
469  durations_forward_mw_ss_real_pole[irepeat]);
470  printf(" error = %40.5e\n\n",
471  errors_mw_ss_real_pole[irepeat]);
472 
474  printf(" - Null, Nan test result (Input Harmonic Coefficients) : %d, %d\n", null_test(flm_orig, L*L), nan_test(flm_orig, L*L));
475  printf(" - Null, Nan test result (Output Harmonic Coefficients) : %d, %d\n", null_test(flm_syn, L*L), nan_test(flm_syn, L*L));
476  }
477 
478  // =========================================================================
479  // DH
480  printf("DH test no. %d\n", irepeat);
481 
482  ssht_test_gen_flm_complex(flm_orig, L, spin, seed);
483  time_start = clock();
484  ssht_core_dh_inverse_sov(f_dh, flm_orig, L, spin, verbosity);
485  time_end = clock();
486  durations_inverse_dh[irepeat] = (time_end - time_start) / (double)CLOCKS_PER_SEC;
487 
488  time_start = clock();
489  ssht_core_dh_forward_sov(flm_syn, f_dh, L, spin, verbosity);
490  time_end = clock();
491  durations_forward_dh[irepeat] = (time_end - time_start) / (double)CLOCKS_PER_SEC;
492 
493  errors_dh[irepeat] = 0.0;
494  for (i = 0; i < L * L; i++)
495  errors_dh[irepeat] = MAX(cabs(flm_orig[i] - flm_syn[i]), errors_dh[irepeat]);
496 
497  printf(" duration_inverse (s) = %40.4f\n",
498  durations_inverse_dh[irepeat]);
499  printf(" duration_forward (s) = %40.4f\n",
500  durations_forward_dh[irepeat]);
501  printf(" error = %40.5e\n\n",
502  errors_dh[irepeat]);
503 
505  printf(" - Null, Nan test result (Input Harmonic Coefficients) : %d, %d\n", null_test(flm_orig, L*L), nan_test(flm_orig, L*L));
506  printf(" - Null, Nan test result (Output Harmonic Coefficients) : %d, %d\n", null_test(flm_syn, L*L), nan_test(flm_syn, L*L));
507 
508  // =========================================================================
509  // GL
510  printf("GL test no. %d\n", irepeat);
511 
512  ssht_test_gen_flm_complex(flm_orig, L, spin, seed);
513  time_start = clock();
514  ssht_core_gl_inverse_sov(f_gl, flm_orig, L, spin, verbosity);
515  time_end = clock();
516  durations_inverse_gl[irepeat] = (time_end - time_start) / (double)CLOCKS_PER_SEC;
517 
518  time_start = clock();
519  ssht_core_gl_forward_sov(flm_syn, f_gl, L, spin, verbosity);
520  time_end = clock();
521  durations_forward_gl[irepeat] = (time_end - time_start) / (double)CLOCKS_PER_SEC;
522 
523  errors_gl[irepeat] = 0.0;
524  for (i = 0; i < L * L; i++)
525  errors_gl[irepeat] = MAX(cabs(flm_orig[i] - flm_syn[i]), errors_gl[irepeat]);
526 
527  printf(" duration_inverse (s) = %40.4f\n",
528  durations_inverse_gl[irepeat]);
529  printf(" duration_forward (s) = %40.4f\n",
530  durations_forward_gl[irepeat]);
531  printf(" error = %40.5e\n\n",
532  errors_gl[irepeat]);
533 
535  printf(" - Null, Nan test result (Input Harmonic Coefficients) : %d, %d\n", null_test(flm_orig, L*L), nan_test(flm_orig, L*L));
536  printf(" - Null, Nan test result (Output Harmonic Coefficients) : %d, %d\n", null_test(flm_syn, L*L), nan_test(flm_syn, L*L));
537 
538  // =========================================================================
539  // MW
540  printf("MW test no. %d\n", irepeat);
541 
542  ssht_test_gen_flm_complex(flm_orig, L, spin, seed);
543  time_start = clock();
544  ssht_core_mw_inverse_sov_sym(f_mw, flm_orig, L, spin, dl_method, verbosity);
545  time_end = clock();
546  durations_inverse_mw[irepeat] = (time_end - time_start) / (double)CLOCKS_PER_SEC;
547 
548  time_start = clock();
549  ssht_core_mw_forward_sov_conv_sym(flm_syn, f_mw, L, spin, dl_method, verbosity);
550  time_end = clock();
551  durations_forward_mw[irepeat] = (time_end - time_start) / (double)CLOCKS_PER_SEC;
552 
553  errors_mw[irepeat] = 0.0;
554  for (i = 0; i < L * L; i++)
555  errors_mw[irepeat] = MAX(cabs(flm_orig[i] - flm_syn[i]), errors_mw[irepeat]);
556 
557  printf(" duration_inverse (s) = %40.4f\n",
558  durations_inverse_mw[irepeat]);
559  printf(" duration_forward (s) = %40.4f\n",
560  durations_forward_mw[irepeat]);
561  printf(" error = %40.5e\n\n",
562  errors_mw[irepeat]);
564  printf(" - Null, Nan test result (Input Harmonic Coefficients) : %d, %d\n", null_test(flm_orig, L*L), nan_test(flm_orig, L*L));
565  printf(" - Null, Nan test result (Output Harmonic Coefficients) : %d, %d\n", null_test(flm_syn, L*L), nan_test(flm_syn, L*L));
566 
567  // =========================================================================
568  // MW with lower band-limit
569  printf("MW lower band-limit test no. %d\n", irepeat);
570 
571  ssht_test_gen_lb_flm_complex(flm_orig, L0, L, spin, seed);
572 
573  time_start = clock();
574  ssht_core_mw_lb_inverse_sov_sym(f_mw_lb, flm_orig, L0, L, spin, dl_method, verbosity);
575  time_end = clock();
576  durations_inverse_mw_lb[irepeat] = (time_end - time_start) / (double)CLOCKS_PER_SEC;
577 
578  time_start = clock();
579  ssht_core_mw_lb_forward_sov_conv_sym(flm_syn, f_mw_lb, L0, L, spin, dl_method, verbosity);
580  time_end = clock();
581  durations_forward_mw_lb[irepeat] = (time_end - time_start) / (double)CLOCKS_PER_SEC;
582 
583  errors_mw_lb[irepeat] = 0.0;
584  for (i = 0; i < L * L; i++)
585  errors_mw_lb[irepeat] = MAX(cabs(flm_orig[i] - flm_syn[i]), errors_mw_lb[irepeat]);
586 
587  printf(" duration_inverse (s) = %40.4f\n",
588  durations_inverse_mw_lb[irepeat]);
589  printf(" duration_forward (s) = %40.4f\n",
590  durations_forward_mw_lb[irepeat]);
591  printf(" error = %40.5e\n\n",
592  errors_mw_lb[irepeat]);
593 
595  printf(" - Null, Nan test result (Input Harmonic Coefficients) : %d, %d\n", null_test(flm_orig, L*L), nan_test(flm_orig, L*L));
596  printf(" - Null, Nan test result (Output Harmonic Coefficients) : %d, %d\n", null_test(flm_syn, L*L), nan_test(flm_syn, L*L));
597 
598  // =========================================================================
599  // MW pole
600  printf("MW pole test no. %d\n", irepeat);
601 
602  ssht_test_gen_flm_complex(flm_orig, L, spin, seed);
603  time_start = clock();
604  ssht_core_mw_inverse_sov_sym_pole(f_mw_pole, &f_mw_sp, &phi_sp,
605  flm_orig, L, spin,
606  dl_method, verbosity);
607  time_end = clock();
608  durations_inverse_mw_pole[irepeat] = (time_end - time_start) / (double)CLOCKS_PER_SEC;
609 
610  time_start = clock();
611  ssht_core_mw_forward_sov_conv_sym_pole(flm_syn, f_mw_pole, f_mw_sp, phi_sp,
612  L, spin,
613  dl_method, verbosity);
614  time_end = clock();
615  durations_forward_mw_pole[irepeat] = (time_end - time_start) / (double)CLOCKS_PER_SEC;
616 
617  errors_mw_pole[irepeat] = 0.0;
618  for (i = 0; i < L * L; i++)
619  errors_mw_pole[irepeat] = MAX(cabs(flm_orig[i] - flm_syn[i]), errors_mw_pole[irepeat]);
620 
621  printf(" duration_inverse (s) = %40.4f\n",
622  durations_inverse_mw_pole[irepeat]);
623  printf(" duration_forward (s) = %40.4f\n",
624  durations_forward_mw_pole[irepeat]);
625  printf(" error = %40.5e\n\n",
626  errors_mw_pole[irepeat]);
627 
629  printf(" - Null, Nan test result (Input Harmonic Coefficients) : %d, %d\n", null_test(flm_orig, L*L), nan_test(flm_orig, L*L));
630  printf(" - Null, Nan test result (Output Harmonic Coefficients) : %d, %d\n", null_test(flm_syn, L*L), nan_test(flm_syn, L*L));
631 
632  // =========================================================================
633  // MW SS
634  printf("MW SS test no. %d\n", irepeat);
635 
636  ssht_test_gen_flm_complex(flm_orig, L, spin, seed);
637  time_start = clock();
638  ssht_core_mw_inverse_sov_sym_ss(f_mw_ss, flm_orig, L, spin,
639  dl_method, verbosity);
640  //ssht_core_mwdirect_inverse_ss(f_mw_ss, flm_orig, L, spin, verbosity);
641  time_end = clock();
642  durations_inverse_mw_ss[irepeat] = (time_end - time_start) / (double)CLOCKS_PER_SEC;
643 
644  time_start = clock();
645  ssht_core_mw_forward_sov_conv_sym_ss(flm_syn, f_mw_ss, L, spin,
646  dl_method, verbosity);
647  time_end = clock();
648  durations_forward_mw_ss[irepeat] = (time_end - time_start) / (double)CLOCKS_PER_SEC;
649 
650  errors_mw_ss[irepeat] = 0.0;
651  for (i = 0; i < L * L; i++)
652  errors_mw_ss[irepeat] = MAX(cabs(flm_orig[i] - flm_syn[i]), errors_mw_ss[irepeat]);
653 
654  printf(" duration_inverse (s) = %40.4f\n",
655  durations_inverse_mw_ss[irepeat]);
656  printf(" duration_forward (s) = %40.4f\n",
657  durations_forward_mw_ss[irepeat]);
658  printf(" error = %40.5e\n\n",
659  errors_mw_ss[irepeat]);
660 
662  printf(" - Null, Nan test result (Input Harmonic Coefficients) : %d, %d\n", null_test(flm_orig, L*L), nan_test(flm_orig, L*L));
663  printf(" - Null, Nan test result (Output Harmonic Coefficients) : %d, %d\n", null_test(flm_syn, L*L), nan_test(flm_syn, L*L));
664 
665  // =========================================================================
666  // MW SS with lower band-limit
667  printf("MW SS (lower band-limit) test no. %d\n", irepeat);
668 
669  ssht_test_gen_lb_flm_complex(flm_orig, L0, L, spin, seed);
670  time_start = clock();
671  ssht_core_mw_lb_inverse_sov_sym_ss(f_mw_lb_ss, flm_orig, L0, L, spin,
672  dl_method, verbosity);
673 
674  time_end = clock();
675  durations_inverse_mw_lb_ss[irepeat] = (time_end - time_start) / (double)CLOCKS_PER_SEC;
676 
677  time_start = clock();
678  ssht_core_mw_lb_forward_sov_conv_sym_ss(flm_syn, f_mw_lb_ss, L0, L, spin,
679  dl_method, verbosity);
680  time_end = clock();
681  durations_forward_mw_lb_ss[irepeat] = (time_end - time_start) / (double)CLOCKS_PER_SEC;
682 
683  errors_mw_lb_ss[irepeat] = 0.0;
684  for (i = 0; i < L * L; i++)
685  errors_mw_lb_ss[irepeat] = MAX(cabs(flm_orig[i] - flm_syn[i]), errors_mw_lb_ss[irepeat]);
686 
687  printf(" duration_inverse (s) = %40.4f\n",
688  durations_inverse_mw_lb_ss[irepeat]);
689  printf(" duration_forward (s) = %40.4f\n",
690  durations_forward_mw_lb_ss[irepeat]);
691  printf(" error = %40.5e\n\n",
692  errors_mw_lb_ss[irepeat]);
693 
695  printf(" - Null, Nan test result (Input Harmonic Coefficients) : %d, %d\n", null_test(flm_orig, L*L), nan_test(flm_orig, L*L));
696  printf(" - Null, Nan test result (Output Harmonic Coefficients) : %d, %d\n", null_test(flm_syn, L*L), nan_test(flm_syn, L*L));
697 
698  // =========================================================================
699  // MW SS pole
700  printf("MW SS pole test no. %d\n", irepeat);
701 
702  ssht_test_gen_flm_complex(flm_orig, L, spin, seed);
703  time_start = clock();
705  &f_mw_ss_np, &phi_np,
706  &f_mw_ss_sp, &phi_sp,
707  flm_orig, L, spin,
708  dl_method, verbosity);
709  time_end = clock();
710  durations_inverse_mw_ss_pole[irepeat] = (time_end - time_start) / (double)CLOCKS_PER_SEC;
711 
712  time_start = clock();
713  ssht_core_mw_forward_sov_conv_sym_ss_pole(flm_syn, f_mw_ss_pole,
714  f_mw_ss_np, phi_np,
715  f_mw_ss_sp, phi_sp,
716  L, spin,
717  dl_method, verbosity);
718  time_end = clock();
719  durations_forward_mw_ss_pole[irepeat] = (time_end - time_start) / (double)CLOCKS_PER_SEC;
720 
721  errors_mw_ss_pole[irepeat] = 0.0;
722  for (i = 0; i < L * L; i++)
723  errors_mw_ss_pole[irepeat] = MAX(cabs(flm_orig[i] - flm_syn[i]), errors_mw_ss_pole[irepeat]);
724 
725  printf(" duration_inverse (s) = %40.4f\n",
726  durations_inverse_mw_ss_pole[irepeat]);
727  printf(" duration_forward (s) = %40.4f\n",
728  durations_forward_mw_ss_pole[irepeat]);
729  printf(" error = %40.5e\n\n",
730  errors_mw_ss_pole[irepeat]);
732  printf(" - Null, Nan test result (Input Harmonic Coefficients) : %d, %d\n", null_test(flm_orig, L*L), nan_test(flm_orig, L*L));
733  printf(" - Null, Nan test result (Output Harmonic Coefficients) : %d, %d\n", null_test(flm_syn, L*L), nan_test(flm_syn, L*L));
734  }
735 
736  // =========================================================================
737  // Summarise results
738 
739  printf("================================================================\n");
740  printf("Summary\n\n");
741  printf("NREPEAT = %40d\n", NREPEAT);
742  printf("L = %40d\n", L);
743  printf("spin = %40d\n\n", spin);
744 
745  if (spin == 0) {
746 
747  printf("DH real\n");
748  tmp = 0.0;
749  for (i = 0; i < NREPEAT; i++)
750  tmp += durations_forward_dh_real[i];
751  printf(" Average forward transform time (s) = %26.4f\n", tmp / (double)NREPEAT);
752  tmp = 0.0;
753  for (i = 0; i < NREPEAT; i++)
754  tmp += durations_inverse_dh_real[i];
755  printf(" Average inverse transform time (s) = %26.4f\n", tmp / (double)NREPEAT);
756  tmp = 0.0;
757  for (i = 0; i < NREPEAT; i++)
758  tmp += errors_dh_real[i];
759  printf(" Average max error = %26.5e\n\n", tmp / (double)NREPEAT);
760 
761  printf("GL real\n");
762  tmp = 0.0;
763  for (i = 0; i < NREPEAT; i++)
764  tmp += durations_forward_gl_real[i];
765  printf(" Average forward transform time (s) = %26.4f\n", tmp / (double)NREPEAT);
766  tmp = 0.0;
767  for (i = 0; i < NREPEAT; i++)
768  tmp += durations_inverse_gl_real[i];
769  printf(" Average inverse transform time (s) = %26.4f\n", tmp / (double)NREPEAT);
770  tmp = 0.0;
771  for (i = 0; i < NREPEAT; i++)
772  tmp += errors_gl_real[i];
773  printf(" Average max error = %26.5e\n\n", tmp / (double)NREPEAT);
774 
775  printf("MW real\n");
776  tmp = 0.0;
777  for (i = 0; i < NREPEAT; i++)
778  tmp += durations_forward_mw_real[i];
779  printf(" Average forward transform time (s) = %26.4f\n", tmp / (double)NREPEAT);
780  tmp = 0.0;
781  for (i = 0; i < NREPEAT; i++)
782  tmp += durations_inverse_mw_real[i];
783  printf(" Average inverse transform time (s) = %26.4f\n", tmp / (double)NREPEAT);
784  tmp = 0.0;
785  for (i = 0; i < NREPEAT; i++)
786  tmp += errors_mw_real[i];
787  printf(" Average max error = %26.5e\n\n", tmp / (double)NREPEAT);
788 
789  printf("MW real with lower band-limit\n");
790  tmp = 0.0;
791  for (i = 0; i < NREPEAT; i++)
792  tmp += durations_forward_mw_lb_real[i];
793  printf(" Average forward transform time (s) = %26.4f\n", tmp / (double)NREPEAT);
794  tmp = 0.0;
795  for (i = 0; i < NREPEAT; i++)
796  tmp += durations_inverse_mw_lb_real[i];
797  printf(" Average inverse transform time (s) = %26.4f\n", tmp / (double)NREPEAT);
798  tmp = 0.0;
799  for (i = 0; i < NREPEAT; i++)
800  tmp += errors_mw_lb_real[i];
801  printf(" Average max error = %26.5e\n\n", tmp / (double)NREPEAT);
802 
803  printf("MW real pole\n");
804  tmp = 0.0;
805  for (i = 0; i < NREPEAT; i++)
806  tmp += durations_forward_mw_real_pole[i];
807  printf(" Average forward transform time (s) = %26.4f\n", tmp / (double)NREPEAT);
808  tmp = 0.0;
809  for (i = 0; i < NREPEAT; i++)
810  tmp += durations_inverse_mw_real_pole[i];
811  printf(" Average inverse transform time (s) = %26.4f\n", tmp / (double)NREPEAT);
812  tmp = 0.0;
813  for (i = 0; i < NREPEAT; i++)
814  tmp += errors_mw_real_pole[i];
815  printf(" Average max error = %26.5e\n\n", tmp / (double)NREPEAT);
816 
817  printf("MW SS real\n");
818  tmp = 0.0;
819  for (i = 0; i < NREPEAT; i++)
820  tmp += durations_forward_mw_lb_ss_real[i];
821  printf(" Average forward transform time (s) = %26.4f\n", tmp / (double)NREPEAT);
822  tmp = 0.0;
823  for (i = 0; i < NREPEAT; i++)
824  tmp += durations_inverse_mw_lb_ss_real[i];
825  printf(" Average inverse transform time (s) = %26.4f\n", tmp / (double)NREPEAT);
826  tmp = 0.0;
827  for (i = 0; i < NREPEAT; i++)
828  tmp += errors_mw_lb_ss_real[i];
829  printf(" Average max error = %26.5e\n\n", tmp / (double)NREPEAT);
830 
831  printf("MW SS real with lower band-limit\n");
832  tmp = 0.0;
833  for (i = 0; i < NREPEAT; i++)
834  tmp += durations_forward_mw_lb_ss_real[i];
835  printf(" Average forward transform time (s) = %26.4f\n", tmp / (double)NREPEAT);
836  tmp = 0.0;
837  for (i = 0; i < NREPEAT; i++)
838  tmp += durations_inverse_mw_lb_ss_real[i];
839  printf(" Average inverse transform time (s) = %26.4f\n", tmp / (double)NREPEAT);
840  tmp = 0.0;
841  for (i = 0; i < NREPEAT; i++)
842  tmp += errors_mw_lb_ss_real[i];
843  printf(" Average max error = %26.5e\n\n", tmp / (double)NREPEAT);
844 
845  printf("MW SS real pole\n");
846  tmp = 0.0;
847  for (i = 0; i < NREPEAT; i++)
848  tmp += durations_forward_mw_ss_real_pole[i];
849  printf(" Average forward transform time (s) = %26.4f\n", tmp / (double)NREPEAT);
850  tmp = 0.0;
851  for (i = 0; i < NREPEAT; i++)
852  tmp += durations_inverse_mw_ss_real_pole[i];
853  printf(" Average inverse transform time (s) = %26.4f\n", tmp / (double)NREPEAT);
854  tmp = 0.0;
855  for (i = 0; i < NREPEAT; i++)
856  tmp += errors_mw_ss_real_pole[i];
857  printf(" Average max error = %26.5e\n\n", tmp / (double)NREPEAT);
858  }
859 
860  printf("DH\n");
861  tmp = 0.0;
862  for (i = 0; i < NREPEAT; i++)
863  tmp += durations_forward_dh[i];
864  printf(" Average forward transform time (s) = %26.4f\n", tmp / (double)NREPEAT);
865  tmp = 0.0;
866  for (i = 0; i < NREPEAT; i++)
867  tmp += durations_inverse_dh[i];
868  printf(" Average inverse transform time (s) = %26.4f\n", tmp / (double)NREPEAT);
869  tmp = 0.0;
870  for (i = 0; i < NREPEAT; i++)
871  tmp += errors_dh[i];
872  printf(" Average max error = %26.5e\n\n", tmp / (double)NREPEAT);
873 
874  printf("GL\n");
875  tmp = 0.0;
876  for (i = 0; i < NREPEAT; i++)
877  tmp += durations_forward_gl[i];
878  printf(" Average forward transform time (s) = %26.4f\n", tmp / (double)NREPEAT);
879  tmp = 0.0;
880  for (i = 0; i < NREPEAT; i++)
881  tmp += durations_inverse_gl[i];
882  printf(" Average inverse transform time (s) = %26.4f\n", tmp / (double)NREPEAT);
883  tmp = 0.0;
884  for (i = 0; i < NREPEAT; i++)
885  tmp += errors_gl[i];
886  printf(" Average max error = %26.5e\n\n", tmp / (double)NREPEAT);
887 
888  printf("MW\n");
889  tmp = 0.0;
890  for (i = 0; i < NREPEAT; i++)
891  tmp += durations_forward_mw[i];
892  printf(" Average forward transform time (s) = %26.4f\n", tmp / (double)NREPEAT);
893  tmp = 0.0;
894  for (i = 0; i < NREPEAT; i++)
895  tmp += durations_inverse_mw[i];
896  printf(" Average inverse transform time (s) = %26.4f\n", tmp / (double)NREPEAT);
897  tmp = 0.0;
898  for (i = 0; i < NREPEAT; i++)
899  tmp += errors_mw[i];
900  printf(" Average max error = %26.5e\n\n", tmp / (double)NREPEAT);
901 
902  printf("MW with lower band-limit\n");
903  tmp = 0.0;
904  for (i = 0; i < NREPEAT; i++)
905  tmp += durations_forward_mw_lb[i];
906  printf(" Average forward transform time (s) = %26.4f\n", tmp / (double)NREPEAT);
907  tmp = 0.0;
908  for (i = 0; i < NREPEAT; i++)
909  tmp += durations_inverse_mw_lb[i];
910  printf(" Average inverse transform time (s) = %26.4f\n", tmp / (double)NREPEAT);
911  tmp = 0.0;
912  for (i = 0; i < NREPEAT; i++)
913  tmp += errors_mw_lb[i];
914  printf(" Average max error = %26.5e\n\n", tmp / (double)NREPEAT);
915 
916  printf("MW pole\n");
917  tmp = 0.0;
918  for (i = 0; i < NREPEAT; i++)
919  tmp += durations_forward_mw_pole[i];
920  printf(" Average forward transform time (s) = %26.4f\n", tmp / (double)NREPEAT);
921  tmp = 0.0;
922  for (i = 0; i < NREPEAT; i++)
923  tmp += durations_inverse_mw_pole[i];
924  printf(" Average inverse transform time (s) = %26.4f\n", tmp / (double)NREPEAT);
925  tmp = 0.0;
926  for (i = 0; i < NREPEAT; i++)
927  tmp += errors_mw_pole[i];
928  printf(" Average max error = %26.5e\n\n", tmp / (double)NREPEAT);
929 
930  printf("MW SS\n");
931  tmp = 0.0;
932  for (i = 0; i < NREPEAT; i++)
933  tmp += durations_forward_mw_ss[i];
934  printf(" Average forward transform time (s) = %26.4f\n", tmp / (double)NREPEAT);
935  tmp = 0.0;
936  for (i = 0; i < NREPEAT; i++)
937  tmp += durations_inverse_mw_ss[i];
938  printf(" Average inverse transform time (s) = %26.4f\n", tmp / (double)NREPEAT);
939  tmp = 0.0;
940  for (i = 0; i < NREPEAT; i++)
941  tmp += errors_mw_ss[i];
942  printf(" Average max error = %26.5e\n\n", tmp / (double)NREPEAT);
943 
944  printf("MW SS with lower band-limit\n");
945  tmp = 0.0;
946  for (i = 0; i < NREPEAT; i++)
947  tmp += durations_forward_mw_lb_ss[i];
948  printf(" Average forward transform time (s) = %26.4f\n", tmp / (double)NREPEAT);
949  tmp = 0.0;
950  for (i = 0; i < NREPEAT; i++)
951  tmp += durations_inverse_mw_lb_ss[i];
952  printf(" Average inverse transform time (s) = %26.4f\n", tmp / (double)NREPEAT);
953  tmp = 0.0;
954  for (i = 0; i < NREPEAT; i++)
955  tmp += errors_mw_lb_ss[i];
956  printf(" Average max error = %26.5e\n\n", tmp / (double)NREPEAT);
957 
958  printf("MW SS pole\n");
959  tmp = 0.0;
960  for (i = 0; i < NREPEAT; i++)
961  tmp += durations_forward_mw_ss_pole[i];
962  printf(" Average forward transform time (s) = %26.4f\n", tmp / (double)NREPEAT);
963  tmp = 0.0;
964  for (i = 0; i < NREPEAT; i++)
965  tmp += durations_inverse_mw_ss_pole[i];
966  printf(" Average inverse transform time (s) = %26.4f\n", tmp / (double)NREPEAT);
967  tmp = 0.0;
968  for (i = 0; i < NREPEAT; i++)
969  tmp += errors_mw_ss_pole[i];
970  printf(" Average max error = %26.5e\n\n", tmp / (double)NREPEAT);
971 
972  // Free memory.
973  free(flm_orig);
974  free(flm_syn);
975  free(f_mw);
976  free(f_mw_lb);
977  free(f_mw_ss);
978  free(f_mw_lb_ss);
979  free(f_mw_real);
980  free(f_mw_lb_real);
981  free(f_mw_ss_real);
982  free(f_mw_lb_ss_real);
983  free(f_mw_pole);
984  free(f_mw_ss_pole);
985  free(f_mw_real_pole);
986  free(f_mw_ss_real_pole);
987  free(f_gl);
988  free(f_gl_real);
989  free(f_dh);
990  free(f_dh_real);
991 
992  return 0;
993 }
994 
1006 void ssht_test_gen_flm_real(complex double* flm, int L, int seed)
1007 {
1008 
1009  int el, m, msign, i, i_op;
1010 
1011  for (el = 0; el < L; el++) {
1012  m = 0;
1013  ssht_sampling_elm2ind(&i, el, m);
1014  flm[i] = (2.0 * ran2_dp(seed) - 1.0);
1015  for (m = 1; m <= el; m++) {
1016  ssht_sampling_elm2ind(&i, el, m);
1017  flm[i] = (2.0 * ran2_dp(seed) - 1.0) + I * (2.0 * ran2_dp(seed) - 1.0);
1018  ssht_sampling_elm2ind(&i_op, el, -m);
1019  msign = m & 1;
1020  msign = 1 - msign - msign; // (-1)^m
1021  flm[i_op] = msign * conj(flm[i]);
1022  }
1023  }
1024 }
1025 
1038 void ssht_test_gen_lb_flm_real(complex double* flm, int L0, int L, int seed)
1039 {
1040 
1041  int el, m, msign, i, i_op;
1042 
1043  for (el = 0; el < L0; el++) {
1044  m = 0;
1045  ssht_sampling_elm2ind(&i, el, m);
1046  flm[i] = 0.0;
1047  for (m = 1; m <= el; m++) {
1048  ssht_sampling_elm2ind(&i, el, m);
1049  flm[i] = 0.0;
1050  ssht_sampling_elm2ind(&i, el, -m);
1051  flm[i] = 0.0;
1052  }
1053  }
1054 
1055  for (el = L0; el < L; el++) {
1056  m = 0;
1057  ssht_sampling_elm2ind(&i, el, m);
1058  flm[i] = (2.0 * ran2_dp(seed) - 1.0);
1059  for (m = 1; m <= el; m++) {
1060  ssht_sampling_elm2ind(&i, el, m);
1061  flm[i] = (2.0 * ran2_dp(seed) - 1.0) + I * (2.0 * ran2_dp(seed) - 1.0);
1062  ssht_sampling_elm2ind(&i_op, el, -m);
1063  msign = m & 1;
1064  msign = 1 - msign - msign; // (-1)^m
1065  flm[i_op] = msign * conj(flm[i]);
1066  }
1067  }
1068 }
1069 
1082 void ssht_test_gen_flm_complex(complex double* flm, int L, int spin, int seed)
1083 {
1084 
1085  int i, i_lo;
1086 
1087  ssht_sampling_elm2ind(&i_lo, abs(spin), 0);
1088  for (i = i_lo; i < L * L; i++)
1089  flm[i] = (2.0 * ran2_dp(seed) - 1.0) + I * (2.0 * ran2_dp(seed) - 1.0);
1090 }
1091 
1105 void ssht_test_gen_lb_flm_complex(complex double* flm, int L0, int L, int spin, int seed)
1106 {
1107 
1108  int i, i_lo;
1109 
1110  ssht_sampling_elm2ind(&i_lo, abs(spin), 0);
1111  for (i = 0; i < MAX(i_lo, L0 * L0); i++){ flm[i] = 0.0; }
1112  for (i = MAX(i_lo, L0 * L0); i < L * L; i++){ flm[i] = (2.0 * ran2_dp(seed) - 1.0) + I * (2.0 * ran2_dp(seed) - 1.0);}
1113 
1114 }
1115 
1128 double ran2_dp(int idum)
1129 {
1130 
1131  int IM1 = 2147483563, IM2 = 2147483399, IMM1 = IM1 - 1,
1132  IA1 = 40014, IA2 = 40692, IQ1 = 53668, IQ2 = 52774, IR1 = 12211, IR2 = 3791,
1133  NTAB = 32, NDIV = 1 + IMM1 / NTAB;
1134 
1135  double AM = 1. / IM1, EPS = 1.2e-7, RNMX = 1. - EPS;
1136  int j, k;
1137  static int iv[32], iy, idum2 = 123456789;
1138  // N.B. in C static variables are initialised to 0 by default.
1139 
1140  if (idum <= 0) {
1141  idum = (-idum > 1 ? -idum : 1); // max(-idum,1);
1142  idum2 = idum;
1143  for (j = NTAB + 8; j >= 1; j--) {
1144  k = idum / IQ1;
1145  idum = IA1 * (idum - k * IQ1) - k * IR1;
1146  if (idum < 0)
1147  idum = idum + IM1;
1148  if (j < NTAB)
1149  iv[j - 1] = idum;
1150  }
1151  iy = iv[0];
1152  }
1153  k = idum / IQ1;
1154  idum = IA1 * (idum - k * IQ1) - k * IR1;
1155  if (idum < 0)
1156  idum = idum + IM1;
1157  k = idum2 / IQ2;
1158  idum2 = IA2 * (idum2 - k * IQ2) - k * IR2;
1159  if (idum2 < 0)
1160  idum2 = idum2 + IM2;
1161  j = 1 + iy / NDIV;
1162  iy = iv[j - 1] - idum2;
1163  iv[j - 1] = idum;
1164  if (iy < 1)
1165  iy = iy + IMM1;
1166  return (AM * iy < RNMX ? AM * iy : RNMX); // min(AM*iy,RNMX);
1167 }
ssht_core_mw_forward_sov_conv_sym
void ssht_core_mw_forward_sov_conv_sym(SSHT_COMPLEX(double) *flm, const SSHT_COMPLEX(double) *f, int L, int spin, ssht_dl_method_t dl_method, int verbosity)
Definition: ssht_core.c:810
ssht_core_dh_inverse_sov
void ssht_core_dh_inverse_sov(SSHT_COMPLEX(double) *f, const SSHT_COMPLEX(double) *flm, int L, int spin, int verbosity)
Definition: ssht_core.c:3996
ssht_core_gl_inverse_sov_real
void ssht_core_gl_inverse_sov_real(double *f, const SSHT_COMPLEX(double) *flm, int L, int verbosity)
Definition: ssht_core.c:3552
ssht_core_gl_forward_sov
void ssht_core_gl_forward_sov(SSHT_COMPLEX(double) *flm, const SSHT_COMPLEX(double) *f, int L, int spin, int verbosity)
Definition: ssht_core.c:3685
ssht_test_gen_flm_real
void ssht_test_gen_flm_real(complex double *flm, int L, int seed)
Definition: ssht_test.c:1006
ssht_core_mw_lb_forward_sov_conv_sym_real
void ssht_core_mw_lb_forward_sov_conv_sym_real(SSHT_COMPLEX(double) *flm, const double *f, int L0, int L, ssht_dl_method_t dl_method, int verbosity)
Definition: ssht_core.c:1230
ssht_core_mw_lb_forward_sov_conv_sym
void ssht_core_mw_lb_forward_sov_conv_sym(SSHT_COMPLEX(double) *flm, const SSHT_COMPLEX(double) *f, int L0, int L, int spin, ssht_dl_method_t dl_method, int verbosity)
Definition: ssht_core.c:837
ssht_core_mw_forward_sov_conv_sym_real_pole
void ssht_core_mw_forward_sov_conv_sym_real_pole(SSHT_COMPLEX(double) *flm, const double *f, double f_sp, int L, ssht_dl_method_t dl_method, int verbosity)
Definition: ssht_core.c:1750
ssht.h
null_test
int null_test(const complex double *X, int n)
Definition: ssht_test.c:45
ssht_core_gl_forward_sov_real
void ssht_core_gl_forward_sov_real(SSHT_COMPLEX(double) *flm, const double *f, int L, int verbosity)
Definition: ssht_core.c:3832
ssht_core_mw_inverse_sov_sym_ss_pole
void ssht_core_mw_inverse_sov_sym_ss_pole(SSHT_COMPLEX(double) *f, SSHT_COMPLEX(double) *f_np, double *phi_np, SSHT_COMPLEX(double) *f_sp, double *phi_sp, const SSHT_COMPLEX(double) *flm, int L, int spin, ssht_dl_method_t dl_method, int verbosity)
Definition: ssht_core.c:3209
ssht_core_mw_inverse_sov_sym_real
void ssht_core_mw_inverse_sov_sym_real(double *f, const SSHT_COMPLEX(double) *flm, int L, ssht_dl_method_t dl_method, int verbosity)
Definition: ssht_core.c:322
ssht_core_mw_inverse_sov_sym_ss_real
void ssht_core_mw_inverse_sov_sym_ss_real(double *f, const SSHT_COMPLEX(double) *flm, int L, ssht_dl_method_t dl_method, int verbosity)
Definition: ssht_core.c:2068
ran2_dp
double ran2_dp(int idum)
Definition: ssht_test.c:1128
ssht_core_gl_inverse_sov
void ssht_core_gl_inverse_sov(SSHT_COMPLEX(double) *f, const SSHT_COMPLEX(double) *flm, int L, int spin, int verbosity)
Definition: ssht_core.c:3421
ssht_test_gen_lb_flm_real
void ssht_test_gen_lb_flm_real(complex double *flm, int L0, int L, int seed)
Definition: ssht_test.c:1038
ssht_test_gen_lb_flm_complex
void ssht_test_gen_lb_flm_complex(complex double *flm, int L_zero, int L, int spin, int seed)
Definition: ssht_test.c:1105
ssht_core_mw_forward_sov_conv_sym_pole
void ssht_core_mw_forward_sov_conv_sym_pole(SSHT_COMPLEX(double) *flm, const SSHT_COMPLEX(double) *f, SSHT_COMPLEX(double) f_sp, double phi_sp, int L, int spin, ssht_dl_method_t dl_method, int verbosity)
Definition: ssht_core.c:1702
main
int main(int argc, char *argv[])
Definition: ssht_test.c:68
ssht_core_mw_inverse_sov_sym
void ssht_core_mw_inverse_sov_sym(SSHT_COMPLEX(double) *f, const SSHT_COMPLEX(double) *flm, int L, int spin, ssht_dl_method_t dl_method, int verbosity)
Definition: ssht_core.c:51
ssht_core_mw_lb_forward_sov_conv_sym_ss
void ssht_core_mw_lb_forward_sov_conv_sym_ss(SSHT_COMPLEX(double) *flm, const SSHT_COMPLEX(double) *f, int L0, int L, int spin, ssht_dl_method_t dl_method, int verbosity)
Definition: ssht_core.c:2452
NREPEAT
#define NREPEAT
Definition: ssht_test.c:27
ssht_core_dh_forward_sov_real
void ssht_core_dh_forward_sov_real(SSHT_COMPLEX(double) *flm, const double *f, int L, int verbosity)
Definition: ssht_core.c:4373
MAX
#define MAX(a, b)
Definition: ssht_test.c:28
ssht_core_mw_forward_sov_conv_sym_ss
void ssht_core_mw_forward_sov_conv_sym_ss(SSHT_COMPLEX(double) *flm, const SSHT_COMPLEX(double) *f, int L, int spin, ssht_dl_method_t dl_method, int verbosity)
Definition: ssht_core.c:2425
ssht_core_mw_inverse_sov_sym_pole
void ssht_core_mw_inverse_sov_sym_pole(SSHT_COMPLEX(double) *f, SSHT_COMPLEX(double) *f_sp, double *phi_sp, const SSHT_COMPLEX(double) *flm, int L, int spin, ssht_dl_method_t dl_method, int verbosity)
Definition: ssht_core.c:1608
SSHT_ERROR_MEM_ALLOC_CHECK
#define SSHT_ERROR_MEM_ALLOC_CHECK(pointer)
Definition: ssht_error.h:28
ssht_core_dh_inverse_sov_real
void ssht_core_dh_inverse_sov_real(double *f, const SSHT_COMPLEX(double) *flm, int L, int verbosity)
Definition: ssht_core.c:4115
ssht_core_mw_forward_sov_conv_sym_ss_real
void ssht_core_mw_forward_sov_conv_sym_ss_real(SSHT_COMPLEX(double) *flm, const double *f, int L, ssht_dl_method_t dl_method, int verbosity)
Definition: ssht_core.c:2820
ssht_core_mw_inverse_sov_sym_real_pole
void ssht_core_mw_inverse_sov_sym_real_pole(double *f, double *f_sp, const SSHT_COMPLEX(double) *flm, int L, ssht_dl_method_t dl_method, int verbosity)
Definition: ssht_core.c:1654
ssht_core_mw_forward_sov_conv_sym_ss_pole
void ssht_core_mw_forward_sov_conv_sym_ss_pole(SSHT_COMPLEX(double) *flm, const SSHT_COMPLEX(double) *f, SSHT_COMPLEX(double) f_np, double phi_np, SSHT_COMPLEX(double) f_sp, double phi_sp, int L, int spin, ssht_dl_method_t dl_method, int verbosity)
Definition: ssht_core.c:3316
nan_test
int nan_test(const complex double *X, int n)
Definition: ssht_test.c:61
ssht_core_mw_inverse_sov_sym_ss_real_pole
void ssht_core_mw_inverse_sov_sym_ss_real_pole(double *f, double *f_np, double *f_sp, const SSHT_COMPLEX(double) *flm, int L, ssht_dl_method_t dl_method, int verbosity)
Definition: ssht_core.c:3261
ssht_core_mw_lb_inverse_sov_sym_real
void ssht_core_mw_lb_inverse_sov_sym_real(double *f, const SSHT_COMPLEX(double) *flm, int L0, int L, ssht_dl_method_t dl_method, int verbosity)
Definition: ssht_core.c:348
ssht_core_mw_lb_inverse_sov_sym
void ssht_core_mw_lb_inverse_sov_sym(SSHT_COMPLEX(double) *f, const SSHT_COMPLEX(double) *flm, int L0, int L, int spin, ssht_dl_method_t dl_method, int verbosity)
Definition: ssht_core.c:77
ssht_core_mw_inverse_sov_sym_ss
void ssht_core_mw_inverse_sov_sym_ss(SSHT_COMPLEX(double) *f, const SSHT_COMPLEX(double) *flm, int L, int spin, ssht_dl_method_t dl_method, int verbosity)
Definition: ssht_core.c:1798
ssht_core_dh_forward_sov
void ssht_core_dh_forward_sov(SSHT_COMPLEX(double) *flm, const SSHT_COMPLEX(double) *f, int L, int spin, int verbosity)
Definition: ssht_core.c:4236
ssht_core_mw_lb_inverse_sov_sym_ss
void ssht_core_mw_lb_inverse_sov_sym_ss(SSHT_COMPLEX(double) *f, const SSHT_COMPLEX(double) *flm, int L0, int L, int spin, ssht_dl_method_t dl_method, int verbosity)
Definition: ssht_core.c:1824
ssht_core_mw_forward_sov_conv_sym_real
void ssht_core_mw_forward_sov_conv_sym_real(SSHT_COMPLEX(double) *flm, const double *f, int L, ssht_dl_method_t dl_method, int verbosity)
Definition: ssht_core.c:1203
ssht_core_mw_forward_sov_conv_sym_ss_real_pole
void ssht_core_mw_forward_sov_conv_sym_ss_real_pole(SSHT_COMPLEX(double) *flm, const double *f, double f_np, double f_sp, int L, ssht_dl_method_t dl_method, int verbosity)
Definition: ssht_core.c:3369
SSHT_DL_RISBO
@ SSHT_DL_RISBO
Definition: ssht_dl.h:21
ssht_dl_method_t
ssht_dl_method_t
Definition: ssht_dl.h:21
ssht_test_gen_flm_complex
void ssht_test_gen_flm_complex(complex double *flm, int L, int spin, int seed)
Definition: ssht_test.c:1082
ssht_core_mw_lb_forward_sov_conv_sym_ss_real
void ssht_core_mw_lb_forward_sov_conv_sym_ss_real(SSHT_COMPLEX(double) *flm, const double *f, int L0, int L, ssht_dl_method_t dl_method, int verbosity)
Definition: ssht_core.c:2847
ssht_core_mw_lb_inverse_sov_sym_ss_real
void ssht_core_mw_lb_inverse_sov_sym_ss_real(double *f, const SSHT_COMPLEX(double) *flm, int L0, int L, ssht_dl_method_t dl_method, int verbosity)
Definition: ssht_core.c:2095