My Project
ssht_adjoint.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 
14 #include <stdio.h>
15 #include <stdlib.h>
16 #include <string.h>
17 #include <math.h>
18 #include <complex.h> // Must be before fftw3.h
19 #include <fftw3.h>
20 
21 #include "ssht_types.h"
22 #include "ssht_error.h"
23 #include "ssht_dl.h"
24 #include "ssht_sampling.h"
25 #include "ssht_adjoint.h"
26 
27 
28 //============================================================================
29 // MW algorithms
30 //============================================================================
31 
32 
49  const SSHT_COMPLEX(double) *f,
50  int L, int spin,
51  ssht_dl_method_t dl_method,
52  int verbosity) {
53 
54  int el, m, mm, ind, t;
55  int eltmp;
56  double *sqrt_tbl, *signs;
57  int el2pel, inds_offset;
58  int *inds;
59  double ssign, elfactor;
60  fftw_plan plan;
61  SSHT_COMPLEX(double) *inout;
62  SSHT_COMPLEX(double) *Fmt, *Fmm;
63  int f_stride, Fmt_stride, Fmt_offset, Fmm_stride, Fmm_offset;
64  double *dl;
65  double *dl8 = NULL;
66  int dl_offset, dl_stride;
67  SSHT_COMPLEX(double) *expsm, *expsmm;
68  int exps_offset;
69  int elmmsign, elssign;
70  int spinneg;
71 
72  // Allocate memory.
73  sqrt_tbl = (double*)calloc(2*(L-1)+2, sizeof(double));
75  signs = (double*)calloc(L+1, sizeof(double));
77  expsm = (SSHT_COMPLEX(double)*)calloc(2*L-1, sizeof(SSHT_COMPLEX(double)));
79  expsmm = (SSHT_COMPLEX(double)*)calloc(2*L-1, sizeof(SSHT_COMPLEX(double)));
81  inds = (int*)calloc(2*L-1, sizeof(int));
83 
84  // Perform precomputations.
85  for (el=0; el<=2*(L-1)+1; el++)
86  sqrt_tbl[el] = sqrt((double)el);
87  for (m=0; m<=L-1; m=m+2) {
88  signs[m] = 1.0;
89  signs[m+1] = -1.0;
90  }
91  ssign = signs[abs(spin)];
92  spinneg = spin <= 0 ? spin : -spin;
93  exps_offset = L-1;
94  for (m=-(L-1); m<=L-1; m++)
95  expsm[m + exps_offset] = cexp(I*SSHT_PION2*(m+spin));
96  for (mm=-(L-1); mm<=L-1; mm++)
97  expsmm[mm + exps_offset] = cexp(-I*mm*SSHT_PI/(2.0*L-1.0));
98 
99  // Print messages depending on verbosity level.
100  if (verbosity > 0) {
101  printf("%s%s\n", SSHT_PROMPT,
102  "Computing adjoint inverse transform using MW sampling with ");
103  printf("%s%s%d%s%d%s\n", SSHT_PROMPT, "parameters (L,spin,reality) = (",
104  L, ",", spin, ", FALSE)");
105  if (verbosity > 1)
106  printf("%s%s\n", SSHT_PROMPT,
107  "Using routine ssht_adjoint_mw_inverse_sov_sym...");
108  }
109 
110  // Compute Fourier transform over phi, i.e. compute Fmt.
111  Fmt = (SSHT_COMPLEX(double)*)calloc((2*L-1)*(2*L-1), sizeof(SSHT_COMPLEX(double)));
113  Fmt_stride = 2*L-1;
114  Fmt_offset = L-1;
115  f_stride = 2*L-1;
116  inout = (SSHT_COMPLEX(double)*)calloc(2*L-1, sizeof(SSHT_COMPLEX(double)));
118  plan = fftw_plan_dft_1d(2*L-1, inout, inout, FFTW_FORWARD, FFTW_MEASURE);
119  for (t=0; t<=L-1; t++) {
120  memcpy(inout, &f[t*f_stride], f_stride*sizeof(SSHT_COMPLEX(double)));
121  fftw_execute_dft(plan, inout, inout);
122  for(m=0; m<=L-1; m++)
123  Fmt[(m+Fmt_offset)*Fmt_stride + t] = inout[m];
124  for(m=-(L-1); m<=-1; m++)
125  Fmt[(m+Fmt_offset)*Fmt_stride + t] = inout[m+2*L-1];
126  }
127 
128  // Apply adjoint of periodic entension.
129  for (m=-(L-1); m<=L-1; m++)
130  for (t=L; t<=2*L-2; t++)
131  Fmt[(m+Fmt_offset)*Fmt_stride + t] = 0.0;
132 
133  // Compute Fourier transform over theta, i.e. compute Fmm.
134  Fmm = (SSHT_COMPLEX(double)*)calloc((2*L-1)*(2*L-1), sizeof(SSHT_COMPLEX(double)));
136  Fmm_stride = 2*L-1;
137  Fmm_offset = L-1;
138  for (m=-(L-1); m<=L-1; m++) {
139  memcpy(inout, &Fmt[(m+Fmt_offset)*Fmt_stride], Fmt_stride*sizeof(SSHT_COMPLEX(double)));
140  fftw_execute_dft(plan, inout, inout);
141  for(mm=0; mm<=L-1; mm++)
142  Fmm[(mm+Fmm_offset)*Fmm_stride + m + Fmm_offset] =
143  inout[mm];
144  for(mm=-(L-1); mm<=-1; mm++)
145  Fmm[(mm+Fmm_offset)*Fmm_stride + m + Fmm_offset] =
146  inout[mm+2*L-1];
147  }
148  fftw_destroy_plan(plan);
149  free(inout);
150 
151  // Apply phase modulation to account for sampling offset.
152  for (mm=-(L-1); mm<=L-1; mm++)
153  for (m=-(L-1); m<=L-1; m++)
154  Fmm[(mm+Fmm_offset)*Fmm_stride + m + Fmm_offset] *=
155  expsmm[mm + exps_offset];
156 
157  // Compute flm.
160  if (dl_method == SSHT_DL_RISBO) {
163  }
164  dl_offset = ssht_dl_get_offset(L, SSHT_DL_QUARTER);
165  dl_stride = ssht_dl_get_stride(L, SSHT_DL_QUARTER);
166  inds_offset = L-1;
167  for (el=0; el<=L-1; el++) {
168  for (m=-el; m<=el; m++) {
169  ssht_sampling_elm2ind(&ind, el, m);
170  flm[ind] = 0.0;
171  }
172  }
173  for (el=abs(spin); el<=L-1; el++) {
174 
175  // Compute Wigner plane.
176  switch (dl_method) {
177 
178  case SSHT_DL_RISBO:
179  if (el!=0 && el==abs(spin)) {
180  for(eltmp=0; eltmp<=abs(spin); eltmp++)
183  eltmp, sqrt_tbl, signs);
185  dl8, L,
188  el,
189  signs);
190  }
191  else {
194  el, sqrt_tbl, signs);
196  dl8, L,
199  el,
200  signs);
201  }
202  break;
203 
204  case SSHT_DL_TRAPANI:
205  if (el!=0 && el==abs(spin)) {
206  for(eltmp=0; eltmp<=abs(spin); eltmp++)
209  eltmp, sqrt_tbl);
212  el, signs);
213  }
214  else {
217  el, sqrt_tbl);
220  el, signs);
221  }
222  break;
223 
224  default:
225  SSHT_ERROR_GENERIC("Invalid dl method")
226  }
227 
228  // Compute flm.
229  elfactor = sqrt((double)(2.0*el+1.0)/(4.0*SSHT_PI));
230  el2pel = el *el + el;
231  for (m=-el; m<=el; m++)
232  inds[m + inds_offset] = el2pel + m;
233  elssign = spin <= 0 ? 1.0 : signs[el];
234 
235  for (m=-el; m<=-1; m++) {
236  // mm = 0
237  ind = inds[m + inds_offset];
238  flm[ind] +=
239  ssign
240  * elfactor
241  * expsm[m + exps_offset]
242  * signs[el] * dl[0*dl_stride - m + dl_offset]
243  * elssign * dl[0*dl_stride - spinneg + dl_offset]
244  * Fmm[(0+Fmm_offset)*Fmm_stride + m + Fmm_offset];
245  }
246  for (m=0; m<=el; m++) {
247  // mm = 0
248  ind = inds[m + inds_offset];
249  flm[ind] +=
250  ssign
251  * elfactor
252  * expsm[m + exps_offset]
253  * dl[0*dl_stride + m + dl_offset]
254  * elssign * dl[0*dl_stride - spinneg + dl_offset]
255  * Fmm[(0+Fmm_offset)*Fmm_stride + m + Fmm_offset];
256  }
257 
258  for (mm=1; mm<=el; mm++) {
259  elmmsign = signs[el] * signs[mm];
260  elssign = spin <= 0 ? 1.0 : elmmsign;
261 
262  for (m=-el; m<=-1; m++) {
263  ind = inds[m + inds_offset];
264  flm[ind] +=
265  ssign
266  * elfactor
267  * expsm[m + exps_offset]
268  * elmmsign * dl[mm*dl_stride - m + dl_offset]
269  * elssign * dl[mm*dl_stride - spinneg + dl_offset]
270  * ( Fmm[(mm+Fmm_offset)*Fmm_stride + m + Fmm_offset]
271  + signs[-m] * ssign
272  * Fmm[(-mm+Fmm_offset)*Fmm_stride + m + Fmm_offset]);
273  }
274  for (m=0; m<=el; m++) {
275  ind = inds[m + inds_offset];
276  flm[ind] +=
277  ssign
278  * elfactor
279  * expsm[m + exps_offset]
280  * dl[mm*dl_stride + m + dl_offset]
281  * elssign * dl[mm*dl_stride - spinneg + dl_offset]
282  * ( Fmm[(mm+Fmm_offset)*Fmm_stride + m + Fmm_offset]
283  + signs[m] * ssign
284  * Fmm[(-mm+Fmm_offset)*Fmm_stride + m + Fmm_offset]);
285  }
286 
287  }
288 
289  }
290 
291  // Free memory.
292  free(dl);
293  if (dl_method == SSHT_DL_RISBO)
294  free(dl8);
295  free(Fmt);
296  free(Fmm);
297  free(sqrt_tbl);
298  free(signs);
299  free(expsm);
300  free(expsmm);
301  free(inds);
302 
303  // Print finished if verbosity set.
304  if (verbosity > 0)
305  printf("%s%s", SSHT_PROMPT, "Adjoint inverse transform computed!");
306 
307 }
308 
309 
326  const double *f,
327  int L,
328  ssht_dl_method_t dl_method,
329  int verbosity) {
330 
331  int el, m, mm, ind, ind_nm, t;
332  int eltmp;
333  double *sqrt_tbl, *signs;
334  int el2pel, inds_offset;
335  int *inds;
336  double ssign, elfactor;
337  fftw_plan plan;
338  double *in_real;
339  SSHT_COMPLEX(double) *inout, *out;
340  SSHT_COMPLEX(double) *Fmt, *Fmm;
341  int f_stride, Fmt_stride, Fmt_offset, Fmm_stride, Fmm_offset;
342  double *dl;
343  double *dl8 = NULL;
344  int dl_offset, dl_stride;
345  SSHT_COMPLEX(double) *expsm, *expsmm;
346  int exps_offset;
347  int elmmsign, elssign;
348  int spinneg;
349  int spin = 0;
350 
351  // Allocate memory.
352  sqrt_tbl = (double*)calloc(2*(L-1)+2, sizeof(double));
354  signs = (double*)calloc(L+1, sizeof(double));
356  expsm = (SSHT_COMPLEX(double)*)calloc(L, sizeof(SSHT_COMPLEX(double)));
358  expsmm = (SSHT_COMPLEX(double)*)calloc(2*L-1, sizeof(SSHT_COMPLEX(double)));
360  inds = (int*)calloc(L, sizeof(int));
362 
363  // Perform precomputations.
364  for (el=0; el<=2*(L-1)+1; el++)
365  sqrt_tbl[el] = sqrt((double)el);
366  for (m=0; m<=L-1; m=m+2) {
367  signs[m] = 1.0;
368  signs[m+1] = -1.0;
369  }
370  ssign = signs[abs(spin)];
371  spinneg = spin <= 0 ? spin : -spin;
372  exps_offset = L-1;
373  for (m=0; m<=L-1; m++)
374  expsm[m] = cexp(I*SSHT_PION2*(m+spin));
375  for (mm=-(L-1); mm<=L-1; mm++)
376  expsmm[mm + exps_offset] = cexp(-I*mm*SSHT_PI/(2.0*L-1.0));
377 
378  // Print messages depending on verbosity level.
379  if (verbosity > 0) {
380  printf("%s%s\n", SSHT_PROMPT,
381  "Computing adjoint inverse transform using MW sampling with ");
382  printf("%s%s%d%s%d%s\n", SSHT_PROMPT, "parameters (L,spin,reality) = (",
383  L, ",", spin, ", TRUE)");
384  if (verbosity > 1)
385  printf("%s%s\n", SSHT_PROMPT,
386  "Using routine ssht_adjoint_mw_inverse_sov_sym_real...");
387  }
388 
389  // Compute Fourier transform over phi, i.e. compute Fmt.
390  Fmt = (SSHT_COMPLEX(double)*)calloc(L*(2*L-1), sizeof(SSHT_COMPLEX(double)));
392  Fmt_stride = 2*L-1;
393  f_stride = 2*L-1;
394  in_real = (double*)calloc(2*L-1, sizeof(double));
396  out = (SSHT_COMPLEX(double)*)calloc(L, sizeof(SSHT_COMPLEX(double)));
398  plan = fftw_plan_dft_r2c_1d(2*L-1, in_real, out, FFTW_MEASURE);
399  for (t=0; t<=L-1; t++) {
400  memcpy(in_real, &f[t*f_stride], f_stride*sizeof(double));
401  fftw_execute_dft_r2c(plan, in_real, out);
402  for(m=0; m<=L-1; m++)
403  Fmt[m*Fmt_stride + t] = out[m];
404  }
405  free(in_real);
406  free(out);
407  fftw_destroy_plan(plan);
408 
409  // Apply adjoint of periodic extension.
410  for (m=0; m<=L-1; m++)
411  for (t=L; t<=2*L-2; t++)
412  Fmt[m*Fmt_stride + t] = 0.0;
413 
414  // Compute Fourier transform over theta, i.e. compute Fmm.
415  Fmm = (SSHT_COMPLEX(double)*)calloc((2*L-1)*L, sizeof(SSHT_COMPLEX(double)));
417  Fmm_stride = L;
418  Fmm_offset = L-1;
419  inout = (SSHT_COMPLEX(double)*)calloc(2*L-1, sizeof(SSHT_COMPLEX(double)));
421  plan = fftw_plan_dft_1d(2*L-1, inout, inout, FFTW_FORWARD, FFTW_MEASURE);
422  for (m=0; m<=L-1; m++) {
423  memcpy(inout, &Fmt[m*Fmt_stride], Fmt_stride*sizeof(SSHT_COMPLEX(double)));
424  fftw_execute_dft(plan, inout, inout);
425  for(mm=0; mm<=L-1; mm++)
426  Fmm[(mm+Fmm_offset)*Fmm_stride + m] =
427  inout[mm];
428  for(mm=-(L-1); mm<=-1; mm++)
429  Fmm[(mm+Fmm_offset)*Fmm_stride + m] =
430  inout[mm+2*L-1];
431  }
432  fftw_destroy_plan(plan);
433  free(inout);
434 
435  // Apply phase modulation to account for sampling offset.
436  for (mm=-(L-1); mm<=L-1; mm++)
437  for (m=0; m<=L-1; m++)
438  Fmm[(mm+Fmm_offset)*Fmm_stride + m] *=
439  expsmm[mm + exps_offset];
440 
441  // Compute flm.
444  if (dl_method == SSHT_DL_RISBO) {
447  }
448  dl_offset = ssht_dl_get_offset(L, SSHT_DL_QUARTER);
449  dl_stride = ssht_dl_get_stride(L, SSHT_DL_QUARTER);
450  inds_offset = 0;
451  for (el=0; el<=L-1; el++) {
452  for (m=0; m<=el; m++) {
453  ssht_sampling_elm2ind(&ind, el, m);
454  flm[ind] = 0.0;
455  }
456  }
457  for (el=abs(spin); el<=L-1; el++) {
458 
459  // Compute Wigner plane.
460  switch (dl_method) {
461 
462  case SSHT_DL_RISBO:
463  if (el!=0 && el==abs(spin)) {
464  for(eltmp=0; eltmp<=abs(spin); eltmp++)
467  eltmp, sqrt_tbl, signs);
469  dl8, L,
472  el,
473  signs);
474  }
475  else {
478  el, sqrt_tbl, signs);
480  dl8, L,
483  el,
484  signs);
485  }
486  break;
487 
488  case SSHT_DL_TRAPANI:
489  if (el!=0 && el==abs(spin)) {
490  for(eltmp=0; eltmp<=abs(spin); eltmp++)
493  eltmp, sqrt_tbl);
496  el, signs);
497  }
498  else {
501  el, sqrt_tbl);
504  el, signs);
505  }
506  break;
507 
508  default:
509  SSHT_ERROR_GENERIC("Invalid dl method")
510  }
511 
512 
513  // Compute flm.
514  elfactor = sqrt((double)(2.0*el+1.0)/(4.0*SSHT_PI));
515  el2pel = el *el + el;
516  for (m=0; m<=el; m++)
517  inds[m + inds_offset] = el2pel + m;
518  elssign = spin <= 0 ? 1.0 : signs[el];
519 
520  for (m=0; m<=el; m++) {
521  // mm = 0
522  ind = inds[m + inds_offset];
523  flm[ind] +=
524  ssign
525  * elfactor
526  * expsm[m]
527  * dl[0*dl_stride + m + dl_offset]
528  * elssign * dl[0*dl_stride - spinneg + dl_offset]
529  * Fmm[(0+Fmm_offset)*Fmm_stride + m];
530  }
531 
532  for (mm=1; mm<=el; mm++) {
533  elmmsign = signs[el] * signs[mm];
534  elssign = spin <= 0 ? 1.0 : elmmsign;
535 
536  for (m=0; m<=el; m++) {
537  ind = inds[m + inds_offset];
538  flm[ind] +=
539  ssign
540  * elfactor
541  * expsm[m]
542  * dl[mm*dl_stride + m + dl_offset]
543  * elssign * dl[mm*dl_stride - spinneg + dl_offset]
544  * ( Fmm[(mm+Fmm_offset)*Fmm_stride + m]
545  + signs[m] * ssign
546  * Fmm[(-mm+Fmm_offset)*Fmm_stride + m]);
547  }
548 
549  }
550 
551  }
552 
553  // Set flm values for negative m using conjugate symmetry.
554  for (el=abs(spin); el<=L-1; el++) {
555  for (m=1; m<=el; m++) {
556  ssht_sampling_elm2ind(&ind, el, m);
557  ssht_sampling_elm2ind(&ind_nm, el, -m);
558  flm[ind_nm] = signs[m] * conj(flm[ind]);
559  }
560  }
561 
562  // Free memory.
563  free(dl);
564  if (dl_method == SSHT_DL_RISBO)
565  free(dl8);
566  free(Fmt);
567  free(Fmm);
568  free(sqrt_tbl);
569  free(signs);
570  free(expsm);
571  free(expsmm);
572  free(inds);
573 
574  // Print finished if verbosity set.
575  if (verbosity > 0)
576  printf("%s%s", SSHT_PROMPT, "Adjoint inverse transform computed!");
577 
578 }
579 
580 
597 void ssht_adjoint_mw_forward_sov_sym(SSHT_COMPLEX(double) *f, const SSHT_COMPLEX(double) *flm,
598  int L, int spin,
599  ssht_dl_method_t dl_method,
600  int verbosity) {
601 
602  int el, m, mm, ind;
603  int eltmp;
604  double *sqrt_tbl, *signs;
605  int el2pel, inds_offset;
606  int *inds;
607  double ssign, elfactor;
608  SSHT_COMPLEX(double) mmfactor;
609  double *dl;
610  double *dl8 = NULL;
611  int dl_offset, dl_stride;
612  SSHT_COMPLEX(double) *exps;
613  int exps_offset;
614  double elmmsign, elssign;
615  int spinneg;
616  SSHT_COMPLEX(double) *Fmm;
617  int Fmm_offset, Fmm_stride;
618 
619  fftw_plan plan, plan_bwd, plan_fwd;
620  SSHT_COMPLEX(double) *Ftm, *Gmm;
621  SSHT_COMPLEX(double) *w, *wr;
622  int w_offset;
623  SSHT_COMPLEX(double) *Fmm_pad, *tmp_pad;
624  int f_stride, Ftm_stride, Ftm_offset;
625  int r, t, p;
626  SSHT_COMPLEX(double) *inout;
627 
628  // Allocate memory.
629  sqrt_tbl = (double*)calloc(2*(L-1)+2, sizeof(double));
631  signs = (double*)calloc(L+1, sizeof(double));
633  exps = (SSHT_COMPLEX(double)*)calloc(2*L-1, sizeof(SSHT_COMPLEX(double)));
635  inds = (int*)calloc(2*L-1, sizeof(int));
637 
638  // Perform precomputations.
639  for (el=0; el<=2*(L-1)+1; el++)
640  sqrt_tbl[el] = sqrt((double)el);
641  for (m=0; m<=L-1; m=m+2) {
642  signs[m] = 1.0;
643  signs[m+1] = -1.0;
644  }
645  ssign = signs[abs(spin)];
646  spinneg = spin <= 0 ? spin : -spin;
647  exps_offset = L-1;
648  for (m=-(L-1); m<=L-1; m++)
649  exps[m + exps_offset] = cexp(-I*SSHT_PION2*(m+spin));
650 
651  // Print messages depending on verbosity level.
652  if (verbosity > 0) {
653  printf("%s%s\n", SSHT_PROMPT,
654  "Computing adjoint forward transform using MW sampling with ");
655  printf("%s%s%d%s%d%s\n", SSHT_PROMPT, "parameters (L,spin,reality) = (",
656  L, ",", spin, ", FALSE)");
657  if (verbosity > 1)
658  printf("%s%s\n", SSHT_PROMPT,
659  "Using routine ssht_adjoint_mw_forward_sov_sym...");
660  }
661 
662  // Compute Fmm.
663  Fmm = (SSHT_COMPLEX(double)*)calloc((2*L-1)*(2*L-1), sizeof(SSHT_COMPLEX(double)));
665  Fmm_offset = L-1;
666  Fmm_stride = 2*L-1;
669  if (dl_method == SSHT_DL_RISBO) {
672  }
673  dl_offset = ssht_dl_get_offset(L, SSHT_DL_QUARTER);
674  dl_stride = ssht_dl_get_stride(L, SSHT_DL_QUARTER);
675  inds_offset = L-1;
676  for (el=abs(spin); el<=L-1; el++) {
677 
678  // Compute Wigner plane.
679  switch (dl_method) {
680 
681  case SSHT_DL_RISBO:
682  if (el!=0 && el==abs(spin)) {
683  for(eltmp=0; eltmp<=abs(spin); eltmp++)
686  eltmp, sqrt_tbl, signs);
688  dl8, L,
691  el,
692  signs);
693  }
694  else {
697  el, sqrt_tbl, signs);
699  dl8, L,
702  el,
703  signs);
704  }
705  break;
706 
707  case SSHT_DL_TRAPANI:
708  if (el!=0 && el==abs(spin)) {
709  for(eltmp=0; eltmp<=abs(spin); eltmp++)
712  eltmp, sqrt_tbl);
715  el, signs);
716  }
717  else {
720  el, sqrt_tbl);
723  el, signs);
724  }
725  break;
726 
727  default:
728  SSHT_ERROR_GENERIC("Invalid dl method")
729  }
730 
731  // Compute Fmm.
732  elfactor = sqrt((double)(2.0*el+1.0)/(4.0*SSHT_PI));
733  el2pel = el *el + el;
734  for (m=-el; m<=el; m++)
735  inds[m + inds_offset] = el2pel + m;
736  for (mm=0; mm<=el; mm++) {
737  elmmsign = signs[el] * signs[mm];
738  elssign = spin <= 0 ? 1.0 : elmmsign;
739 
740  for (m=-el; m<=-1; m++) {
741  ind = inds[m + inds_offset];
742  Fmm[(m + Fmm_offset)*Fmm_stride + mm + Fmm_offset] +=
743  ssign
744  * elfactor
745  * exps[m + exps_offset]
746  * elmmsign * dl[mm*dl_stride - m + dl_offset]
747  * elssign * dl[mm*dl_stride - spinneg + dl_offset]
748  * flm[ind];
749  }
750  for (m=0; m<=el; m++) {
751  ind = inds[m + inds_offset];
752  Fmm[(m + Fmm_offset)*Fmm_stride + mm + Fmm_offset] +=
753  ssign
754  * elfactor
755  * exps[m + exps_offset]
756  * dl[mm*dl_stride + m + dl_offset]
757  * elssign * dl[mm*dl_stride - spinneg + dl_offset]
758  * flm[ind];
759  }
760 
761  }
762 
763  }
764 
765  // Free dl memory.
766  free(dl);
767  if (dl_method == SSHT_DL_RISBO)
768  free(dl8);
769 
770  // Use symmetry to compute Fmm for negative mm.
771  for (m=-(L-1); m<=L-1; m++)
772  for (mm=-(L-1); mm<=-1; mm++)
773  Fmm[(m + Fmm_offset)*Fmm_stride + mm + Fmm_offset] =
774  signs[abs(m)] * ssign
775  * Fmm[(m + Fmm_offset)*Fmm_stride - mm + Fmm_offset];
776 
777  // Compute weights.
778  w = (SSHT_COMPLEX(double)*)calloc(4*L-3, sizeof(SSHT_COMPLEX(double)));
780  w_offset = 2*(L-1);
781  for (mm=-2*(L-1); mm<=2*(L-1); mm++)
782  w[mm+w_offset] = ssht_sampling_weight_mw(mm);
783 
784  // Compute IFFT of w to give wr.
785  wr = (SSHT_COMPLEX(double)*)calloc(4*L-3, sizeof(SSHT_COMPLEX(double)));
787  inout = (SSHT_COMPLEX(double)*)calloc(4*L-3, sizeof(SSHT_COMPLEX(double)));
789  plan_bwd = fftw_plan_dft_1d(4*L-3, inout, inout, FFTW_BACKWARD, FFTW_MEASURE);
790  plan_fwd = fftw_plan_dft_1d(4*L-3, inout, inout, FFTW_FORWARD, FFTW_MEASURE);
791  for (mm=1; mm<=2*L-2; mm++)
792  inout[mm + w_offset] = w[mm - 2*(L-1) - 1 + w_offset];
793  for (mm=-2*(L-1); mm<=0; mm++)
794  inout[mm + w_offset] = w[mm + 2*(L-1) + w_offset];
795  fftw_execute_dft(plan_bwd, inout, inout);
796  for (mm=0; mm<=2*L-2; mm++)
797  wr[mm + w_offset] = inout[mm - 2*(L-1) + w_offset];
798  for (mm=-2*(L-1); mm<=-1; mm++)
799  wr[mm + w_offset] = inout[mm + 2*(L-1) + 1 + w_offset];
800 
801  // Compute Gmm by convolution implemented as product in real space.
802  Fmm_pad = (SSHT_COMPLEX(double)*)calloc(4*L-3, sizeof(SSHT_COMPLEX(double)));
804  tmp_pad = (SSHT_COMPLEX(double)*)calloc(4*L-3, sizeof(SSHT_COMPLEX(double)));
806  Gmm = (SSHT_COMPLEX(double)*)calloc((2*L-1)*(2*L-1), sizeof(SSHT_COMPLEX(double)));
808  for (m=-(L-1); m<=L-1; m++) {
809 
810  // Zero-pad Fmm.
811  for (mm=-2*(L-1); mm<=-L; mm++)
812  Fmm_pad[mm+w_offset] = 0.0;
813  for (mm=L; mm<=2*(L-1); mm++)
814  Fmm_pad[mm+w_offset] = 0.0;
815  for (mm=-(L-1); mm<=L-1; mm++)
816  Fmm_pad[mm + w_offset] =
817  Fmm[(m+Fmm_offset)*Fmm_stride + mm + Fmm_offset];
818 
819  // Compute IFFT of Fmm.
820  for (mm=1; mm<=2*L-2; mm++)
821  inout[mm + w_offset] = Fmm_pad[mm - 2*(L-1) - 1 + w_offset];
822  for (mm=-2*(L-1); mm<=0; mm++)
823  inout[mm + w_offset] = Fmm_pad[mm + 2*(L-1) + w_offset];
824  fftw_execute_dft(plan_bwd, inout, inout);
825  for (mm=0; mm<=2*L-2; mm++)
826  Fmm_pad[mm + w_offset] = inout[mm - 2*(L-1) + w_offset];
827  for (mm=-2*(L-1); mm<=-1; mm++)
828  Fmm_pad[mm + w_offset] = inout[mm + 2*(L-1) + 1 + w_offset];
829 
830  // Compute product of Fmm and weight in real space.
831  for (r=-2*(L-1); r<=2*(L-1); r++)
832  Fmm_pad[r + w_offset] *= wr[-r + w_offset];
833 
834  // Compute Gmm by FFT.
835  for (mm=1; mm<=2*L-2; mm++)
836  inout[mm + w_offset] = Fmm_pad[mm - 2*(L-1) - 1 + w_offset];
837  for (mm=-2*(L-1); mm<=0; mm++)
838  inout[mm + w_offset] = Fmm_pad[mm + 2*(L-1) + w_offset];
839  fftw_execute_dft(plan_fwd, inout, inout);
840  for (mm=0; mm<=2*L-2; mm++)
841  Fmm_pad[mm + w_offset] = inout[mm - 2*(L-1) + w_offset];
842  for (mm=-2*(L-1); mm<=-1; mm++)
843  Fmm_pad[mm + w_offset] = inout[mm + 2*(L-1) + 1 + w_offset];
844 
845  // Extract section of Gmm of interest.
846  for (mm=-(L-1); mm<=L-1; mm++)
847  Gmm[(m+Fmm_offset)*Fmm_stride + mm + Fmm_offset] =
848  Fmm_pad[mm + w_offset] * 2.0 * SSHT_PI / (4.0*L-3.0);
849 
850  }
851  fftw_destroy_plan(plan_bwd);
852  fftw_destroy_plan(plan_fwd);
853  free(inout);
854 
855  // Apply phase modulation to account for sampling offset.
856  for (mm=-(L-1); mm<=L-1; mm++) {
857  mmfactor = cexp(I*mm*SSHT_PI/(2.0*L-1.0));
858  for (m=-(L-1); m<=L-1; m++)
859  Gmm[(m + Fmm_offset)*Fmm_stride + mm + Fmm_offset] *=
860  mmfactor;
861  }
862 
863  // Compute Fourier transform over theta.
864  inout = (SSHT_COMPLEX(double)*)calloc(2*L-1, sizeof(SSHT_COMPLEX(double)));
866  plan = fftw_plan_dft_1d(2*L-1, inout, inout, FFTW_BACKWARD, FFTW_MEASURE);
867  Ftm = (SSHT_COMPLEX(double)*)calloc((2*L-1)*(2*L-1), sizeof(SSHT_COMPLEX(double)));
869  Ftm_stride = 2*L-1;
870  Ftm_offset = L-1;
871  for (m=-(L-1); m<=L-1; m++) {
872 
873  for(mm=0; mm<=L-1; mm++)
874  inout[mm] = Gmm[(m+Fmm_offset)*Fmm_stride + mm + Fmm_offset];
875  for(mm=-(L-1); mm<=-1; mm++)
876  inout[mm+2*L-1] = Gmm[(m+Fmm_offset)*Fmm_stride + mm + Fmm_offset];
877  fftw_execute_dft(plan, inout, inout);
878 
879  for(t=0; t<=2*L-2; t++)
880  Ftm[t*Ftm_stride + m + Ftm_offset] = inout[t] / (2.0*L-1.0);
881 
882  }
883 
884  // Adjoint of periodic extension of Ftm.
885  for(t=0; t<=L-2; t++)
886  for (m=-(L-1); m<=L-1; m++)
887  Ftm[t*Ftm_stride + m + Ftm_offset] =
888  Ftm[t*Ftm_stride + m + Ftm_offset]
889  + signs[abs(m)] * ssign * Ftm[(2*L-2-t)*Ftm_stride + m + Ftm_offset];
890 
891  // Compute Fourier transform over phi.
892  f_stride = 2*L-1;
893  for(t=0; t<=L-1; t++) {
894 
895  for(m=0; m<=L-1; m++)
896  inout[m] = Ftm[t*Ftm_stride + m + Ftm_offset];
897  for(m=-(L-1); m<=-1; m++)
898  inout[m+2*L-1] = Ftm[t*Ftm_stride + m + Ftm_offset];
899  fftw_execute_dft(plan, inout, inout);
900 
901  for(p=0; p<=2*L-2; p++)
902  f[t*f_stride + p] = inout[p] / (2.0*L-1.0);
903 
904  }
905  fftw_destroy_plan(plan);
906  free(inout);
907 
908  // Free memory.
909  free(Fmm);
910  free(Ftm);
911  free(w);
912  free(wr);
913  free(Fmm_pad);
914  free(tmp_pad);
915  free(Gmm);
916 
917  // Free precomputation memory.
918  free(sqrt_tbl);
919  free(signs);
920  free(exps);
921  free(inds);
922 
923  // Print finished if verbosity set.
924  if (verbosity > 0)
925  printf("%s%s", SSHT_PROMPT, "Adjoint forward transform computed!");
926 
927 }
928 
929 
947  const SSHT_COMPLEX(double) *flm,
948  int L,
949  ssht_dl_method_t dl_method,
950  int verbosity) {
951 
952  int el, m, mm, ind;
953  int eltmp;
954  double *sqrt_tbl, *signs;
955  int el2pel, inds_offset;
956  int *inds;
957  double ssign, elfactor;
958  SSHT_COMPLEX(double) mmfactor;
959  double *dl;
960  double *dl8 = NULL;
961  int dl_offset, dl_stride;
962  SSHT_COMPLEX(double) *exps;
963  int exps_offset;
964  double elmmsign, elssign;
965  int spinneg;
966  SSHT_COMPLEX(double) *Fmm;
967  int Fmm_offset, Fmm_stride;
968 
969  fftw_plan plan, plan_bwd, plan_fwd;
970  SSHT_COMPLEX(double) *Ftm, *Gmm;
971  SSHT_COMPLEX(double) *w, *wr;
972  int w_offset;
973  SSHT_COMPLEX(double) *Fmm_pad, *tmp_pad;
974  int f_stride, Ftm_stride;
975  int r, t, p;
976  SSHT_COMPLEX(double) *inout;
977  SSHT_COMPLEX(double) *in;
978  double *out_real;
979  int Gmm_stride;
980  int spin = 0;
981 
982  // Allocate memory.
983  sqrt_tbl = (double*)calloc(2*(L-1)+2, sizeof(double));
985  signs = (double*)calloc(L+1, sizeof(double));
987  exps = (SSHT_COMPLEX(double)*)calloc(2*L-1, sizeof(SSHT_COMPLEX(double)));
989  inds = (int*)calloc(2*L-1, sizeof(int));
991 
992  // Perform precomputations.
993  for (el=0; el<=2*(L-1)+1; el++)
994  sqrt_tbl[el] = sqrt((double)el);
995  for (m=0; m<=L-1; m=m+2) {
996  signs[m] = 1.0;
997  signs[m+1] = -1.0;
998  }
999  ssign = signs[abs(spin)];
1000  spinneg = spin <= 0 ? spin : -spin;
1001  exps_offset = L-1;
1002  for (m=-(L-1); m<=L-1; m++)
1003  exps[m + exps_offset] = cexp(-I*SSHT_PION2*(m+spin));
1004 
1005  // Print messages depending on verbosity level.
1006  if (verbosity > 0) {
1007  printf("%s%s\n", SSHT_PROMPT,
1008  "Computing adjoint forward transform using MW sampling with ");
1009  printf("%s%s%d%s%d%s\n", SSHT_PROMPT, "parameters (L,spin,reality) = (",
1010  L, ",", spin, ", FALSE)");
1011  if (verbosity > 1)
1012  printf("%s%s\n", SSHT_PROMPT,
1013  "Using routine ssht_adjoint_mw_forward_sov_sym_real...");
1014  }
1015 
1016  // Compute Fmm.
1017  Fmm = (SSHT_COMPLEX(double)*)calloc(L*(2*L-1), sizeof(SSHT_COMPLEX(double)));
1019  Fmm_offset = L-1;
1020  Fmm_stride = 2*L-1;
1023  if (dl_method == SSHT_DL_RISBO) {
1026  }
1027  dl_offset = ssht_dl_get_offset(L, SSHT_DL_QUARTER);
1028  dl_stride = ssht_dl_get_stride(L, SSHT_DL_QUARTER);
1029  inds_offset = L-1;
1030  for (el=abs(spin); el<=L-1; el++) {
1031 
1032  // Compute Wigner plane.
1033  switch (dl_method) {
1034 
1035  case SSHT_DL_RISBO:
1036  if (el!=0 && el==abs(spin)) {
1037  for(eltmp=0; eltmp<=abs(spin); eltmp++)
1040  eltmp, sqrt_tbl, signs);
1042  dl8, L,
1045  el,
1046  signs);
1047  }
1048  else {
1051  el, sqrt_tbl, signs);
1053  dl8, L,
1056  el,
1057  signs);
1058  }
1059  break;
1060 
1061  case SSHT_DL_TRAPANI:
1062  if (el!=0 && el==abs(spin)) {
1063  for(eltmp=0; eltmp<=abs(spin); eltmp++)
1066  eltmp, sqrt_tbl);
1069  el, signs);
1070  }
1071  else {
1074  el, sqrt_tbl);
1077  el, signs);
1078  }
1079  break;
1080 
1081  default:
1082  SSHT_ERROR_GENERIC("Invalid dl method")
1083  }
1084 
1085  // Compute Fmm.
1086  elfactor = sqrt((double)(2.0*el+1.0)/(4.0*SSHT_PI));
1087  el2pel = el *el + el;
1088  for (m=0; m<=el; m++)
1089  inds[m + inds_offset] = el2pel + m;
1090  for (mm=0; mm<=el; mm++) {
1091  elmmsign = signs[el] * signs[mm];
1092  elssign = spin <= 0 ? 1.0 : elmmsign;
1093 
1094  for (m=0; m<=el; m++) {
1095  ind = inds[m + inds_offset];
1096  Fmm[m*Fmm_stride + mm + Fmm_offset] +=
1097  ssign
1098  * elfactor
1099  * exps[m + exps_offset]
1100  * dl[mm*dl_stride + m + dl_offset]
1101  * elssign * dl[mm*dl_stride - spinneg + dl_offset]
1102  * flm[ind];
1103  }
1104 
1105  }
1106 
1107  }
1108 
1109  // Free dl memory.
1110  free(dl);
1111  if (dl_method == SSHT_DL_RISBO)
1112  free(dl8);
1113 
1114  // Use symmetry to compute Fmm for negative mm.
1115  for (m=0; m<=L-1; m++)
1116  for (mm=-(L-1); mm<=-1; mm++)
1117  Fmm[m*Fmm_stride + mm + Fmm_offset] =
1118  signs[abs(m)] * ssign
1119  * Fmm[m*Fmm_stride - mm + Fmm_offset];
1120 
1121  // Compute weights.
1122  w = (SSHT_COMPLEX(double)*)calloc(4*L-3, sizeof(SSHT_COMPLEX(double)));
1124  w_offset = 2*(L-1);
1125  for (mm=-2*(L-1); mm<=2*(L-1); mm++)
1126  w[mm+w_offset] = ssht_sampling_weight_mw(mm);
1127 
1128  // Compute IFFT of w to give wr.
1129  wr = (SSHT_COMPLEX(double)*)calloc(4*L-3, sizeof(SSHT_COMPLEX(double)));
1131  inout = (SSHT_COMPLEX(double)*)calloc(4*L-3, sizeof(SSHT_COMPLEX(double)));
1133  plan_bwd = fftw_plan_dft_1d(4*L-3, inout, inout, FFTW_BACKWARD, FFTW_MEASURE);
1134  plan_fwd = fftw_plan_dft_1d(4*L-3, inout, inout, FFTW_FORWARD, FFTW_MEASURE);
1135  for (mm=1; mm<=2*L-2; mm++)
1136  inout[mm + w_offset] = w[mm - 2*(L-1) - 1 + w_offset];
1137  for (mm=-2*(L-1); mm<=0; mm++)
1138  inout[mm + w_offset] = w[mm + 2*(L-1) + w_offset];
1139  fftw_execute_dft(plan_bwd, inout, inout);
1140  for (mm=0; mm<=2*L-2; mm++)
1141  wr[mm + w_offset] = inout[mm - 2*(L-1) + w_offset];
1142  for (mm=-2*(L-1); mm<=-1; mm++)
1143  wr[mm + w_offset] = inout[mm + 2*(L-1) + 1 + w_offset];
1144 
1145  // Compute Gmm by convolution implemented as product in real space.
1146  Fmm_pad = (SSHT_COMPLEX(double)*)calloc(4*L-3, sizeof(SSHT_COMPLEX(double)));
1148  tmp_pad = (SSHT_COMPLEX(double)*)calloc(4*L-3, sizeof(SSHT_COMPLEX(double)));
1150  Gmm = (SSHT_COMPLEX(double)*)calloc(L*(2*L-1), sizeof(SSHT_COMPLEX(double)));
1151  Gmm_stride = 2*L-1;
1153  for (m=0; m<=L-1; m++) {
1154 
1155  // Zero-pad Fmm.
1156  for (mm=-2*(L-1); mm<=-L; mm++)
1157  Fmm_pad[mm+w_offset] = 0.0;
1158  for (mm=L; mm<=2*(L-1); mm++)
1159  Fmm_pad[mm+w_offset] = 0.0;
1160  for (mm=-(L-1); mm<=L-1; mm++)
1161  Fmm_pad[mm + w_offset] =
1162  Fmm[m*Fmm_stride + mm + Fmm_offset];
1163 
1164  // Compute IFFT of Fmm.
1165  for (mm=1; mm<=2*L-2; mm++)
1166  inout[mm + w_offset] = Fmm_pad[mm - 2*(L-1) - 1 + w_offset];
1167  for (mm=-2*(L-1); mm<=0; mm++)
1168  inout[mm + w_offset] = Fmm_pad[mm + 2*(L-1) + w_offset];
1169  fftw_execute_dft(plan_bwd, inout, inout);
1170  for (mm=0; mm<=2*L-2; mm++)
1171  Fmm_pad[mm + w_offset] = inout[mm - 2*(L-1) + w_offset];
1172  for (mm=-2*(L-1); mm<=-1; mm++)
1173  Fmm_pad[mm + w_offset] = inout[mm + 2*(L-1) + 1 + w_offset];
1174 
1175  // Compute product of Fmm and weight in real space.
1176  for (r=-2*(L-1); r<=2*(L-1); r++)
1177  Fmm_pad[r + w_offset] *= wr[-r + w_offset];
1178 
1179  // Compute Gmm by FFT.
1180  for (mm=1; mm<=2*L-2; mm++)
1181  inout[mm + w_offset] = Fmm_pad[mm - 2*(L-1) - 1 + w_offset];
1182  for (mm=-2*(L-1); mm<=0; mm++)
1183  inout[mm + w_offset] = Fmm_pad[mm + 2*(L-1) + w_offset];
1184  fftw_execute_dft(plan_fwd, inout, inout);
1185  for (mm=0; mm<=2*L-2; mm++)
1186  Fmm_pad[mm + w_offset] = inout[mm - 2*(L-1) + w_offset];
1187  for (mm=-2*(L-1); mm<=-1; mm++)
1188  Fmm_pad[mm + w_offset] = inout[mm + 2*(L-1) + 1 + w_offset];
1189 
1190  // Extract section of Gmm of interest.
1191  for (mm=-(L-1); mm<=L-1; mm++)
1192  Gmm[m*Gmm_stride + mm + Fmm_offset] =
1193  Fmm_pad[mm + w_offset] * 2.0 * SSHT_PI / (4.0*L-3.0);
1194 
1195  }
1196  fftw_destroy_plan(plan_bwd);
1197  fftw_destroy_plan(plan_fwd);
1198  free(inout);
1199 
1200  // Apply phase modulation to account for sampling offset.
1201  for (mm=-(L-1); mm<=L-1; mm++) {
1202  mmfactor = cexp(I*mm*SSHT_PI/(2.0*L-1.0));
1203  for (m=0; m<=L-1; m++)
1204  Gmm[m*Gmm_stride + mm + Fmm_offset] *=
1205  mmfactor;
1206  }
1207 
1208  // Compute Fourier transform over theta.
1209  inout = (SSHT_COMPLEX(double)*)calloc(2*L-1, sizeof(SSHT_COMPLEX(double)));
1211  plan = fftw_plan_dft_1d(2*L-1, inout, inout, FFTW_BACKWARD, FFTW_MEASURE);
1212  Ftm = (SSHT_COMPLEX(double)*)calloc((2*L-1)*L, sizeof(SSHT_COMPLEX(double)));
1214  Ftm_stride = L;
1215  for (m=0; m<=L-1; m++) {
1216 
1217  for(mm=0; mm<=L-1; mm++)
1218  inout[mm] = Gmm[m*Gmm_stride + mm + Fmm_offset];
1219  for(mm=-(L-1); mm<=-1; mm++)
1220  inout[mm+2*L-1] = Gmm[m*Gmm_stride + mm + Fmm_offset];
1221  fftw_execute_dft(plan, inout, inout);
1222 
1223  for(t=0; t<=2*L-2; t++)
1224  Ftm[t*Ftm_stride + m] = inout[t] / (2.0*L-1.0);
1225 
1226  }
1227  fftw_destroy_plan(plan);
1228  free(inout);
1229 
1230  // Adjoint of periodic extension of Ftm.
1231  for(t=0; t<=L-2; t++)
1232  for (m=0; m<=L-1; m++)
1233  Ftm[t*Ftm_stride + m] =
1234  Ftm[t*Ftm_stride + m]
1235  + signs[abs(m)] * ssign * Ftm[(2*L-2-t)*Ftm_stride + m];
1236 
1237  // Compute Fourier transform over phi.
1238  out_real = (double*)calloc(2*L-1, sizeof(double));
1239  SSHT_ERROR_MEM_ALLOC_CHECK(out_real)
1240  in = (SSHT_COMPLEX(double)*)calloc(L, sizeof(SSHT_COMPLEX(double)));
1242  plan = fftw_plan_dft_c2r_1d(2*L-1, in, out_real, FFTW_MEASURE);
1243  f_stride = 2*L-1;
1244  for(t=0; t<=L-1; t++) {
1245 
1246  memcpy(in, &Ftm[t*Ftm_stride], Ftm_stride*sizeof(SSHT_COMPLEX(double)));
1247  fftw_execute_dft_c2r(plan, in, out_real);
1248 
1249  for(p=0; p<=2*L-2; p++)
1250  f[t*f_stride + p] = out_real[p] / (2.0*L-1.0);
1251 
1252  }
1253  fftw_destroy_plan(plan);
1254  free(out_real);
1255  free(in);
1256 
1257  // Free memory.
1258  free(Fmm);
1259  free(Ftm);
1260  free(w);
1261  free(wr);
1262  free(Fmm_pad);
1263  free(tmp_pad);
1264  free(Gmm);
1265 
1266  // Free precomputation memory.
1267  free(sqrt_tbl);
1268  free(signs);
1269  free(exps);
1270  free(inds);
1271 
1272  // Print finished if verbosity set.
1273  if (verbosity > 0)
1274  printf("%s%s", SSHT_PROMPT, "Adjoint forward transform computed!");
1275 
1276 }
1277 
1278 
1279 //============================================================================
1280 // MW South pole interfaces
1281 //============================================================================
1282 
1283 
1305  SSHT_COMPLEX(double) *f_sp, double *phi_sp,
1306  SSHT_COMPLEX(double) *flm,
1307  int L, int spin,
1308  ssht_dl_method_t dl_method,
1309  int verbosity) {
1310 
1311  SSHT_COMPLEX(double)* f_full;
1312  int f_stride = 2*L-1;
1313 
1314  // Allocate full array.
1315  f_full = (SSHT_COMPLEX(double)*)calloc(L*(2*L-1), sizeof(SSHT_COMPLEX(double)));
1317 
1318  // Perform inverse transform.
1319  ssht_adjoint_mw_forward_sov_sym(f_full, flm, L, spin,
1320  dl_method, verbosity);
1321 
1322  // Copy output function values, including separate point for South pole.
1323  memcpy(f, f_full, (L-1)*(2*L-1)*sizeof(SSHT_COMPLEX(double)));
1324  *f_sp = f_full[(L-1)*f_stride + 0];
1325  *phi_sp = ssht_sampling_mw_p2phi(0, L);
1326 
1327  // Free memory.
1328  free(f_full);
1329 
1330 }
1331 
1332 
1351  double *f_sp,
1352  SSHT_COMPLEX(double) *flm,
1353  int L,
1354  ssht_dl_method_t dl_method,
1355  int verbosity) {
1356 
1357  double *f_full;
1358  int f_stride = 2*L-1;
1359 
1360  // Allocate full array.
1361  f_full = (double*)calloc(L*(2*L-1), sizeof(double));
1363 
1364  // Perform inverse transform.
1365  ssht_adjoint_mw_forward_sov_sym_real(f_full, flm, L,
1366  dl_method, verbosity);
1367 
1368  // Copy output function values, including separate point for South pole.
1369  memcpy(f, f_full, (L-1)*(2*L-1)*sizeof(double));
1370  *f_sp = f_full[(L-1)*f_stride + 0];
1371 
1372  // Free memory.
1373  free(f_full);
1374 
1375 }
1376 
1377 
1399  SSHT_COMPLEX(double) f_sp, double phi_sp,
1400  int L, int spin,
1401  ssht_dl_method_t dl_method,
1402  int verbosity) {
1403 
1404  SSHT_COMPLEX(double) *f_full;
1405  int p, f_stride = 2*L-1;
1406  double phi;
1407 
1408  // Copy function values to full array.
1409  f_full = (SSHT_COMPLEX(double)*)calloc(L*(2*L-1), sizeof(SSHT_COMPLEX(double)));
1411  memcpy(f_full, f, (L-1)*(2*L-1)*sizeof(SSHT_COMPLEX(double)));
1412 
1413  // Define South pole for all phi.
1414  for (p=0; p<=2*L-2; p++) {
1415  phi = ssht_sampling_mw_p2phi(p, L);
1416  f_full[(L-1)*f_stride + p] = f_sp * cexp(I*spin*(phi-phi_sp));
1417  }
1418 
1419  // Perform forward transform.
1420  ssht_adjoint_mw_inverse_sov_sym(flm, f_full, L, spin,
1421  dl_method, verbosity);
1422 
1423  // Free memory.
1424  free(f_full);
1425 
1426 }
1427 
1428 
1447  double *f,
1448  double f_sp,
1449  int L,
1450  ssht_dl_method_t dl_method,
1451  int verbosity) {
1452 
1453  double *f_full;
1454  int p, f_stride = 2*L-1;
1455 
1456  // Copy function values to full array.
1457  f_full = (double*)calloc(L*(2*L-1), sizeof(double));
1459  memcpy(f_full, f, (L-1)*(2*L-1)*sizeof(double));
1460 
1461  // Define South pole for all phi.
1462  for (p=0; p<=2*L-2; p++)
1463  f_full[(L-1)*f_stride + p] = f_sp;
1464 
1465  // Perform forward transform.
1466  ssht_adjoint_mw_inverse_sov_sym_real(flm, f_full, L,
1467  dl_method, verbosity);
1468 
1469  // Free memory.
1470  free(f_full);
1471 
1472 }
1473 
1474 
1475 //============================================================================
1476 // MW SS algorithms
1477 //============================================================================
1478 
1479 
1496  int L, int spin,
1497  ssht_dl_method_t dl_method,
1498  int verbosity) {
1499 
1500  int el, m, mm, ind, t;
1501  int eltmp;
1502  double *sqrt_tbl, *signs;
1503  int el2pel, inds_offset;
1504  int *inds;
1505  double ssign, elfactor;
1506  fftw_plan plan;
1507  SSHT_COMPLEX(double) *inout;
1508  SSHT_COMPLEX(double) *Fmt, *Fmm;
1509  int f_stride, Fmt_stride, Fmt_offset, Fmm_stride, Fmm_offset;
1510  double *dl;
1511  double *dl8 = NULL;
1512  int dl_offset, dl_stride;
1513  SSHT_COMPLEX(double) *expsm, *expsmm;
1514  int exps_offset;
1515  int elmmsign, elssign;
1516  int spinneg;
1517 
1518  // Allocate memory.
1519  sqrt_tbl = (double*)calloc(2*(L-1)+2, sizeof(double));
1520  SSHT_ERROR_MEM_ALLOC_CHECK(sqrt_tbl)
1521  signs = (double*)calloc(L+1, sizeof(double));
1523  expsm = (SSHT_COMPLEX(double)*)calloc(2*L-1, sizeof(SSHT_COMPLEX(double)));
1525  expsmm = (SSHT_COMPLEX(double)*)calloc(2*L-1, sizeof(SSHT_COMPLEX(double)));
1527  inds = (int*)calloc(2*L-1, sizeof(int));
1529 
1530  // Perform precomputations.
1531  for (el=0; el<=2*(L-1)+1; el++)
1532  sqrt_tbl[el] = sqrt((double)el);
1533  for (m=0; m<=L-1; m=m+2) {
1534  signs[m] = 1.0;
1535  signs[m+1] = -1.0;
1536  }
1537  ssign = signs[abs(spin)];
1538  spinneg = spin <= 0 ? spin : -spin;
1539  exps_offset = L-1;
1540  for (m=-(L-1); m<=L-1; m++)
1541  expsm[m + exps_offset] = cexp(I*SSHT_PION2*(m+spin));
1542  for (mm=-(L-1); mm<=L-1; mm++)
1543  expsmm[mm + exps_offset] = cexp(-I*mm*SSHT_PI/(2.0*L-1.0));
1544 
1545  // Print messages depending on verbosity level.
1546  if (verbosity > 0) {
1547  printf("%s%s\n", SSHT_PROMPT,
1548  "Computing adjoint inverse transform using MW symmetric sampling with ");
1549  printf("%s%s%d%s%d%s\n", SSHT_PROMPT, "parameters (L,spin,reality) = (",
1550  L, ",", spin, ", FALSE)");
1551  if (verbosity > 1)
1552  printf("%s%s\n", SSHT_PROMPT,
1553  "Using routine ssht_adjoint_mw_inverse_sov_sym_ss...");
1554  }
1555 
1556  // Compute Fourier transform over phi, i.e. compute Fmt.
1557  // Note that t and p indices of fext are increased in size by
1558  // one compared to usual sampling.
1559  Fmt = (SSHT_COMPLEX(double)*)calloc((2*L)*(2*L), sizeof(SSHT_COMPLEX(double)));
1561  Fmt_stride = 2*L;
1562  Fmt_offset = L-1;
1563  f_stride = 2*L;
1564  inout = (SSHT_COMPLEX(double)*)calloc(2*L, sizeof(SSHT_COMPLEX(double)));
1566  plan = fftw_plan_dft_1d(2*L, inout, inout, FFTW_FORWARD, FFTW_MEASURE);
1567  for (t=0; t<=L; t++) {
1568  memcpy(inout, &f[t*f_stride], f_stride*sizeof(SSHT_COMPLEX(double)));
1569  fftw_execute_dft(plan, inout, inout);
1570  for(m=0; m<=L; m++)
1571  Fmt[(m+Fmt_offset)*Fmt_stride + t] = inout[m];
1572  for(m=-(L-1); m<=-1; m++)
1573  Fmt[(m+Fmt_offset)*Fmt_stride + t] = inout[m+2*L-1+1];
1574  }
1575  fftw_destroy_plan(plan);
1576  free(inout);
1577 
1578  // Apply adjoint of periodic extension.
1579  for (m=-(L-1); m<=L; m++)
1580  for (t=L+1; t<=2*L-1; t++)
1581  Fmt[(m+Fmt_offset)*Fmt_stride + t] = 0.0;
1582 
1583  // Compute Fourier transform over theta, i.e. compute Fmm.
1584  // Note that m and mm indices are increased in size by one.
1585  Fmm = (SSHT_COMPLEX(double)*)calloc((2*L)*(2*L), sizeof(SSHT_COMPLEX(double)));
1587  Fmm_stride = 2*L;
1588  Fmm_offset = L-1;
1589  inout = (SSHT_COMPLEX(double)*)calloc(2*L, sizeof(SSHT_COMPLEX(double)));
1591  plan = fftw_plan_dft_1d(2*L, inout, inout, FFTW_FORWARD, FFTW_MEASURE);
1592  for (m=-(L-1); m<=L; m++) {
1593  memcpy(inout, &Fmt[(m+Fmt_offset)*Fmt_stride], Fmt_stride*sizeof(SSHT_COMPLEX(double)));
1594  fftw_execute_dft(plan, inout, inout);
1595  for(mm=0; mm<=L; mm++)
1596  Fmm[(mm+Fmm_offset)*Fmm_stride + m + Fmm_offset] =
1597  inout[mm];
1598  for(mm=-(L-1); mm<=-1; mm++)
1599  Fmm[(mm+Fmm_offset)*Fmm_stride + m + Fmm_offset] =
1600  inout[mm+2*L-1+1];
1601  }
1602  fftw_destroy_plan(plan);
1603  free(inout);
1604 
1605  // Compute flm.
1608  if (dl_method == SSHT_DL_RISBO) {
1611  }
1612  dl_offset = ssht_dl_get_offset(L, SSHT_DL_QUARTER);
1613  dl_stride = ssht_dl_get_stride(L, SSHT_DL_QUARTER);
1614  inds_offset = L-1;
1615  for (el=0; el<=L-1; el++) {
1616  for (m=-el; m<=el; m++) {
1617  ssht_sampling_elm2ind(&ind, el, m);
1618  flm[ind] = 0.0;
1619  }
1620  }
1621  for (el=abs(spin); el<=L-1; el++) {
1622 
1623  // Compute Wigner plane.
1624  switch (dl_method) {
1625 
1626  case SSHT_DL_RISBO:
1627  if (el!=0 && el==abs(spin)) {
1628  for(eltmp=0; eltmp<=abs(spin); eltmp++)
1631  eltmp, sqrt_tbl, signs);
1633  dl8, L,
1636  el,
1637  signs);
1638  }
1639  else {
1642  el, sqrt_tbl, signs);
1644  dl8, L,
1647  el,
1648  signs);
1649  }
1650  break;
1651 
1652  case SSHT_DL_TRAPANI:
1653  if (el!=0 && el==abs(spin)) {
1654  for(eltmp=0; eltmp<=abs(spin); eltmp++)
1657  eltmp, sqrt_tbl);
1660  el, signs);
1661  }
1662  else {
1665  el, sqrt_tbl);
1668  el, signs);
1669  }
1670  break;
1671 
1672  default:
1673  SSHT_ERROR_GENERIC("Invalid dl method")
1674  }
1675 
1676  // Compute flm.
1677  elfactor = sqrt((double)(2.0*el+1.0)/(4.0*SSHT_PI));
1678  el2pel = el *el + el;
1679  for (m=-el; m<=el; m++)
1680  inds[m + inds_offset] = el2pel + m;
1681  elssign = spin <= 0 ? 1.0 : signs[el];
1682 
1683  for (m=-el; m<=-1; m++) {
1684  // mm = 0
1685  ind = inds[m + inds_offset];
1686  flm[ind] +=
1687  ssign
1688  * elfactor
1689  * expsm[m + exps_offset]
1690  * signs[el] * dl[0*dl_stride - m + dl_offset]
1691  * elssign * dl[0*dl_stride - spinneg + dl_offset]
1692  * Fmm[(0+Fmm_offset)*Fmm_stride + m + Fmm_offset];
1693  }
1694  for (m=0; m<=el; m++) {
1695  // mm = 0
1696  ind = inds[m + inds_offset];
1697  flm[ind] +=
1698  ssign
1699  * elfactor
1700  * expsm[m + exps_offset]
1701  * dl[0*dl_stride + m + dl_offset]
1702  * elssign * dl[0*dl_stride - spinneg + dl_offset]
1703  * Fmm[(0+Fmm_offset)*Fmm_stride + m + Fmm_offset];
1704  }
1705 
1706  for (mm=1; mm<=el; mm++) {
1707  elmmsign = signs[el] * signs[mm];
1708  elssign = spin <= 0 ? 1.0 : elmmsign;
1709 
1710  for (m=-el; m<=-1; m++) {
1711  ind = inds[m + inds_offset];
1712  flm[ind] +=
1713  ssign
1714  * elfactor
1715  * expsm[m + exps_offset]
1716  * elmmsign * dl[mm*dl_stride - m + dl_offset]
1717  * elssign * dl[mm*dl_stride - spinneg + dl_offset]
1718  * ( Fmm[(mm+Fmm_offset)*Fmm_stride + m + Fmm_offset]
1719  + signs[-m] * ssign
1720  * Fmm[(-mm+Fmm_offset)*Fmm_stride + m + Fmm_offset]);
1721  }
1722  for (m=0; m<=el; m++) {
1723  ind = inds[m + inds_offset];
1724  flm[ind] +=
1725  ssign
1726  * elfactor
1727  * expsm[m + exps_offset]
1728  * dl[mm*dl_stride + m + dl_offset]
1729  * elssign * dl[mm*dl_stride - spinneg + dl_offset]
1730  * ( Fmm[(mm+Fmm_offset)*Fmm_stride + m + Fmm_offset]
1731  + signs[m] * ssign
1732  * Fmm[(-mm+Fmm_offset)*Fmm_stride + m + Fmm_offset]);
1733  }
1734 
1735  }
1736 
1737  }
1738 
1739  // Free memory.
1740  free(dl);
1741  if (dl_method == SSHT_DL_RISBO)
1742  free(dl8);
1743  free(Fmt);
1744  free(Fmm);
1745  free(sqrt_tbl);
1746  free(signs);
1747  free(expsm);
1748  free(expsmm);
1749  free(inds);
1750 
1751  // Print finished if verbosity set.
1752  if (verbosity > 0)
1753  printf("%s%s", SSHT_PROMPT, "Adjoint inverse transform computed!");
1754 
1755 }
1756 
1757 
1774  int L,
1775  ssht_dl_method_t dl_method,
1776  int verbosity) {
1777 
1778  int el, m, mm, ind, ind_nm, t;
1779  int eltmp;
1780  double *sqrt_tbl, *signs;
1781  int el2pel, inds_offset;
1782  int *inds;
1783  double ssign, elfactor;
1784  fftw_plan plan;
1785  double *in_real;
1786  SSHT_COMPLEX(double) *inout, *out;
1787  SSHT_COMPLEX(double) *Fmt, *Fmm;
1788  int f_stride, Fmt_stride, Fmt_offset, Fmm_stride, Fmm_offset;
1789  double *dl;
1790  double *dl8 = NULL;
1791  int dl_offset, dl_stride;
1792  SSHT_COMPLEX(double) *expsm;
1793  int exps_offset;
1794  int elmmsign, elssign;
1795  int spinneg;
1796  int spin = 0;
1797 
1798  // Allocate memory.
1799  sqrt_tbl = (double*)calloc(2*(L-1)+2, sizeof(double));
1800  SSHT_ERROR_MEM_ALLOC_CHECK(sqrt_tbl)
1801  signs = (double*)calloc(L+1, sizeof(double));
1803  expsm = (SSHT_COMPLEX(double)*)calloc(L, sizeof(SSHT_COMPLEX(double)));
1805  inds = (int*)calloc(L, sizeof(int));
1807 
1808  // Perform precomputations.
1809  for (el=0; el<=2*(L-1)+1; el++)
1810  sqrt_tbl[el] = sqrt((double)el);
1811  for (m=0; m<=L-1; m=m+2) {
1812  signs[m] = 1.0;
1813  signs[m+1] = -1.0;
1814  }
1815  ssign = signs[abs(spin)];
1816  spinneg = spin <= 0 ? spin : -spin;
1817  for (m=0; m<=L-1; m++)
1818  expsm[m] = cexp(I*SSHT_PION2*(m+spin));
1819 
1820  // Print messages depending on verbosity level.
1821  if (verbosity > 0) {
1822  printf("%s%s\n", SSHT_PROMPT,
1823  "Computing adjoint inverse transform using MW symmetric sampling with ");
1824  printf("%s%s%d%s%d%s\n", SSHT_PROMPT, "parameters (L,spin,reality) = (",
1825  L, ",", spin, ", TRUE)");
1826  if (verbosity > 1)
1827  printf("%s%s\n", SSHT_PROMPT,
1828  "Using routine ssht_adjoint_mw_inverse_sov_sym_ss_real...");
1829  }
1830 
1831  // Compute Fourier transform over phi, i.e. compute Fmt.
1832  // Note that t and p indices of fext are increased in size by
1833  // one compared to usual sampling.
1834  Fmt = (SSHT_COMPLEX(double)*)calloc((L+1)*(2*L), sizeof(SSHT_COMPLEX(double)));
1836  Fmt_stride = 2*L;
1837  f_stride = 2*L;
1838  in_real = (double*)calloc(2*L, sizeof(double));
1840  out = (SSHT_COMPLEX(double)*)calloc(L+1, sizeof(SSHT_COMPLEX(double)));
1842  plan = fftw_plan_dft_r2c_1d(2*L, in_real, out, FFTW_MEASURE);
1843  for (t=0; t<=L; t++) {
1844  memcpy(in_real, &f[t*f_stride], f_stride*sizeof(double));
1845  fftw_execute_dft_r2c(plan, in_real, out);
1846  for(m=0; m<=L; m++)
1847  Fmt[m*Fmt_stride + t] = out[m];
1848  }
1849  free(in_real);
1850  free(out);
1851  fftw_destroy_plan(plan);
1852 
1853  // Apply adjoint of periodic extension.
1854  for (m=0; m<=L; m++)
1855  for (t=L+1; t<=2*L-1; t++)
1856  Fmt[m*Fmt_stride + t] = 0.0;
1857 
1858  // Compute Fourier transform over theta, i.e. compute Fmm.
1859  // Note that m and mm indices are increased in size by one.
1860  Fmm = (SSHT_COMPLEX(double)*)calloc((2*L)*(L+1), sizeof(SSHT_COMPLEX(double)));
1862  Fmm_stride = L+1;
1863  Fmm_offset = L-1;
1864  inout = (SSHT_COMPLEX(double)*)calloc(2*L, sizeof(SSHT_COMPLEX(double)));
1866  plan = fftw_plan_dft_1d(2*L, inout, inout, FFTW_FORWARD, FFTW_MEASURE);
1867  for (m=0; m<=L; m++) {
1868  memcpy(inout, &Fmt[m*Fmt_stride], Fmt_stride*sizeof(SSHT_COMPLEX(double)));
1869  fftw_execute_dft(plan, inout, inout);
1870  for(mm=0; mm<=L; mm++)
1871  Fmm[(mm+Fmm_offset)*Fmm_stride + m] =
1872  inout[mm];
1873  for(mm=-(L-1); mm<=-1; mm++)
1874  Fmm[(mm+Fmm_offset)*Fmm_stride + m] =
1875  inout[mm+2*L-1+1];
1876  }
1877  fftw_destroy_plan(plan);
1878  free(inout);
1879 
1880  // Compute flm.
1883  if (dl_method == SSHT_DL_RISBO) {
1886  }
1887  dl_offset = ssht_dl_get_offset(L, SSHT_DL_QUARTER);
1888  dl_stride = ssht_dl_get_stride(L, SSHT_DL_QUARTER);
1889  inds_offset = 0;
1890  for (el=0; el<=L-1; el++) {
1891  for (m=0; m<=el; m++) {
1892  ssht_sampling_elm2ind(&ind, el, m);
1893  flm[ind] = 0.0;
1894  }
1895  }
1896  for (el=abs(spin); el<=L-1; el++) {
1897 
1898  // Compute Wigner plane.
1899  switch (dl_method) {
1900 
1901  case SSHT_DL_RISBO:
1902  if (el!=0 && el==abs(spin)) {
1903  for(eltmp=0; eltmp<=abs(spin); eltmp++)
1906  eltmp, sqrt_tbl, signs);
1908  dl8, L,
1911  el,
1912  signs);
1913  }
1914  else {
1917  el, sqrt_tbl, signs);
1919  dl8, L,
1922  el,
1923  signs);
1924  }
1925  break;
1926 
1927  case SSHT_DL_TRAPANI:
1928  if (el!=0 && el==abs(spin)) {
1929  for(eltmp=0; eltmp<=abs(spin); eltmp++)
1932  eltmp, sqrt_tbl);
1935  el, signs);
1936  }
1937  else {
1940  el, sqrt_tbl);
1943  el, signs);
1944  }
1945  break;
1946 
1947  default:
1948  SSHT_ERROR_GENERIC("Invalid dl method")
1949  }
1950 
1951  // Compute flm.
1952  elfactor = sqrt((double)(2.0*el+1.0)/(4.0*SSHT_PI));
1953  el2pel = el *el + el;
1954  for (m=0; m<=el; m++)
1955  inds[m + inds_offset] = el2pel + m;
1956  elssign = spin <= 0 ? 1.0 : signs[el];
1957 
1958  for (m=0; m<=el; m++) {
1959  // mm = 0
1960  ind = inds[m + inds_offset];
1961  flm[ind] +=
1962  ssign
1963  * elfactor
1964  * expsm[m]
1965  * dl[0*dl_stride + m + dl_offset]
1966  * elssign * dl[0*dl_stride - spinneg + dl_offset]
1967  * Fmm[(0+Fmm_offset)*Fmm_stride + m];
1968  }
1969 
1970  for (mm=1; mm<=el; mm++) {
1971  elmmsign = signs[el] * signs[mm];
1972  elssign = spin <= 0 ? 1.0 : elmmsign;
1973 
1974  for (m=0; m<=el; m++) {
1975  ind = inds[m + inds_offset];
1976  flm[ind] +=
1977  ssign
1978  * elfactor
1979  * expsm[m]
1980  * dl[mm*dl_stride + m + dl_offset]
1981  * elssign * dl[mm*dl_stride - spinneg + dl_offset]
1982  * ( Fmm[(mm+Fmm_offset)*Fmm_stride + m]
1983  + signs[m] * ssign
1984  * Fmm[(-mm+Fmm_offset)*Fmm_stride + m]);
1985  }
1986 
1987  }
1988 
1989  }
1990 
1991  // Set flm values for negative m using conjugate symmetry.
1992  for (el=abs(spin); el<=L-1; el++) {
1993  for (m=1; m<=el; m++) {
1994  ssht_sampling_elm2ind(&ind, el, m);
1995  ssht_sampling_elm2ind(&ind_nm, el, -m);
1996  flm[ind_nm] = signs[m] * conj(flm[ind]);
1997  }
1998  }
1999 
2000  // Free memory.
2001  free(dl);
2002  if (dl_method == SSHT_DL_RISBO)
2003  free(dl8);
2004  free(Fmt);
2005  free(Fmm);
2006  free(sqrt_tbl);
2007  free(signs);
2008  free(expsm);
2009  free(inds);
2010 
2011  // Print finished if verbosity set.
2012  if (verbosity > 0)
2013  printf("%s%s", SSHT_PROMPT, "Adjoint inverse transform computed!");
2014 
2015 }
2016 
2017 
2035  int L, int spin,
2036  ssht_dl_method_t dl_method,
2037  int verbosity) {
2038 
2039  int el, m, mm, ind;
2040  int eltmp;
2041  double *sqrt_tbl, *signs;
2042  int el2pel, inds_offset;
2043  int *inds;
2044  double ssign, elfactor;
2045  double *dl;
2046  double *dl8 = NULL;
2047  int dl_offset, dl_stride;
2048  SSHT_COMPLEX(double) *exps;
2049  int exps_offset;
2050  double elmmsign, elssign;
2051  int spinneg;
2052  SSHT_COMPLEX(double) *Fmm;
2053  int Fmm_offset, Fmm_stride;
2054 
2055  fftw_plan plan, plan_bwd, plan_fwd;
2056  SSHT_COMPLEX(double) *Ftm, *Gmm;
2057  SSHT_COMPLEX(double) *w, *wr;
2058  int w_offset;
2059  SSHT_COMPLEX(double) *Fmm_pad, *tmp_pad;
2060  int f_stride, Ftm_stride, Ftm_offset;
2061  int r, t, p;
2062  SSHT_COMPLEX(double) *inout;
2063 
2064  int Gmm_offset, Gmm_stride;
2065 
2066  // Allocate memory.
2067  sqrt_tbl = (double*)calloc(2*(L-1)+2, sizeof(double));
2068  SSHT_ERROR_MEM_ALLOC_CHECK(sqrt_tbl)
2069  signs = (double*)calloc(L+1, sizeof(double));
2071  exps = (SSHT_COMPLEX(double)*)calloc(2*L-1, sizeof(SSHT_COMPLEX(double)));
2073  inds = (int*)calloc(2*L-1, sizeof(int));
2075 
2076  // Perform precomputations.
2077  for (el=0; el<=2*(L-1)+1; el++)
2078  sqrt_tbl[el] = sqrt((double)el);
2079  for (m=0; m<=L-1; m=m+2) {
2080  signs[m] = 1.0;
2081  signs[m+1] = -1.0;
2082  }
2083  ssign = signs[abs(spin)];
2084  spinneg = spin <= 0 ? spin : -spin;
2085  exps_offset = L-1;
2086  for (m=-(L-1); m<=L-1; m++)
2087  exps[m + exps_offset] = cexp(-I*SSHT_PION2*(m+spin));
2088 
2089  // Print messages depending on verbosity level.
2090  if (verbosity > 0) {
2091  printf("%s%s\n", SSHT_PROMPT,
2092  "Computing adjoint forward transform using MW symmetric sampling with ");
2093  printf("%s%s%d%s%d%s\n", SSHT_PROMPT, "parameters (L,spin,reality) = (",
2094  L, ",", spin, ", FALSE)");
2095  if (verbosity > 1)
2096  printf("%s%s\n", SSHT_PROMPT,
2097  "Using routine ssht_adjoint_mw_forward_sov_sym_ss...");
2098  }
2099 
2100  // Compute Fmm.
2101  Fmm = (SSHT_COMPLEX(double)*)calloc((2*L-1)*(2*L-1), sizeof(SSHT_COMPLEX(double)));
2103  Fmm_offset = L-1;
2104  Fmm_stride = 2*L-1;
2107  if (dl_method == SSHT_DL_RISBO) {
2110  }
2111  dl_offset = ssht_dl_get_offset(L, SSHT_DL_QUARTER);
2112  dl_stride = ssht_dl_get_stride(L, SSHT_DL_QUARTER);
2113  inds_offset = L-1;
2114  for (el=abs(spin); el<=L-1; el++) {
2115 
2116  // Compute Wigner plane.
2117  switch (dl_method) {
2118 
2119  case SSHT_DL_RISBO:
2120  if (el!=0 && el==abs(spin)) {
2121  for(eltmp=0; eltmp<=abs(spin); eltmp++)
2124  eltmp, sqrt_tbl, signs);
2126  dl8, L,
2129  el,
2130  signs);
2131  }
2132  else {
2135  el, sqrt_tbl, signs);
2137  dl8, L,
2140  el,
2141  signs);
2142  }
2143  break;
2144 
2145  case SSHT_DL_TRAPANI:
2146  if (el!=0 && el==abs(spin)) {
2147  for(eltmp=0; eltmp<=abs(spin); eltmp++)
2150  eltmp, sqrt_tbl);
2153  el, signs);
2154  }
2155  else {
2158  el, sqrt_tbl);
2161  el, signs);
2162  }
2163  break;
2164 
2165  default:
2166  SSHT_ERROR_GENERIC("Invalid dl method")
2167  }
2168 
2169  // Compute Fmm.
2170  elfactor = sqrt((double)(2.0*el+1.0)/(4.0*SSHT_PI));
2171  el2pel = el *el + el;
2172  for (m=-el; m<=el; m++)
2173  inds[m + inds_offset] = el2pel + m;
2174  for (mm=0; mm<=el; mm++) {
2175  elmmsign = signs[el] * signs[mm];
2176  elssign = spin <= 0 ? 1.0 : elmmsign;
2177 
2178  for (m=-el; m<=-1; m++) {
2179  ind = inds[m + inds_offset];
2180  Fmm[(m + Fmm_offset)*Fmm_stride + mm + Fmm_offset] +=
2181  ssign
2182  * elfactor
2183  * exps[m + exps_offset]
2184  * elmmsign * dl[mm*dl_stride - m + dl_offset]
2185  * elssign * dl[mm*dl_stride - spinneg + dl_offset]
2186  * flm[ind];
2187  }
2188  for (m=0; m<=el; m++) {
2189  ind = inds[m + inds_offset];
2190  Fmm[(m + Fmm_offset)*Fmm_stride + mm + Fmm_offset] +=
2191  ssign
2192  * elfactor
2193  * exps[m + exps_offset]
2194  * dl[mm*dl_stride + m + dl_offset]
2195  * elssign * dl[mm*dl_stride - spinneg + dl_offset]
2196  * flm[ind];
2197  }
2198 
2199  }
2200 
2201  }
2202 
2203  // Free dl memory.
2204  free(dl);
2205  if (dl_method == SSHT_DL_RISBO)
2206  free(dl8);
2207 
2208  // Use symmetry to compute Fmm for negative mm.
2209  for (m=-(L-1); m<=L-1; m++)
2210  for (mm=-(L-1); mm<=-1; mm++)
2211  Fmm[(m + Fmm_offset)*Fmm_stride + mm + Fmm_offset] =
2212  signs[abs(m)] * ssign
2213  * Fmm[(m + Fmm_offset)*Fmm_stride - mm + Fmm_offset];
2214 
2215  // Compute weights.
2216  w = (SSHT_COMPLEX(double)*)calloc(4*L-3, sizeof(SSHT_COMPLEX(double)));
2218  w_offset = 2*(L-1);
2219  for (mm=-2*(L-1); mm<=2*(L-1); mm++)
2220  w[mm+w_offset] = ssht_sampling_weight_mw(mm);
2221 
2222  // Compute IFFT of w to give wr.
2223  wr = (SSHT_COMPLEX(double)*)calloc(4*L-3, sizeof(SSHT_COMPLEX(double)));
2225  inout = (SSHT_COMPLEX(double)*)calloc(4*L-3, sizeof(SSHT_COMPLEX(double)));
2227  plan_bwd = fftw_plan_dft_1d(4*L-3, inout, inout, FFTW_BACKWARD, FFTW_MEASURE);
2228  plan_fwd = fftw_plan_dft_1d(4*L-3, inout, inout, FFTW_FORWARD, FFTW_MEASURE);
2229  for (mm=1; mm<=2*L-2; mm++)
2230  inout[mm + w_offset] = w[mm - 2*(L-1) - 1 + w_offset];
2231  for (mm=-2*(L-1); mm<=0; mm++)
2232  inout[mm + w_offset] = w[mm + 2*(L-1) + w_offset];
2233  fftw_execute_dft(plan_bwd, inout, inout);
2234  for (mm=0; mm<=2*L-2; mm++)
2235  wr[mm + w_offset] = inout[mm - 2*(L-1) + w_offset];
2236  for (mm=-2*(L-1); mm<=-1; mm++)
2237  wr[mm + w_offset] = inout[mm + 2*(L-1) + 1 + w_offset];
2238 
2239  // Compute Gmm by convolution implemented as product in real space.
2240  Fmm_pad = (SSHT_COMPLEX(double)*)calloc(4*L-3, sizeof(SSHT_COMPLEX(double)));
2242  tmp_pad = (SSHT_COMPLEX(double)*)calloc(4*L-3, sizeof(SSHT_COMPLEX(double)));
2244  Gmm = (SSHT_COMPLEX(double)*)calloc((2*L)*(2*L), sizeof(SSHT_COMPLEX(double)));
2246  Gmm_stride = 2*L;
2247  Gmm_offset = L-1;
2248  for (m=-(L-1); m<=L-1; m++) {
2249 
2250  // Zero-pad Fmm.
2251  for (mm=-2*(L-1); mm<=-L; mm++)
2252  Fmm_pad[mm+w_offset] = 0.0;
2253  for (mm=L; mm<=2*(L-1); mm++)
2254  Fmm_pad[mm+w_offset] = 0.0;
2255  for (mm=-(L-1); mm<=L-1; mm++)
2256  Fmm_pad[mm + w_offset] =
2257  Fmm[(m+Fmm_offset)*Fmm_stride + mm + Fmm_offset];
2258 
2259  // Compute IFFT of Fmm.
2260  for (mm=1; mm<=2*L-2; mm++)
2261  inout[mm + w_offset] = Fmm_pad[mm - 2*(L-1) - 1 + w_offset];
2262  for (mm=-2*(L-1); mm<=0; mm++)
2263  inout[mm + w_offset] = Fmm_pad[mm + 2*(L-1) + w_offset];
2264  fftw_execute_dft(plan_bwd, inout, inout);
2265  for (mm=0; mm<=2*L-2; mm++)
2266  Fmm_pad[mm + w_offset] = inout[mm - 2*(L-1) + w_offset];
2267  for (mm=-2*(L-1); mm<=-1; mm++)
2268  Fmm_pad[mm + w_offset] = inout[mm + 2*(L-1) + 1 + w_offset];
2269 
2270  // Compute product of Fmm and weight in real space.
2271  for (r=-2*(L-1); r<=2*(L-1); r++)
2272  Fmm_pad[r + w_offset] *= wr[-r + w_offset];
2273 
2274  // Compute Gmm by FFT.
2275  for (mm=1; mm<=2*L-2; mm++)
2276  inout[mm + w_offset] = Fmm_pad[mm - 2*(L-1) - 1 + w_offset];
2277  for (mm=-2*(L-1); mm<=0; mm++)
2278  inout[mm + w_offset] = Fmm_pad[mm + 2*(L-1) + w_offset];
2279  fftw_execute_dft(plan_fwd, inout, inout);
2280  for (mm=0; mm<=2*L-2; mm++)
2281  Fmm_pad[mm + w_offset] = inout[mm - 2*(L-1) + w_offset];
2282  for (mm=-2*(L-1); mm<=-1; mm++)
2283  Fmm_pad[mm + w_offset] = inout[mm + 2*(L-1) + 1 + w_offset];
2284 
2285  // Extract section of Gmm of interest.
2286  for (mm=-(L-1); mm<=L-1; mm++)
2287  Gmm[(m+Gmm_offset)*Gmm_stride + mm + Gmm_offset] =
2288  Fmm_pad[mm + w_offset] * 2.0 * SSHT_PI / (4.0*L-3.0);
2289 
2290  }
2291  fftw_destroy_plan(plan_bwd);
2292  fftw_destroy_plan(plan_fwd);
2293  free(inout);
2294 
2295  // Compute Fourier transform over theta.
2296  inout = (SSHT_COMPLEX(double)*)calloc(2*L, sizeof(SSHT_COMPLEX(double)));
2298  plan = fftw_plan_dft_1d(2*L, inout, inout, FFTW_BACKWARD, FFTW_MEASURE);
2299  Ftm = (SSHT_COMPLEX(double)*)calloc((2*L)*(2*L), sizeof(SSHT_COMPLEX(double)));
2301  Ftm_stride = 2*L;
2302  Ftm_offset = L-1;
2303  for (m=-(L-1); m<=L; m++) {
2304 
2305  for(mm=0; mm<=L; mm++)
2306  inout[mm] = Gmm[(m+Gmm_offset)*Gmm_stride + mm + Gmm_offset];
2307  for(mm=-(L-1); mm<=-1; mm++)
2308  inout[mm+2*L-1+1] = Gmm[(m+Gmm_offset)*Gmm_stride + mm + Gmm_offset];
2309  fftw_execute_dft(plan, inout, inout);
2310 
2311  for(t=0; t<=2*L-1; t++)
2312  Ftm[t*Ftm_stride + m + Ftm_offset] = inout[t] / (2.0*L);
2313 
2314  }
2315 
2316  // Adjoint of periodic extension of Ftm.
2317  for(t=1; t<=L-1; t++)
2318  for (m=-(L-1); m<=L; m++)
2319  Ftm[t*Ftm_stride + m + Ftm_offset] =
2320  Ftm[t*Ftm_stride + m + Ftm_offset]
2321  + signs[abs(m)] * ssign * Ftm[(2*L-t)*Ftm_stride + m + Ftm_offset];
2322 
2323  // Compute Fourier transform over phi.
2324  f_stride = 2*L;
2325  for(t=0; t<=L; t++) {
2326 
2327  for(m=0; m<=L; m++)
2328  inout[m] = Ftm[t*Ftm_stride + m + Ftm_offset];
2329  for(m=-(L-1); m<=-1; m++)
2330  inout[m+2*L-1+1] = Ftm[t*Ftm_stride + m + Ftm_offset];
2331  fftw_execute_dft(plan, inout, inout);
2332 
2333  for(p=0; p<=2*L-1; p++)
2334  f[t*f_stride + p] = inout[p] / (2.0*L);
2335 
2336  }
2337  fftw_destroy_plan(plan);
2338  free(inout);
2339 
2340  // Free memory.
2341  free(Fmm);
2342  free(Ftm);
2343  free(w);
2344  free(wr);
2345  free(Fmm_pad);
2346  free(tmp_pad);
2347  free(Gmm);
2348 
2349  // Free precomputation memory.
2350  free(sqrt_tbl);
2351  free(signs);
2352  free(exps);
2353  free(inds);
2354 
2355  // Print finished if verbosity set.
2356  if (verbosity > 0)
2357  printf("%s%s", SSHT_PROMPT, "Adjoint forward transform computed!");
2358 
2359 }
2360 
2361 
2379  SSHT_COMPLEX(double) *flm,
2380  int L,
2381  ssht_dl_method_t dl_method,
2382  int verbosity) {
2383 
2384  int el, m, mm, ind;
2385  int eltmp;
2386  double *sqrt_tbl, *signs;
2387  int el2pel, inds_offset;
2388  int *inds;
2389  double ssign, elfactor;
2390  double *dl;
2391  double *dl8 = NULL;
2392  int dl_offset, dl_stride;
2393  SSHT_COMPLEX(double) *exps;
2394  int exps_offset;
2395  double elmmsign, elssign;
2396  int spinneg;
2397  SSHT_COMPLEX(double) *Fmm;
2398  int Fmm_offset, Fmm_stride;
2399 
2400  fftw_plan plan, plan_bwd, plan_fwd;
2401  SSHT_COMPLEX(double) *Ftm, *Gmm;
2402  SSHT_COMPLEX(double) *w, *wr;
2403  int w_offset;
2404  SSHT_COMPLEX(double) *Fmm_pad, *tmp_pad;
2405  int f_stride, Ftm_stride;
2406  int r, t, p;
2407  SSHT_COMPLEX(double) *inout;
2408  SSHT_COMPLEX(double) *in;
2409  double *out_real;
2410  int Gmm_offset, Gmm_stride;
2411  int spin = 0;
2412 
2413  // Allocate memory.
2414  sqrt_tbl = (double*)calloc(2*(L-1)+2, sizeof(double));
2415  SSHT_ERROR_MEM_ALLOC_CHECK(sqrt_tbl)
2416  signs = (double*)calloc(L+1, sizeof(double));
2418  exps = (SSHT_COMPLEX(double)*)calloc(2*L-1, sizeof(SSHT_COMPLEX(double)));
2420  inds = (int*)calloc(2*L-1, sizeof(int));
2422 
2423  // Perform precomputations.
2424  for (el=0; el<=2*(L-1)+1; el++)
2425  sqrt_tbl[el] = sqrt((double)el);
2426  for (m=0; m<=L-1; m=m+2) {
2427  signs[m] = 1.0;
2428  signs[m+1] = -1.0;
2429  }
2430  ssign = signs[abs(spin)];
2431  spinneg = spin <= 0 ? spin : -spin;
2432  exps_offset = L-1;
2433  for (m=-(L-1); m<=L-1; m++)
2434  exps[m + exps_offset] = cexp(-I*SSHT_PION2*(m+spin));
2435 
2436  // Print messages depending on verbosity level.
2437  if (verbosity > 0) {
2438  printf("%s%s\n", SSHT_PROMPT,
2439  "Computing adjoint forward transform using MW symmetric sampling with ");
2440  printf("%s%s%d%s%d%s\n", SSHT_PROMPT, "parameters (L,spin,reality) = (",
2441  L, ",", spin, ", FALSE)");
2442  if (verbosity > 1)
2443  printf("%s%s\n", SSHT_PROMPT,
2444  "Using routine ssht_adjoint_mw_forward_sov_sym_ss_real...");
2445  }
2446 
2447  // Compute Fmm.
2448  Fmm = (SSHT_COMPLEX(double)*)calloc(L*(2*L-1), sizeof(SSHT_COMPLEX(double)));
2450  Fmm_offset = L-1;
2451  Fmm_stride = 2*L-1;
2454  if (dl_method == SSHT_DL_RISBO) {
2457  }
2458  dl_offset = ssht_dl_get_offset(L, SSHT_DL_QUARTER);
2459  dl_stride = ssht_dl_get_stride(L, SSHT_DL_QUARTER);
2460  inds_offset = L-1;
2461  for (el=abs(spin); el<=L-1; el++) {
2462 
2463  // Compute Wigner plane.
2464  switch (dl_method) {
2465 
2466  case SSHT_DL_RISBO:
2467  if (el!=0 && el==abs(spin)) {
2468  for(eltmp=0; eltmp<=abs(spin); eltmp++)
2471  eltmp, sqrt_tbl, signs);
2473  dl8, L,
2476  el,
2477  signs);
2478  }
2479  else {
2482  el, sqrt_tbl, signs);
2484  dl8, L,
2487  el,
2488  signs);
2489  }
2490  break;
2491 
2492  case SSHT_DL_TRAPANI:
2493  if (el!=0 && el==abs(spin)) {
2494  for(eltmp=0; eltmp<=abs(spin); eltmp++)
2497  eltmp, sqrt_tbl);
2500  el, signs);
2501  }
2502  else {
2505  el, sqrt_tbl);
2508  el, signs);
2509  }
2510  break;
2511 
2512  default:
2513  SSHT_ERROR_GENERIC("Invalid dl method")
2514  }
2515 
2516  // Compute Fmm.
2517  elfactor = sqrt((double)(2.0*el+1.0)/(4.0*SSHT_PI));
2518  el2pel = el *el + el;
2519  for (m=0; m<=el; m++)
2520  inds[m + inds_offset] = el2pel + m;
2521  for (mm=0; mm<=el; mm++) {
2522  elmmsign = signs[el] * signs[mm];
2523  elssign = spin <= 0 ? 1.0 : elmmsign;
2524 
2525  for (m=0; m<=el; m++) {
2526  ind = inds[m + inds_offset];
2527  Fmm[m*Fmm_stride + mm + Fmm_offset] +=
2528  ssign
2529  * elfactor
2530  * exps[m + exps_offset]
2531  * dl[mm*dl_stride + m + dl_offset]
2532  * elssign * dl[mm*dl_stride - spinneg + dl_offset]
2533  * flm[ind];
2534  }
2535 
2536  }
2537 
2538  }
2539 
2540  // Free dl memory.
2541  free(dl);
2542  if (dl_method == SSHT_DL_RISBO)
2543  free(dl8);
2544 
2545  // Use symmetry to compute Fmm for negative mm.
2546  for (m=0; m<=L-1; m++)
2547  for (mm=-(L-1); mm<=-1; mm++)
2548  Fmm[m*Fmm_stride + mm + Fmm_offset] =
2549  signs[abs(m)] * ssign
2550  * Fmm[m*Fmm_stride - mm + Fmm_offset];
2551 
2552  // Compute weights.
2553  w = (SSHT_COMPLEX(double)*)calloc(4*L-3, sizeof(SSHT_COMPLEX(double)));
2555  w_offset = 2*(L-1);
2556  for (mm=-2*(L-1); mm<=2*(L-1); mm++)
2557  w[mm+w_offset] = ssht_sampling_weight_mw(mm);
2558 
2559  // Compute IFFT of w to give wr.
2560  wr = (SSHT_COMPLEX(double)*)calloc(4*L-3, sizeof(SSHT_COMPLEX(double)));
2562  inout = (SSHT_COMPLEX(double)*)calloc(4*L-3, sizeof(SSHT_COMPLEX(double)));
2564  plan_bwd = fftw_plan_dft_1d(4*L-3, inout, inout, FFTW_BACKWARD, FFTW_MEASURE);
2565  plan_fwd = fftw_plan_dft_1d(4*L-3, inout, inout, FFTW_FORWARD, FFTW_MEASURE);
2566  for (mm=1; mm<=2*L-2; mm++)
2567  inout[mm + w_offset] = w[mm - 2*(L-1) - 1 + w_offset];
2568  for (mm=-2*(L-1); mm<=0; mm++)
2569  inout[mm + w_offset] = w[mm + 2*(L-1) + w_offset];
2570  fftw_execute_dft(plan_bwd, inout, inout);
2571  for (mm=0; mm<=2*L-2; mm++)
2572  wr[mm + w_offset] = inout[mm - 2*(L-1) + w_offset];
2573  for (mm=-2*(L-1); mm<=-1; mm++)
2574  wr[mm + w_offset] = inout[mm + 2*(L-1) + 1 + w_offset];
2575 
2576  // Compute Gmm by convolution implemented as product in real space.
2577  Fmm_pad = (SSHT_COMPLEX(double)*)calloc(4*L-3, sizeof(SSHT_COMPLEX(double)));
2579  tmp_pad = (SSHT_COMPLEX(double)*)calloc(4*L-3, sizeof(SSHT_COMPLEX(double)));
2581  Gmm = (SSHT_COMPLEX(double)*)calloc((L+1)*(2*L), sizeof(SSHT_COMPLEX(double)));
2582  Gmm_stride = 2*L;
2583  Gmm_offset = L-1;
2585  for (m=0; m<=L-1; m++) {
2586 
2587  // Zero-pad Fmm.
2588  for (mm=-2*(L-1); mm<=-L; mm++)
2589  Fmm_pad[mm+w_offset] = 0.0;
2590  for (mm=L; mm<=2*(L-1); mm++)
2591  Fmm_pad[mm+w_offset] = 0.0;
2592  for (mm=-(L-1); mm<=L-1; mm++)
2593  Fmm_pad[mm + w_offset] =
2594  Fmm[m*Fmm_stride + mm + Fmm_offset];
2595 
2596  // Compute IFFT of Fmm.
2597  for (mm=1; mm<=2*L-2; mm++)
2598  inout[mm + w_offset] = Fmm_pad[mm - 2*(L-1) - 1 + w_offset];
2599  for (mm=-2*(L-1); mm<=0; mm++)
2600  inout[mm + w_offset] = Fmm_pad[mm + 2*(L-1) + w_offset];
2601  fftw_execute_dft(plan_bwd, inout, inout);
2602  for (mm=0; mm<=2*L-2; mm++)
2603  Fmm_pad[mm + w_offset] = inout[mm - 2*(L-1) + w_offset];
2604  for (mm=-2*(L-1); mm<=-1; mm++)
2605  Fmm_pad[mm + w_offset] = inout[mm + 2*(L-1) + 1 + w_offset];
2606 
2607  // Compute product of Fmm and weight in real space.
2608  for (r=-2*(L-1); r<=2*(L-1); r++)
2609  Fmm_pad[r + w_offset] *= wr[-r + w_offset];
2610 
2611  // Compute Gmm by FFT.
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_fwd, 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  // Extract section of Gmm of interest.
2623  for (mm=-(L-1); mm<=L-1; mm++)
2624  Gmm[m*Gmm_stride + mm + Gmm_offset] =
2625  Fmm_pad[mm + w_offset] * 2.0 * SSHT_PI / (4.0*L-3.0);
2626 
2627  }
2628  fftw_destroy_plan(plan_bwd);
2629  fftw_destroy_plan(plan_fwd);
2630  free(inout);
2631 
2632  // Compute Fourier transform over theta.
2633  inout = (SSHT_COMPLEX(double)*)calloc(2*L, sizeof(SSHT_COMPLEX(double)));
2635  plan = fftw_plan_dft_1d(2*L, inout, inout, FFTW_BACKWARD, FFTW_MEASURE);
2636  Ftm = (SSHT_COMPLEX(double)*)calloc((2*L)*(L+1), sizeof(SSHT_COMPLEX(double)));
2638  Ftm_stride = L+1;
2639  for (m=0; m<=L; m++) {
2640 
2641  for(mm=0; mm<=L; mm++)
2642  inout[mm] = Gmm[m*Gmm_stride + mm + Gmm_offset];
2643  for(mm=-(L-1); mm<=-1; mm++)
2644  inout[mm+2*L-1+1] = Gmm[m*Gmm_stride + mm + Gmm_offset];
2645  fftw_execute_dft(plan, inout, inout);
2646 
2647  for(t=0; t<=2*L-1; t++)
2648  Ftm[t*Ftm_stride + m] = inout[t] / (2.0*L);
2649 
2650  }
2651  fftw_destroy_plan(plan);
2652  free(inout);
2653 
2654  // Adjoint of periodic extension of Ftm.
2655  for(t=1; t<=L-1; t++)
2656  for (m=0; m<=L; m++)
2657  Ftm[t*Ftm_stride + m] =
2658  Ftm[t*Ftm_stride + m]
2659  + signs[abs(m)] * ssign * Ftm[(2*L-t)*Ftm_stride + m];
2660 
2661  // Compute Fourier transform over phi.
2662  out_real = (double*)calloc(2*L, sizeof(double));
2663  SSHT_ERROR_MEM_ALLOC_CHECK(out_real)
2664  in = (SSHT_COMPLEX(double)*)calloc(L+1, sizeof(SSHT_COMPLEX(double)));
2666  plan = fftw_plan_dft_c2r_1d(2*L, in, out_real, FFTW_MEASURE);
2667  f_stride = 2*L;
2668  for(t=0; t<=L; t++) {
2669 
2670  memcpy(in, &Ftm[t*Ftm_stride], Ftm_stride*sizeof(SSHT_COMPLEX(double)));
2671  fftw_execute_dft_c2r(plan, in, out_real);
2672 
2673  for(p=0; p<=2*L-1; p++)
2674  f[t*f_stride + p] = out_real[p] / (2.0*L);
2675 
2676  }
2677  fftw_destroy_plan(plan);
2678  free(out_real);
2679  free(in);
2680 
2681  // Free memory.
2682  free(Fmm);
2683  free(Ftm);
2684  free(w);
2685  free(wr);
2686  free(Fmm_pad);
2687  free(tmp_pad);
2688  free(Gmm);
2689 
2690  // Free precomputation memory.
2691  free(sqrt_tbl);
2692  free(signs);
2693  free(exps);
2694  free(inds);
2695 
2696  // Print finished if verbosity set.
2697  if (verbosity > 0)
2698  printf("%s%s", SSHT_PROMPT, "Adjoint forward transform computed!");
2699 
2700 }
2701 
2702 
2703 //============================================================================
2704 // MW SS Noth-South pole interfaces
2705 //============================================================================
2706 
2707 
2732  SSHT_COMPLEX(double) *f_np, double *phi_np,
2733  SSHT_COMPLEX(double) *f_sp, double *phi_sp,
2734  SSHT_COMPLEX(double) *flm,
2735  int L, int spin,
2736  ssht_dl_method_t dl_method,
2737  int verbosity) {
2738 
2739  SSHT_COMPLEX(double)* f_full;
2740  int t, f_stride = 2*L;
2741 
2742  // Allocate full array.
2743  f_full = (SSHT_COMPLEX(double)*)calloc((L+1)*(2*L), sizeof(SSHT_COMPLEX(double)));
2745 
2746  // Perform inverse transform.
2747  ssht_adjoint_mw_forward_sov_sym_ss(f_full, flm, L, spin,
2748  dl_method, verbosity);
2749 
2750  // Copy output function values, including separate points for poles.
2751  for (t=1; t<=L-1; t++)
2752  memcpy(&f[(t-1)*f_stride], &f_full[t*f_stride],
2753  (2*L)*sizeof(SSHT_COMPLEX(double)));
2754  *f_np = f_full[0];
2755  *phi_np = ssht_sampling_mw_ss_p2phi(0, L);
2756  *f_sp = f_full[L*f_stride + 0];
2757  *phi_sp = ssht_sampling_mw_ss_p2phi(0, L);
2758 
2759  // Free memory.
2760  free(f_full);
2761 
2762 }
2763 
2764 
2784  double *f_np,
2785  double *f_sp,
2786  SSHT_COMPLEX(double) *flm,
2787  int L,
2788  ssht_dl_method_t dl_method,
2789  int verbosity) {
2790 
2791  double *f_full;
2792  int t, f_stride = 2*L;
2793 
2794  // Allocate full array.
2795  f_full = (double*)calloc((L+1)*(2*L), sizeof(double));
2797 
2798  // Perform inverse transform.
2800  dl_method, verbosity);
2801 
2802  // Copy output function values, including separate points for poles.
2803  for (t=1; t<=L-1; t++)
2804  memcpy(&f[(t-1)*f_stride], &f_full[t*f_stride],
2805  (2*L)*sizeof(double));
2806  *f_np = f_full[0];
2807  *f_sp = f_full[L*f_stride + 0];
2808 
2809  // Free memory.
2810  free(f_full);
2811 
2812 }
2813 
2814 
2839  SSHT_COMPLEX(double) f_np, double phi_np,
2840  SSHT_COMPLEX(double) f_sp, double phi_sp,
2841  int L, int spin,
2842  ssht_dl_method_t dl_method,
2843  int verbosity) {
2844 
2845  SSHT_COMPLEX(double) *f_full;
2846  int t, p, f_stride = 2*L;
2847  double phi;
2848 
2849  // Copy function values to full array.
2850  f_full = (SSHT_COMPLEX(double)*)calloc((L+1)*(2*L), sizeof(SSHT_COMPLEX(double)));
2852  for (t=1; t<=L-1; t++)
2853  memcpy(&f_full[t*f_stride], &f[(t-1)*f_stride],
2854  (2*L)*sizeof(SSHT_COMPLEX(double)));
2855 
2856  // Define poles for all phi.
2857  for (p=0; p<=2*L-1; p++) {
2858  phi = ssht_sampling_mw_ss_p2phi(p, L);
2859  f_full[0*f_stride + p] = f_np * cexp(-I*spin*(phi-phi_np));
2860  f_full[L*f_stride + p] = f_sp * cexp(I*spin*(phi-phi_sp));
2861  }
2862 
2863  // Perform forward transform.
2864  ssht_adjoint_mw_inverse_sov_sym_ss(flm, f_full, L, spin,
2865  dl_method, verbosity);
2866 
2867  // Free memory.
2868  free(f_full);
2869 
2870 }
2871 
2872 
2892  double *f,
2893  double f_np,
2894  double f_sp,
2895  int L,
2896  ssht_dl_method_t dl_method,
2897  int verbosity) {
2898 
2899  double *f_full;
2900  int t, p, f_stride = 2*L;
2901 
2902  // Copy function values to full array.
2903  f_full = (double*)calloc((L+1)*(2*L), sizeof(double));
2905  for (t=1; t<=L-1; t++)
2906  memcpy(&f_full[t*f_stride], &f[(t-1)*f_stride],
2907  (2*L)*sizeof(double));
2908 
2909  // Define poles for all phi.
2910  for (p=0; p<=2*L-1; p++) {
2911  f_full[0*f_stride + p] = f_np;
2912  f_full[L*f_stride + p] = f_sp;
2913  }
2914 
2915  // Perform forward transform.
2917  dl_method, verbosity);
2918 
2919  // Free memory.
2920  free(f_full);
2921 
2922 }
ssht_adjoint_mw_forward_sov_sym_real
void ssht_adjoint_mw_forward_sov_sym_real(double *f, const SSHT_COMPLEX(double) *flm, int L, ssht_dl_method_t dl_method, int verbosity)
Definition: ssht_adjoint.c:946
SSHT_PION2
#define SSHT_PION2
Definition: ssht_types.h:39
SSHT_PI
#define SSHT_PI
Definition: ssht_types.h:38
ssht_adjoint_mw_forward_sov_sym_real_pole
void ssht_adjoint_mw_forward_sov_sym_real_pole(double *f, double *f_sp, SSHT_COMPLEX(double) *flm, int L, ssht_dl_method_t dl_method, int verbosity)
Definition: ssht_adjoint.c:1350
ssht_adjoint_mw_forward_sov_sym
void ssht_adjoint_mw_forward_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_adjoint.c:597
ssht_error.h
SSHT_DL_QUARTER_EXTENDED
@ SSHT_DL_QUARTER_EXTENDED
Definition: ssht_dl.h:16
ssht_adjoint_mw_forward_sov_sym_ss_real_pole
void ssht_adjoint_mw_forward_sov_sym_ss_real_pole(double *f, double *f_np, double *f_sp, SSHT_COMPLEX(double) *flm, int L, ssht_dl_method_t dl_method, int verbosity)
Definition: ssht_adjoint.c:2783
SSHT_ERROR_GENERIC
#define SSHT_ERROR_GENERIC(comment)
Definition: ssht_error.h:19
ssht_adjoint_mw_inverse_sov_sym_real_pole
void ssht_adjoint_mw_inverse_sov_sym_real_pole(SSHT_COMPLEX(double) *flm, double *f, double f_sp, int L, ssht_dl_method_t dl_method, int verbosity)
Definition: ssht_adjoint.c:1446
SSHT_DL_TRAPANI
@ SSHT_DL_TRAPANI
Definition: ssht_dl.h:22
ssht_dl.h
ssht_dl_get_stride
int ssht_dl_get_stride(int L, ssht_dl_size_t dl_size)
Definition: ssht_dl.c:146
ssht_adjoint_mw_inverse_sov_sym_ss_pole
void ssht_adjoint_mw_inverse_sov_sym_ss_pole(SSHT_COMPLEX(double) *flm, 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_adjoint.c:2838
ssht_dl_calloc
double * ssht_dl_calloc(int L, ssht_dl_size_t dl_size)
Definition: ssht_dl.c:64
ssht_sampling.h
ssht_adjoint_mw_inverse_sov_sym_real
void ssht_adjoint_mw_inverse_sov_sym_real(SSHT_COMPLEX(double) *flm, const double *f, int L, ssht_dl_method_t dl_method, int verbosity)
Definition: ssht_adjoint.c:325
ssht_adjoint_mw_forward_sov_sym_pole
void ssht_adjoint_mw_forward_sov_sym_pole(SSHT_COMPLEX(double) *f, SSHT_COMPLEX(double) *f_sp, double *phi_sp, SSHT_COMPLEX(double) *flm, int L, int spin, ssht_dl_method_t dl_method, int verbosity)
Definition: ssht_adjoint.c:1304
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_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_adjoint_mw_inverse_sov_sym_ss_real_pole
void ssht_adjoint_mw_inverse_sov_sym_ss_real_pole(SSHT_COMPLEX(double) *flm, double *f, double f_np, double f_sp, int L, ssht_dl_method_t dl_method, int verbosity)
Definition: ssht_adjoint.c:2891
ssht_adjoint_mw_inverse_sov_sym_ss_real
void ssht_adjoint_mw_inverse_sov_sym_ss_real(SSHT_COMPLEX(double) *flm, double *f, int L, ssht_dl_method_t dl_method, int verbosity)
Definition: ssht_adjoint.c:1773
SSHT_ERROR_MEM_ALLOC_CHECK
#define SSHT_ERROR_MEM_ALLOC_CHECK(pointer)
Definition: ssht_error.h:28
ssht_adjoint_mw_inverse_sov_sym_pole
void ssht_adjoint_mw_inverse_sov_sym_pole(SSHT_COMPLEX(double) *flm, 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_adjoint.c:1398
ssht_adjoint_mw_forward_sov_sym_ss_pole
void ssht_adjoint_mw_forward_sov_sym_ss_pole(SSHT_COMPLEX(double) *f, SSHT_COMPLEX(double) *f_np, double *phi_np, SSHT_COMPLEX(double) *f_sp, double *phi_sp, SSHT_COMPLEX(double) *flm, int L, int spin, ssht_dl_method_t dl_method, int verbosity)
Definition: ssht_adjoint.c:2731
ssht_adjoint.h
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_adjoint_mw_inverse_sov_sym_ss
void ssht_adjoint_mw_inverse_sov_sym_ss(SSHT_COMPLEX(double) *flm, SSHT_COMPLEX(double) *f, int L, int spin, ssht_dl_method_t dl_method, int verbosity)
Definition: ssht_adjoint.c:1495
SSHT_PROMPT
#define SSHT_PROMPT
Definition: ssht_types.h:45
ssht_dl_get_offset
int ssht_dl_get_offset(int L, ssht_dl_size_t dl_size)
Definition: ssht_dl.c:111
ssht_adjoint_mw_forward_sov_sym_ss
void ssht_adjoint_mw_forward_sov_sym_ss(SSHT_COMPLEX(double) *f, SSHT_COMPLEX(double) *flm, int L, int spin, ssht_dl_method_t dl_method, int verbosity)
Definition: ssht_adjoint.c:2034
SSHT_COMPLEX
#define SSHT_COMPLEX(TYPE)
Definition: ssht_types.h:52
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_types.h
SSHT_DL_QUARTER
@ SSHT_DL_QUARTER
Definition: ssht_dl.h:15
ssht_adjoint_mw_forward_sov_sym_ss_real
void ssht_adjoint_mw_forward_sov_sym_ss_real(double *f, SSHT_COMPLEX(double) *flm, int L, ssht_dl_method_t dl_method, int verbosity)
Definition: ssht_adjoint.c:2378
ssht_adjoint_mw_inverse_sov_sym
void ssht_adjoint_mw_inverse_sov_sym(SSHT_COMPLEX(double) *flm, const SSHT_COMPLEX(double) *f, int L, int spin, ssht_dl_method_t dl_method, int verbosity)
Definition: ssht_adjoint.c:48