My Project
ssht_core.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 
5 
13 #include <stdio.h>
14 #include <stdlib.h>
15 #include <string.h>
16 #include <math.h>
17 
18 #include <complex.h>
19 
20 #include <fftw3.h>
21 
22 #include "ssht_types.h"
23 #include "ssht_error.h"
24 #include "ssht_dl.h"
25 #include "ssht_sampling.h"
26 #include "ssht_core.h"
27 
28 #define MAX(a,b) ((a) > (b) ? (a) : (b))
29 
30 
31 //============================================================================
32 // MW algorithms
33 //============================================================================
34 
35 
51 void ssht_core_mw_inverse_sov_sym(SSHT_COMPLEX(double) *f, const SSHT_COMPLEX(double) *flm,
52  int L, int spin,
53  ssht_dl_method_t dl_method,
54  int verbosity) {
56  0, L, spin,
57  dl_method,
58  verbosity);
59 }
60 
77 void ssht_core_mw_lb_inverse_sov_sym(SSHT_COMPLEX(double) *f, const SSHT_COMPLEX(double) *flm,
78  int L0, int L, int spin,
79  ssht_dl_method_t dl_method,
80  int verbosity) {
81 
82  int el, m, mm;
83  //int t, p;
84  int eltmp;
85  double *sqrt_tbl, *signs;
86  int el2pel, m_offset;
87  double ssign, elfactor;
88  SSHT_COMPLEX(double) mmfactor;
89  double *dl;
90  double *dl8 = NULL;
91  int dl_offset, dl_stride;
92  SSHT_COMPLEX(double) *exps, *m_factors;
93  int exps_offset;
94  double elmmsign, elssign;
95  int spinneg;
96  SSHT_COMPLEX(double) *Fmm, *fext;
97  int Fmm_offset, Fmm_stride, fext_stride;
98  fftw_plan plan;
99 
100  // Allocate memory.
101  sqrt_tbl = (double*)calloc(2*(L-1)+2, sizeof(double));
103  signs = (double*)calloc(L+1, sizeof(double));
105  exps = (SSHT_COMPLEX(double)*)calloc(2*L-1, sizeof(SSHT_COMPLEX(double)));
107  m_factors = calloc(2*L-1, sizeof *m_factors);
108  SSHT_ERROR_MEM_ALLOC_CHECK(m_factors)
109 
110  // Perform precomputations.
111  for (el=0; el<=2*(L-1)+1; el++)
112  sqrt_tbl[el] = sqrt((double)el);
113  for (m=0; m<=L-1; m=m+2) {
114  signs[m] = 1.0;
115  signs[m+1] = -1.0;
116  }
117  ssign = signs[abs(spin)];
118  spinneg = spin <= 0 ? spin : -spin;
119  exps_offset = L-1;
120  for (m=-(L-1); m<=L-1; m++)
121  exps[m + exps_offset] = cexp(-I*SSHT_PION2*(m+spin));
122 
123  // Print messages depending on verbosity level.
124  if (verbosity > 0) {
125  printf("%s%s\n", SSHT_PROMPT,
126  "Computing inverse transform using MW sampling with ");
127  printf("%s%s%d%s%d%s\n", SSHT_PROMPT, "parameters (L,spin,reality) = (",
128  L, ",", spin, ", FALSE)");
129  if (verbosity > 1)
130  printf("%s%s\n", SSHT_PROMPT,
131  "Using routine ssht_core_mw_inverse_sov_sym...");
132  }
133 
134  // Compute Fmm.
135  Fmm = (SSHT_COMPLEX(double)*)calloc((2*L-1)*(2*L-1), sizeof(SSHT_COMPLEX(double)));
137  Fmm_offset = L-1;
138  Fmm_stride = 2*L-1;
141  if (dl_method == SSHT_DL_RISBO) {
144  }
145  dl_offset = ssht_dl_get_offset(L, SSHT_DL_QUARTER);
146  dl_stride = ssht_dl_get_stride(L, SSHT_DL_QUARTER);
147  m_offset = L-1;
148  for (el=MAX(L0, abs(spin)); el<=L-1; el++) {
149 
150  // Compute Wigner plane.
151  switch (dl_method) {
152 
153  case SSHT_DL_RISBO:
154  if (el!=0 && el==MAX(L0, abs(spin))) {
155  for(eltmp=0; eltmp<=MAX(L0, abs(spin)); eltmp++)
158  eltmp, sqrt_tbl, signs);
160  dl8, L,
163  el,
164  signs);
165  }
166  else {
169  el, sqrt_tbl, signs);
171  dl8, L,
174  el,
175  signs);
176  }
177  break;
178 
179  case SSHT_DL_TRAPANI:
180  if (el!=0 && el==MAX(L0, abs(spin))) {
181  for(eltmp=0; eltmp<=MAX(L0, abs(spin)); eltmp++)
184  eltmp, sqrt_tbl);
187  el, signs);
188  }
189  else {
192  el, sqrt_tbl);
195  el, signs);
196  }
197  break;
198 
199  default:
200  SSHT_ERROR_GENERIC("Invalid dl method")
201  }
202 
203  // Compute Fmm.
204  elfactor = ssign * sqrt((double)(2.0*el+1.0)/(4.0*SSHT_PI));
205  el2pel = el*el + el;
206  for (m=-el; m<=el; m++)
207  m_factors[m + m_offset] = flm[el2pel + m] * exps[m + exps_offset];
208 
209  for (mm=0; mm<=el; mm++) {
210  double mm_factor;
211  int mm_offset = mm*dl_stride;
212  elmmsign = signs[el] * signs[mm];
213  elssign = spin <= 0 ? 1.0 : elmmsign;
214  mm_factor = elfactor * elssign * dl[mm_offset - spinneg + dl_offset];
215  for (m=-el; m<=-1; m++) {
216  Fmm[(mm + Fmm_offset)*Fmm_stride + m + Fmm_offset] +=
217  mm_factor
218  * m_factors[m + m_offset]
219  * elmmsign * dl[mm_offset - m + dl_offset];
220  }
221  for (m=0; m<=el; m++) {
222  Fmm[(mm + Fmm_offset)*Fmm_stride + m + Fmm_offset] +=
223  mm_factor
224  * m_factors[m + m_offset]
225  * dl[mm_offset + m + dl_offset];
226  }
227 
228  }
229 
230  }
231 
232  // Free dl memory.
233  free(dl);
234  if (dl_method == SSHT_DL_RISBO)
235  free(dl8);
236 
237  // Use symmetry to compute Fmm for negative mm.
238  for (mm=-(L-1); mm<=-1; mm++)
239  for (m=-(L-1); m<=L-1; m++)
240  Fmm[(mm + Fmm_offset)*Fmm_stride + m + Fmm_offset] =
241  signs[abs(m)] * ssign
242  * Fmm[(-mm + Fmm_offset)*Fmm_stride + m + Fmm_offset];
243 
244  // Apply phase modulation to account for sampling offset.
245  for (mm=-(L-1); mm<=L-1; mm++) {
246  mmfactor = cexp(I*mm*SSHT_PI/(2.0*L-1.0));
247  for (m=-(L-1); m<=L-1; m++)
248  Fmm[(mm + Fmm_offset)*Fmm_stride + m + Fmm_offset] *=
249  mmfactor;
250  }
251 
252  // Allocate space for function values.
253  fext = (SSHT_COMPLEX(double)*)calloc((2*L-1)*(2*L-1), sizeof(SSHT_COMPLEX(double)));
255  fext_stride = 2*L-1;
256 
257  // Apply spatial shift.
258  for (mm=0; mm<=L-1; mm++)
259  for (m=0; m<=L-1; m++)
260  fext[mm*fext_stride + m] =
261  Fmm[(mm + Fmm_offset)*Fmm_stride + m + Fmm_offset];
262  for (mm=0; mm<=L-1; mm++)
263  for (m=-(L-1); m<=-1; m++)
264  fext[mm*fext_stride + (m+2*L-1)] =
265  Fmm[(mm + Fmm_offset)*Fmm_stride + m + Fmm_offset];
266  for (mm=-(L-1); mm<=-1; mm++)
267  for (m=0; m<=L-1; m++)
268  fext[(mm + 2*L-1)*fext_stride + m] =
269  Fmm[(mm + Fmm_offset)*Fmm_stride + m + Fmm_offset];
270  for (mm=-(L-1); mm<=-1; mm++)
271  for (m=-(L-1); m<=-1; m++)
272  fext[(mm+2*L-1)*fext_stride + m + 2*L-1] =
273  Fmm[(mm + Fmm_offset)*Fmm_stride + m + Fmm_offset];
274 
275  // Perform 2D FFT.
276  plan = fftw_plan_dft_2d(2*L-1, 2*L-1, Fmm, Fmm,
277  FFTW_BACKWARD, FFTW_ESTIMATE);
278  fftw_execute_dft(plan, fext, fext);
279  fftw_destroy_plan(plan);
280 
281  // Free Fmm memory.
282  free(Fmm);
283 
284  // Extract f from version of f extended to the torus (fext).
285  memcpy(f, fext, L*(2*L-1)*sizeof(SSHT_COMPLEX(double)));
286  /* Memcpy equivalent to:
287  for (t=0; t<=L-1; t++)
288  for (p=0; p<=2*L-2; p++)
289  f[t*fext_stride + p] = fext[t*fext_stride + p];
290  */
291 
292  // Free fext memory.
293  free(fext);
294 
295  // Print finished if verbosity set.
296  if (verbosity > 0)
297  printf("%s%s", SSHT_PROMPT, "Inverse transform computed!");
298 
299  // Free precomputation memory.
300  free(sqrt_tbl);
301  free(signs);
302  free(exps);
303  free(m_factors);
304 
305 }
306 
322 void ssht_core_mw_inverse_sov_sym_real(double *f, const SSHT_COMPLEX(double) *flm,
323  int L,
324  ssht_dl_method_t dl_method,
325  int verbosity) {
327  0, L,
328  dl_method,
329  verbosity);
330 }
331 
348 void ssht_core_mw_lb_inverse_sov_sym_real(double *f, const SSHT_COMPLEX(double) *flm,
349  int L0, int L,
350  ssht_dl_method_t dl_method,
351  int verbosity) {
352 
353  int el, m, mm;
354  //int t, p;
355  int eltmp;
356  double *sqrt_tbl, *signs;
357  int el2pel, m_offset;
358  double ssign, elfactor;
359  SSHT_COMPLEX(double) *m_factors;
360  SSHT_COMPLEX(double) mmfactor;
361  double *dl;
362  double *dl8 = NULL;
363  int dl_offset, dl_stride;
364  SSHT_COMPLEX(double) *exps;
365  int exps_offset;
366  double elmmsign, elssign;
367  int spinneg;
368  SSHT_COMPLEX(double) *Fmm, *Fmm_shift;
369  double *fext_real;
370  int Fmm_offset, Fmm_stride;
371  // int fext_stride;
372  fftw_plan plan;
373  int spin = 0;
374 
375  // Allocate memory.
376  sqrt_tbl = (double*)calloc(2*(L-1)+2, sizeof(double));
378  signs = (double*)calloc(L+1, sizeof(double));
380  exps = (SSHT_COMPLEX(double)*)calloc(2*L-1, sizeof(SSHT_COMPLEX(double)));
382  m_factors = calloc(2*L-1, sizeof *m_factors);
383  SSHT_ERROR_MEM_ALLOC_CHECK(m_factors)
384 
385  // Perform precomputations.
386  for (el=0; el<=2*(L-1)+1; el++)
387  sqrt_tbl[el] = sqrt((double)el);
388  for (m=0; m<=L-1; m=m+2) {
389  signs[m] = 1.0;
390  signs[m+1] = -1.0;
391  }
392  ssign = signs[abs(spin)];
393  spinneg = spin <= 0 ? spin : -spin;
394  exps_offset = L-1;
395  for (m=-(L-1); m<=L-1; m++)
396  exps[m + exps_offset] = cexp(-I*SSHT_PION2*(m+spin));
397 
398  // Print messages depending on verbosity level.
399  if (verbosity > 0) {
400  printf("%s%s\n", SSHT_PROMPT,
401  "Computing inverse transform using MW sampling with ");
402  printf("%s%s%d%s%d%s\n", SSHT_PROMPT, "parameters (L,spin,reality) = (",
403  L, ",", spin, ", TRUE)");
404  if (verbosity > 1)
405  printf("%s%s\n", SSHT_PROMPT,
406  "Using routine ssht_core_mw_inverse_sov_sym_real...");
407  }
408 
409  // Compute Fmm.
410  Fmm = (SSHT_COMPLEX(double)*)calloc((2*L-1)*L, sizeof(SSHT_COMPLEX(double)));
412  Fmm_offset = L-1;
413  Fmm_stride = L;
416  if (dl_method == SSHT_DL_RISBO) {
419  }
420  dl_offset = ssht_dl_get_offset(L, SSHT_DL_QUARTER);
421  dl_stride = ssht_dl_get_stride(L, SSHT_DL_QUARTER);
422  m_offset = L-1;
423  for (el=MAX(L0, abs(spin)); el<=L-1; el++) {
424 
425  // Compute Wigner plane.
426  switch (dl_method) {
427 
428  case SSHT_DL_RISBO:
429  if (el!=0 && el==MAX(L0, abs(spin))) {
430  for(eltmp=0; eltmp<=MAX(L0, abs(spin)); eltmp++)
433  eltmp, sqrt_tbl, signs);
435  dl8, L,
438  el,
439  signs);
440  }
441  else {
444  el, sqrt_tbl, signs);
446  dl8, L,
449  el,
450  signs);
451  }
452  break;
453 
454  case SSHT_DL_TRAPANI:
455  if (el!=0 && el==MAX(L0, abs(spin))) {
456  for(eltmp=0; eltmp<=MAX(L0, abs(spin)); eltmp++)
459  eltmp, sqrt_tbl);
462  el, signs);
463  }
464  else {
467  el, sqrt_tbl);
470  el, signs);
471  }
472  break;
473 
474  default:
475  SSHT_ERROR_GENERIC("Invalid dl method")
476  }
477 
478 
479  // Compute Fmm.
480  elfactor = ssign * sqrt((double)(2.0*el+1.0)/(4.0*SSHT_PI));
481  el2pel = el*el + el;
482  for (m=-el; m<=el; m++)
483  m_factors[m + m_offset] = flm[el2pel + m] * exps[m + exps_offset];
484 
485  for (mm=0; mm<=el; mm++) {
486  double mm_factor;
487  int mm_offset = mm*dl_stride;
488  elmmsign = signs[el] * signs[mm];
489  elssign = spin <= 0 ? 1.0 : elmmsign;
490 
491  mm_factor = elfactor * elssign * dl[mm_offset - spinneg + dl_offset];
492 
493  for (m=0; m<=el; m++) {
494  Fmm[(mm + Fmm_offset)*Fmm_stride + m] +=
495  mm_factor
496  * m_factors[m + m_offset]
497  * dl[mm_offset + m + dl_offset];
498  }
499 
500  }
501 
502  }
503 
504  // Free dl memory.
505  free(dl);
506  if (dl_method == SSHT_DL_RISBO)
507  free(dl8);
508 
509  // Use symmetry to compute Fmm for negative mm.
510  for (mm=-(L-1); mm<=-1; mm++)
511  for (m=0; m<=L-1; m++)
512  Fmm[(mm + Fmm_offset)*Fmm_stride + m] =
513  signs[abs(m)] * ssign
514  * Fmm[(-mm + Fmm_offset)*Fmm_stride + m];
515 
516  // Apply phase modulation to account for sampling offset.
517  for (mm=-(L-1); mm<=L-1; mm++) {
518  mmfactor = cexp(I*mm*SSHT_PI/(2.0*L-1.0));
519  for (m=0; m<=L-1; m++)
520  Fmm[(mm + Fmm_offset)*Fmm_stride + m] *=
521  mmfactor;
522  }
523 
524  // Apply spatial shift.
525  Fmm_shift = (SSHT_COMPLEX(double)*)calloc((2*L-1)*L, sizeof(SSHT_COMPLEX(double)));
526  SSHT_ERROR_MEM_ALLOC_CHECK(Fmm_shift)
527  for (mm=0; mm<=L-1; mm++)
528  for (m=0; m<=L-1; m++)
529  Fmm_shift[mm*Fmm_stride + m] =
530  Fmm[(mm + Fmm_offset)*Fmm_stride + m];
531  for (mm=-(L-1); mm<=-1; mm++)
532  for (m=0; m<=L-1; m++)
533  Fmm_shift[(mm + 2*L-1)*Fmm_stride + m] =
534  Fmm[(mm + Fmm_offset)*Fmm_stride + m];
535 
536  // Allocate space for function values.
537  fext_real = (double*)calloc((2*L-1)*(2*L-1), sizeof(double));
538  SSHT_ERROR_MEM_ALLOC_CHECK(fext_real)
539  // fext_stride = 2*L-1;
540 
541  // Perform 2D FFT.
542  plan = fftw_plan_dft_c2r_2d(2*L-1, 2*L-1, Fmm_shift, fext_real,
543  FFTW_ESTIMATE);
544  fftw_execute_dft_c2r(plan, Fmm_shift, fext_real);
545  fftw_destroy_plan(plan);
546 
547  // Free Fmm memory.
548  free(Fmm);
549  free(Fmm_shift);
550 
551  // Extract f from version of f extended to the torus (fext).
552  memcpy(f, fext_real, L*(2*L-1)*sizeof(double));
553  /* Memcpy equivalent to:
554  for (t=0; t<=L-1; t++)
555  for (p=0; p<=2*L-2; p++)
556  f[t*fext_stride + p] = fext[t*fext_stride + p];
557  */
558 
559  // Free fext memory.
560  free(fext_real);
561 
562  // Print finished if verbosity set.
563  if (verbosity > 0)
564  printf("%s%s", SSHT_PROMPT, "Inverse transform computed!");
565 
566  // Free precomputation memory.
567  free(sqrt_tbl);
568  free(signs);
569  free(exps);
570  free(m_factors);
571 
572 }
573 
574 
590 void ssht_core_mwdirect_inverse(SSHT_COMPLEX(double) *f, const SSHT_COMPLEX(double) *flm,
591  int L, int spin, int verbosity) {
592 
593  int t, p, m, el, ind, eltmp;
594  double *dl;
595  double *sqrt_tbl;
596  double theta, phi, elfactor;
597  int ssign;
598  int dl_offset, dl_stride, f_stride;
599 
600  // Allocate memory.
601  sqrt_tbl = (double*)calloc(2*(L-1)+2, sizeof(double));
603 
604  // Perform precomputations.
605  for (el=0; el<=2*(L-1)+1; el++)
606  sqrt_tbl[el] = sqrt((double)el);
607  ssign = spin & 1;
608  ssign = 1 - ssign - ssign; // (-1)^spin
609 
610  // Print messages depending on verbosity level.
611  if (verbosity > 0) {
612  printf("%s%s\n", SSHT_PROMPT,
613  "Computing inverse transform using MW sampling with ");
614  printf("%s%s%d%s%d%s\n", SSHT_PROMPT, "parameters (L,spin,reality) = (",
615  L, ",", spin, ", FALSE)");
616  if (verbosity > 1)
617  printf("%s%s\n", SSHT_PROMPT,
618  "Using routine ssht_core_mwdirect_inverse...");
619  }
620 
621  // Initialise f with zeros.
622  f_stride = 2*L-1;
623  for (t=0; t<=L-1; t++)
624  for (p=0; p<=2*L-2; p++)
625  f[t*f_stride + p] = 0.0;
626 
627  // Compute inverse transform.
628  dl = ssht_dl_calloc(L, SSHT_DL_FULL);
630  dl_offset = ssht_dl_get_offset(L, SSHT_DL_FULL);
631  dl_stride = ssht_dl_get_stride(L, SSHT_DL_FULL);
632  for (t=0; t<=L-1; t++) {
633  theta = ssht_sampling_mw_t2theta(t, L);
634  for (el=abs(spin); el<=L-1; el++) {
635  elfactor = sqrt((double)(2.0*el+1.0)/(4.0*SSHT_PI));
636  if (el!=0 && el==abs(spin)) {
637  for(eltmp=0; eltmp<=abs(spin); eltmp++)
638  ssht_dl_beta_risbo_full_table(dl, theta, L,
639  SSHT_DL_FULL,
640  eltmp, sqrt_tbl);
641  }
642  else {
643  ssht_dl_beta_risbo_full_table(dl, theta, L,
644  SSHT_DL_FULL,
645  el, sqrt_tbl);
646  }
647 
648  for (m=-el; m<=el; m++) {
649  ssht_sampling_elm2ind(&ind, el, m);
650  for (p=0; p<=2*L-2; p++) {
651  phi = ssht_sampling_mw_p2phi(p, L);
652  f[t*f_stride + p] +=
653  ssign
654  * elfactor
655  * cexp(I*m*phi)
656  * dl[(m+dl_offset)*dl_stride - spin + dl_offset]
657  * flm[ind];
658  }
659  }
660 
661  }
662  }
663 
664  free(sqrt_tbl);
665  free(dl);
666 
667  // Print finished if verbosity set.
668  if (verbosity > 0)
669  printf("%s%s", SSHT_PROMPT, "Inverse transform computed!");
670 
671 }
672 
673 
689 void ssht_core_mwdirect_inverse_sov(SSHT_COMPLEX(double) *f, const SSHT_COMPLEX(double) *flm,
690  int L, int spin, int verbosity) {
691 
692  int t, p, m, el, ind;
693  int ftm_stride, ftm_offset, f_stride;
694  double *dlm1p1_line, *dl_line;
695  double *dl_ptr;
696  double *sqrt_tbl, *signs;
697  SSHT_COMPLEX(double) *ftm, *inout;
698  double theta, ssign, elfactor;
699  fftw_plan plan;
700 
701  // Allocate memory.
702  sqrt_tbl = (double*)calloc(2*(L-1)+2, sizeof(double));
704  signs = (double*)calloc(L+1, sizeof(double));
706 
707  // Perform precomputations.
708  for (el=0; el<=2*(L-1)+1; el++)
709  sqrt_tbl[el] = sqrt((double)el);
710  for (m=0; m<=L-1; m=m+2) {
711  signs[m] = 1.0;
712  signs[m+1] = -1.0;
713  }
714  ssign = signs[abs(spin)];
715 
716  // Print messages depending on verbosity level.
717  if (verbosity > 0) {
718  printf("%s%s\n", SSHT_PROMPT,
719  "Computing inverse transform using MW sampling with ");
720  printf("%s%s%d%s%d%s\n", SSHT_PROMPT, "parameters (L,spin,reality) = (",
721  L, ",", spin, ", FALSE)");
722  if (verbosity > 1)
723  printf("%s%s\n", SSHT_PROMPT,
724  "Using routine ssht_core_mwdirect_inverse_sov...");
725  }
726 
727  // Compute ftm.
728  ftm = (SSHT_COMPLEX(double)*)calloc(L*(2*L-1), sizeof(SSHT_COMPLEX(double)));
730  ftm_stride = 2*L-1;
731  ftm_offset = L-1;
732  dlm1p1_line = (double*)calloc(2*L-1, sizeof(double));
733  SSHT_ERROR_MEM_ALLOC_CHECK(dlm1p1_line)
734  dl_line = (double*)calloc(2*L-1, sizeof(double));
736  for (t=0; t<=L-1; t++) {
737  theta = ssht_sampling_mw_t2theta(t, L);
738  for (el=abs(spin); el<=L-1; el++) {
739  elfactor = sqrt((double)(2.0*el+1.0)/(4.0*SSHT_PI));
740 
741  // Compute dl line for given spin.
742  ssht_dl_beta_kostelec_line_table(dlm1p1_line, dl_line,
743  theta, L, -spin, el,
744  sqrt_tbl, signs);
745  // Switch current and previous dls.
746  dl_ptr = dl_line;
747  dl_line = dlm1p1_line;
748  dlm1p1_line = dl_ptr;
749 
750  for (m=-el; m<=el; m++) {
751  ssht_sampling_elm2ind(&ind, el, m);
752  ftm[t*ftm_stride + m + ftm_offset] +=
753  ssign
754  * elfactor
755  * dl_line[m + L-1]
756  * flm[ind];
757  }
758  }
759  }
760 
761  // Free dl memory.
762  free(dlm1p1_line);
763  free(dl_line);
764 
765  // Compute f.
766  inout = (SSHT_COMPLEX(double)*)calloc(2*L-1, sizeof(SSHT_COMPLEX(double)));
768  f_stride = 2*L-1;
769  plan = fftw_plan_dft_1d(2*L-1, inout, inout, FFTW_BACKWARD, FFTW_MEASURE);
770  for (t=0; t<=L-1; t++) {
771  for (m=0; m<=L-1; m++)
772  inout[m] = ftm[t*ftm_stride + m + ftm_offset];
773  for (m=-(L-1); m<=-1; m++)
774  inout[m+2*L-1] = ftm[t*ftm_stride + m + ftm_offset];
775  fftw_execute_dft(plan, inout, inout);
776  for (p=0; p<=2*L-2; p++)
777  f[t*f_stride + p] = inout[p];
778  }
779  fftw_destroy_plan(plan);
780 
781  // Free memory.
782  free(ftm);
783  free(inout);
784  free(signs);
785  free(sqrt_tbl);
786 
787  // Print finished if verbosity set.
788  if (verbosity > 0)
789  printf("%s%s", SSHT_PROMPT, "Inverse transform computed!");
790 
791 }
792 
793 
810 void ssht_core_mw_forward_sov_conv_sym(SSHT_COMPLEX(double) *flm, const SSHT_COMPLEX(double) *f,
811  int L, int spin,
812  ssht_dl_method_t dl_method,
813  int verbosity) {
815  0, L, spin,
816  dl_method,
817  verbosity);
818 }
819 
838  int L0, int L, int spin,
839  ssht_dl_method_t dl_method,
840  int verbosity) {
841 
842  int el, m, mm, ind, t, r;
843  int eltmp;
844  double *sqrt_tbl, *signs;
845  int el2pel;
846  double ssign, elfactor;
847  fftw_plan plan, plan_bwd, plan_fwd;
848  SSHT_COMPLEX(double) *inout;
849  SSHT_COMPLEX(double) *Fmt, *Fmm, *Gmm, *m_mm_factor;
850  SSHT_COMPLEX(double) *w, *wr;
851  SSHT_COMPLEX(double) *Fmm_pad, *tmp_pad;
852  int f_stride, Fmt_stride, Fmt_offset, Fmm_stride, Fmm_offset;
853  double *dl;
854  double *dl8 = NULL;
855  int dl_offset, dl_stride;
856  int w_offset;
857  SSHT_COMPLEX(double) *expsm, *expsmm;
858  int exps_offset;
859  int elmmsign, elssign;
860  int spinneg;
861 
862  // Allocate memory.
863  sqrt_tbl = (double*)calloc(2*(L-1)+2, sizeof(double));
865  signs = (double*)calloc(L+1, sizeof(double));
867  expsm = (SSHT_COMPLEX(double)*)calloc(2*L-1, sizeof(SSHT_COMPLEX(double)));
869  expsmm = (SSHT_COMPLEX(double)*)calloc(2*L-1, sizeof(SSHT_COMPLEX(double)));
871 
872  // Perform precomputations.
873  for (el=0; el<=2*(L-1)+1; el++)
874  sqrt_tbl[el] = sqrt((double)el);
875  for (m=0; m<=L-1; m=m+2) {
876  signs[m] = 1.0;
877  signs[m+1] = -1.0;
878  }
879  ssign = signs[abs(spin)];
880  spinneg = spin <= 0 ? spin : -spin;
881  exps_offset = L-1;
882  for (m=-(L-1); m<=L-1; m++)
883  expsm[m + exps_offset] = cexp(I*SSHT_PION2*(m+spin));
884  for (mm=-(L-1); mm<=L-1; mm++)
885  expsmm[mm + exps_offset] = cexp(-I*mm*SSHT_PI/(2.0*L-1.0));
886 
887  // Print messages depending on verbosity level.
888  if (verbosity > 0) {
889  printf("%s%s\n", SSHT_PROMPT,
890  "Computing forward transform using MW sampling with ");
891  printf("%s%s%d%s%d%s\n", SSHT_PROMPT, "parameters (L,spin,reality) = (",
892  L, ",", spin, ", FALSE)");
893  if (verbosity > 1)
894  printf("%s%s\n", SSHT_PROMPT,
895  "Using routine ssht_core_mw_forward_sov_conv_sym...");
896  }
897 
898  // Compute Fourier transform over phi, i.e. compute Fmt.
899  Fmt = (SSHT_COMPLEX(double)*)calloc((2*L-1)*(2*L-1), sizeof(SSHT_COMPLEX(double)));
901  Fmt_stride = 2*L-1;
902  Fmt_offset = L-1;
903  f_stride = 2*L-1;
904  inout = (SSHT_COMPLEX(double)*)calloc(2*L-1, sizeof(SSHT_COMPLEX(double)));
906  plan = fftw_plan_dft_1d(2*L-1, inout, inout, FFTW_FORWARD, FFTW_ESTIMATE);
907  for (t=0; t<=L-1; t++) {
908  memcpy(inout, &f[t*f_stride], f_stride*sizeof(SSHT_COMPLEX(double)));
909  fftw_execute_dft(plan, inout, inout);
910  for(m=0; m<=L-1; m++)
911  Fmt[(m+Fmt_offset)*Fmt_stride + t] = inout[m] / (2.0*L-1.0);
912  for(m=-(L-1); m<=-1; m++)
913  Fmt[(m+Fmt_offset)*Fmt_stride + t] = inout[m+2*L-1] / (2.0*L-1.0);
914  }
915 
916  // Extend Fmt periodically.
917  for (m=-(L-1); m<=L-1; m++)
918  for (t=L; t<=2*L-2; t++)
919  Fmt[(m+Fmt_offset)*Fmt_stride + t] =
920  signs[abs(m)] * ssign * Fmt[(m+Fmt_offset)*Fmt_stride + (2*L-2-t)];
921 
922  // Compute Fourier transform over theta, i.e. compute Fmm.
923  Fmm = (SSHT_COMPLEX(double)*)calloc((2*L-1)*(2*L-1), sizeof(SSHT_COMPLEX(double)));
925  Fmm_stride = 2*L-1;
926  Fmm_offset = L-1;
927  for (m=-(L-1); m<=L-1; m++) {
928  memcpy(inout, &Fmt[(m+Fmt_offset)*Fmt_stride], Fmt_stride*sizeof(SSHT_COMPLEX(double)));
929  fftw_execute_dft(plan, inout, inout);
930  for(mm=0; mm<=L-1; mm++)
931  Fmm[(m+Fmm_offset)*Fmm_stride + mm + Fmm_offset] =
932  inout[mm] / (2.0*L-1.0);
933  for(mm=-(L-1); mm<=-1; mm++)
934  Fmm[(m+Fmm_offset)*Fmm_stride + mm + Fmm_offset] =
935  inout[mm+2*L-1] / (2.0*L-1.0);
936  }
937  fftw_destroy_plan(plan);
938  free(inout);
939 
940  // Apply phase modulation to account for sampling offset.
941  for (m=-(L-1); m<=L-1; m++)
942  for (mm=-(L-1); mm<=L-1; mm++)
943  Fmm[(m+Fmm_offset)*Fmm_stride + mm + Fmm_offset] *=
944  expsmm[mm + exps_offset];
945 
946  // Compute weights.
947  w = (SSHT_COMPLEX(double)*)calloc(4*L-3, sizeof(SSHT_COMPLEX(double)));
949  w_offset = 2*(L-1);
950  for (mm=-2*(L-1); mm<=2*(L-1); mm++)
951  w[mm+w_offset] = ssht_sampling_weight_mw(mm);
952 
953  // Compute IFFT of w to give wr.
954  wr = (SSHT_COMPLEX(double)*)calloc(4*L-3, sizeof(SSHT_COMPLEX(double)));
956  inout = (SSHT_COMPLEX(double)*)calloc(4*L-3, sizeof(SSHT_COMPLEX(double)));
958  plan_bwd = fftw_plan_dft_1d(4*L-3, inout, inout, FFTW_BACKWARD, FFTW_MEASURE);
959  plan_fwd = fftw_plan_dft_1d(4*L-3, inout, inout, FFTW_FORWARD, FFTW_MEASURE);
960  for (mm=1; mm<=2*L-2; mm++)
961  inout[mm + w_offset] = w[mm - 2*(L-1) - 1 + w_offset];
962  for (mm=-2*(L-1); mm<=0; mm++)
963  inout[mm + w_offset] = w[mm + 2*(L-1) + w_offset];
964  fftw_execute_dft(plan_bwd, inout, inout);
965  for (mm=0; mm<=2*L-2; mm++)
966  wr[mm + w_offset] = inout[mm - 2*(L-1) + w_offset];
967  for (mm=-2*(L-1); mm<=-1; mm++)
968  wr[mm + w_offset] = inout[mm + 2*(L-1) + 1 + w_offset];
969 
970  // Compute Gmm by convolution implemented as product in real space.
971  Fmm_pad = (SSHT_COMPLEX(double)*)calloc(4*L-3, sizeof(SSHT_COMPLEX(double)));
973  tmp_pad = (SSHT_COMPLEX(double)*)calloc(4*L-3, sizeof(SSHT_COMPLEX(double)));
975  Gmm = (SSHT_COMPLEX(double)*)calloc((2*L-1)*(2*L-1), sizeof(SSHT_COMPLEX(double)));
977  for (m=-(L-1); m<=L-1; m++) {
978 
979  // Zero-pad Fmm.
980  for (mm=-2*(L-1); mm<=-L; mm++)
981  Fmm_pad[mm+w_offset] = 0.0;
982  for (mm=L; mm<=2*(L-1); mm++)
983  Fmm_pad[mm+w_offset] = 0.0;
984  for (mm=-(L-1); mm<=L-1; mm++)
985  Fmm_pad[mm + w_offset] =
986  Fmm[(m+Fmm_offset)*Fmm_stride + mm + Fmm_offset];
987 
988  // Compute IFFT of Fmm.
989  for (mm=1; mm<=2*L-2; mm++)
990  inout[mm + w_offset] = Fmm_pad[mm - 2*(L-1) - 1 + w_offset];
991  for (mm=-2*(L-1); mm<=0; mm++)
992  inout[mm + w_offset] = Fmm_pad[mm + 2*(L-1) + w_offset];
993  fftw_execute_dft(plan_bwd, inout, inout);
994  for (mm=0; mm<=2*L-2; mm++)
995  Fmm_pad[mm + w_offset] = inout[mm - 2*(L-1) + w_offset];
996  for (mm=-2*(L-1); mm<=-1; mm++)
997  Fmm_pad[mm + w_offset] = inout[mm + 2*(L-1) + 1 + w_offset];
998 
999  // Compute product of Fmm and weight in real space.
1000  for (r=-2*(L-1); r<=2*(L-1); r++)
1001  Fmm_pad[r + w_offset] *= wr[-r + w_offset];
1002 
1003  // Compute Gmm by FFT.
1004  for (mm=1; mm<=2*L-2; mm++)
1005  inout[mm + w_offset] = Fmm_pad[mm - 2*(L-1) - 1 + w_offset];
1006  for (mm=-2*(L-1); mm<=0; mm++)
1007  inout[mm + w_offset] = Fmm_pad[mm + 2*(L-1) + w_offset];
1008  fftw_execute_dft(plan_fwd, inout, inout);
1009  for (mm=0; mm<=2*L-2; mm++)
1010  Fmm_pad[mm + w_offset] = inout[mm - 2*(L-1) + w_offset];
1011  for (mm=-2*(L-1); mm<=-1; mm++)
1012  Fmm_pad[mm + w_offset] = inout[mm + 2*(L-1) + 1 + w_offset];
1013 
1014  // Extract section of Gmm of interest.
1015  for (mm=-(L-1); mm<=L-1; mm++)
1016  Gmm[(mm+Fmm_offset)*Fmm_stride + m + Fmm_offset] =
1017  Fmm_pad[mm + w_offset] * 2.0 * SSHT_PI / (4.0*L-3.0);
1018 
1019  }
1020  fftw_destroy_plan(plan_bwd);
1021  fftw_destroy_plan(plan_fwd);
1022 
1023  // Precompute factors depending on particular Gmm, to be used later
1024  // in computing the flm.
1025  m_mm_factor = calloc(L*(2*L-1), sizeof *m_mm_factor);
1026  SSHT_ERROR_MEM_ALLOC_CHECK(m_mm_factor)
1027 
1028  for (mm = 1; mm < L; mm++) {
1029  for (m = -L+1; m < L; m++) {
1030  m_mm_factor[mm*Fmm_stride + m + Fmm_offset] =
1031  expsm[m + exps_offset]
1032  * ( Gmm[(mm+Fmm_offset)*Fmm_stride + m + Fmm_offset]
1033  + signs[abs(m)] * ssign
1034  * Gmm[(-mm+Fmm_offset)*Fmm_stride + m + Fmm_offset] );
1035  }
1036  }
1037 
1038  // Compute flm.
1041  if (dl_method == SSHT_DL_RISBO) {
1044  }
1045  dl_offset = ssht_dl_get_offset(L, SSHT_DL_QUARTER);
1046  dl_stride = ssht_dl_get_stride(L, SSHT_DL_QUARTER);
1047  for (el=0; el<=L-1; el++) {
1048  for (m=-el; m<=el; m++) {
1049  ssht_sampling_elm2ind(&ind, el, m);
1050  flm[ind] = 0.0;
1051  }
1052  }
1053  for (el=MAX(L0, abs(spin)); el<=L-1; el++) {
1054  double mm_factor;
1055 
1056  // Compute Wigner plane.
1057  switch (dl_method) {
1058 
1059  case SSHT_DL_RISBO:
1060  if (el!=0 && el==MAX(L0, abs(spin))) {
1061  for(eltmp=0; eltmp<=MAX(L0, abs(spin)); eltmp++)
1064  eltmp, sqrt_tbl, signs);
1066  dl8, L,
1069  el,
1070  signs);
1071  }
1072  else {
1075  el, sqrt_tbl, signs);
1077  dl8, L,
1080  el,
1081  signs);
1082  }
1083  break;
1084 
1085  case SSHT_DL_TRAPANI:
1086  if (el!=0 && el==MAX(L0, abs(spin))) {
1087  for(eltmp=0; eltmp<=MAX(L0, abs(spin)); eltmp++)
1090  eltmp, sqrt_tbl);
1093  el, signs);
1094  }
1095  else {
1098  el, sqrt_tbl);
1101  el, signs);
1102  }
1103  break;
1104 
1105  default:
1106  SSHT_ERROR_GENERIC("Invalid dl method")
1107  }
1108 
1109  // Compute flm.
1110  elfactor = sqrt((double)(2.0*el+1.0)/(4.0*SSHT_PI));
1111  el2pel = el *el + el;
1112  elssign = spin <= 0 ? 1.0 : signs[el];
1113 
1114  mm_factor = ssign * elfactor * elssign * dl[- spinneg + dl_offset];
1115 
1116  for (m=-el; m<=-1; m++) {
1117  // mm = 0
1118  flm[m + el2pel] +=
1119  mm_factor
1120  * signs[el]
1121  * expsm[m + exps_offset]
1122  * dl[0*dl_stride - m + dl_offset]
1123  * Gmm[(0+Fmm_offset)*Fmm_stride + m + Fmm_offset];
1124  }
1125  for (m=0; m<=el; m++) {
1126  // mm = 0
1127  flm[m + el2pel] +=
1128  mm_factor
1129  * expsm[m + exps_offset]
1130  * dl[0*dl_stride + m + dl_offset]
1131  * Gmm[(0+Fmm_offset)*Fmm_stride + m + Fmm_offset];
1132  }
1133 
1134  for (mm=1; mm<=el; mm++) {
1135  int mm_offset = mm * dl_stride;
1136  elmmsign = signs[el] * signs[mm];
1137  elssign = spin <= 0 ? 1.0 : elmmsign;
1138 
1139  mm_factor =
1140  ssign
1141  * elfactor
1142  * elssign
1143  * dl[mm_offset - spinneg + dl_offset];
1144 
1145  for (m=-el; m<=-1; m++) {
1146  flm[m + el2pel] +=
1147  mm_factor
1148  * elmmsign * dl[mm_offset - m + dl_offset]
1149  * m_mm_factor[mm*Fmm_stride + m + Fmm_offset];
1150  }
1151  for (m=0; m<=el; m++) {
1152  flm[m + el2pel] +=
1153  mm_factor
1154  * dl[mm_offset + m + dl_offset]
1155  * m_mm_factor[mm*Fmm_stride + m + Fmm_offset];
1156  }
1157 
1158  }
1159 
1160  }
1161 
1162  // Free memory.
1163  free(dl);
1164  if (dl_method == SSHT_DL_RISBO)
1165  free(dl8);
1166  free(Fmt);
1167  free(Fmm);
1168  free(inout);
1169  free(w);
1170  free(wr);
1171  free(Fmm_pad);
1172  free(tmp_pad);
1173  free(Gmm);
1174  free(m_mm_factor);
1175  free(sqrt_tbl);
1176  free(signs);
1177  free(expsm);
1178  free(expsmm);
1179 
1180  // Print finished if verbosity set.
1181  if (verbosity > 0)
1182  printf("%s%s", SSHT_PROMPT, "Forward transform computed!");
1183 
1184 }
1185 
1186 
1203 void ssht_core_mw_forward_sov_conv_sym_real(SSHT_COMPLEX(double) *flm, const double *f,
1204  int L,
1205  ssht_dl_method_t dl_method,
1206  int verbosity) {
1208  0, L,
1209  dl_method,
1210  verbosity);
1211 }
1212 
1230 void ssht_core_mw_lb_forward_sov_conv_sym_real(SSHT_COMPLEX(double) *flm, const double *f,
1231  int L0, int L,
1232  ssht_dl_method_t dl_method,
1233  int verbosity) {
1234 
1235  int el, m, mm, ind, ind_nm, t, r;
1236  int eltmp;
1237  double *sqrt_tbl, *signs;
1238  int el2pel;
1239  double ssign, elfactor;
1240  fftw_plan plan, plan_bwd, plan_fwd;
1241  double *in_real;
1242  SSHT_COMPLEX(double) *inout, *out;
1243  SSHT_COMPLEX(double) *Fmt, *Fmm, *Gmm, *m_mm_factor;
1244  SSHT_COMPLEX(double) *w, *wr;
1245  SSHT_COMPLEX(double) *Fmm_pad, *tmp_pad;
1246  int f_stride, Fmt_stride, Fmt_offset, Fmm_stride, Fmm_offset, Gmm_stride;
1247  double *dl;
1248  double *dl8 = NULL;
1249  int dl_offset, dl_stride;
1250  int w_offset;
1251  SSHT_COMPLEX(double) *expsm, *expsmm;
1252  int exps_offset;
1253  int elmmsign, elssign;
1254  int spinneg;
1255  int spin = 0;
1256 
1257  // Allocate memory.
1258  sqrt_tbl = (double*)calloc(2*(L-1)+2, sizeof(double));
1259  SSHT_ERROR_MEM_ALLOC_CHECK(sqrt_tbl)
1260  signs = (double*)calloc(L+1, sizeof(double));
1262  expsm = (SSHT_COMPLEX(double)*)calloc(L, sizeof(SSHT_COMPLEX(double)));
1264  expsmm = (SSHT_COMPLEX(double)*)calloc(2*L-1, sizeof(SSHT_COMPLEX(double)));
1266 
1267  // Perform precomputations.
1268  for (el=0; el<=2*(L-1)+1; el++)
1269  sqrt_tbl[el] = sqrt((double)el);
1270  for (m=0; m<=L-1; m=m+2) {
1271  signs[m] = 1.0;
1272  signs[m+1] = -1.0;
1273  }
1274  ssign = signs[abs(spin)];
1275  spinneg = spin <= 0 ? spin : -spin;
1276  exps_offset = L-1;
1277  for (m=0; m<=L-1; m++)
1278  expsm[m] = cexp(I*SSHT_PION2*(m+spin));
1279  for (mm=-(L-1); mm<=L-1; mm++)
1280  expsmm[mm + exps_offset] = cexp(-I*mm*SSHT_PI/(2.0*L-1.0));
1281 
1282  // Print messages depending on verbosity level.
1283  if (verbosity > 0) {
1284  printf("%s%s\n", SSHT_PROMPT,
1285  "Computing forward transform using MW sampling with ");
1286  printf("%s%s%d%s%d%s\n", SSHT_PROMPT, "parameters (L,spin,reality) = (",
1287  L, ",", spin, ", TRUE)");
1288  if (verbosity > 1)
1289  printf("%s%s\n", SSHT_PROMPT,
1290  "Using routine ssht_core_mw_forward_sov_conv_sym_real...");
1291  }
1292 
1293  // Compute Fourier transform over phi, i.e. compute Fmt.
1294  Fmt = (SSHT_COMPLEX(double)*)calloc(L*(2*L-1), sizeof(SSHT_COMPLEX(double)));
1296  Fmt_stride = 2*L-1;
1297  f_stride = 2*L-1;
1298  in_real = (double*)calloc(2*L-1, sizeof(double));
1300  out = (SSHT_COMPLEX(double)*)calloc(L, sizeof(SSHT_COMPLEX(double)));
1302  plan = fftw_plan_dft_r2c_1d(2*L-1, in_real, out, FFTW_MEASURE);
1303  for (t=0; t<=L-1; t++) {
1304  memcpy(in_real, &f[t*f_stride], f_stride*sizeof(double));
1305  fftw_execute_dft_r2c(plan, in_real, out);
1306  for(m=0; m<=L-1; m++)
1307  Fmt[m*Fmt_stride + t] = out[m] / (2.0*L-1.0);
1308  }
1309  free(in_real);
1310  free(out);
1311  fftw_destroy_plan(plan);
1312 
1313  // Extend Fmt periodically.
1314  for (m=0; m<=L-1; m++)
1315  for (t=L; t<=2*L-2; t++)
1316  Fmt[m*Fmt_stride + t] =
1317  signs[abs(m)] * ssign * Fmt[m*Fmt_stride + (2*L-2-t)];
1318 
1319  // Compute Fourier transform over theta, i.e. compute Fmm.
1320  Fmm = (SSHT_COMPLEX(double)*)calloc(L*(2*L-1), sizeof(SSHT_COMPLEX(double)));
1322  Fmm_stride = 2*L-1;
1323  Fmm_offset = L-1;
1324  inout = (SSHT_COMPLEX(double)*)calloc(2*L-1, sizeof(SSHT_COMPLEX(double)));
1326  plan = fftw_plan_dft_1d(2*L-1, inout, inout, FFTW_FORWARD, FFTW_MEASURE);
1327  for (m=0; m<=L-1; m++) {
1328  memcpy(inout, &Fmt[m*Fmt_stride], Fmt_stride*sizeof(SSHT_COMPLEX(double)));
1329  fftw_execute_dft(plan, inout, inout);
1330  for(mm=0; mm<=L-1; mm++)
1331  Fmm[m*Fmm_stride + mm + Fmm_offset] =
1332  inout[mm] / (2.0*L-1.0);
1333  for(mm=-(L-1); mm<=-1; mm++)
1334  Fmm[m*Fmm_stride + mm + Fmm_offset] =
1335  inout[mm+2*L-1] / (2.0*L-1.0);
1336  }
1337  fftw_destroy_plan(plan);
1338  free(inout);
1339 
1340  // Apply phase modulation to account for sampling offset.
1341  for (m=0; m<=L-1; m++)
1342  for (mm=-(L-1); mm<=L-1; mm++)
1343  Fmm[m*Fmm_stride + mm + Fmm_offset] *=
1344  expsmm[mm + exps_offset];
1345 
1346  // Compute weights.
1347  w = (SSHT_COMPLEX(double)*)calloc(4*L-3, sizeof(SSHT_COMPLEX(double)));
1349  w_offset = 2*(L-1);
1350  for (mm=-2*(L-1); mm<=2*(L-1); mm++)
1351  w[mm+w_offset] = ssht_sampling_weight_mw(mm);
1352 
1353  // Compute IFFT of w to give wr.
1354  wr = (SSHT_COMPLEX(double)*)calloc(4*L-3, sizeof(SSHT_COMPLEX(double)));
1356  inout = (SSHT_COMPLEX(double)*)calloc(4*L-3, sizeof(SSHT_COMPLEX(double)));
1358  plan_bwd = fftw_plan_dft_1d(4*L-3, inout, inout, FFTW_BACKWARD, FFTW_MEASURE);
1359  plan_fwd = fftw_plan_dft_1d(4*L-3, inout, inout, FFTW_FORWARD, FFTW_MEASURE);
1360  for (mm=1; mm<=2*L-2; mm++)
1361  inout[mm + w_offset] = w[mm - 2*(L-1) - 1 + w_offset];
1362  for (mm=-2*(L-1); mm<=0; mm++)
1363  inout[mm + w_offset] = w[mm + 2*(L-1) + w_offset];
1364  fftw_execute_dft(plan_bwd, inout, inout);
1365  for (mm=0; mm<=2*L-2; mm++)
1366  wr[mm + w_offset] = inout[mm - 2*(L-1) + w_offset];
1367  for (mm=-2*(L-1); mm<=-1; mm++)
1368  wr[mm + w_offset] = inout[mm + 2*(L-1) + 1 + w_offset];
1369 
1370  // Compute Gmm by convolution implemented as product in real space.
1371  Fmm_pad = (SSHT_COMPLEX(double)*)calloc(4*L-3, sizeof(SSHT_COMPLEX(double)));
1373  tmp_pad = (SSHT_COMPLEX(double)*)calloc(4*L-3, sizeof(SSHT_COMPLEX(double)));
1375  Gmm = (SSHT_COMPLEX(double)*)calloc((2*L-1)*L, sizeof(SSHT_COMPLEX(double)));
1376  Gmm_stride = L;
1378  for (m=0; m<=L-1; m++) {
1379 
1380  // Zero-pad Fmm.
1381  for (mm=-2*(L-1); mm<=-L; mm++)
1382  Fmm_pad[mm+w_offset] = 0.0;
1383  for (mm=L; mm<=2*(L-1); mm++)
1384  Fmm_pad[mm+w_offset] = 0.0;
1385  for (mm=-(L-1); mm<=L-1; mm++)
1386  Fmm_pad[mm + w_offset] =
1387  Fmm[m*Fmm_stride + mm + Fmm_offset];
1388 
1389  // Compute IFFT of Fmm.
1390  for (mm=1; mm<=2*L-2; mm++)
1391  inout[mm + w_offset] = Fmm_pad[mm - 2*(L-1) - 1 + w_offset];
1392  for (mm=-2*(L-1); mm<=0; mm++)
1393  inout[mm + w_offset] = Fmm_pad[mm + 2*(L-1) + w_offset];
1394  fftw_execute_dft(plan_bwd, inout, inout);
1395  for (mm=0; mm<=2*L-2; mm++)
1396  Fmm_pad[mm + w_offset] = inout[mm - 2*(L-1) + w_offset];
1397  for (mm=-2*(L-1); mm<=-1; mm++)
1398  Fmm_pad[mm + w_offset] = inout[mm + 2*(L-1) + 1 + w_offset];
1399 
1400  // Compute product of Fmm and weight in real space.
1401  for (r=-2*(L-1); r<=2*(L-1); r++)
1402  Fmm_pad[r + w_offset] *= wr[-r + w_offset];
1403 
1404  // Compute Gmm by FFT.
1405  for (mm=1; mm<=2*L-2; mm++)
1406  inout[mm + w_offset] = Fmm_pad[mm - 2*(L-1) - 1 + w_offset];
1407  for (mm=-2*(L-1); mm<=0; mm++)
1408  inout[mm + w_offset] = Fmm_pad[mm + 2*(L-1) + w_offset];
1409  fftw_execute_dft(plan_fwd, inout, inout);
1410  for (mm=0; mm<=2*L-2; mm++)
1411  Fmm_pad[mm + w_offset] = inout[mm - 2*(L-1) + w_offset];
1412  for (mm=-2*(L-1); mm<=-1; mm++)
1413  Fmm_pad[mm + w_offset] = inout[mm + 2*(L-1) + 1 + w_offset];
1414 
1415  // Extract section of Gmm of interest.
1416  for (mm=-(L-1); mm<=L-1; mm++)
1417  Gmm[(mm+Fmm_offset)*Gmm_stride + m] =
1418  Fmm_pad[mm + w_offset] * 2.0 * SSHT_PI / (4.0*L-3.0);
1419 
1420  }
1421  fftw_destroy_plan(plan_bwd);
1422  fftw_destroy_plan(plan_fwd);
1423 
1424  // Precompute factors depending on particular Gmm, to be used later
1425  // in computing the flm.
1426  m_mm_factor = calloc(L*L, sizeof *m_mm_factor);
1427  SSHT_ERROR_MEM_ALLOC_CHECK(m_mm_factor)
1428 
1429  for (mm = 1; mm < L; mm++) {
1430  for (m = 0; m < L; m++) {
1431  m_mm_factor[mm*Gmm_stride + m] =
1432  expsm[m]
1433  * ( Gmm[(mm+Fmm_offset)*Gmm_stride + m]
1434  + signs[m] * ssign
1435  * Gmm[(-mm+Fmm_offset)*Gmm_stride + m] );
1436  }
1437  }
1438 
1439  // Compute flm.
1442  if (dl_method == SSHT_DL_RISBO) {
1445  }
1446  dl_offset = ssht_dl_get_offset(L, SSHT_DL_QUARTER);
1447  dl_stride = ssht_dl_get_stride(L, SSHT_DL_QUARTER);
1448  for (el=0; el<=L-1; el++) {
1449  for (m=0; m<=el; m++) {
1450  ssht_sampling_elm2ind(&ind, el, m);
1451  flm[ind] = 0.0;
1452  }
1453  }
1454  for (el=MAX(L0, abs(spin)); el<=L-1; el++) {
1455  double mm_factor;
1456 
1457  // Compute Wigner plane.
1458  switch (dl_method) {
1459 
1460  case SSHT_DL_RISBO:
1461  if (el!=0 && el==MAX(L0, abs(spin))) {
1462  for(eltmp=0; eltmp<=MAX(L0, abs(spin)); eltmp++)
1465  eltmp, sqrt_tbl, signs);
1467  dl8, L,
1470  el,
1471  signs);
1472  }
1473  else {
1476  el, sqrt_tbl, signs);
1478  dl8, L,
1481  el,
1482  signs);
1483  }
1484  break;
1485 
1486  case SSHT_DL_TRAPANI:
1487  if (el!=0 && el==MAX(L0, abs(spin))) {
1488  for(eltmp=0; eltmp<=MAX(L0, abs(spin)); eltmp++)
1491  eltmp, sqrt_tbl);
1494  el, signs);
1495  }
1496  else {
1499  el, sqrt_tbl);
1502  el, signs);
1503  }
1504  break;
1505 
1506  default:
1507  SSHT_ERROR_GENERIC("Invalid dl method")
1508  }
1509 
1510 
1511  // Compute flm.
1512  elfactor = sqrt((double)(2.0*el+1.0)/(4.0*SSHT_PI));
1513  el2pel = el *el + el;
1514  elssign = spin <= 0 ? 1.0 : signs[el];
1515 
1516  mm_factor = ssign * elfactor * elssign * dl[- spinneg + dl_offset];
1517 
1518  for (m=0; m<=el; m++) {
1519  // mm = 0
1520  flm[el2pel + m] +=
1521  mm_factor
1522  * expsm[m]
1523  * dl[0*dl_stride + m + dl_offset]
1524  * Gmm[(0+Fmm_offset)*Gmm_stride + m];
1525  }
1526 
1527  for (mm=1; mm<=el; mm++) {
1528  int mm_offset = mm * dl_stride;
1529  elmmsign = signs[el] * signs[mm];
1530  elssign = spin <= 0 ? 1.0 : elmmsign;
1531 
1532  mm_factor =
1533  ssign
1534  * elfactor
1535  * elssign
1536  * dl[mm_offset - spinneg + dl_offset];
1537 
1538  for (m=0; m<=el; m++) {
1539  flm[el2pel + m] +=
1540  mm_factor
1541  * dl[mm_offset + m + dl_offset]
1542  * m_mm_factor[mm*Gmm_stride + m];
1543  }
1544 
1545  }
1546 
1547  }
1548 
1549  // Set flm values for negative m using conjugate symmetry.
1550  for (el=abs(spin); el<=L-1; el++) {
1551  for (m=1; m<=el; m++) {
1552  ssht_sampling_elm2ind(&ind, el, m);
1553  ssht_sampling_elm2ind(&ind_nm, el, -m);
1554  flm[ind_nm] = signs[m] * conj(flm[ind]);
1555  }
1556  }
1557 
1558  // Free memory.
1559  free(dl);
1560  if (dl_method == SSHT_DL_RISBO)
1561  free(dl8);
1562  free(Fmt);
1563  free(Fmm);
1564  free(inout);
1565  free(w);
1566  free(wr);
1567  free(Fmm_pad);
1568  free(tmp_pad);
1569  free(Gmm);
1570  free(sqrt_tbl);
1571  free(signs);
1572  free(expsm);
1573  free(expsmm);
1574  free(m_mm_factor);
1575 
1576  // Print finished if verbosity set.
1577  if (verbosity > 0)
1578  printf("%s%s", SSHT_PROMPT, "Forward transform computed!");
1579 
1580 }
1581 
1582 
1583 //============================================================================
1584 // MW South pole interfaces
1585 //============================================================================
1586 
1587 
1609  SSHT_COMPLEX(double) *f_sp, double *phi_sp,
1610  const SSHT_COMPLEX(double) *flm,
1611  int L, int spin,
1612  ssht_dl_method_t dl_method,
1613  int verbosity) {
1614 
1615  SSHT_COMPLEX(double)* f_full;
1616  int f_stride = 2*L-1;
1617 
1618  // Allocate full array.
1619  f_full = (SSHT_COMPLEX(double)*)calloc(L*(2*L-1), sizeof(SSHT_COMPLEX(double)));
1621 
1622  // Perform inverse transform.
1623  ssht_core_mw_inverse_sov_sym(f_full, flm, L, spin,
1624  dl_method, verbosity);
1625 
1626  // Copy output function values, including separate point for South pole.
1627  memcpy(f, f_full, (L-1)*(2*L-1)*sizeof(SSHT_COMPLEX(double)));
1628  *f_sp = f_full[(L-1)*f_stride + 0];
1629  *phi_sp = ssht_sampling_mw_p2phi(0, L);
1630 
1631  // Free memory.
1632  free(f_full);
1633 
1634 }
1635 
1636 
1655  double *f_sp,
1656  const SSHT_COMPLEX(double) *flm,
1657  int L,
1658  ssht_dl_method_t dl_method,
1659  int verbosity) {
1660 
1661  double *f_full;
1662  int f_stride = 2*L-1;
1663 
1664  // Allocate full array.
1665  f_full = (double*)calloc(L*(2*L-1), sizeof(double));
1667 
1668  // Perform inverse transform.
1669  ssht_core_mw_inverse_sov_sym_real(f_full, flm, L,
1670  dl_method, verbosity);
1671 
1672  // Copy output function values, including separate point for South pole.
1673  memcpy(f, f_full, (L-1)*(2*L-1)*sizeof(double));
1674  *f_sp = f_full[(L-1)*f_stride + 0];
1675 
1676  // Free memory.
1677  free(f_full);
1678 
1679 }
1680 
1681 
1703  SSHT_COMPLEX(double) f_sp, double phi_sp,
1704  int L, int spin,
1705  ssht_dl_method_t dl_method,
1706  int verbosity) {
1707 
1708  SSHT_COMPLEX(double) *f_full;
1709  int p, f_stride = 2*L-1;
1710  double phi;
1711 
1712  // Copy function values to full array.
1713  f_full = (SSHT_COMPLEX(double)*)calloc(L*(2*L-1), sizeof(SSHT_COMPLEX(double)));
1715  memcpy(f_full, f, (L-1)*(2*L-1)*sizeof(SSHT_COMPLEX(double)));
1716 
1717  // Define South pole for all phi.
1718  for (p=0; p<=2*L-2; p++) {
1719  phi = ssht_sampling_mw_p2phi(p, L);
1720  f_full[(L-1)*f_stride + p] = f_sp * cexp(I*spin*(phi-phi_sp));
1721  }
1722 
1723  // Perform forward transform.
1724  ssht_core_mw_forward_sov_conv_sym(flm, f_full, L, spin,
1725  dl_method, verbosity);
1726 
1727  // Free memory.
1728  free(f_full);
1729 
1730 }
1731 
1732 
1751  const double *f,
1752  double f_sp,
1753  int L,
1754  ssht_dl_method_t dl_method,
1755  int verbosity) {
1756 
1757  double *f_full;
1758  int p, f_stride = 2*L-1;
1759 
1760  // Copy function values to full array.
1761  f_full = (double*)calloc(L*(2*L-1), sizeof(double));
1763  memcpy(f_full, f, (L-1)*(2*L-1)*sizeof(double));
1764 
1765  // Define South pole for all phi.
1766  for (p=0; p<=2*L-2; p++)
1767  f_full[(L-1)*f_stride + p] = f_sp;
1768 
1769  // Perform forward transform.
1771  dl_method, verbosity);
1772 
1773  // Free memory.
1774  free(f_full);
1775 
1776 }
1777 
1778 
1779 //============================================================================
1780 // MW SS algorithms
1781 //============================================================================
1782 
1798 void ssht_core_mw_inverse_sov_sym_ss(SSHT_COMPLEX(double) *f, const SSHT_COMPLEX(double) *flm,
1799  int L, int spin,
1800  ssht_dl_method_t dl_method,
1801  int verbosity) {
1803  0, L, spin,
1804  dl_method,
1805  verbosity);
1806 }
1807 
1825  int L0, int L, int spin,
1826  ssht_dl_method_t dl_method,
1827  int verbosity) {
1828 
1829  int el, m, mm, ind;
1830  //int t, p;
1831  int eltmp;
1832  double *sqrt_tbl, *signs;
1833  int el2pel, inds_offset;
1834  int *inds;
1835  double ssign, elfactor;
1836  double *dl;
1837  double *dl8 = NULL;
1838  int dl_offset, dl_stride;
1839  SSHT_COMPLEX(double) *exps;
1840  int exps_offset;
1841  double elmmsign, elssign;
1842  int spinneg;
1843  SSHT_COMPLEX(double) *Fmm, *fext;
1844  int Fmm_offset, Fmm_stride, fext_stride;
1845  fftw_plan plan;
1846 
1847  // Allocate memory.
1848  sqrt_tbl = (double*)calloc(2*(L-1)+2, sizeof(double));
1849  SSHT_ERROR_MEM_ALLOC_CHECK(sqrt_tbl)
1850  signs = (double*)calloc(L+1, sizeof(double));
1852  exps = (SSHT_COMPLEX(double)*)calloc(2*L-1, sizeof(SSHT_COMPLEX(double)));
1854  inds = (int*)calloc(2*L-1, sizeof(int));
1856 
1857  // Perform precomputations.
1858  for (el=0; el<=2*(L-1)+1; el++)
1859  sqrt_tbl[el] = sqrt((double)el);
1860  for (m=0; m<=L-1; m=m+2) {
1861  signs[m] = 1.0;
1862  signs[m+1] = -1.0;
1863  }
1864  ssign = signs[abs(spin)];
1865  spinneg = spin <= 0 ? spin : -spin;
1866  exps_offset = L-1;
1867  for (m=-(L-1); m<=L-1; m++)
1868  exps[m + exps_offset] = cexp(-I*SSHT_PION2*(m+spin));
1869 
1870  // Print messages depending on verbosity level.
1871  if (verbosity > 0) {
1872  printf("%s%s\n", SSHT_PROMPT,
1873  "Computing inverse transform using MW symmetric sampling with ");
1874  printf("%s%s%d%s%d%s\n", SSHT_PROMPT, "parameters (L,spin,reality) = (",
1875  L, ",", spin, ", FALSE)");
1876  if (verbosity > 1)
1877  printf("%s%s\n", SSHT_PROMPT,
1878  "Using routine ssht_core_mw_inverse_sov_sym_ss...");
1879  }
1880 
1881  // Compute Fmm.
1882  // Note that m and mm indices are increased in size by one and
1883  // will be filled with zeros by calloc.
1884  Fmm = (SSHT_COMPLEX(double)*)calloc((2*L)*(2*L), sizeof(SSHT_COMPLEX(double)));
1886  Fmm_offset = L-1;
1887  Fmm_stride = 2*L;
1890  if (dl_method == SSHT_DL_RISBO) {
1893  }
1894  dl_offset = ssht_dl_get_offset(L, SSHT_DL_QUARTER);
1895  dl_stride = ssht_dl_get_stride(L, SSHT_DL_QUARTER);
1896  inds_offset = L-1;
1897  for (el=MAX(L0, abs(spin)); el<=L-1; el++) {
1898 
1899  // Compute Wigner plane.
1900  switch (dl_method) {
1901 
1902  case SSHT_DL_RISBO:
1903  if (el!=0 && el==MAX(L0, abs(spin))) {
1904  for(eltmp=0; eltmp<=MAX(L0, abs(spin)); eltmp++)
1907  eltmp, sqrt_tbl, signs);
1909  dl8, L,
1912  el,
1913  signs);
1914  }
1915  else {
1918  el, sqrt_tbl, signs);
1920  dl8, L,
1923  el,
1924  signs);
1925  }
1926  break;
1927 
1928  case SSHT_DL_TRAPANI:
1929  if (el!=0 && el==MAX(L0, abs(spin))) {
1930  for(eltmp=0; eltmp<=MAX(L0, abs(spin)); eltmp++)
1933  eltmp, sqrt_tbl);
1936  el, signs);
1937  }
1938  else {
1941  el, sqrt_tbl);
1944  el, signs);
1945  }
1946  break;
1947 
1948  default:
1949  SSHT_ERROR_GENERIC("Invalid dl method")
1950  }
1951 
1952  // Compute Fmm.
1953  elfactor = sqrt((double)(2.0*el+1.0)/(4.0*SSHT_PI));
1954  el2pel = el *el + el;
1955  for (m=-el; m<=el; m++)
1956  inds[m + inds_offset] = el2pel + m;
1957  for (mm=0; mm<=el; mm++) {
1958  elmmsign = signs[el] * signs[mm];
1959  elssign = spin <= 0 ? 1.0 : elmmsign;
1960 
1961  for (m=-el; m<=-1; m++) {
1962  ind = inds[m + inds_offset];
1963  Fmm[(mm + Fmm_offset)*Fmm_stride + m + Fmm_offset] +=
1964  ssign
1965  * elfactor
1966  * exps[m + exps_offset]
1967  * elmmsign * dl[mm*dl_stride - m + dl_offset]
1968  * elssign * dl[mm*dl_stride - spinneg + dl_offset]
1969  * flm[ind];
1970  }
1971  for (m=0; m<=el; m++) {
1972  ind = inds[m + inds_offset];
1973  Fmm[(mm + Fmm_offset)*Fmm_stride + m + Fmm_offset] +=
1974  ssign
1975  * elfactor
1976  * exps[m + exps_offset]
1977  * dl[mm*dl_stride + m + dl_offset]
1978  * elssign * dl[mm*dl_stride - spinneg + dl_offset]
1979  * flm[ind];
1980  }
1981 
1982  }
1983 
1984  }
1985 
1986  // Free dl memory.
1987  free(dl);
1988  if (dl_method == SSHT_DL_RISBO)
1989  free(dl8);
1990 
1991  // Use symmetry to compute Fmm for negative mm.
1992  for (mm=-(L-1); mm<=-1; mm++)
1993  for (m=-(L-1); m<=L; m++)
1994  Fmm[(mm + Fmm_offset)*Fmm_stride + m + Fmm_offset] =
1995  signs[abs(m)] * ssign
1996  * Fmm[(-mm + Fmm_offset)*Fmm_stride + m + Fmm_offset];
1997 
1998  // Allocate space for function values.
1999  // Note that t and p indices of fext are increased in size by
2000  // one compared to usual sampling.
2001  fext = (SSHT_COMPLEX(double)*)calloc((2*L)*(2*L), sizeof(SSHT_COMPLEX(double)));
2003  fext_stride = 2*L;
2004 
2005  // Apply spatial shift.
2006  for (mm=0; mm<=L; mm++)
2007  for (m=0; m<=L; m++)
2008  fext[mm*fext_stride + m] =
2009  Fmm[(mm + Fmm_offset)*Fmm_stride + m + Fmm_offset];
2010  for (mm=0; mm<=L; mm++)
2011  for (m=-(L-1); m<=-1; m++)
2012  fext[mm*fext_stride + (m+2*L-1+1)] =
2013  Fmm[(mm + Fmm_offset)*Fmm_stride + m + Fmm_offset];
2014  for (mm=-(L-1); mm<=-1; mm++)
2015  for (m=0; m<=L; m++)
2016  fext[(mm + 2*L-1+1)*fext_stride + m] =
2017  Fmm[(mm + Fmm_offset)*Fmm_stride + m + Fmm_offset];
2018  for (mm=-(L-1); mm<=-1; mm++)
2019  for (m=-(L-1); m<=-1; m++)
2020  fext[(mm+2*L-1+1)*fext_stride + m + 2*L-1+1] =
2021  Fmm[(mm + Fmm_offset)*Fmm_stride + m + Fmm_offset];
2022 
2023  // Perform 2D FFT.
2024  plan = fftw_plan_dft_2d(2*L, 2*L, Fmm, Fmm,
2025  FFTW_BACKWARD, FFTW_ESTIMATE);
2026  fftw_execute_dft(plan, fext, fext);
2027  fftw_destroy_plan(plan);
2028 
2029  // Free Fmm memory.
2030  free(Fmm);
2031 
2032  // Extract f from version of f extended to the torus (fext).
2033  // Note that t and p indices of fext are increased in size by
2034  // one compared to usual sampling.
2035  memcpy(f, fext, (L+1)*(2*L)*sizeof(SSHT_COMPLEX(double)));
2036 
2037  // Free fext memory.
2038  free(fext);
2039 
2040  // Print finished if verbosity set.
2041  if (verbosity > 0)
2042  printf("%s%s", SSHT_PROMPT, "Inverse transform computed!");
2043 
2044  // Free precomputation memory.
2045  free(sqrt_tbl);
2046  free(signs);
2047  free(exps);
2048  free(inds);
2049 
2050 }
2051 
2068 void ssht_core_mw_inverse_sov_sym_ss_real(double *f, const SSHT_COMPLEX(double) *flm,
2069  int L,
2070  ssht_dl_method_t dl_method,
2071  int verbosity) {
2073  0, L,
2074  dl_method,
2075  verbosity);
2076 }
2077 
2095 void ssht_core_mw_lb_inverse_sov_sym_ss_real(double *f, const SSHT_COMPLEX(double) *flm,
2096  int L0, int L,
2097  ssht_dl_method_t dl_method,
2098  int verbosity) {
2099 
2100  int el, m, mm, ind;
2101  //int t, p;
2102  int eltmp;
2103  double *sqrt_tbl, *signs;
2104  int el2pel, inds_offset;
2105  int *inds;
2106  double ssign, elfactor;
2107  double *dl;
2108  double *dl8 = NULL;
2109  int dl_offset, dl_stride;
2110  SSHT_COMPLEX(double) *exps;
2111  int exps_offset;
2112  double elmmsign, elssign;
2113  int spinneg;
2114  SSHT_COMPLEX(double) *Fmm, *Fmm_shift;
2115  double *fext_real;
2116  int Fmm_offset, Fmm_stride;
2117  fftw_plan plan;
2118  int spin = 0;
2119 
2120  // Allocate memory.
2121  sqrt_tbl = (double*)calloc(2*(L-1)+2, sizeof(double));
2122  SSHT_ERROR_MEM_ALLOC_CHECK(sqrt_tbl)
2123  signs = (double*)calloc(L+1, sizeof(double));
2125  exps = (SSHT_COMPLEX(double)*)calloc(2*L-1, sizeof(SSHT_COMPLEX(double)));
2127  inds = (int*)calloc(2*L-1, sizeof(int));
2129 
2130  // Perform precomputations.
2131  for (el=0; el<=2*(L-1)+1; el++)
2132  sqrt_tbl[el] = sqrt((double)el);
2133  for (m=0; m<=L-1; m=m+2) {
2134  signs[m] = 1.0;
2135  signs[m+1] = -1.0;
2136  }
2137  ssign = signs[abs(spin)];
2138  spinneg = spin <= 0 ? spin : -spin;
2139  exps_offset = L-1;
2140  for (m=-(L-1); m<=L-1; m++)
2141  exps[m + exps_offset] = cexp(-I*SSHT_PION2*(m+spin));
2142 
2143  // Print messages depending on verbosity level.
2144  if (verbosity > 0) {
2145  printf("%s%s\n", SSHT_PROMPT,
2146  "Computing inverse transform using MW symmetric sampling with ");
2147  printf("%s%s%d%s%d%s\n", SSHT_PROMPT, "parameters (L,spin,reality) = (",
2148  L, ",", spin, ", TRUE)");
2149  if (verbosity > 1)
2150  printf("%s%s\n", SSHT_PROMPT,
2151  "Using routine ssht_core_mw_inverse_sov_sym_ss_real...");
2152  }
2153 
2154  // Compute Fmm.
2155  // Note that m and mm indices are increased in size by one and
2156  // will be filled with zeros by calloc.
2157  Fmm = (SSHT_COMPLEX(double)*)calloc((2*L)*(L+1), sizeof(SSHT_COMPLEX(double)));
2159  Fmm_offset = L-1;
2160  Fmm_stride = L+1;
2163  if (dl_method == SSHT_DL_RISBO) {
2166  }
2167  dl_offset = ssht_dl_get_offset(L, SSHT_DL_QUARTER);
2168  dl_stride = ssht_dl_get_stride(L, SSHT_DL_QUARTER);
2169  inds_offset = L-1;
2170  for (el=MAX(L0, abs(spin)); el<=L-1; el++) {
2171 
2172  // Compute Wigner plane.
2173  switch (dl_method) {
2174 
2175  case SSHT_DL_RISBO:
2176  if (el!=0 && el==MAX(L0, abs(spin))) {
2177  for(eltmp=0; eltmp<=MAX(L0, abs(spin)); eltmp++)
2180  eltmp, sqrt_tbl, signs);
2182  dl8, L,
2185  el,
2186  signs);
2187  }
2188  else {
2191  el, sqrt_tbl, signs);
2193  dl8, L,
2196  el,
2197  signs);
2198  }
2199  break;
2200 
2201  case SSHT_DL_TRAPANI:
2202  if (el!=0 && el==MAX(L0, abs(spin))) {
2203  for(eltmp=0; eltmp<=MAX(L0, abs(spin)); eltmp++)
2206  eltmp, sqrt_tbl);
2209  el, signs);
2210  }
2211  else {
2214  el, sqrt_tbl);
2217  el, signs);
2218  }
2219  break;
2220 
2221  default:
2222  SSHT_ERROR_GENERIC("Invalid dl method")
2223  }
2224 
2225  // Compute Fmm.
2226  elfactor = sqrt((double)(2.0*el+1.0)/(4.0*SSHT_PI));
2227  el2pel = el *el + el;
2228  for (m=0; m<=el; m++)
2229  inds[m + inds_offset] = el2pel + m;
2230  for (mm=0; mm<=el; mm++) {
2231  elmmsign = signs[el] * signs[mm];
2232  elssign = spin <= 0 ? 1.0 : elmmsign;
2233 
2234  for (m=0; m<=el; m++) {
2235  ind = inds[m + inds_offset];
2236  Fmm[(mm + Fmm_offset)*Fmm_stride + m] +=
2237  ssign
2238  * elfactor
2239  * exps[m + exps_offset]
2240  * dl[mm*dl_stride + m + dl_offset]
2241  * elssign * dl[mm*dl_stride - spinneg + dl_offset]
2242  * flm[ind];
2243  }
2244 
2245  }
2246 
2247  }
2248 
2249  // Free dl memory.
2250  free(dl);
2251  if (dl_method == SSHT_DL_RISBO)
2252  free(dl8);
2253 
2254  // Use symmetry to compute Fmm for negative mm.
2255  for (mm=-(L-1); mm<=-1; mm++)
2256  for (m=0; m<=L; m++)
2257  Fmm[(mm + Fmm_offset)*Fmm_stride + m] =
2258  signs[abs(m)] * ssign
2259  * Fmm[(-mm + Fmm_offset)*Fmm_stride + m];
2260 
2261  // Apply spatial shift.
2262  Fmm_shift = (SSHT_COMPLEX(double)*)calloc((2*L)*(L+1), sizeof(SSHT_COMPLEX(double)));
2263  SSHT_ERROR_MEM_ALLOC_CHECK(Fmm_shift)
2264  for (mm=0; mm<=L-1; mm++)
2265  for (m=0; m<=L-1; m++)
2266  Fmm_shift[mm*Fmm_stride + m] =
2267  Fmm[(mm + Fmm_offset)*Fmm_stride + m];
2268  for (mm=-(L-1); mm<=-1; mm++)
2269  for (m=0; m<=L-1; m++)
2270  Fmm_shift[(mm + 2*L-1+1)*Fmm_stride + m] =
2271  Fmm[(mm + Fmm_offset)*Fmm_stride + m];
2272 
2273  // Allocate space for function values.
2274  // Note that t and p indices of fext are increased in size by
2275  // one compared to usual sampling.
2276  fext_real = (double*)calloc((2*L)*(2*L), sizeof(double));
2277  SSHT_ERROR_MEM_ALLOC_CHECK(fext_real)
2278 
2279  // Perform 2D FFT.
2280  plan = fftw_plan_dft_c2r_2d(2*L, 2*L, Fmm_shift, fext_real,
2281  FFTW_ESTIMATE);
2282  fftw_execute_dft_c2r(plan, Fmm_shift, fext_real);
2283  fftw_destroy_plan(plan);
2284 
2285  // Free Fmm memory.
2286  free(Fmm);
2287  free(Fmm_shift);
2288 
2289  // Extract f from version of f extended to the torus (fext).
2290  // Note that t and p indices of fext are increased in size by
2291  // one compared to usual sampling.
2292  memcpy(f, fext_real, (L+1)*(2*L)*sizeof(double));
2293 
2294  // Free fext memory.
2295  free(fext_real);
2296 
2297  // Print finished if verbosity set.
2298  if (verbosity > 0)
2299  printf("%s%s", SSHT_PROMPT, "Inverse transform computed!");
2300 
2301  // Free precomputation memory.
2302  free(sqrt_tbl);
2303  free(signs);
2304  free(exps);
2305  free(inds);
2306 
2307 }
2308 
2309 
2326 void ssht_core_mwdirect_inverse_ss(SSHT_COMPLEX(double) *f, const SSHT_COMPLEX(double) *flm,
2327  int L, int spin, int verbosity) {
2328 
2329  int t, p, m, el, ind, eltmp;
2330  double *dl;
2331  double *sqrt_tbl;
2332  double theta, phi, elfactor;
2333  int ssign;
2334  int dl_offset, dl_stride, f_stride;
2335 
2336  // Allocate memory.
2337  sqrt_tbl = (double*)calloc(2*(L-1)+2, sizeof(double));
2338  SSHT_ERROR_MEM_ALLOC_CHECK(sqrt_tbl)
2339 
2340  // Perform precomputations.
2341  for (el=0; el<=2*(L-1)+1; el++)
2342  sqrt_tbl[el] = sqrt((double)el);
2343  ssign = spin & 1;
2344  ssign = 1 - ssign - ssign; // (-1)^spin
2345 
2346  // Print messages depending on verbosity level.
2347  if (verbosity > 0) {
2348  printf("%s%s\n", SSHT_PROMPT,
2349  "Computing inverse transform using MW sampling with ");
2350  printf("%s%s%d%s%d%s\n", SSHT_PROMPT, "parameters (L,spin,reality) = (",
2351  L, ",", spin, ", FALSE)");
2352  if (verbosity > 1)
2353  printf("%s%s\n", SSHT_PROMPT,
2354  "Using routine ssht_core_mwdirect_inverse_ss...");
2355  }
2356 
2357  // Initialise f with zeros.
2358  f_stride = 2*L;
2359  for (t=0; t<=L; t++)
2360  for (p=0; p<=2*L-1; p++)
2361  f[t*f_stride + p] = 0.0;
2362 
2363  // Compute inverse transform.
2364  dl = ssht_dl_calloc(L, SSHT_DL_FULL);
2366  dl_offset = ssht_dl_get_offset(L, SSHT_DL_FULL);
2367  dl_stride = ssht_dl_get_stride(L, SSHT_DL_FULL);
2368  for (t=0; t<=L; t++) {
2369  theta = ssht_sampling_mw_ss_t2theta(t, L);
2370  for (el=abs(spin); el<=L-1; el++) {
2371  elfactor = sqrt((double)(2.0*el+1.0)/(4.0*SSHT_PI));
2372  if (el!=0 && el==abs(spin)) {
2373  for(eltmp=0; eltmp<=abs(spin); eltmp++)
2374  ssht_dl_beta_risbo_full_table(dl, theta, L,
2375  SSHT_DL_FULL,
2376  eltmp, sqrt_tbl);
2377  }
2378  else {
2379  ssht_dl_beta_risbo_full_table(dl, theta, L,
2380  SSHT_DL_FULL,
2381  el, sqrt_tbl);
2382  }
2383 
2384  for (m=-el; m<=el; m++) {
2385  ssht_sampling_elm2ind(&ind, el, m);
2386  for (p=0; p<=2*L-1; p++) {
2387  phi = ssht_sampling_mw_ss_p2phi(p, L);
2388  f[t*f_stride + p] +=
2389  ssign
2390  * elfactor
2391  * cexp(I*m*phi)
2392  * dl[(m+dl_offset)*dl_stride - spin + dl_offset]
2393  * flm[ind];
2394  }
2395  }
2396 
2397  }
2398  }
2399 
2400  free(sqrt_tbl);
2401  free(dl);
2402 
2403  // Print finished if verbosity set.
2404  if (verbosity > 0)
2405  printf("%s%s", SSHT_PROMPT, "Inverse transform computed!");
2406 
2407 }
2408 
2426  int L, int spin,
2427  ssht_dl_method_t dl_method,
2428  int verbosity) {
2430  0, L, spin,
2431  dl_method,
2432  verbosity);
2433 }
2434 
2453  int L0, int L, int spin,
2454  ssht_dl_method_t dl_method,
2455  int verbosity) {
2456 
2457  int el, m, mm, ind, t, r;
2458  int eltmp;
2459  double *sqrt_tbl, *signs;
2460  int el2pel, inds_offset;
2461  int *inds;
2462  double ssign, elfactor;
2463  fftw_plan plan, plan_bwd, plan_fwd;
2464  SSHT_COMPLEX(double) *inout;
2465  SSHT_COMPLEX(double) *Fmt, *Fmm, *Gmm;
2466  SSHT_COMPLEX(double) *w, *wr;
2467  SSHT_COMPLEX(double) *Fmm_pad, *tmp_pad;
2468  int f_stride, Fmt_stride, Fmt_offset, Fmm_stride, Fmm_offset;
2469  double *dl;
2470  double *dl8 = NULL;
2471  int dl_offset, dl_stride;
2472  int w_offset;
2473  SSHT_COMPLEX(double) *expsm, *expsmm;
2474  int exps_offset;
2475  int elmmsign, elssign;
2476  int spinneg;
2477  int Gmm_stride, Gmm_offset;
2478 
2479  // Allocate memory.
2480  sqrt_tbl = (double*)calloc(2*(L-1)+2, sizeof(double));
2481  SSHT_ERROR_MEM_ALLOC_CHECK(sqrt_tbl)
2482  signs = (double*)calloc(L+1, sizeof(double));
2484  expsm = (SSHT_COMPLEX(double)*)calloc(2*L-1, sizeof(SSHT_COMPLEX(double)));
2486  expsmm = (SSHT_COMPLEX(double)*)calloc(2*L-1, sizeof(SSHT_COMPLEX(double)));
2488  inds = (int*)calloc(2*L-1, sizeof(int));
2490 
2491  // Perform precomputations.
2492  for (el=0; el<=2*(L-1)+1; el++)
2493  sqrt_tbl[el] = sqrt((double)el);
2494  for (m=0; m<=L-1; m=m+2) {
2495  signs[m] = 1.0;
2496  signs[m+1] = -1.0;
2497  }
2498  ssign = signs[abs(spin)];
2499  spinneg = spin <= 0 ? spin : -spin;
2500  exps_offset = L-1;
2501  for (m=-(L-1); m<=L-1; m++)
2502  expsm[m + exps_offset] = cexp(I*SSHT_PION2*(m+spin));
2503  for (mm=-(L-1); mm<=L-1; mm++)
2504  expsmm[mm + exps_offset] = cexp(-I*mm*SSHT_PI/(2.0*L-1.0));
2505 
2506  // Print messages depending on verbosity level.
2507  if (verbosity > 0) {
2508  printf("%s%s\n", SSHT_PROMPT,
2509  "Computing forward transform using MW symmetric sampling with ");
2510  printf("%s%s%d%s%d%s\n", SSHT_PROMPT, "parameters (L,spin,reality) = (",
2511  L, ",", spin, ", FALSE)");
2512  if (verbosity > 1)
2513  printf("%s%s\n", SSHT_PROMPT,
2514  "Using routine ssht_core_mw_forward_sov_conv_sym_ss...");
2515  }
2516 
2517  // Compute Fourier transform over phi, i.e. compute Fmt.
2518  // Note that t and p indices of fext are increased in size by
2519  // one compared to usual sampling.
2520  Fmt = (SSHT_COMPLEX(double)*)calloc((2*L)*(2*L), sizeof(SSHT_COMPLEX(double)));
2522  Fmt_stride = 2*L;
2523  Fmt_offset = L-1;
2524  f_stride = 2*L;
2525  inout = (SSHT_COMPLEX(double)*)calloc(2*L, sizeof(SSHT_COMPLEX(double)));
2527  plan = fftw_plan_dft_1d(2*L, inout, inout, FFTW_FORWARD, FFTW_MEASURE);
2528  for (t=0; t<=L; t++) {
2529  memcpy(inout, &f[t*f_stride], f_stride*sizeof(SSHT_COMPLEX(double)));
2530  fftw_execute_dft(plan, inout, inout);
2531  for(m=0; m<=L; m++)
2532  Fmt[(m+Fmt_offset)*Fmt_stride + t] = inout[m] / (2.0*L);
2533  for(m=-(L-1); m<=-1; m++)
2534  Fmt[(m+Fmt_offset)*Fmt_stride + t] = inout[m+2*L-1+1] / (2.0*L);
2535  }
2536  fftw_destroy_plan(plan);
2537  free(inout);
2538 
2539  // Extend Fmt periodically.
2540  for (m=-(L-1); m<=L; m++)
2541  for (t=L+1; t<=2*L-1; t++)
2542  Fmt[(m+Fmt_offset)*Fmt_stride + t] =
2543  signs[abs(m)] * ssign * Fmt[(m+Fmt_offset)*Fmt_stride + (2*L-t)];
2544 
2545  // Compute Fourier transform over theta, i.e. compute Fmm.
2546  // Note that m and mm indices are increased in size by one.
2547  Fmm = (SSHT_COMPLEX(double)*)calloc((2*L)*(2*L), sizeof(SSHT_COMPLEX(double)));
2549  Fmm_stride = 2*L;
2550  Fmm_offset = L-1;
2551  inout = (SSHT_COMPLEX(double)*)calloc(2*L, sizeof(SSHT_COMPLEX(double)));
2553  plan = fftw_plan_dft_1d(2*L, inout, inout, FFTW_FORWARD, FFTW_MEASURE);
2554  for (m=-(L-1); m<=L; m++) {
2555  memcpy(inout, &Fmt[(m+Fmt_offset)*Fmt_stride], Fmt_stride*sizeof(SSHT_COMPLEX(double)));
2556  fftw_execute_dft(plan, inout, inout);
2557  for(mm=0; mm<=L; mm++)
2558  Fmm[(m+Fmm_offset)*Fmm_stride + mm + Fmm_offset] =
2559  inout[mm] / (2.0*L);
2560  for(mm=-(L-1); mm<=-1; mm++)
2561  Fmm[(m+Fmm_offset)*Fmm_stride + mm + Fmm_offset] =
2562  inout[mm+2*L-1+1] / (2.0*L);
2563  }
2564  fftw_destroy_plan(plan);
2565  free(inout);
2566 
2567  // Compute weights.
2568  w = (SSHT_COMPLEX(double)*)calloc(4*L-3, sizeof(SSHT_COMPLEX(double)));
2570  w_offset = 2*(L-1);
2571  for (mm=-2*(L-1); mm<=2*(L-1); mm++)
2572  w[mm+w_offset] = ssht_sampling_weight_mw(mm);
2573 
2574  // Compute IFFT of w to give wr.
2575  wr = (SSHT_COMPLEX(double)*)calloc(4*L-3, sizeof(SSHT_COMPLEX(double)));
2577  inout = (SSHT_COMPLEX(double)*)calloc(4*L-3, sizeof(SSHT_COMPLEX(double)));
2579  plan_bwd = fftw_plan_dft_1d(4*L-3, inout, inout, FFTW_BACKWARD, FFTW_MEASURE);
2580  plan_fwd = fftw_plan_dft_1d(4*L-3, inout, inout, FFTW_FORWARD, FFTW_MEASURE);
2581  for (mm=1; mm<=2*L-2; mm++)
2582  inout[mm + w_offset] = w[mm - 2*(L-1) - 1 + w_offset];
2583  for (mm=-2*(L-1); mm<=0; mm++)
2584  inout[mm + w_offset] = w[mm + 2*(L-1) + w_offset];
2585  fftw_execute_dft(plan_bwd, inout, inout);
2586  for (mm=0; mm<=2*L-2; mm++)
2587  wr[mm + w_offset] = inout[mm - 2*(L-1) + w_offset];
2588  for (mm=-2*(L-1); mm<=-1; mm++)
2589  wr[mm + w_offset] = inout[mm + 2*(L-1) + 1 + w_offset];
2590 
2591  // Compute Gmm by convolution implemented as product in real space.
2592  Fmm_pad = (SSHT_COMPLEX(double)*)calloc(4*L-3, sizeof(SSHT_COMPLEX(double)));
2594  tmp_pad = (SSHT_COMPLEX(double)*)calloc(4*L-3, sizeof(SSHT_COMPLEX(double)));
2596  Gmm = (SSHT_COMPLEX(double)*)calloc((2*L-1)*(2*L-1), sizeof(SSHT_COMPLEX(double)));
2598  Gmm_stride = 2*L-1;
2599  Gmm_offset = L-1;
2600  for (m=-(L-1); m<=L-1; m++) {
2601 
2602  // Zero-pad Fmm.
2603  for (mm=-2*(L-1); mm<=-L; mm++)
2604  Fmm_pad[mm+w_offset] = 0.0;
2605  for (mm=L; mm<=2*(L-1); mm++)
2606  Fmm_pad[mm+w_offset] = 0.0;
2607  for (mm=-(L-1); mm<=L-1; mm++)
2608  Fmm_pad[mm + w_offset] =
2609  Fmm[(m+Fmm_offset)*Fmm_stride + mm + Fmm_offset];
2610 
2611  // Compute IFFT of Fmm.
2612  for (mm=1; mm<=2*L-2; mm++)
2613  inout[mm + w_offset] = Fmm_pad[mm - 2*(L-1) - 1 + w_offset];
2614  for (mm=-2*(L-1); mm<=0; mm++)
2615  inout[mm + w_offset] = Fmm_pad[mm + 2*(L-1) + w_offset];
2616  fftw_execute_dft(plan_bwd, inout, inout);
2617  for (mm=0; mm<=2*L-2; mm++)
2618  Fmm_pad[mm + w_offset] = inout[mm - 2*(L-1) + w_offset];
2619  for (mm=-2*(L-1); mm<=-1; mm++)
2620  Fmm_pad[mm + w_offset] = inout[mm + 2*(L-1) + 1 + w_offset];
2621 
2622  // Compute product of Fmm and weight in real space.
2623  for (r=-2*(L-1); r<=2*(L-1); r++)
2624  Fmm_pad[r + w_offset] *= wr[-r + w_offset];
2625 
2626  // Compute Gmm by FFT.
2627  for (mm=1; mm<=2*L-2; mm++)
2628  inout[mm + w_offset] = Fmm_pad[mm - 2*(L-1) - 1 + w_offset];
2629  for (mm=-2*(L-1); mm<=0; mm++)
2630  inout[mm + w_offset] = Fmm_pad[mm + 2*(L-1) + w_offset];
2631  fftw_execute_dft(plan_fwd, inout, inout);
2632  for (mm=0; mm<=2*L-2; mm++)
2633  Fmm_pad[mm + w_offset] = inout[mm - 2*(L-1) + w_offset];
2634  for (mm=-2*(L-1); mm<=-1; mm++)
2635  Fmm_pad[mm + w_offset] = inout[mm + 2*(L-1) + 1 + w_offset];
2636 
2637  // Extract section of Gmm of interest.
2638  for (mm=-(L-1); mm<=L-1; mm++)
2639  Gmm[(mm+Gmm_offset)*Gmm_stride + m + Gmm_offset] =
2640  Fmm_pad[mm + w_offset] * 2.0 * SSHT_PI / (4.0*L-3.0);
2641 
2642  }
2643  fftw_destroy_plan(plan_bwd);
2644  fftw_destroy_plan(plan_fwd);
2645 
2646  // Compute flm.
2649  if (dl_method == SSHT_DL_RISBO) {
2652  }
2653  dl_offset = ssht_dl_get_offset(L, SSHT_DL_QUARTER);
2654  dl_stride = ssht_dl_get_stride(L, SSHT_DL_QUARTER);
2655  inds_offset = L-1;
2656  for (el=0; el<=L-1; el++) {
2657  for (m=-el; m<=el; m++) {
2658  ssht_sampling_elm2ind(&ind, el, m);
2659  flm[ind] = 0.0;
2660  }
2661  }
2662  for (el=MAX(L0, abs(spin)); el<=L-1; el++) {
2663 
2664  // Compute Wigner plane.
2665  switch (dl_method) {
2666 
2667  case SSHT_DL_RISBO:
2668  if (el!=0 && el==MAX(L0, abs(spin))) {
2669  for(eltmp=0; eltmp<=MAX(L0, abs(spin)); eltmp++)
2672  eltmp, sqrt_tbl, signs);
2674  dl8, L,
2677  el,
2678  signs);
2679  }
2680  else {
2683  el, sqrt_tbl, signs);
2685  dl8, L,
2688  el,
2689  signs);
2690  }
2691  break;
2692 
2693  case SSHT_DL_TRAPANI:
2694  if (el!=0 && el==MAX(L0, abs(spin))) {
2695  for(eltmp=0; eltmp<=MAX(L0, abs(spin)); eltmp++)
2698  eltmp, sqrt_tbl);
2701  el, signs);
2702  }
2703  else {
2706  el, sqrt_tbl);
2709  el, signs);
2710  }
2711  break;
2712 
2713  default:
2714  SSHT_ERROR_GENERIC("Invalid dl method")
2715  }
2716 
2717  // Compute flm.
2718  elfactor = sqrt((double)(2.0*el+1.0)/(4.0*SSHT_PI));
2719  el2pel = el *el + el;
2720  for (m=-el; m<=el; m++)
2721  inds[m + inds_offset] = el2pel + m;
2722  elssign = spin <= 0 ? 1.0 : signs[el];
2723 
2724  for (m=-el; m<=-1; m++) {
2725  // mm = 0
2726  ind = inds[m + inds_offset];
2727  flm[ind] +=
2728  ssign
2729  * elfactor
2730  * expsm[m + exps_offset]
2731  * signs[el] * dl[0*dl_stride - m + dl_offset]
2732  * elssign * dl[0*dl_stride - spinneg + dl_offset]
2733  * Gmm[(0+Gmm_offset)*Gmm_stride + m + Gmm_offset];
2734  }
2735  for (m=0; m<=el; m++) {
2736  // mm = 0
2737  ind = inds[m + inds_offset];
2738  flm[ind] +=
2739  ssign
2740  * elfactor
2741  * expsm[m + exps_offset]
2742  * dl[0*dl_stride + m + dl_offset]
2743  * elssign * dl[0*dl_stride - spinneg + dl_offset]
2744  * Gmm[(0+Gmm_offset)*Gmm_stride + m + Gmm_offset];
2745  }
2746 
2747  for (mm=1; mm<=el; mm++) {
2748  elmmsign = signs[el] * signs[mm];
2749  elssign = spin <= 0 ? 1.0 : elmmsign;
2750 
2751  for (m=-el; m<=-1; m++) {
2752  ind = inds[m + inds_offset];
2753  flm[ind] +=
2754  ssign
2755  * elfactor
2756  * expsm[m + exps_offset]
2757  * elmmsign * dl[mm*dl_stride - m + dl_offset]
2758  * elssign * dl[mm*dl_stride - spinneg + dl_offset]
2759  * ( Gmm[(mm+Gmm_offset)*Gmm_stride + m + Gmm_offset]
2760  + signs[-m] * ssign
2761  * Gmm[(-mm+Gmm_offset)*Gmm_stride + m + Gmm_offset]);
2762  }
2763  for (m=0; m<=el; m++) {
2764  ind = inds[m + inds_offset];
2765  flm[ind] +=
2766  ssign
2767  * elfactor
2768  * expsm[m + exps_offset]
2769  * dl[mm*dl_stride + m + dl_offset]
2770  * elssign * dl[mm*dl_stride - spinneg + dl_offset]
2771  * ( Gmm[(mm+Gmm_offset)*Gmm_stride + m + Gmm_offset]
2772  + signs[m] * ssign
2773  * Gmm[(-mm+Gmm_offset)*Gmm_stride + m + Gmm_offset]);
2774  }
2775 
2776  }
2777 
2778  }
2779 
2780  // Free memory.
2781  free(dl);
2782  if (dl_method == SSHT_DL_RISBO)
2783  free(dl8);
2784  free(Fmt);
2785  free(Fmm);
2786  free(inout);
2787  free(w);
2788  free(wr);
2789  free(Fmm_pad);
2790  free(tmp_pad);
2791  free(Gmm);
2792  free(sqrt_tbl);
2793  free(signs);
2794  free(expsm);
2795  free(expsmm);
2796  free(inds);
2797 
2798  // Print finished if verbosity set.
2799  if (verbosity > 0)
2800  printf("%s%s", SSHT_PROMPT, "Forward transform computed!");
2801 
2802 }
2803 
2820 void ssht_core_mw_forward_sov_conv_sym_ss_real(SSHT_COMPLEX(double) *flm, const double *f,
2821  int L,
2822  ssht_dl_method_t dl_method,
2823  int verbosity) {
2825  0, L,
2826  dl_method,
2827  verbosity);
2828 }
2829 
2848  int L0, int L,
2849  ssht_dl_method_t dl_method,
2850  int verbosity) {
2851 
2852  int el, m, mm, ind, ind_nm, t, r;
2853  int eltmp;
2854  double *sqrt_tbl, *signs;
2855  int el2pel, inds_offset;
2856  int *inds;
2857  double ssign, elfactor;
2858  fftw_plan plan, plan_bwd, plan_fwd;
2859  double *in_real;
2860  SSHT_COMPLEX(double) *inout, *out;
2861  SSHT_COMPLEX(double) *Fmt, *Fmm, *Gmm;
2862  SSHT_COMPLEX(double) *w, *wr;
2863  SSHT_COMPLEX(double) *Fmm_pad, *tmp_pad;
2864  int f_stride, Fmt_stride, Fmt_offset, Fmm_stride, Fmm_offset, Gmm_stride;
2865  double *dl;
2866  double *dl8 = NULL;
2867  int dl_offset, dl_stride;
2868  int w_offset;
2869  SSHT_COMPLEX(double) *expsm;
2870  int exps_offset;
2871  int elmmsign, elssign;
2872  int spinneg;
2873  int spin = 0;
2874 
2875  // Allocate memory.
2876  sqrt_tbl = (double*)calloc(2*(L-1)+2, sizeof(double));
2877  SSHT_ERROR_MEM_ALLOC_CHECK(sqrt_tbl)
2878  signs = (double*)calloc(L+1, sizeof(double));
2880  expsm = (SSHT_COMPLEX(double)*)calloc(L, sizeof(SSHT_COMPLEX(double)));
2882  inds = (int*)calloc(L, sizeof(int));
2884 
2885  // Perform precomputations.
2886  for (el=0; el<=2*(L-1)+1; el++)
2887  sqrt_tbl[el] = sqrt((double)el);
2888  for (m=0; m<=L-1; m=m+2) {
2889  signs[m] = 1.0;
2890  signs[m+1] = -1.0;
2891  }
2892  ssign = signs[abs(spin)];
2893  spinneg = spin <= 0 ? spin : -spin;
2894  for (m=0; m<=L-1; m++)
2895  expsm[m] = cexp(I*SSHT_PION2*(m+spin));
2896 
2897  // Print messages depending on verbosity level.
2898  if (verbosity > 0) {
2899  printf("%s%s\n", SSHT_PROMPT,
2900  "Computing forward transform using MW symmetric sampling with ");
2901  printf("%s%s%d%s%d%s\n", SSHT_PROMPT, "parameters (L,spin,reality) = (",
2902  L, ",", spin, ", TRUE)");
2903  if (verbosity > 1)
2904  printf("%s%s\n", SSHT_PROMPT,
2905  "Using routine ssht_core_mw_forward_sov_conv_sym_ss_real...");
2906  }
2907 
2908  // Compute Fourier transform over phi, i.e. compute Fmt.
2909  // Note that t and p indices of fext are increased in size by
2910  // one compared to usual sampling.
2911  Fmt = (SSHT_COMPLEX(double)*)calloc((L+1)*(2*L), sizeof(SSHT_COMPLEX(double)));
2913  Fmt_stride = 2*L;
2914  f_stride = 2*L;
2915  in_real = (double*)calloc(2*L, sizeof(double));
2917  out = (SSHT_COMPLEX(double)*)calloc(L+1, sizeof(SSHT_COMPLEX(double)));
2919  plan = fftw_plan_dft_r2c_1d(2*L, in_real, out, FFTW_MEASURE);
2920  for (t=0; t<=L; t++) {
2921  memcpy(in_real, &f[t*f_stride], f_stride*sizeof(double));
2922  fftw_execute_dft_r2c(plan, in_real, out);
2923  for(m=0; m<=L; m++)
2924  Fmt[m*Fmt_stride + t] = out[m] / (2.0*L);
2925 
2926  }
2927  free(in_real);
2928  free(out);
2929  fftw_destroy_plan(plan);
2930 
2931  // Extend Fmt periodically.
2932  for (m=0; m<=L; m++)
2933  for (t=L+1; t<=2*L-1; t++)
2934  Fmt[m*Fmt_stride + t] =
2935  signs[abs(m)] * ssign * Fmt[m*Fmt_stride + (2*L-t)];
2936 
2937  // Compute Fourier transform over theta, i.e. compute Fmm.
2938  // Note that m and mm indices are increased in size by one.
2939  Fmm = (SSHT_COMPLEX(double)*)calloc((L+1)*(2*L), sizeof(SSHT_COMPLEX(double)));
2941  Fmm_stride = 2*L;
2942  Fmm_offset = L-1;
2943  inout = (SSHT_COMPLEX(double)*)calloc(2*L, sizeof(SSHT_COMPLEX(double)));
2945  plan = fftw_plan_dft_1d(2*L, inout, inout, FFTW_FORWARD, FFTW_MEASURE);
2946  for (m=0; m<=L; m++) {
2947  memcpy(inout, &Fmt[m*Fmt_stride], Fmt_stride*sizeof(SSHT_COMPLEX(double)));
2948  fftw_execute_dft(plan, inout, inout);
2949  for(mm=0; mm<=L; mm++)
2950  Fmm[m*Fmm_stride + mm + Fmm_offset] =
2951  inout[mm] / (2.0*L);
2952  for(mm=-(L-1); mm<=-1; mm++)
2953  Fmm[m*Fmm_stride + mm + Fmm_offset] =
2954  inout[mm+2*L-1+1] / (2.0*L);
2955  }
2956  fftw_destroy_plan(plan);
2957  free(inout);
2958 
2959  // Compute weights.
2960  w = (SSHT_COMPLEX(double)*)calloc(4*L-3, sizeof(SSHT_COMPLEX(double)));
2962  w_offset = 2*(L-1);
2963  for (mm=-2*(L-1); mm<=2*(L-1); mm++)
2964  w[mm+w_offset] = ssht_sampling_weight_mw(mm);
2965 
2966  // Compute IFFT of w to give wr.
2967  wr = (SSHT_COMPLEX(double)*)calloc(4*L-3, sizeof(SSHT_COMPLEX(double)));
2969  inout = (SSHT_COMPLEX(double)*)calloc(4*L-3, sizeof(SSHT_COMPLEX(double)));
2971  plan_bwd = fftw_plan_dft_1d(4*L-3, inout, inout, FFTW_BACKWARD, FFTW_MEASURE);
2972  plan_fwd = fftw_plan_dft_1d(4*L-3, inout, inout, FFTW_FORWARD, FFTW_MEASURE);
2973  for (mm=1; mm<=2*L-2; mm++)
2974  inout[mm + w_offset] = w[mm - 2*(L-1) - 1 + w_offset];
2975  for (mm=-2*(L-1); mm<=0; mm++)
2976  inout[mm + w_offset] = w[mm + 2*(L-1) + w_offset];
2977  fftw_execute_dft(plan_bwd, inout, inout);
2978  for (mm=0; mm<=2*L-2; mm++)
2979  wr[mm + w_offset] = inout[mm - 2*(L-1) + w_offset];
2980  for (mm=-2*(L-1); mm<=-1; mm++)
2981  wr[mm + w_offset] = inout[mm + 2*(L-1) + 1 + w_offset];
2982 
2983  // Compute Gmm by convolution implemented as product in real space.
2984  Fmm_pad = (SSHT_COMPLEX(double)*)calloc(4*L-3, sizeof(SSHT_COMPLEX(double)));
2986  tmp_pad = (SSHT_COMPLEX(double)*)calloc(4*L-3, sizeof(SSHT_COMPLEX(double)));
2988  Gmm = (SSHT_COMPLEX(double)*)calloc((2*L-1)*L, sizeof(SSHT_COMPLEX(double)));
2989  Gmm_stride = L;
2991  for (m=0; m<=L-1; m++) {
2992 
2993  // Zero-pad Fmm.
2994  for (mm=-2*(L-1); mm<=-L; mm++)
2995  Fmm_pad[mm+w_offset] = 0.0;
2996  for (mm=L; mm<=2*(L-1); mm++)
2997  Fmm_pad[mm+w_offset] = 0.0;
2998  for (mm=-(L-1); mm<=L-1; mm++)
2999  Fmm_pad[mm + w_offset] =
3000  Fmm[m*Fmm_stride + mm + Fmm_offset];
3001 
3002  // Compute IFFT of Fmm.
3003  for (mm=1; mm<=2*L-2; mm++)
3004  inout[mm + w_offset] = Fmm_pad[mm - 2*(L-1) - 1 + w_offset];
3005  for (mm=-2*(L-1); mm<=0; mm++)
3006  inout[mm + w_offset] = Fmm_pad[mm + 2*(L-1) + w_offset];
3007  fftw_execute_dft(plan_bwd, inout, inout);
3008  for (mm=0; mm<=2*L-2; mm++)
3009  Fmm_pad[mm + w_offset] = inout[mm - 2*(L-1) + w_offset];
3010  for (mm=-2*(L-1); mm<=-1; mm++)
3011  Fmm_pad[mm + w_offset] = inout[mm + 2*(L-1) + 1 + w_offset];
3012 
3013  // Compute product of Fmm and weight in real space.
3014  for (r=-2*(L-1); r<=2*(L-1); r++)
3015  Fmm_pad[r + w_offset] *= wr[-r + w_offset];
3016 
3017  // Compute Gmm by FFT.
3018  for (mm=1; mm<=2*L-2; mm++)
3019  inout[mm + w_offset] = Fmm_pad[mm - 2*(L-1) - 1 + w_offset];
3020  for (mm=-2*(L-1); mm<=0; mm++)
3021  inout[mm + w_offset] = Fmm_pad[mm + 2*(L-1) + w_offset];
3022  fftw_execute_dft(plan_fwd, inout, inout);
3023  for (mm=0; mm<=2*L-2; mm++)
3024  Fmm_pad[mm + w_offset] = inout[mm - 2*(L-1) + w_offset];
3025  for (mm=-2*(L-1); mm<=-1; mm++)
3026  Fmm_pad[mm + w_offset] = inout[mm + 2*(L-1) + 1 + w_offset];
3027 
3028  // Extract section of Gmm of interest.
3029  for (mm=-(L-1); mm<=L-1; mm++)
3030  Gmm[(mm+Fmm_offset)*Gmm_stride + m] =
3031  Fmm_pad[mm + w_offset] * 2.0 * SSHT_PI / (4.0*L-3.0);
3032 
3033  }
3034  fftw_destroy_plan(plan_bwd);
3035  fftw_destroy_plan(plan_fwd);
3036 
3037  // Compute flm.
3040  if (dl_method == SSHT_DL_RISBO) {
3043  }
3044  dl_offset = ssht_dl_get_offset(L, SSHT_DL_QUARTER);
3045  dl_stride = ssht_dl_get_stride(L, SSHT_DL_QUARTER);
3046  inds_offset = 0;
3047  for (el=0; el<=L-1; el++) {
3048  for (m=0; m<=el; m++) {
3049  ssht_sampling_elm2ind(&ind, el, m);
3050  flm[ind] = 0.0;
3051  }
3052  }
3053  for (el=MAX(L0, abs(spin)); el<=L-1; el++) {
3054 
3055  // Compute Wigner plane.
3056  switch (dl_method) {
3057 
3058  case SSHT_DL_RISBO:
3059  if (el!=0 && el==MAX(L0, abs(spin))) {
3060  for(eltmp=0; eltmp<=MAX(L0, abs(spin)); eltmp++)
3063  eltmp, sqrt_tbl, signs);
3065  dl8, L,
3068  el,
3069  signs);
3070  }
3071  else {
3074  el, sqrt_tbl, signs);
3076  dl8, L,
3079  el,
3080  signs);
3081  }
3082  break;
3083 
3084  case SSHT_DL_TRAPANI:
3085  if (el!=0 && el==MAX(L0, abs(spin))) {
3086  for(eltmp=0; eltmp<=MAX(L0, abs(spin)); eltmp++)
3089  eltmp, sqrt_tbl);
3092  el, signs);
3093  }
3094  else {
3097  el, sqrt_tbl);
3100  el, signs);
3101  }
3102  break;
3103 
3104  default:
3105  SSHT_ERROR_GENERIC("Invalid dl method")
3106  }
3107 
3108  // Compute flm.
3109  elfactor = sqrt((double)(2.0*el+1.0)/(4.0*SSHT_PI));
3110  el2pel = el *el + el;
3111  for (m=0; m<=el; m++)
3112  inds[m + inds_offset] = el2pel + m;
3113  elssign = spin <= 0 ? 1.0 : signs[el];
3114 
3115  for (m=0; m<=el; m++) {
3116  // mm = 0
3117  ind = inds[m + inds_offset];
3118  flm[ind] +=
3119  ssign
3120  * elfactor
3121  * expsm[m]
3122  * dl[0*dl_stride + m + dl_offset]
3123  * elssign * dl[0*dl_stride - spinneg + dl_offset]
3124  * Gmm[(0+Fmm_offset)*Gmm_stride + m];
3125  }
3126 
3127  for (mm=1; mm<=el; mm++) {
3128  elmmsign = signs[el] * signs[mm];
3129  elssign = spin <= 0 ? 1.0 : elmmsign;
3130 
3131  for (m=0; m<=el; m++) {
3132  ind = inds[m + inds_offset];
3133  flm[ind] +=
3134  ssign
3135  * elfactor
3136  * expsm[m]
3137  * dl[mm*dl_stride + m + dl_offset]
3138  * elssign * dl[mm*dl_stride - spinneg + dl_offset]
3139  * ( Gmm[(mm+Fmm_offset)*Gmm_stride + m]
3140  + signs[m] * ssign
3141  * Gmm[(-mm+Fmm_offset)*Gmm_stride + m]);
3142  }
3143 
3144  }
3145 
3146  }
3147 
3148  // Set flm values for negative m using conjugate symmetry.
3149  for (el=abs(spin); el<=L-1; el++) {
3150  for (m=1; m<=el; m++) {
3151  ssht_sampling_elm2ind(&ind, el, m);
3152  ssht_sampling_elm2ind(&ind_nm, el, -m);
3153  flm[ind_nm] = signs[m] * conj(flm[ind]);
3154  }
3155  }
3156 
3157  // Free memory.
3158  free(dl);
3159  if (dl_method == SSHT_DL_RISBO)
3160  free(dl8);
3161  free(Fmt);
3162  free(Fmm);
3163  free(inout);
3164  free(w);
3165  free(wr);
3166  free(Fmm_pad);
3167  free(tmp_pad);
3168  free(Gmm);
3169  free(sqrt_tbl);
3170  free(signs);
3171  free(expsm);
3172  free(inds);
3173 
3174  // Print finished if verbosity set.
3175  if (verbosity > 0)
3176  printf("%s%s", SSHT_PROMPT, "Forward transform computed!");
3177 
3178 }
3179 
3180 
3181 //============================================================================
3182 // MW SS Noth-South pole interfaces
3183 //============================================================================
3184 
3185 
3210  SSHT_COMPLEX(double) *f_np, double *phi_np,
3211  SSHT_COMPLEX(double) *f_sp, double *phi_sp,
3212  const SSHT_COMPLEX(double) *flm,
3213  int L, int spin,
3214  ssht_dl_method_t dl_method,
3215  int verbosity) {
3216 
3217  SSHT_COMPLEX(double)* f_full;
3218  int t, f_stride = 2*L;
3219 
3220  // Allocate full array.
3221  f_full = (SSHT_COMPLEX(double)*)calloc((L+1)*(2*L), sizeof(SSHT_COMPLEX(double)));
3223 
3224  // Perform inverse transform.
3225  ssht_core_mw_inverse_sov_sym_ss(f_full, flm, L, spin,
3226  dl_method, verbosity);
3227 
3228  // Copy output function values, including separate points for poles.
3229  for (t=1; t<=L-1; t++)
3230  memcpy(&f[(t-1)*f_stride], &f_full[t*f_stride],
3231  (2*L)*sizeof(SSHT_COMPLEX(double)));
3232  *f_np = f_full[0];
3233  *phi_np = ssht_sampling_mw_ss_p2phi(0, L);
3234  *f_sp = f_full[L*f_stride + 0];
3235  *phi_sp = ssht_sampling_mw_ss_p2phi(0, L);
3236 
3237  // Free memory.
3238  free(f_full);
3239 
3240 }
3241 
3242 
3262  double *f_np,
3263  double *f_sp,
3264  const SSHT_COMPLEX(double) *flm,
3265  int L,
3266  ssht_dl_method_t dl_method,
3267  int verbosity) {
3268 
3269  double *f_full;
3270  int t, f_stride = 2*L;
3271 
3272  // Allocate full array.
3273  f_full = (double*)calloc((L+1)*(2*L), sizeof(double));
3275 
3276  // Perform inverse transform.
3277  ssht_core_mw_inverse_sov_sym_ss_real(f_full, flm, L,
3278  dl_method, verbosity);
3279 
3280  // Copy output function values, including separate points for poles.
3281  for (t=1; t<=L-1; t++)
3282  memcpy(&f[(t-1)*f_stride], &f_full[t*f_stride],
3283  (2*L)*sizeof(double));
3284  *f_np = f_full[0];
3285  *f_sp = f_full[L*f_stride + 0];
3286 
3287  // Free memory.
3288  free(f_full);
3289 
3290 }
3291 
3292 
3317  SSHT_COMPLEX(double) f_np, double phi_np,
3318  SSHT_COMPLEX(double) f_sp, double phi_sp,
3319  int L, int spin,
3320  ssht_dl_method_t dl_method,
3321  int verbosity) {
3322 
3323  SSHT_COMPLEX(double) *f_full;
3324  int t, p, f_stride = 2*L;
3325  double phi;
3326 
3327  // Copy function values to full array.
3328  f_full = (SSHT_COMPLEX(double)*)calloc((L+1)*(2*L), sizeof(SSHT_COMPLEX(double)));
3330  for (t=1; t<=L-1; t++)
3331  memcpy(&f_full[t*f_stride], &f[(t-1)*f_stride],
3332  (2*L)*sizeof(SSHT_COMPLEX(double)));
3333 
3334  // Define poles for all phi.
3335  for (p=0; p<=2*L-1; p++) {
3336  phi = ssht_sampling_mw_ss_p2phi(p, L);
3337  f_full[0*f_stride + p] = f_np * cexp(-I*spin*(phi-phi_np));
3338  f_full[L*f_stride + p] = f_sp * cexp(I*spin*(phi-phi_sp));
3339  }
3340 
3341  // Perform forward transform.
3342  ssht_core_mw_forward_sov_conv_sym_ss(flm, f_full, L, spin,
3343  dl_method, verbosity);
3344 
3345  // Free memory.
3346  free(f_full);
3347 
3348 }
3349 
3350 
3370  const double *f,
3371  double f_np,
3372  double f_sp,
3373  int L,
3374  ssht_dl_method_t dl_method,
3375  int verbosity) {
3376 
3377  double *f_full;
3378  int t, p, f_stride = 2*L;
3379 
3380  // Copy function values to full array.
3381  f_full = (double*)calloc((L+1)*(2*L), sizeof(double));
3383  for (t=1; t<=L-1; t++)
3384  memcpy(&f_full[t*f_stride], &f[(t-1)*f_stride],
3385  (2*L)*sizeof(double));
3386 
3387  // Define poles for all phi.
3388  for (p=0; p<=2*L-1; p++) {
3389  f_full[0*f_stride + p] = f_np;
3390  f_full[L*f_stride + p] = f_sp;
3391  }
3392 
3393  // Perform forward transform.
3395  dl_method, verbosity);
3396 
3397  // Free memory.
3398  free(f_full);
3399 
3400 }
3401 
3402 
3403 //============================================================================
3404 // GL algorithms
3405 //============================================================================
3406 
3407 
3421 void ssht_core_gl_inverse_sov(SSHT_COMPLEX(double) *f, const SSHT_COMPLEX(double) *flm,
3422  int L, int spin, int verbosity) {
3423 
3424  int t, p, m, el, ind;
3425  int ftm_stride, ftm_offset, f_stride;
3426  double *dlm1p1_line, *dl_line;
3427  double *dl_ptr;
3428  double *sqrt_tbl, *signs;
3429  SSHT_COMPLEX(double) *ftm, *inout;
3430  double theta, ssign, elfactor;
3431  fftw_plan plan;
3432  double *thetas, *weights;
3433 
3434  // Allocate memory.
3435  sqrt_tbl = (double*)calloc(2*(L-1)+2, sizeof(double));
3436  SSHT_ERROR_MEM_ALLOC_CHECK(sqrt_tbl)
3437  signs = (double*)calloc(L+1, sizeof(double));
3439 
3440  // Perform precomputations.
3441  for (el=0; el<=2*(L-1)+1; el++)
3442  sqrt_tbl[el] = sqrt((double)el);
3443  for (m=0; m<=L-1; m=m+2) {
3444  signs[m] = 1.0;
3445  signs[m+1] = -1.0;
3446  }
3447  ssign = signs[abs(spin)];
3448 
3449  // Print messages depending on verbosity level.
3450  if (verbosity > 0) {
3451  printf("%s%s\n", SSHT_PROMPT,
3452  "Computing inverse transform using GL sampling with ");
3453  printf("%s%s%d%s%d%s\n", SSHT_PROMPT, "parameters (L,spin,reality) = (",
3454  L, ",", spin, ", FALSE)");
3455  if (verbosity > 1)
3456  printf("%s%s\n", SSHT_PROMPT,
3457  "Using routine ssht_core_gl_inverse_sov...");
3458  }
3459 
3460  // Compute weights and theta positions.
3461  thetas = (double*)calloc(L, sizeof(double));
3463  weights = (double*)calloc(L, sizeof(double));
3465  ssht_sampling_gl_thetas_weights(thetas, weights, L);
3466 
3467  // Compute ftm.
3468  ftm = (SSHT_COMPLEX(double)*)calloc(L*(2*L-1), sizeof(SSHT_COMPLEX(double)));
3470  ftm_stride = 2*L-1;
3471  ftm_offset = L-1;
3472  dlm1p1_line = (double*)calloc(2*L-1, sizeof(double));
3473  SSHT_ERROR_MEM_ALLOC_CHECK(dlm1p1_line)
3474  dl_line = (double*)calloc(2*L-1, sizeof(double));
3476  for (t=0; t<=L-1; t++) {
3477  theta = thetas[t];
3478  for (el=abs(spin); el<=L-1; el++) {
3479  elfactor = sqrt((double)(2.0*el+1.0)/(4.0*SSHT_PI));
3480 
3481  // Compute dl line for given spin.
3482  ssht_dl_beta_kostelec_line_table(dlm1p1_line, dl_line,
3483  theta, L, -spin, el,
3484  sqrt_tbl, signs);
3485  // Switch current and previous dls.
3486  dl_ptr = dl_line;
3487  dl_line = dlm1p1_line;
3488  dlm1p1_line = dl_ptr;
3489 
3490  for (m=-el; m<=el; m++) {
3491  ssht_sampling_elm2ind(&ind, el, m);
3492  ftm[t*ftm_stride + m + ftm_offset] +=
3493  ssign
3494  * elfactor
3495  * dl_line[m + L-1]
3496  * flm[ind];
3497  }
3498  }
3499  }
3500 
3501  // Free dl memory.
3502  free(dlm1p1_line);
3503  free(dl_line);
3504 
3505  // Free memory for thetas and weights.
3506  free(thetas);
3507  free(weights);
3508 
3509  // Compute f.
3510  inout = (SSHT_COMPLEX(double)*)calloc(2*L-1, sizeof(SSHT_COMPLEX(double)));
3512  f_stride = 2*L-1;
3513  plan = fftw_plan_dft_1d(2*L-1, inout, inout, FFTW_BACKWARD, FFTW_MEASURE);
3514  for (t=0; t<=L-1; t++) {
3515  for (m=0; m<=L-1; m++)
3516  inout[m] = ftm[t*ftm_stride + m + ftm_offset];
3517  for (m=-(L-1); m<=-1; m++)
3518  inout[m+2*L-1] = ftm[t*ftm_stride + m + ftm_offset];
3519  fftw_execute_dft(plan, inout, inout);
3520  for (p=0; p<=2*L-2; p++)
3521  f[t*f_stride + p] = inout[p];
3522  }
3523  fftw_destroy_plan(plan);
3524 
3525  // Free memory.
3526  free(ftm);
3527  free(inout);
3528  free(signs);
3529  free(sqrt_tbl);
3530 
3531  // Print finished if verbosity set.
3532  if (verbosity > 0)
3533  printf("%s%s", SSHT_PROMPT, "Inverse transform computed!");
3534 
3535 }
3536 
3537 
3552 void ssht_core_gl_inverse_sov_real(double *f, const SSHT_COMPLEX(double) *flm,
3553  int L, int verbosity) {
3554 
3555  int t, p, m, el, ind;
3556  int ftm_stride, ftm_offset, f_stride;
3557  double *dlm1p1_line, *dl_line;
3558  double *dl_ptr;
3559  double *sqrt_tbl, *signs;
3560  SSHT_COMPLEX(double) *ftm;
3561  SSHT_COMPLEX(double) *in;
3562  double *out;
3563  double theta, ssign, elfactor;
3564  fftw_plan plan;
3565  double *thetas, *weights;
3566  int spin = 0;
3567 
3568  // Allocate memory.
3569  sqrt_tbl = (double*)calloc(2*(L-1)+2, sizeof(double));
3570  SSHT_ERROR_MEM_ALLOC_CHECK(sqrt_tbl)
3571  signs = (double*)calloc(L+1, sizeof(double));
3573 
3574  // Perform precomputations.
3575  for (el=0; el<=2*(L-1)+1; el++)
3576  sqrt_tbl[el] = sqrt((double)el);
3577  for (m=0; m<=L-1; m=m+2) {
3578  signs[m] = 1.0;
3579  signs[m+1] = -1.0;
3580  }
3581  ssign = signs[abs(spin)];
3582 
3583  // Print messages depending on verbosity level.
3584  if (verbosity > 0) {
3585  printf("%s%s\n", SSHT_PROMPT,
3586  "Computing inverse transform using GL sampling with ");
3587  printf("%s%s%d%s%d%s\n", SSHT_PROMPT, "parameters (L,spin,reality) = (",
3588  L, ",", spin, ", TRUE)");
3589  if (verbosity > 1)
3590  printf("%s%s\n", SSHT_PROMPT,
3591  "Using routine ssht_core_gl_inverse_sov_real...");
3592  }
3593 
3594  // Compute weights and theta positions.
3595  thetas = (double*)calloc(L, sizeof(double));
3597  weights = (double*)calloc(L, sizeof(double));
3599  ssht_sampling_gl_thetas_weights(thetas, weights, L);
3600 
3601  // Compute ftm.
3602  ftm = (SSHT_COMPLEX(double)*)calloc(L*L, sizeof(SSHT_COMPLEX(double)));
3604  ftm_stride = L;
3605  ftm_offset = 0;
3606  dlm1p1_line = (double*)calloc(L, sizeof(double));
3607  SSHT_ERROR_MEM_ALLOC_CHECK(dlm1p1_line)
3608  dl_line = (double*)calloc(L, sizeof(double));
3610  for (t=0; t<=L-1; t++) {
3611  theta = thetas[t];
3612  for (el=abs(spin); el<=L-1; el++) {
3613  elfactor = sqrt((double)(2.0*el+1.0)/(4.0*SSHT_PI));
3614 
3615  // Compute half dl line for given spin.
3616  ssht_dl_beta_kostelec_halfline_table(dlm1p1_line, dl_line,
3617  theta, L, -spin, el,
3618  sqrt_tbl, signs);
3619  // Switch current and previous dls.
3620  dl_ptr = dl_line;
3621  dl_line = dlm1p1_line;
3622  dlm1p1_line = dl_ptr;
3623 
3624  for (m=0; m<=el; m++) {
3625  ssht_sampling_elm2ind(&ind, el, m);
3626  ftm[t*ftm_stride + m + ftm_offset] +=
3627  ssign
3628  * elfactor
3629  * dl_line[m]
3630  * flm[ind];
3631  }
3632  }
3633  }
3634 
3635  // Free dl memory.
3636  free(dlm1p1_line);
3637  free(dl_line);
3638 
3639  // Free memory for thetas and weights.
3640  free(thetas);
3641  free(weights);
3642 
3643  // Compute f.
3644  in = (SSHT_COMPLEX(double)*)calloc(L, sizeof(SSHT_COMPLEX(double)));
3646  out = (double*)calloc(2*L-1, sizeof(double));
3648  plan = fftw_plan_dft_c2r_1d(2*L-1, in, out, FFTW_MEASURE);
3649  f_stride = 2*L-1;
3650  for (t=0; t<=L-1; t++) {
3651  memcpy(in, &ftm[t*ftm_stride], L*sizeof(SSHT_COMPLEX(double)));
3652  fftw_execute_dft_c2r(plan, in, out);
3653  for (p=0; p<=2*L-2; p++)
3654  f[t*f_stride + p] = out[p];
3655  }
3656  fftw_destroy_plan(plan);
3657 
3658  // Free memory.
3659  free(ftm);
3660  free(in);
3661  free(out);
3662  free(signs);
3663  free(sqrt_tbl);
3664 
3665  // Print finished if verbosity set.
3666  if (verbosity > 0)
3667  printf("%s%s", SSHT_PROMPT, "Inverse transform computed!");
3668 
3669 }
3670 
3671 
3685 void ssht_core_gl_forward_sov(SSHT_COMPLEX(double) *flm, const SSHT_COMPLEX(double) *f,
3686  int L, int spin, int verbosity) {
3687 
3688  int t, m, el, ind;
3689  int f_stride;
3690  double *dlm1p1_line, *dl_line;
3691  double *dl_ptr;
3692  int el2pel, inds_offset;
3693  int *inds;
3694  double *sqrt_tbl, *signs;
3695  int Ftm_stride, Ftm_offset;
3696  SSHT_COMPLEX(double) *Ftm, *inout;
3697  double theta, ssign, elfactor;
3698  fftw_plan plan;
3699  double *thetas, *weights;
3700  double w;
3701 
3702  // Allocate memory.
3703  sqrt_tbl = (double*)calloc(2*(L-1)+2, sizeof(double));
3704  SSHT_ERROR_MEM_ALLOC_CHECK(sqrt_tbl)
3705  signs = (double*)calloc(L+1, sizeof(double));
3707  inds = (int*)calloc(2*L-1, sizeof(int));
3709 
3710  // Perform precomputations.
3711  for (el=0; el<=2*(L-1)+1; el++)
3712  sqrt_tbl[el] = sqrt((double)el);
3713  for (m=0; m<=L-1; m=m+2) {
3714  signs[m] = 1.0;
3715  signs[m+1] = -1.0;
3716  }
3717  ssign = signs[abs(spin)];
3718 
3719  // Print messages depending on verbosity level.
3720  if (verbosity > 0) {
3721  printf("%s%s\n", SSHT_PROMPT,
3722  "Computing forward transform using GL sampling with ");
3723  printf("%s%s%d%s%d%s\n", SSHT_PROMPT, "parameters (L,spin,reality) = (",
3724  L, ",", spin, ", FALSE)");
3725  if (verbosity > 1)
3726  printf("%s%s\n", SSHT_PROMPT,
3727  "Using routine ssht_core_gl_forward_sov...");
3728  }
3729 
3730  // Compute weights and theta positions.
3731  thetas = (double*)calloc(L, sizeof(double));
3733  weights = (double*)calloc(L, sizeof(double));
3735  ssht_sampling_gl_thetas_weights(thetas, weights, L);
3736 
3737  // Compute Fourier transform over phi, i.e. compute Ftm.
3738  Ftm = (SSHT_COMPLEX(double)*)calloc(L*(2*L-1), sizeof(SSHT_COMPLEX(double)));
3740  Ftm_stride = 2*L-1;
3741  Ftm_offset = L-1;
3742  f_stride = 2*L-1;
3743  inout = (SSHT_COMPLEX(double)*)calloc(2*L-1, sizeof(SSHT_COMPLEX(double)));
3745  plan = fftw_plan_dft_1d(2*L-1, inout, inout, FFTW_FORWARD, FFTW_MEASURE);
3746  for (t=0; t<=L-1; t++) {
3747  memcpy(inout, &f[t*f_stride], f_stride*sizeof(SSHT_COMPLEX(double)));
3748  fftw_execute_dft(plan, inout, inout);
3749  for(m=0; m<=L-1; m++)
3750  Ftm[t*Ftm_stride + m + Ftm_offset] =
3751  inout[m] * 2.0 * SSHT_PI / (2.0*L-1.0) ;
3752  for(m=-(L-1); m<=-1; m++)
3753  Ftm[t*Ftm_stride + m + Ftm_offset] =
3754  inout[m+2*L-1] * 2.0 * SSHT_PI / (2.0*L-1.0);
3755  }
3756  fftw_destroy_plan(plan);
3757 
3758  // Compute flm.
3759  dlm1p1_line = (double*)calloc(2*L-1, sizeof(double));
3760  SSHT_ERROR_MEM_ALLOC_CHECK(dlm1p1_line)
3761  dl_line = (double*)calloc(2*L-1, sizeof(double));
3763  inds_offset = L-1;
3764  for (el=0; el<=L-1; el++) {
3765  for (m=-el; m<=el; m++) {
3766  ssht_sampling_elm2ind(&ind, el, m);
3767  flm[ind] = 0.0;
3768  }
3769  }
3770  for (t=0; t<=L-1; t++) {
3771  theta = thetas[t];
3772  w = weights[t];
3773  for (el=abs(spin); el<=L-1; el++) {
3774  elfactor = sqrt((double)(2.0*el+1.0)/(4.0*SSHT_PI));
3775  el2pel = el *el + el;
3776  for (m=-el; m<=el; m++)
3777  inds[m + inds_offset] = el2pel + m;
3778 
3779  // Compute dl line for given spin.
3780  ssht_dl_beta_kostelec_line_table(dlm1p1_line, dl_line,
3781  theta, L, -spin, el,
3782  sqrt_tbl, signs);
3783  // Switch current and previous dls.
3784  dl_ptr = dl_line;
3785  dl_line = dlm1p1_line;
3786  dlm1p1_line = dl_ptr;
3787 
3788  for (m=-el; m<=el; m++) {
3789  ind = inds[m + inds_offset];
3790  flm[ind] +=
3791  ssign
3792  * elfactor
3793  * w
3794  * dl_line[m+L-1]
3795  * Ftm[t*Ftm_stride + m + Ftm_offset];
3796  }
3797  }
3798  }
3799 
3800  // Free memory.
3801  free(dlm1p1_line);
3802  free(dl_line);
3803  free(thetas);
3804  free(weights);
3805  free(Ftm);
3806  free(inout);
3807  free(signs);
3808  free(sqrt_tbl);
3809  free(inds);
3810 
3811  // Print finished if verbosity set.
3812  if (verbosity > 0)
3813  printf("%s%s", SSHT_PROMPT, "Forward transform computed!");
3814 
3815 }
3816 
3817 
3832 void ssht_core_gl_forward_sov_real(SSHT_COMPLEX(double) *flm, const double *f,
3833  int L, int verbosity) {
3834 
3835  int t, m, el, ind, ind_nm;
3836  int f_stride;
3837  double *dlm1p1_line, *dl_line;
3838  double *dl_ptr;
3839  int el2pel, inds_offset;
3840  int *inds;
3841  double *sqrt_tbl, *signs;
3842  int Ftm_stride, Ftm_offset;
3843  SSHT_COMPLEX(double) *Ftm;
3844  double *in_real;
3845  SSHT_COMPLEX(double) *out;
3846  double theta, ssign, elfactor;
3847  fftw_plan plan;
3848  double *thetas, *weights;
3849  double w;
3850  int spin = 0;
3851 
3852  // Allocate memory.
3853  sqrt_tbl = (double*)calloc(2*(L-1)+2, sizeof(double));
3854  SSHT_ERROR_MEM_ALLOC_CHECK(sqrt_tbl)
3855  signs = (double*)calloc(L+1, sizeof(double));
3857  inds = (int*)calloc(L, sizeof(int));
3859 
3860  // Perform precomputations.
3861  for (el=0; el<=2*(L-1)+1; el++)
3862  sqrt_tbl[el] = sqrt((double)el);
3863  for (m=0; m<=L-1; m=m+2) {
3864  signs[m] = 1.0;
3865  signs[m+1] = -1.0;
3866  }
3867  ssign = signs[abs(spin)];
3868 
3869  // Print messages depending on verbosity level.
3870  if (verbosity > 0) {
3871  printf("%s%s\n", SSHT_PROMPT,
3872  "Computing forward transform using GL sampling with ");
3873  printf("%s%s%d%s%d%s\n", SSHT_PROMPT, "parameters (L,spin,reality) = (",
3874  L, ",", spin, ", TRUE)");
3875  if (verbosity > 1)
3876  printf("%s%s\n", SSHT_PROMPT,
3877  "Using routine ssht_core_gl_forward_sov_real...");
3878  }
3879 
3880  // Compute weights and theta positions.
3881  thetas = (double*)calloc(L, sizeof(double));
3883  weights = (double*)calloc(L, sizeof(double));
3885  ssht_sampling_gl_thetas_weights(thetas, weights, L);
3886 
3887  // Compute Fourier transform over phi, i.e. compute Ftm.
3888  Ftm = (SSHT_COMPLEX(double)*)calloc(L*L, sizeof(SSHT_COMPLEX(double)));
3890  Ftm_stride = L;
3891  Ftm_offset = 0;
3892  f_stride = 2*L-1;
3893  in_real = (double*)calloc(2*L-1, sizeof(double));
3895  out = (SSHT_COMPLEX(double)*)calloc(L, sizeof(SSHT_COMPLEX(double)));
3897  plan = fftw_plan_dft_r2c_1d(2*L-1, in_real, out, FFTW_MEASURE);
3898  for (t=0; t<=L-1; t++) {
3899  memcpy(in_real, &f[t*f_stride], f_stride*sizeof(double));
3900  fftw_execute_dft_r2c(plan, in_real, out);
3901  for(m=0; m<=L-1; m++)
3902  Ftm[t*Ftm_stride + m + Ftm_offset] =
3903  out[m] * 2.0 * SSHT_PI / (2.0*L-1.0);
3904  }
3905  free(in_real);
3906  free(out);
3907  fftw_destroy_plan(plan);
3908 
3909  // Compute flm.
3910  dlm1p1_line = (double*)calloc(L, sizeof(double));
3911  SSHT_ERROR_MEM_ALLOC_CHECK(dlm1p1_line)
3912  dl_line = (double*)calloc(L, sizeof(double));
3914  inds_offset = 0;
3915  for (el=0; el<=L-1; el++) {
3916  for (m=0; m<=el; m++) {
3917  ssht_sampling_elm2ind(&ind, el, m);
3918  flm[ind] = 0.0;
3919  }
3920  }
3921  for (t=0; t<=L-1; t++) {
3922  theta = thetas[t];
3923  w = weights[t];
3924  for (el=abs(spin); el<=L-1; el++) {
3925  elfactor = sqrt((double)(2.0*el+1.0)/(4.0*SSHT_PI));
3926  el2pel = el *el + el;
3927  for (m=0; m<=el; m++)
3928  inds[m + inds_offset] = el2pel + m;
3929 
3930  // Compute half dl line for given spin.
3931  ssht_dl_beta_kostelec_halfline_table(dlm1p1_line, dl_line,
3932  theta, L, -spin, el,
3933  sqrt_tbl, signs);
3934 
3935  // Switch current and previous dls.
3936  dl_ptr = dl_line;
3937  dl_line = dlm1p1_line;
3938  dlm1p1_line = dl_ptr;
3939 
3940  for (m=0; m<=el; m++) {
3941  ind = inds[m + inds_offset];
3942  flm[ind] +=
3943  ssign
3944  * elfactor
3945  * w
3946  * dl_line[m]
3947  * Ftm[t*Ftm_stride + m + Ftm_offset];
3948  }
3949  }
3950  }
3951 
3952  // Set flm values for negative m using conjugate symmetry.
3953  for (el=abs(spin); el<=L-1; el++) {
3954  for (m=1; m<=el; m++) {
3955  ssht_sampling_elm2ind(&ind, el, m);
3956  ssht_sampling_elm2ind(&ind_nm, el, -m);
3957  flm[ind_nm] = signs[m] * conj(flm[ind]);
3958  }
3959  }
3960 
3961  // Free memory.
3962  free(dlm1p1_line);
3963  free(dl_line);
3964  free(thetas);
3965  free(weights);
3966  free(Ftm);
3967  free(signs);
3968  free(sqrt_tbl);
3969  free(inds);
3970 
3971  // Print finished if verbosity set.
3972  if (verbosity > 0)
3973  printf("%s%s", SSHT_PROMPT, "Forward transform computed!");
3974 
3975 }
3976 
3977 
3978 //============================================================================
3979 // DH algorithms
3980 //============================================================================
3981 
3982 
3996 void ssht_core_dh_inverse_sov(SSHT_COMPLEX(double) *f, const SSHT_COMPLEX(double) *flm,
3997  int L, int spin, int verbosity) {
3998 
3999  int t, p, m, el, ind;
4000  int ftm_stride, ftm_offset, f_stride;
4001  double *dlm1p1_line, *dl_line;
4002  double *dl_ptr;
4003  double *sqrt_tbl, *signs;
4004  SSHT_COMPLEX(double) *ftm, *inout;
4005  double theta, ssign, elfactor;
4006  fftw_plan plan;
4007 
4008  // Allocate memory.
4009  sqrt_tbl = (double*)calloc(2*(L-1)+2, sizeof(double));
4010  SSHT_ERROR_MEM_ALLOC_CHECK(sqrt_tbl)
4011  signs = (double*)calloc(L+1, sizeof(double));
4013 
4014  // Perform precomputations.
4015  for (el=0; el<=2*(L-1)+1; el++)
4016  sqrt_tbl[el] = sqrt((double)el);
4017  for (m=0; m<=L-1; m=m+2) {
4018  signs[m] = 1.0;
4019  signs[m+1] = -1.0;
4020  }
4021  ssign = signs[abs(spin)];
4022 
4023  // Print messages depending on verbosity level.
4024  if (verbosity > 0) {
4025  printf("%s%s\n", SSHT_PROMPT,
4026  "Computing inverse transform using DH sampling with ");
4027  printf("%s%s%d%s%d%s\n", SSHT_PROMPT, "parameters (L,spin,reality) = (",
4028  L, ",", spin, ", FALSE)");
4029  if (verbosity > 1)
4030  printf("%s%s\n", SSHT_PROMPT,
4031  "Using routine ssht_core_dh_inverse_sov...");
4032  }
4033 
4034  // Compute ftm.
4035  ftm = (SSHT_COMPLEX(double)*)calloc((2*L)*(2*L-1), sizeof(SSHT_COMPLEX(double)));
4037  ftm_stride = 2*L-1;
4038  ftm_offset = L-1;
4039  dlm1p1_line = (double*)calloc(2*L-1, sizeof(double));
4040  SSHT_ERROR_MEM_ALLOC_CHECK(dlm1p1_line)
4041  dl_line = (double*)calloc(2*L-1, sizeof(double));
4043  for (t=0; t<=2*L-1; t++) {
4044  theta = ssht_sampling_dh_t2theta(t, L);
4045  for (el=abs(spin); el<=L-1; el++) {
4046  elfactor = sqrt((double)(2.0*el+1.0)/(4.0*SSHT_PI));
4047 
4048  // Compute dl line for given spin.
4049  ssht_dl_beta_kostelec_line_table(dlm1p1_line, dl_line,
4050  theta, L, -spin, el,
4051  sqrt_tbl, signs);
4052  // Switch current and previous dls.
4053  dl_ptr = dl_line;
4054  dl_line = dlm1p1_line;
4055  dlm1p1_line = dl_ptr;
4056 
4057  for (m=-el; m<=el; m++) {
4058  ssht_sampling_elm2ind(&ind, el, m);
4059  ftm[t*ftm_stride + m + ftm_offset] +=
4060  ssign
4061  * elfactor
4062  * dl_line[m + L-1]
4063  * flm[ind];
4064  }
4065  }
4066  }
4067 
4068  // Free dl memory.
4069  free(dlm1p1_line);
4070  free(dl_line);
4071 
4072  // Compute f.
4073  inout = (SSHT_COMPLEX(double)*)calloc(2*L-1, sizeof(SSHT_COMPLEX(double)));
4075  f_stride = 2*L-1;
4076  plan = fftw_plan_dft_1d(2*L-1, inout, inout, FFTW_BACKWARD, FFTW_MEASURE);
4077  for (t=0; t<=2*L-1; t++) {
4078  for (m=0; m<=L-1; m++)
4079  inout[m] = ftm[t*ftm_stride + m + ftm_offset];
4080  for (m=-(L-1); m<=-1; m++)
4081  inout[m+2*L-1] = ftm[t*ftm_stride + m + ftm_offset];
4082  fftw_execute_dft(plan, inout, inout);
4083  for (p=0; p<=2*L-2; p++)
4084  f[t*f_stride + p] = inout[p];
4085  }
4086  fftw_destroy_plan(plan);
4087 
4088  // Free memory.
4089  free(ftm);
4090  free(inout);
4091  free(signs);
4092  free(sqrt_tbl);
4093 
4094  // Print finished if verbosity set.
4095  if (verbosity > 0)
4096  printf("%s%s", SSHT_PROMPT, "Inverse transform computed!");
4097 
4098 }
4099 
4100 
4115 void ssht_core_dh_inverse_sov_real(double *f, const SSHT_COMPLEX(double) *flm,
4116  int L, int verbosity) {
4117 
4118  int t, p, m, el, ind;
4119  int ftm_stride, ftm_offset, f_stride;
4120  double *dlm1p1_line, *dl_line;
4121  double *dl_ptr;
4122  double *sqrt_tbl, *signs;
4123  SSHT_COMPLEX(double) *ftm;
4124  SSHT_COMPLEX(double) *in;
4125  double *out;
4126  double theta, ssign, elfactor;
4127  fftw_plan plan;
4128  int spin = 0;
4129 
4130  // Allocate memory.
4131  sqrt_tbl = (double*)calloc(2*(L-1)+2, sizeof(double));
4132  SSHT_ERROR_MEM_ALLOC_CHECK(sqrt_tbl)
4133  signs = (double*)calloc(L+1, sizeof(double));
4135 
4136  // Perform precomputations.
4137  for (el=0; el<=2*(L-1)+1; el++)
4138  sqrt_tbl[el] = sqrt((double)el);
4139  for (m=0; m<=L-1; m=m+2) {
4140  signs[m] = 1.0;
4141  signs[m+1] = -1.0;
4142  }
4143  ssign = signs[abs(spin)];
4144 
4145  // Print messages depending on verbosity level.
4146  if (verbosity > 0) {
4147  printf("%s%s\n", SSHT_PROMPT,
4148  "Computing inverse transform using DH sampling with ");
4149  printf("%s%s%d%s%d%s\n", SSHT_PROMPT, "parameters (L,spin,reality) = (",
4150  L, ",", spin, ", TRUE)");
4151  if (verbosity > 1)
4152  printf("%s%s\n", SSHT_PROMPT,
4153  "Using routine ssht_core_dh_inverse_sov_real...");
4154  }
4155 
4156  // Compute ftm.
4157  ftm = (SSHT_COMPLEX(double)*)calloc(2*L*L, sizeof(SSHT_COMPLEX(double)));
4159  ftm_stride = L;
4160  ftm_offset = 0;
4161  dlm1p1_line = (double*)calloc(L, sizeof(double));
4162  SSHT_ERROR_MEM_ALLOC_CHECK(dlm1p1_line)
4163  dl_line = (double*)calloc(L, sizeof(double));
4165  for (t=0; t<=2*L-1; t++) {
4166  theta = ssht_sampling_dh_t2theta(t, L);
4167  for (el=abs(spin); el<=L-1; el++) {
4168  elfactor = sqrt((double)(2.0*el+1.0)/(4.0*SSHT_PI));
4169 
4170  // Compute half dl line for given spin.
4171  ssht_dl_beta_kostelec_halfline_table(dlm1p1_line, dl_line,
4172  theta, L, -spin, el,
4173  sqrt_tbl, signs);
4174  // Switch current and previous dls.
4175  dl_ptr = dl_line;
4176  dl_line = dlm1p1_line;
4177  dlm1p1_line = dl_ptr;
4178 
4179  for (m=0; m<=el; m++) {
4180  ssht_sampling_elm2ind(&ind, el, m);
4181  ftm[t*ftm_stride + m + ftm_offset] +=
4182  ssign
4183  * elfactor
4184  * dl_line[m]
4185  * flm[ind];
4186  }
4187  }
4188  }
4189 
4190  // Free dl memory.
4191  free(dlm1p1_line);
4192  free(dl_line);
4193 
4194  // Compute f.
4195  in = (SSHT_COMPLEX(double)*)calloc(L, sizeof(SSHT_COMPLEX(double)));
4197  out = (double*)calloc(2*L-1, sizeof(double));
4199  plan = fftw_plan_dft_c2r_1d(2*L-1, in, out, FFTW_MEASURE);
4200  f_stride = 2*L-1;
4201  for (t=0; t<=2*L-1; t++) {
4202  memcpy(in, &ftm[t*ftm_stride], L*sizeof(SSHT_COMPLEX(double)));
4203  fftw_execute_dft_c2r(plan, in, out);
4204  for (p=0; p<=2*L-2; p++)
4205  f[t*f_stride + p] = out[p];
4206  }
4207  fftw_destroy_plan(plan);
4208 
4209  // Free memory.
4210  free(ftm);
4211  free(in);
4212  free(out);
4213  free(signs);
4214  free(sqrt_tbl);
4215 
4216  // Print finished if verbosity set.
4217  if (verbosity > 0)
4218  printf("%s%s", SSHT_PROMPT, "Inverse transform computed!");
4219 
4220 }
4221 
4222 
4236 void ssht_core_dh_forward_sov(SSHT_COMPLEX(double) *flm, const SSHT_COMPLEX(double) *f,
4237  int L, int spin, int verbosity) {
4238 
4239  int t, m, el, ind;
4240  int f_stride;
4241  double *dlm1p1_line, *dl_line;
4242  double *dl_ptr;
4243  int el2pel, inds_offset;
4244  int *inds;
4245  double *sqrt_tbl, *signs;
4246  int Ftm_stride, Ftm_offset;
4247  SSHT_COMPLEX(double) *Ftm, *inout;
4248  double theta, ssign, elfactor;
4249  fftw_plan plan;
4250  double w;
4251 
4252  // Allocate memory.
4253  sqrt_tbl = (double*)calloc(2*(L-1)+2, sizeof(double));
4254  SSHT_ERROR_MEM_ALLOC_CHECK(sqrt_tbl)
4255  signs = (double*)calloc(L+1, sizeof(double));
4257  inds = (int*)calloc(2*L-1, sizeof(int));
4259 
4260  // Perform precomputations.
4261  for (el=0; el<=2*(L-1)+1; el++)
4262  sqrt_tbl[el] = sqrt((double)el);
4263  for (m=0; m<=L-1; m=m+2) {
4264  signs[m] = 1.0;
4265  signs[m+1] = -1.0;
4266  }
4267  ssign = signs[abs(spin)];
4268 
4269  // Print messages depending on verbosity level.
4270  if (verbosity > 0) {
4271  printf("%s%s\n", SSHT_PROMPT,
4272  "Computing forward transform using DH sampling with ");
4273  printf("%s%s%d%s%d%s\n", SSHT_PROMPT, "parameters (L,spin,reality) = (",
4274  L, ",", spin, ", FALSE)");
4275  if (verbosity > 1)
4276  printf("%s%s\n", SSHT_PROMPT,
4277  "Using routine ssht_core_dh_forward_sov...");
4278  }
4279 
4280  // Compute Fourier transform over phi, i.e. compute Ftm.
4281  Ftm = (SSHT_COMPLEX(double)*)calloc((2*L)*(2*L-1), sizeof(SSHT_COMPLEX(double)));
4283  Ftm_stride = 2*L-1;
4284  Ftm_offset = L-1;
4285  f_stride = 2*L-1;
4286  inout = (SSHT_COMPLEX(double)*)calloc(2*L-1, sizeof(SSHT_COMPLEX(double)));
4288  plan = fftw_plan_dft_1d(2*L-1, inout, inout, FFTW_FORWARD, FFTW_MEASURE);
4289  for (t=0; t<=2*L-1; t++) {
4290  memcpy(inout, &f[t*f_stride], f_stride*sizeof(SSHT_COMPLEX(double)));
4291  fftw_execute_dft(plan, inout, inout);
4292  for(m=0; m<=L-1; m++)
4293  Ftm[t*Ftm_stride + m + Ftm_offset] =
4294  inout[m] * 2.0 * SSHT_PI / (2.0*L-1.0) ;
4295  for(m=-(L-1); m<=-1; m++)
4296  Ftm[t*Ftm_stride + m + Ftm_offset] =
4297  inout[m+2*L-1] * 2.0 * SSHT_PI / (2.0*L-1.0);
4298  }
4299  fftw_destroy_plan(plan);
4300 
4301  // Compute flm.
4302  dlm1p1_line = (double*)calloc(2*L-1, sizeof(double));
4303  SSHT_ERROR_MEM_ALLOC_CHECK(dlm1p1_line)
4304  dl_line = (double*)calloc(2*L-1, sizeof(double));
4306  inds_offset = L-1;
4307  for (el=0; el<=L-1; el++) {
4308  for (m=-el; m<=el; m++) {
4309  ssht_sampling_elm2ind(&ind, el, m);
4310  flm[ind] = 0.0;
4311  }
4312  }
4313  for (t=0; t<=2*L-1; t++) {
4314  theta = ssht_sampling_dh_t2theta(t, L);
4315  w = ssht_sampling_weight_dh(theta, L);
4316  for (el=abs(spin); el<=L-1; el++) {
4317  elfactor = sqrt((double)(2.0*el+1.0)/(4.0*SSHT_PI));
4318  el2pel = el *el + el;
4319  for (m=-el; m<=el; m++)
4320  inds[m + inds_offset] = el2pel + m;
4321 
4322  // Compute dl line for given spin.
4323  ssht_dl_beta_kostelec_line_table(dlm1p1_line, dl_line,
4324  theta, L, -spin, el,
4325  sqrt_tbl, signs);
4326  // Switch current and previous dls.
4327  dl_ptr = dl_line;
4328  dl_line = dlm1p1_line;
4329  dlm1p1_line = dl_ptr;
4330 
4331  for (m=-el; m<=el; m++) {
4332  ind = inds[m + inds_offset];
4333  flm[ind] +=
4334  ssign
4335  * elfactor
4336  * w
4337  * dl_line[m+L-1]
4338  * Ftm[t*Ftm_stride + m + Ftm_offset];
4339  }
4340  }
4341  }
4342 
4343  // Free memory.
4344  free(dlm1p1_line);
4345  free(dl_line);
4346  free(Ftm);
4347  free(inout);
4348  free(signs);
4349  free(sqrt_tbl);
4350  free(inds);
4351 
4352  // Print finished if verbosity set.
4353  if (verbosity > 0)
4354  printf("%s%s", SSHT_PROMPT, "Forward transform computed!");
4355 
4356 }
4357 
4358 
4373 void ssht_core_dh_forward_sov_real(SSHT_COMPLEX(double) *flm, const double *f,
4374  int L, int verbosity) {
4375 
4376  int t, m, el, ind, ind_nm;
4377  int f_stride;
4378  double *dlm1p1_line, *dl_line;
4379  double *dl_ptr;
4380  int el2pel, inds_offset;
4381  int *inds;
4382  double *sqrt_tbl, *signs;
4383  int Ftm_stride, Ftm_offset;
4384  SSHT_COMPLEX(double) *Ftm;
4385  double *in_real;
4386  SSHT_COMPLEX(double) *out;
4387  double theta, ssign, elfactor;
4388  fftw_plan plan;
4389  double w;
4390  int spin = 0;
4391 
4392  // Allocate memory.
4393  sqrt_tbl = (double*)calloc(2*(L-1)+2, sizeof(double));
4394  SSHT_ERROR_MEM_ALLOC_CHECK(sqrt_tbl)
4395  signs = (double*)calloc(L+1, sizeof(double));
4397  inds = (int*)calloc(L, sizeof(int));
4399 
4400  // Perform precomputations.
4401  for (el=0; el<=2*(L-1)+1; el++)
4402  sqrt_tbl[el] = sqrt((double)el);
4403  for (m=0; m<=L-1; m=m+2) {
4404  signs[m] = 1.0;
4405  signs[m+1] = -1.0;
4406  }
4407  ssign = signs[abs(spin)];
4408 
4409  // Print messages depending on verbosity level.
4410  if (verbosity > 0) {
4411  printf("%s%s\n", SSHT_PROMPT,
4412  "Computing forward transform using GL sampling with ");
4413  printf("%s%s%d%s%d%s\n", SSHT_PROMPT, "parameters (L,spin,reality) = (",
4414  L, ",", spin, ", TRUE)");
4415  if (verbosity > 1)
4416  printf("%s%s\n", SSHT_PROMPT,
4417  "Using routine ssht_core_gl_forward_sov_real...");
4418  }
4419 
4420  // Compute Fourier transform over phi, i.e. compute Ftm.
4421  Ftm = (SSHT_COMPLEX(double)*)calloc(2*L*L, sizeof(SSHT_COMPLEX(double)));
4423  Ftm_stride = L;
4424  Ftm_offset = 0;
4425  f_stride = 2*L-1;
4426  in_real = (double*)calloc(2*L-1, sizeof(double));
4428  out = (SSHT_COMPLEX(double)*)calloc(L, sizeof(SSHT_COMPLEX(double)));
4430  plan = fftw_plan_dft_r2c_1d(2*L-1, in_real, out, FFTW_MEASURE);
4431  for (t=0; t<=2*L-1; t++) {
4432  memcpy(in_real, &f[t*f_stride], f_stride*sizeof(double));
4433  fftw_execute_dft_r2c(plan, in_real, out);
4434  for(m=0; m<=L-1; m++)
4435  Ftm[t*Ftm_stride + m + Ftm_offset] =
4436  out[m] * 2.0 * SSHT_PI / (2.0*L-1.0);
4437  }
4438  free(in_real);
4439  free(out);
4440  fftw_destroy_plan(plan);
4441 
4442  // Compute flm.
4443  dlm1p1_line = (double*)calloc(L, sizeof(double));
4444  SSHT_ERROR_MEM_ALLOC_CHECK(dlm1p1_line)
4445  dl_line = (double*)calloc(L, sizeof(double));
4447  inds_offset = 0;
4448  for (el=0; el<=L-1; el++) {
4449  for (m=0; m<=el; m++) {
4450  ssht_sampling_elm2ind(&ind, el, m);
4451  flm[ind] = 0.0;
4452  }
4453  }
4454  for (t=0; t<=2*L-1; t++) {
4455  theta = ssht_sampling_dh_t2theta(t, L);
4456  w = ssht_sampling_weight_dh(theta, L);
4457  for (el=abs(spin); el<=L-1; el++) {
4458  elfactor = sqrt((double)(2.0*el+1.0)/(4.0*SSHT_PI));
4459  el2pel = el *el + el;
4460  for (m=0; m<=el; m++)
4461  inds[m + inds_offset] = el2pel + m;
4462 
4463  // Compute half dl line for given spin.
4464  ssht_dl_beta_kostelec_halfline_table(dlm1p1_line, dl_line,
4465  theta, L, -spin, el,
4466  sqrt_tbl, signs);
4467 
4468  // Switch current and previous dls.
4469  dl_ptr = dl_line;
4470  dl_line = dlm1p1_line;
4471  dlm1p1_line = dl_ptr;
4472 
4473  for (m=0; m<=el; m++) {
4474  ind = inds[m + inds_offset];
4475  flm[ind] +=
4476  ssign
4477  * elfactor
4478  * w
4479  * dl_line[m]
4480  * Ftm[t*Ftm_stride + m + Ftm_offset];
4481  }
4482  }
4483  }
4484 
4485  // Set flm values for negative m using conjugate symmetry.
4486  for (el=abs(spin); el<=L-1; el++) {
4487  for (m=1; m<=el; m++) {
4488  ssht_sampling_elm2ind(&ind, el, m);
4489  ssht_sampling_elm2ind(&ind_nm, el, -m);
4490  flm[ind_nm] = signs[m] * conj(flm[ind]);
4491  }
4492  }
4493 
4494  // Free memory.
4495  free(dlm1p1_line);
4496  free(dl_line);
4497  free(Ftm);
4498  free(signs);
4499  free(sqrt_tbl);
4500  free(inds);
4501 
4502  // Print finished if verbosity set.
4503  if (verbosity > 0)
4504  printf("%s%s", SSHT_PROMPT, "Forward transform computed!");
4505 
4506 }
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_dl_beta_risbo_full_table
void ssht_dl_beta_risbo_full_table(double *dl, double beta, int L, ssht_dl_size_t dl_size, int el, double *sqrt_tbl)
Definition: ssht_dl.c:190
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_PION2
#define SSHT_PION2
Definition: ssht_types.h:39
SSHT_PI
#define SSHT_PI
Definition: ssht_types.h:38
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_error.h
SSHT_DL_QUARTER_EXTENDED
@ SSHT_DL_QUARTER_EXTENDED
Definition: ssht_dl.h:16
SSHT_ERROR_GENERIC
#define SSHT_ERROR_GENERIC(comment)
Definition: ssht_error.h:19
MAX
#define MAX(a, b)
Definition: ssht_core.c:28
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_DL_TRAPANI
@ SSHT_DL_TRAPANI
Definition: ssht_dl.h:22
ssht_dl.h
ssht_sampling_mw_t2theta
double ssht_sampling_mw_t2theta(int t, int L)
Definition: ssht_sampling.c:197
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_dl_get_stride
int ssht_dl_get_stride(int L, ssht_dl_size_t dl_size)
Definition: ssht_dl.c:146
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_DL_FULL
@ SSHT_DL_FULL
Definition: ssht_dl.h:18
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
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_dl_calloc
double * ssht_dl_calloc(int L, ssht_dl_size_t dl_size)
Definition: ssht_dl.c:64
ssht_sampling.h
ssht_dl_beta_kostelec_halfline_table
void ssht_dl_beta_kostelec_halfline_table(double *dlm1p1_line, double *dl_line, double beta, int L, int mm, int el, double *sqrt_tbl, double *signs)
Definition: ssht_dl.c:943
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
ssht_core_mwdirect_inverse_ss
void ssht_core_mwdirect_inverse_ss(SSHT_COMPLEX(double) *f, const SSHT_COMPLEX(double) *flm, int L, int spin, int verbosity)
Definition: ssht_core.c:2326
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_sampling_dh_t2theta
double ssht_sampling_dh_t2theta(int t, int L)
Definition: ssht_sampling.c:361
ssht_dl_halfpi_trapani_eighth_table
void ssht_dl_halfpi_trapani_eighth_table(double *dl, int L, ssht_dl_size_t dl_size, int el, double *sqrt_tbl)
Definition: ssht_dl.c:1077
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
ssht_dl_beta_risbo_eighth_table
void ssht_dl_beta_risbo_eighth_table(double *dl, double beta, int L, ssht_dl_size_t dl_size, int el, double *sqrt_tbl, double *signs)
Definition: ssht_dl.c:449
ssht_dl_halfpi_trapani_fill_eighth2quarter_table
void ssht_dl_halfpi_trapani_fill_eighth2quarter_table(double *dl, int L, ssht_dl_size_t dl_size, int el, double *signs)
Definition: ssht_dl.c:1352
ssht_sampling_mw_p2phi
double ssht_sampling_mw_p2phi(int p, int L)
Definition: ssht_sampling.c:231
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
ssht_core.h
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_mwdirect_inverse
void ssht_core_mwdirect_inverse(SSHT_COMPLEX(double) *f, const SSHT_COMPLEX(double) *flm, int L, int spin, int verbosity)
Definition: ssht_core.c:590
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_sampling_weight_dh
double ssht_sampling_weight_dh(double theta_t, int L)
Definition: ssht_sampling.c:65
ssht_sampling_gl_thetas_weights
void ssht_sampling_gl_thetas_weights(double *thetas, double *weights, int L)
Definition: ssht_sampling.c:94
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
ssht_core_mwdirect_inverse_sov
void ssht_core_mwdirect_inverse_sov(SSHT_COMPLEX(double) *f, const SSHT_COMPLEX(double) *flm, int L, int spin, int verbosity)
Definition: ssht_core.c:689
ssht_dl_beta_risbo_fill_eighth2quarter_table
void ssht_dl_beta_risbo_fill_eighth2quarter_table(double *dl, double *dl8, int L, ssht_dl_size_t dl_size, ssht_dl_size_t dl8_size, int el, double *signs)
Definition: ssht_dl.c:599
ssht_sampling_mw_ss_p2phi
double ssht_sampling_mw_ss_p2phi(int p, int L)
Definition: ssht_sampling.c:307
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_PROMPT
#define SSHT_PROMPT
Definition: ssht_types.h:45
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_sampling_mw_ss_t2theta
double ssht_sampling_mw_ss_t2theta(int t, int L)
Definition: ssht_sampling.c:269
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_dl_get_offset
int ssht_dl_get_offset(int L, ssht_dl_size_t dl_size)
Definition: ssht_dl.c:111
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_COMPLEX
#define SSHT_COMPLEX(TYPE)
Definition: ssht_types.h:52
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_dl_beta_kostelec_line_table
void ssht_dl_beta_kostelec_line_table(double *dlm1p1_line, double *dl_line, double beta, int L, int mm, int el, double *sqrt_tbl, double *signs)
Definition: ssht_dl.c:795
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_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_types.h
SSHT_DL_QUARTER
@ SSHT_DL_QUARTER
Definition: ssht_dl.h:15
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