My Project
so3_adjoint.c
Go to the documentation of this file.
1 // S03 package to perform Wigner transform on the rotation group SO(3)
2 // Copyright (C) 2013 Martin Büttner and Jason McEwen
3 // See LICENSE.txt for license details
4 
13 #include <stdio.h>
14 #include <stdlib.h>
15 #include <string.h>
16 #include <math.h>
17 #include <complex.h> // Must be before fftw3.h
18 #include <fftw3.h>
19 
20 #include "ssht/ssht.h"
21 
22 #include "so3_types.h"
23 #include "so3_error.h"
24 #include "so3_sampling.h"
25 
26 #define MIN(a,b) ((a < b) ? (a) : (b))
27 #define MAX(a,b) ((a > b) ? (a) : (b))
28 
29 typedef void (*inverse_complex_ssht)(complex double *, const complex double *, int, int, int, ssht_dl_method_t, int);
30 typedef void (*inverse_real_ssht)(double *, const complex double *, int, int, ssht_dl_method_t, int);
31 typedef void (*forward_complex_ssht)(complex double *, const complex double *, int, int, int, ssht_dl_method_t, int);
32 typedef void (*forward_real_ssht)(complex double *, const double *, int, int, ssht_dl_method_t, int);
33 
35  complex double *flmn, const complex double *f,
36  const so3_parameters_t *parameters
37 ) {
38  int L0, L, N;
39  so3_sampling_t sampling;
40  so3_storage_t storage;
41  so3_n_mode_t n_mode;
42  ssht_dl_method_t dl_method;
43  int steerable;
44  int verbosity;
45 
46  L0 = parameters->L0;
47  L = parameters->L;
48  N = parameters->N;
49  sampling = parameters->sampling_scheme;
50  storage = parameters->storage;
51  // TODO: Add optimisations for all n-modes.
52  n_mode = parameters->n_mode;
53  dl_method = parameters->dl_method;
54  verbosity = parameters->verbosity;
55  steerable = parameters->steerable;
56 
57  // Print messages depending on verbosity level.
58  if (verbosity > 0)
59  {
60  printf("%sComputing adjoint inverse transform using MW sampling with\n", SO3_PROMPT);
61  printf("%sparameters (L, N, reality) = (%d, %d, FALSE)\n", SO3_PROMPT, L, N);
62  if (verbosity > 1)
63  printf("%sUsing routine so3_adjoint_inverse_direct with storage method %d...\n"
64  , SO3_PROMPT
65  , storage);
66  }
67 
68  int m_stride = 2*L-1;
69  int m_offset = L-1;
70  // unused: int n_stride = 2*N-1;
71  int n_offset = N-1;
72  int mm_stride = 2*L-1;
73  int mm_offset = L-1;
74  int a_stride = 2*L-1;
75  int b_stride = L;
76  int bext_stride = 2*L-1;
77  // unused: int g_stride = 2*N-1;
78 
79  int n_start, n_stop, n_inc;
80 
81  switch (n_mode)
82  {
83  case SO3_N_MODE_ALL:
84  case SO3_N_MODE_L:
85  n_start = -N+1;
86  n_stop = N-1;
87  n_inc = 1;
88  break;
89  case SO3_N_MODE_EVEN:
90  n_start = ((N-1) % 2 == 0) ? -N+1 : -N+2;
91  n_stop = ((N-1) % 2 == 0) ? N-1 : N-2;
92  n_inc = 2;
93  break;
94  case SO3_N_MODE_ODD:
95  n_start = ((N-1) % 2 != 0) ? -N+1 : -N+2;
96  n_stop = ((N-1) % 2 != 0) ? N-1 : N-2;
97  n_inc = 2;
98  break;
99  case SO3_N_MODE_MAXIMUM:
100  n_start = -N+1;
101  n_stop = N-1;
102  n_inc = MAX(1,2*N - 2);
103  break;
104  default:
105  SO3_ERROR_GENERIC("Invalid n-mode.");
106  }
107 
108  double *sqrt_tbl = calloc(2*(L-1)+2, sizeof(*sqrt_tbl));
109  SO3_ERROR_MEM_ALLOC_CHECK(sqrt_tbl);
110  double *signs = calloc(L+1, sizeof(*signs));
111  SO3_ERROR_MEM_ALLOC_CHECK(signs);
112  complex double *exps = calloc(4, sizeof(*exps));
113  SO3_ERROR_MEM_ALLOC_CHECK(exps);
114  complex double *expsmm = calloc(2*L-1, sizeof(*expsmm));
115  SO3_ERROR_MEM_ALLOC_CHECK(expsmm);
116 
117  int el, m, n, mm; // mm is for m'
118  // Perform precomputations.
119  for (el = 0; el <= 2*(L-1)+1; ++el)
120  sqrt_tbl[el] = sqrt((double)el);
121  for (m = 0; m <= L-1; m += 2)
122  {
123  signs[m] = 1.0;
124  signs[m+1] = -1.0;
125  }
126  int i;
127  for (i = 0; i < 4; ++i)
128  exps[i] = cexp(I*SO3_PION2*i);
129  for (mm = -L+1; mm <= L-1; ++mm)
130  expsmm[mm + mm_offset] = cexp(-I*mm*SSHT_PI/(2.0*L-1.0));
131 
132  // Compute Fourier transform over alpha and gamma, i.e. compute Fmn(b).
133  complex double *Fmnb = calloc((2*L-1)*(2*L-1)*(2*N-1), sizeof(*Fmnb));
134  SO3_ERROR_MEM_ALLOC_CHECK(Fmnb);
135  complex double *inout = calloc((2*L-1)*(2*N-1), sizeof(*inout));
136  SO3_ERROR_MEM_ALLOC_CHECK(inout);
137  fftw_plan plan = fftw_plan_dft_2d(
138  2*N-1, 2*L-1,
139  inout, inout,
140  FFTW_FORWARD,
141  FFTW_ESTIMATE);
142 
143  int b, g;
144  for (b = 0; b < L; ++b)
145  {
146  // TODO: This memcpy loop could probably be avoided by using
147  // a more elaborate FFTW plan which performs the FFT directly
148  // over the 1st and 3rd dimensions of f.
149  for (g = 0; g < 2*N-1; ++g)
150  memcpy(
151  inout + g*a_stride,
152  f + 0 + a_stride*(
153  b + b_stride*(
154  g)),
155  a_stride*sizeof(*f));
156  fftw_execute_dft(plan, inout, inout);
157 
158  // Apply spatial shift
159  for (n = n_start; n <= n_stop; n += n_inc)
160  {
161  int n_shift = n < 0 ? 2*N-1 : 0;
162  for (m = -L+1; m <= L-1; ++m)
163  {
164  int m_shift = m < 0 ? 2*L-1 : 0;
165  Fmnb[b + bext_stride*(
166  m + m_offset + m_stride*(
167  n + n_offset))] =
168  inout[m + m_shift + m_stride*(
169  n + n_shift)];
170  }
171  }
172  }
173  fftw_destroy_plan(plan);
174 
175  // Extend Fmnb by filling it with zeroes.
176  for (n = n_start; n <= n_stop; n += n_inc)
177  for (m = -L+1; m <= L-1; ++m)
178  for (b = L; b < 2*L-1; ++b)
179  Fmnb[b + bext_stride*(
180  m + m_offset + m_stride*(
181  n + n_offset))] = 0.0;
182 
183 
184  // Compute Fourier transform over beta, i.e. compute Fmnm'.
185  complex double *Fmnm = calloc((2*L-1)*(2*L-1)*(2*N-1), sizeof(*Fmnm));
186  SO3_ERROR_MEM_ALLOC_CHECK(Fmnm);
187 
188  plan = fftw_plan_dft_1d(
189  2*L-1,
190  inout, inout,
191  FFTW_FORWARD,
192  FFTW_ESTIMATE);
193  for (n = n_start; n <= n_stop; n += n_inc)
194  for (m = -L+1; m <= L-1; ++m)
195  {
196  memcpy(inout,
197  Fmnb + 0 + bext_stride*(
198  m + m_offset + m_stride*(
199  n + n_offset)),
200  bext_stride*sizeof(*Fmnb));
201  fftw_execute_dft(plan, inout, inout);
202 
203  // Apply spatial shift
204  for (mm = -L+1; mm <= L-1; ++mm)
205  {
206  int mm_shift = mm < 0 ? 2*L-1 : 0;
207  Fmnm[m + m_offset + m_stride*(
208  mm + mm_offset + mm_stride*(
209  n + n_offset))] =
210  inout[mm + mm_shift];
211  }
212  }
213  fftw_destroy_plan(plan);
214  free(inout);
215 
216  // Apply phase modulation to account for sampling offset.
217  for (n = n_start; n <= n_stop; n += n_inc)
218  for (mm = -L+1; mm <= L-1; ++mm)
219  for (m = -L+1; m <= L-1; ++m)
220  Fmnm[m + m_offset + m_stride*(
221  mm + mm_offset + mm_stride*(
222  n + n_offset))] *=
223  expsmm[mm + mm_offset];
224 
225  // Compute flmn.
226  double *dl, *dl8 = NULL;
227  dl = ssht_dl_calloc(L, SSHT_DL_QUARTER);
228  SO3_ERROR_MEM_ALLOC_CHECK(dl);
229  if (dl_method == SSHT_DL_RISBO)
230  {
231  dl8 = ssht_dl_calloc(L, SSHT_DL_QUARTER_EXTENDED);
232  SO3_ERROR_MEM_ALLOC_CHECK(dl8);
233  }
234  int dl_offset = ssht_dl_get_offset(L, SSHT_DL_QUARTER);
235  int dl_stride = ssht_dl_get_stride(L, SSHT_DL_QUARTER);
236  for (n = -N+1; n <= N-1; ++n)
237  for (el = abs(n); el < L; ++el)
238  for (m = -el; m <= el; ++m)
239  {
240  int ind;
241  so3_sampling_elmn2ind(&ind, el, m, n, parameters);
242  flmn[ind] = 0.0;
243  }
244 
245  for (el = L0; el < L; ++el)
246  {
247  int eltmp;
248 
249  // Compute Wigner plane.
250  switch (dl_method)
251  {
252  case SSHT_DL_RISBO:
253  if (el != 0 && el == L0)
254  {
255  for(eltmp = 0; eltmp <= L0; ++eltmp)
256  ssht_dl_beta_risbo_eighth_table(dl8, SO3_PION2, L,
257  SSHT_DL_QUARTER_EXTENDED,
258  eltmp, sqrt_tbl, signs);
259  ssht_dl_beta_risbo_fill_eighth2quarter_table(dl,
260  dl8, L,
261  SSHT_DL_QUARTER,
262  SSHT_DL_QUARTER_EXTENDED,
263  el,
264  signs);
265  }
266  else
267  {
268  ssht_dl_beta_risbo_eighth_table(dl8, SO3_PION2, L,
269  SSHT_DL_QUARTER_EXTENDED,
270  el, sqrt_tbl, signs);
271  ssht_dl_beta_risbo_fill_eighth2quarter_table(dl,
272  dl8, L,
273  SSHT_DL_QUARTER,
274  SSHT_DL_QUARTER_EXTENDED,
275  el,
276  signs);
277  }
278  break;
279 
280  case SSHT_DL_TRAPANI:
281  if (el != 0 && el == L0)
282  {
283  for(eltmp = 0; eltmp <= L0; ++eltmp)
284  ssht_dl_halfpi_trapani_eighth_table(dl, L,
285  SSHT_DL_QUARTER,
286  eltmp, sqrt_tbl);
287  ssht_dl_halfpi_trapani_fill_eighth2quarter_table(dl, L,
288  SSHT_DL_QUARTER,
289  el, signs);
290  }
291  else
292  {
293  ssht_dl_halfpi_trapani_eighth_table(dl, L,
294  SSHT_DL_QUARTER,
295  el, sqrt_tbl);
296  ssht_dl_halfpi_trapani_fill_eighth2quarter_table(dl, L,
297  SSHT_DL_QUARTER,
298  el, signs);
299  }
300  break;
301 
302  default:
303  SO3_ERROR_GENERIC("Invalid dl method");
304  }
305 
306  // Compute flmn for current el.
307 
308  switch (n_mode)
309  {
310  case SO3_N_MODE_ALL:
311  n_start = MAX(-N+1,-el);
312  n_stop = MIN( N-1, el);
313  n_inc = 1;
314  break;
315  case SO3_N_MODE_EVEN:
316  n_start = MAX(-N+1,-el);
317  n_start += (-n_start)%2;
318  n_stop = MIN( N-1, el);
319  n_stop -= n_stop%2;
320  n_inc = 2;
321  break;
322  case SO3_N_MODE_ODD:
323  n_start = MAX(-N+1,-el);
324  n_start += 1+n_start%2;
325  n_stop = MIN( N-1, el);
326  n_stop -= 1-n_stop%2;
327  n_inc = 2;
328  break;
329  case SO3_N_MODE_MAXIMUM:
330  if (el < N-1)
331  continue;
332  n_start = -N+1;
333  n_stop = N-1;
334  n_inc = MAX(1,2*N-2);
335  break;
336  case SO3_N_MODE_L:
337  if (el >= N)
338  continue;
339  n_start = -el;
340  n_stop = el;
341  n_inc = MAX(1,2*el);
342  break;
343  default:
344  SO3_ERROR_GENERIC("Invalid n-mode.");
345  }
346 
347  // Factor which depends only on el.
348  double elfactor = (2.0*el+1.0)/(8.0*SO3_PI*SO3_PI);
349 
350  // TODO: Pull out a few multiplications into precomputations
351  // or split up loops to avoid conditionals to check signs.
352  for (mm = -el; mm <= el; ++mm)
353  {
354  // These signs are needed for the symmetry relations of
355  // Wigner symbols.
356  double elmmsign = signs[el] * signs[abs(mm)];
357 
358  for (n = n_start; n <= n_stop; n += n_inc)
359  {
360  double mmsign = mm >= 0 ? 1.0 : signs[el] * signs[abs(n)];
361  double elnsign = n >= 0 ? 1.0 : elmmsign;
362 
363  // Factor which does not depend on m.
364  double elnmm_factor = mmsign * elnsign * elfactor
365  * dl[abs(n) + dl_offset + abs(mm)*dl_stride];
366 
367  for (m = -el; m <= el; ++m)
368  {
369  mmsign = mm >= 0 ? 1.0 : signs[el] * signs[abs(m)];
370  double elmsign = m >= 0 ? 1.0 : elmmsign;
371  int ind;
372  so3_sampling_elmn2ind(&ind, el, m, n, parameters);
373  int mod = ((m-n)%4 + 4)%4;
374  flmn[ind] +=
375  exps[mod]
376  * elnmm_factor
377  * mmsign * elmsign
378  * dl[abs(m) + dl_offset + abs(mm)*dl_stride]
379  * Fmnm[m + m_offset + m_stride*(
380  mm + mm_offset + mm_stride*(
381  n + n_offset))];
382 
383  }
384  }
385  }
386  }
387 
388  free(dl);
389  if (dl_method == SSHT_DL_RISBO)
390  free(dl8);
391  free(Fmnb);
392  free(Fmnm);
393  free(sqrt_tbl);
394  free(signs);
395  free(exps);
396  free(expsmm);
397 
398  if (verbosity > 0)
399  printf("%sAdjoint inverse transform computed!\n", SO3_PROMPT);
400 }
401 
403  complex double *f, const complex double *flmn,
404  const so3_parameters_t *parameters
405 ) {
406  int L0, L, N;
407  so3_sampling_t sampling;
408  so3_storage_t storage;
409  so3_n_mode_t n_mode;
410  ssht_dl_method_t dl_method;
411  int steerable;
412  int verbosity;
413 
414  L0 = parameters->L0;
415  L = parameters->L;
416  N = parameters->N;
417  sampling = parameters->sampling_scheme;
418  storage = parameters->storage;
419  // TODO: Add optimisations for all n-modes.
420  n_mode = parameters->n_mode;
421  dl_method = parameters->dl_method;
422  verbosity = parameters->verbosity;
423  steerable = parameters->steerable;
424 
425  // Print messages depending on verbosity level.
426  if (verbosity > 0)
427  {
428  printf("%sComputing forward adjoint transform using MW sampling with\n", SO3_PROMPT);
429  printf("%sparameters (L, N, reality) = (%d, %d, FALSE)\n", SO3_PROMPT, L, N);
430  if (verbosity > 1)
431  printf("%sUsing routine so3_adjoint_forward_direct with storage method %d...\n"
432  , SO3_PROMPT
433  , storage);
434  }
435 
436  // Iterators
437  int el, m, n, mm, b, g; // mm for m'
438 
439  int m_offset = L-1;
440  int m_stride = 2*L-1;
441  int n_offset = N-1;
442  // unused: int n_stride = 2*N-1;
443  int mm_offset = L-1;
444  int mm_stride = 2*L-1;
445  int a_stride = 2*L-1;
446  int b_stride = L;
447  int bext_stride = 2*L-1;
448  // unused: int g_stride = 2*N-1;
449 
450  complex double* inout; // Used as temporary storage for various FFTs.
451 
452  // Allocate memory.
453  double *sqrt_tbl = calloc(2*(L-1)+2, sizeof(*sqrt_tbl));
454  SO3_ERROR_MEM_ALLOC_CHECK(sqrt_tbl);
455  double *signs = calloc(L+1, sizeof(*signs));
456  SO3_ERROR_MEM_ALLOC_CHECK(signs);
457  complex double *exps = calloc(4, sizeof(*exps));
458  SO3_ERROR_MEM_ALLOC_CHECK(exps);
459  complex double *expsmm = calloc(2*L-1, sizeof(*expsmm));
460  SO3_ERROR_MEM_ALLOC_CHECK(expsmm);
461 
462  // Perform precomputations.
463  for (el = 0; el <= 2*L-1; ++el)
464  sqrt_tbl[el] = sqrt((double)el);
465  for (m = 0; m <= L-1; m += 2)
466  {
467  signs[m] = 1.0;
468  signs[m+1] = -1.0;
469  }
470  int i;
471  for (i = 0; i < 4; ++i)
472  exps[i] = cexp(I*SO3_PION2*i);
473  for (mm = -L+1; mm <= L-1; ++mm)
474  expsmm[mm + mm_offset] = cexp(I*mm*SSHT_PI/(2.0*L-1.0));
475 
476  // Compute Gmnm'
477  // TODO: Currently m is fastest-varying, then n, then m'.
478  // Should this order be changed to m-m'-n?
479  complex double *Gmnm = calloc((2*L-1)*(2*L-1)*(2*N-1), sizeof(*Gmnm));
480  SO3_ERROR_MEM_ALLOC_CHECK(Gmnm);
481 
482  int n_start, n_stop, n_inc;
483 
484 
485  double *dl = ssht_dl_calloc(L, SSHT_DL_QUARTER);
486  SO3_ERROR_MEM_ALLOC_CHECK(dl);
487  double *dl8 = NULL;
488  if (dl_method == SSHT_DL_RISBO)
489  {
490  dl8 = ssht_dl_calloc(L, SSHT_DL_QUARTER_EXTENDED);
491  SO3_ERROR_MEM_ALLOC_CHECK(dl8);
492  }
493  int dl_offset = ssht_dl_get_offset(L, SSHT_DL_QUARTER);
494  int dl_stride = ssht_dl_get_stride(L, SSHT_DL_QUARTER);
495 
496  complex double *mn_factors = calloc((2*L-1)*(2*N-1), sizeof *mn_factors);
497  SO3_ERROR_MEM_ALLOC_CHECK(mn_factors);
498 
499  // TODO: SSHT starts this loop from MAX(L0, abs(spin)).
500  // Can we use a similar optimisation? el can probably
501  // be limited by n, but then we'd need to switch the
502  // loop order, which means we'd have to recompute the
503  // Wigner plane for each n. That seems wrong?
504  for (el = L0; el <= L-1; ++el)
505  {
506  int eltmp;
507  // Compute Wigner plane.
508  switch (dl_method)
509  {
510  case SSHT_DL_RISBO:
511  if (el != 0 && el == L0)
512  {
513  for(eltmp = 0; eltmp <= L0; ++eltmp)
514  ssht_dl_beta_risbo_eighth_table(
515  dl8,
516  SO3_PION2,
517  L,
518  SSHT_DL_QUARTER_EXTENDED,
519  eltmp,
520  sqrt_tbl,
521  signs);
522  ssht_dl_beta_risbo_fill_eighth2quarter_table(dl,
523  dl8, L,
524  SSHT_DL_QUARTER,
525  SSHT_DL_QUARTER_EXTENDED,
526  el,
527  signs);
528  }
529  else
530  {
531  ssht_dl_beta_risbo_eighth_table(
532  dl8,
533  SO3_PION2,
534  L,
535  SSHT_DL_QUARTER_EXTENDED,
536  el,
537  sqrt_tbl,
538  signs);
539  ssht_dl_beta_risbo_fill_eighth2quarter_table(
540  dl,
541  dl8, L,
542  SSHT_DL_QUARTER,
543  SSHT_DL_QUARTER_EXTENDED,
544  el,
545  signs);
546  }
547  break;
548 
549  case SSHT_DL_TRAPANI:
550  if (el != 0 && el == L0)
551  {
552  for(eltmp = 0; eltmp <= L0; ++eltmp)
553  ssht_dl_halfpi_trapani_eighth_table(
554  dl,
555  L,
556  SSHT_DL_QUARTER,
557  eltmp,
558  sqrt_tbl);
559  ssht_dl_halfpi_trapani_fill_eighth2quarter_table(
560  dl,
561  L,
562  SSHT_DL_QUARTER,
563  el,
564  signs);
565  }
566  else
567  {
568  ssht_dl_halfpi_trapani_eighth_table(
569  dl,
570  L,
571  SSHT_DL_QUARTER,
572  el,
573  sqrt_tbl);
574  ssht_dl_halfpi_trapani_fill_eighth2quarter_table(
575  dl,
576  L,
577  SSHT_DL_QUARTER,
578  el,
579  signs);
580  }
581  break;
582 
583  default:
584  SO3_ERROR_GENERIC("Invalid dl method");
585  }
586 
587  // Compute Gmnm' contribution for current el.
588 
589  switch (n_mode)
590  {
591  case SO3_N_MODE_ALL:
592  n_start = MAX(-N+1,-el);
593  n_stop = MIN( N-1, el);
594  n_inc = 1;
595  break;
596  case SO3_N_MODE_EVEN:
597  n_start = MAX(-N+1,-el);
598  n_start += (-n_start)%2;
599  n_stop = MIN( N-1, el);
600  n_stop -= n_stop%2;
601  n_inc = 2;
602  break;
603  case SO3_N_MODE_ODD:
604  n_start = MAX(-N+1,-el);
605  n_start += 1+n_start%2;
606  n_stop = MIN( N-1, el);
607  n_stop -= 1-n_stop%2;
608  n_inc = 2;
609  break;
610  case SO3_N_MODE_MAXIMUM:
611  if (el < N-1)
612  continue;
613  n_start = -N+1;
614  n_stop = N-1;
615  n_inc = MAX(1,2*N-2);
616  break;
617  case SO3_N_MODE_L:
618  if (el >= N)
619  continue;
620  n_start = -el;
621  n_stop = el;
622  n_inc = MAX(1,2*el);
623  break;
624  default:
625  SO3_ERROR_GENERIC("Invalid n-mode.");
626  }
627 
628  // Factors which do not depend on m'.
629  for (n = n_start; n <= n_stop; n += n_inc)
630  for (m = -el; m <= el; ++m)
631  {
632  int ind;
633  so3_sampling_elmn2ind(&ind, el, m, n, parameters);
634  int mod = ((n-m)%4 + 4)%4;
635  mn_factors[m + m_offset + m_stride*(
636  n + n_offset)] =
637  flmn[ind] * exps[mod];
638  }
639 
640  for (mm = 0; mm <= el; ++mm)
641  {
642  // These signs are needed for the symmetry relations of
643  // Wigner symbols.
644  double elmmsign = signs[el] * signs[mm];
645 
646  // TODO: If the conditional for elnsign is a bottleneck
647  // this loop can be split up just like the inner loop.
648  for (n = n_start; n <= n_stop; n += n_inc)
649  {
650  double elnsign = n >= 0 ? 1.0 : elmmsign;
651  // Factor which does not depend on m.
652  double elnmm_factor = elnsign
653  * dl[abs(n) + dl_offset + mm*dl_stride];
654  for (m = -el; m < 0; ++m)
655  Gmnm[mm + mm_offset + mm_stride*(
656  m + m_offset + m_stride*(
657  n + n_offset))] +=
658  elnmm_factor
659  * mn_factors[m + m_offset + m_stride*(
660  n + n_offset)]
661  * elmmsign * dl[-m + dl_offset + mm*dl_stride];
662  for (m = 0; m <= el; ++m)
663  Gmnm[mm + mm_offset + mm_stride*(
664  m + m_offset + m_stride*(
665  n + n_offset))] +=
666  elnmm_factor
667  * mn_factors[m + m_offset + m_stride*(
668  n + n_offset)]
669  * dl[m + dl_offset + mm*dl_stride];
670  }
671  }
672  }
673 
674  // Free dl memory.
675  free(dl);
676  if (dl_method == SSHT_DL_RISBO)
677  free(dl8);
678  free(mn_factors);
679 
680  switch (n_mode)
681  {
682  case SO3_N_MODE_ALL:
683  case SO3_N_MODE_L:
684  n_start = -N+1;
685  n_stop = N-1;
686  n_inc = 1;
687  break;
688  case SO3_N_MODE_EVEN:
689  n_start = ((N-1) % 2 == 0) ? -N+1 : -N+2;
690  n_stop = ((N-1) % 2 == 0) ? N-1 : N-2;
691  n_inc = 2;
692  break;
693  case SO3_N_MODE_ODD:
694  n_start = ((N-1) % 2 != 0) ? -N+1 : -N+2;
695  n_stop = ((N-1) % 2 != 0) ? N-1 : N-2;
696  n_inc = 2;
697  break;
698  case SO3_N_MODE_MAXIMUM:
699  n_start = -N+1;
700  n_stop = N-1;
701  n_inc = MAX(1,2*N - 2);
702  break;
703  default:
704  SO3_ERROR_GENERIC("Invalid n-mode.");
705  }
706 
707  // Use symmetry to compute Gmnm' for negative m'.
708  for (n = n_start; n <= n_stop; n += n_inc)
709  for (m = -L+1; m <= L-1; ++m)
710  for (mm = -L+1; mm < 0; ++mm)
711  Gmnm[mm + mm_offset + mm_stride*(
712  m + m_offset + m_stride*(
713  n + n_offset))] =
714  signs[abs(m+n)%2]
715  * Gmnm[-mm + mm_offset + mm_stride*(
716  m + m_offset + m_stride*(
717  n + n_offset))];
718 
719 
720  // Compute weights.
721  complex double *w = calloc(4*L-3, sizeof(*w));
722  SO3_ERROR_MEM_ALLOC_CHECK(w);
723  int w_offset = 2*(L-1);
724  for (mm = -2*(L-1); mm <= 2*(L-1); ++mm)
725  w[mm+w_offset] = so3_sampling_weight(parameters, mm);
726 
727  // Compute IFFT of w to give wr.
728  complex double *wr = calloc(4*L-3, sizeof(*w));
729  SO3_ERROR_MEM_ALLOC_CHECK(wr);
730  inout = calloc(4*L-3, sizeof(*inout));
731  SO3_ERROR_MEM_ALLOC_CHECK(inout);
732  fftw_plan plan_bwd = fftw_plan_dft_1d(
733  4*L-3,
734  inout, inout,
735  FFTW_BACKWARD,
736  FFTW_MEASURE);
737  fftw_plan plan_fwd = fftw_plan_dft_1d(
738  4*L-3,
739  inout,
740  inout,
741  FFTW_FORWARD,
742  FFTW_MEASURE);
743 
744  // Apply spatial shift.
745  for (mm = 1; mm <= 2*L-2; ++mm)
746  inout[mm + w_offset] = w[mm - 2*(L-1) - 1 + w_offset];
747  for (mm = -2*(L-1); mm <= 0; ++mm)
748  inout[mm + w_offset] = w[mm + 2*(L-1) + w_offset];
749 
750  fftw_execute_dft(plan_bwd, inout, inout);
751 
752  // Apply spatial shift.
753  for (mm = 0; mm <= 2*L-2; ++mm)
754  wr[mm + w_offset] = inout[mm - 2*(L-1) + w_offset];
755  for (mm = -2*(L-1); mm <= -1; ++mm)
756  wr[mm + w_offset] = inout[mm + 2*(L-1) + 1 + w_offset];
757 
758  // Compute Fmnm'' by convolution implemented as product in real space.
759  complex double *Gmnm_pad = calloc(4*L-3, sizeof(*Gmnm_pad));
760  SO3_ERROR_MEM_ALLOC_CHECK(Gmnm_pad);
761  complex double *Fmnm = calloc((2*L-1)*(2*L-1)*(2*N-1), sizeof(*Fmnm));
762  SO3_ERROR_MEM_ALLOC_CHECK(Fmnm);
763  for (n = n_start; n <= n_stop; n += n_inc)
764  for (m = -L+1; m <= L-1; ++m)
765  {
766 
767  // Zero-pad Gmnm'.
768  for (mm = -2*(L-1); mm <= -L; ++mm)
769  Gmnm_pad[mm+w_offset] = 0.0;
770  for (mm = L; mm <= 2*(L-1); ++mm)
771  Gmnm_pad[mm+w_offset] = 0.0;
772  for (mm = -(L-1); mm <= L-1; ++mm)
773  Gmnm_pad[mm + w_offset] =
774  Gmnm[mm + mm_offset + mm_stride*(
775  m + m_offset + m_stride*(
776  n + n_offset))];
777 
778  // Apply spatial shift.
779  for (mm = 1; mm <= 2*L-2; ++mm)
780  inout[mm + w_offset] = Gmnm_pad[mm - 2*(L-1) - 1 + w_offset];
781  for (mm = -2*(L-1); mm <= 0; ++mm)
782  inout[mm + w_offset] = Gmnm_pad[mm + 2*(L-1) + w_offset];
783  // Compute IFFT of Gmnm'.
784  fftw_execute_dft(plan_bwd, inout, inout);
785  // Apply spatial shift.
786  for (mm = 0; mm <= 2*L-2; ++mm)
787  Gmnm_pad[mm + w_offset] = inout[mm - 2*(L-1) + w_offset];
788  for (mm = -2*(L-1); mm <= -1; ++mm)
789  Gmnm_pad[mm + w_offset] = inout[mm + 2*(L-1) + 1 + w_offset];
790 
791  // Compute product of Gmnm' and weight in real space.
792  int r;
793  for (r = -2*(L-1); r <= 2*(L-1); ++r)
794  Gmnm_pad[r + w_offset] *= wr[r + w_offset];
795 
796  // Apply spatial shift.
797  for (mm = 1; mm <= 2*L-2; ++mm)
798  inout[mm + w_offset] = Gmnm_pad[mm - 2*(L-1) - 1 + w_offset];
799  for (mm = -2*(L-1); mm <= 0; ++mm)
800  inout[mm + w_offset] = Gmnm_pad[mm + 2*(L-1) + w_offset];
801  // Compute Fmnm'' by FFT.
802  fftw_execute_dft(plan_fwd, inout, inout);
803  // Apply spatial shift.
804  for (mm = 0; mm <= 2*L-2; ++mm)
805  Gmnm_pad[mm + w_offset] = inout[mm - 2*(L-1) + w_offset];
806  for (mm = -2*(L-1); mm <= -1; ++mm)
807  Gmnm_pad[mm + w_offset] = inout[mm + 2*(L-1) + 1 + w_offset];
808 
809  // Extract section of Fmnm'' of interest.
810  for (mm = -(L-1); mm <= L-1; ++mm)
811  Fmnm[mm + mm_offset + mm_stride*(
812  m + m_offset + m_stride*(
813  n + n_offset))] =
814  Gmnm_pad[mm + w_offset]
815  * 4.0 * SSHT_PI * SSHT_PI / (4.0*L-3.0);
816 
817  }
818  fftw_destroy_plan(plan_bwd);
819  fftw_destroy_plan(plan_fwd);
820 
821  // Apply phase modulation to account for sampling offset.
822  for (n = n_start; n <= n_stop; n += n_inc)
823  for (m = -L+1; m <= L-1; ++m)
824  for (mm = -L+1; mm <= L-1; ++mm)
825  Fmnm[mm + mm_offset + mm_stride*(
826  m + m_offset + m_stride*(
827  n + n_offset))] *=
828  expsmm[mm + mm_offset];
829 
830  // Compute Fourier transform over mm, i.e. compute Fmnb.
831  complex double *Fmnb = calloc((2*L-1)*(2*L-1)*(2*N-1), sizeof(*Fmnb));
832  SO3_ERROR_MEM_ALLOC_CHECK(Fmnm);
833 
834  fftw_plan plan = fftw_plan_dft_1d(
835  2*L-1,
836  inout, inout,
837  FFTW_BACKWARD,
838  FFTW_ESTIMATE);
839  for (n = n_start; n <= n_stop; n += n_inc)
840  for (m = -L+1; m <= L-1; ++m)
841  {
842  // Apply spatial shift and normalisation factor
843  for (mm = -L+1; mm <= L-1; ++mm)
844  {
845  int mm_shift = mm < 0 ? 2*L-1 : 0;
846  inout[mm + mm_shift] =
847  Fmnm[mm + mm_offset + mm_stride*(
848  m + m_offset + m_stride*(
849  n + n_offset))] / (2.0*L-1.0);
850  }
851  fftw_execute_dft(plan, inout, inout);
852  memcpy(Fmnb + 0 + bext_stride*(
853  m + m_offset + m_stride*(
854  n + n_offset)),
855  inout,
856  bext_stride*sizeof(*Fmnb));
857 
858  }
859 
860  fftw_destroy_plan(plan);
861  free(inout);
862 
863 
864 
865  // Adjoint of periodic extension of Ftm.
866  for (n = n_start; n <= n_stop; n += n_inc)
867  for (m = -L+1; m <= L-1; ++m)
868  {
869  int signmn = signs[abs(m+n)%2];
870  for (b = 0; b <= L-2; ++b)
871  Fmnb[b + bext_stride*(
872  m + m_offset + m_stride*(
873  n + n_offset))] +=
874  signmn
875  * Fmnb[(2*L-2-b) + bext_stride*(
876  m + m_offset + m_stride*(
877  n + n_offset))];
878  }
879 
880  // Compute Fourier transform over alpha and gamma, i.e. compute f.
881  inout = calloc((2*L-1)*(2*N-1), sizeof(*inout));
882  SO3_ERROR_MEM_ALLOC_CHECK(inout);
883  plan = fftw_plan_dft_2d(
884  2*N-1, 2*L-1,
885  inout, inout,
886  FFTW_BACKWARD,
887  FFTW_ESTIMATE);
888 
889  double norm_factor = 1.0/(2.0*L-1.0)/(2.0*N-1.0);
890 
891  for (b = 0; b < L; ++b)
892  {
893  for (int i_count=0; i_count<(2*L-1)*(2*N-1); i_count++) {inout[i_count] = 0.0;}
894 
895  // Apply spatial shift and normalisation factor
896  for (n = n_start; n <= n_stop; n += n_inc)
897  {
898  int n_shift = n < 0 ? 2*N-1 : 0;
899  for (m = -L+1; m <= L-1; ++m)
900  {
901  int m_shift = m < 0 ? 2*L-1 : 0;
902  inout[m + m_shift + m_stride*(
903  n + n_shift)] =
904  Fmnb[b + bext_stride*(
905  m + m_offset + m_stride*(
906  n + n_offset))] * norm_factor;
907  }
908  }
909  fftw_execute_dft(plan, inout, inout);
910 
911  // TODO: This memcpy loop could probably be avoided by using
912  // a more elaborate FFTW plan which performs the FFT directly
913  // over the 1st and 3rd dimensions of f.
914  for (g = 0; g < 2*N-1; ++g)
915  memcpy(
916  f + 0 + a_stride*(
917  b + b_stride*(
918  g)),
919  inout + g*a_stride,
920  a_stride*sizeof(*f));
921 
922 
923  }
924  fftw_destroy_plan(plan);
925 
926  free(Fmnb);
927  free(Fmnm);
928  free(inout);
929  free(w);
930  free(wr);
931  free(Gmnm_pad);
932  free(Gmnm);
933  free(sqrt_tbl);
934  free(signs);
935  free(exps);
936  free(expsmm);
937 }
938 
940  complex double *flmn, const double *f,
941  const so3_parameters_t *parameters
942 ) {
943  int L0, L, N;
944  so3_sampling_t sampling;
945  so3_storage_t storage;
946  so3_n_mode_t n_mode;
947  ssht_dl_method_t dl_method;
948  int steerable;
949  int verbosity;
950 
951  L0 = parameters->L0;
952  L = parameters->L;
953  N = parameters->N;
954  sampling = parameters->sampling_scheme;
955  storage = parameters->storage;
956  // TODO: Add optimisations for all n-modes.
957  n_mode = parameters->n_mode;
958  dl_method = parameters->dl_method;
959  verbosity = parameters->verbosity;
960  steerable = parameters->steerable;
961 
962  // Print messages depending on verbosity level.
963  if (verbosity > 0) {
964  printf("%sComputing adjoint inverse transform using MW sampling with\n", SO3_PROMPT);
965  printf("%sparameters (L, N, reality) = (%d, %d, FALSE)\n", SO3_PROMPT, L, N);
966  if (verbosity > 1)
967  printf("%sUsing routine so3_adjoint_inverse_direct_real with storage method %d...\n"
968  , SO3_PROMPT
969  , storage);
970  }
971 
972  int m_stride = 2*L-1;
973  int m_offset = L-1;
974  int n_offset = 0;
975  int n_stride = N;
976  int mm_stride = 2*L-1;
977  int mm_offset = L-1;
978  int a_stride = 2*L-1;
979  int b_stride = L;
980  int bext_stride = 2*L-1;
981  int g_stride = 2*N-1;
982 
983  int n_start, n_stop, n_inc;
984 
985  switch (n_mode)
986  {
987  case SO3_N_MODE_ALL:
988  case SO3_N_MODE_L:
989  n_start = 0;
990  n_stop = N-1;
991  n_inc = 1;
992  break;
993  case SO3_N_MODE_EVEN:
994  n_start = 0;
995  n_stop = ((N-1) % 2 == 0) ? N-1 : N-2;
996  n_inc = 2;
997  break;
998  case SO3_N_MODE_ODD:
999  n_start = 1;
1000  n_stop = ((N-1) % 2 != 0) ? N-1 : N-2;
1001  n_inc = 2;
1002  break;
1003  case SO3_N_MODE_MAXIMUM:
1004  n_start = N-1;
1005  n_stop = N-1;
1006  n_inc = 1;
1007  break;
1008  default:
1009  SO3_ERROR_GENERIC("Invalid n-mode.");
1010  }
1011 
1012  double *sqrt_tbl = calloc(2*(L-1)+2, sizeof(*sqrt_tbl));
1013  SO3_ERROR_MEM_ALLOC_CHECK(sqrt_tbl);
1014  double *signs = calloc(L+1, sizeof(*signs));
1015  SO3_ERROR_MEM_ALLOC_CHECK(signs);
1016  complex double *exps = calloc(4, sizeof(*exps));
1017  SO3_ERROR_MEM_ALLOC_CHECK(exps);
1018  complex double *expsmm = calloc(2*L-1, sizeof(*expsmm));
1019  SO3_ERROR_MEM_ALLOC_CHECK(expsmm);
1020 
1021  int el, m, n, mm; // mm is for m'
1022  // Perform precomputations.
1023  for (el = 0; el <= 2*(L-1)+1; ++el)
1024  sqrt_tbl[el] = sqrt((double)el);
1025  for (m = 0; m <= L-1; m += 2)
1026  {
1027  signs[m] = 1.0;
1028  signs[m+1] = -1.0;
1029  }
1030  int i;
1031  for (i = 0; i < 4; ++i)
1032  exps[i] = cexp(I*SO3_PION2*i);
1033  for (mm = -L+1; mm <= L-1; ++mm)
1034  expsmm[mm + mm_offset] = cexp(-I*mm*SSHT_PI/(2.0*L-1.0));
1035 
1036  // Compute Fourier transform over alpha and gamma, i.e. compute Fmn(b).
1037  complex double *Fmnb = calloc((2*L-1)*(2*L-1)*N, sizeof(*Fmnb));
1038  SO3_ERROR_MEM_ALLOC_CHECK(Fmnb);
1039  double *fft_in = calloc((2*L-1)*(2*N-1), sizeof(*fft_in));
1040  SO3_ERROR_MEM_ALLOC_CHECK(fft_in);
1041  complex double *fft_out = calloc((2*L-1)*N, sizeof(*fft_out));
1042  SO3_ERROR_MEM_ALLOC_CHECK(fft_out);
1043  // Redundant dimension needs to be last
1044  fftw_plan plan = fftw_plan_dft_r2c_2d(
1045  2*L-1, 2*N-1,
1046  fft_in, fft_out,
1047  FFTW_ESTIMATE);
1048 
1049  int a, b, g;
1050  for (b = 0; b < L; ++b)
1051  {
1052  // TODO: This loop could probably be avoided by using
1053  // a more elaborate FFTW plan which performs the FFT directly
1054  // over the 1st and 3rd dimensions of f.
1055  // Instead, for each index in the 2nd dimension, we copy the
1056  // corresponding values in the 1st and 3rd dimension into a
1057  // new 2D array, to perform a standard 2D FFT there. While
1058  // we're at it, we also reshape that array such that gamma
1059  // is the inner dimension, as required by FFTW.
1060  for (a = 0; a < 2*L-1; ++a)
1061  for (g = 0; g < 2*N-1; ++g)
1062  fft_in[g + g_stride*(
1063  a)] =
1064  f[a + a_stride*(
1065  b + b_stride*(
1066  g))];
1067 
1068  fftw_execute(plan);
1069 
1070  // Apply spatial shift, while
1071  // reshaping the dimensions once more.
1072  for (n = n_start; n <= n_stop; n += n_inc)
1073  {
1074  for (m = -L+1; m <= L-1; ++m)
1075  {
1076  int m_shift = m < 0 ? 2*L-1 : 0;
1077  Fmnb[b + bext_stride*(
1078  m + m_offset + m_stride*(
1079  n + n_offset))] =
1080  fft_out[n + n_stride*(
1081  m + m_shift)];
1082  }
1083  }
1084  }
1085  fftw_destroy_plan(plan);
1086 
1087  // Extend Fmnb by filling it with zeroes.
1088  for (n = n_start; n <= n_stop; n += n_inc)
1089  for (m = -L+1; m <= L-1; ++m)
1090  for (b = L; b < 2*L-1; ++b)
1091  Fmnb[b + bext_stride*(
1092  m + m_offset + m_stride*(
1093  n + n_offset))] = 0.0;
1094 
1095  // Compute Fourier transform over beta, i.e. compute Fmnm'.
1096  complex double *Fmnm = calloc((2*L-1)*(2*L-1)*(2*N-1), sizeof(*Fmnm));
1097  SO3_ERROR_MEM_ALLOC_CHECK(Fmnm);
1098  complex double *inout = calloc(2*L-1, sizeof(*inout));
1099  SO3_ERROR_MEM_ALLOC_CHECK(inout);
1100 
1101  plan = fftw_plan_dft_1d(
1102  2*L-1,
1103  inout, inout,
1104  FFTW_FORWARD,
1105  FFTW_ESTIMATE);
1106  for (n = n_start; n <= n_stop; n += n_inc)
1107  for (m = -L+1; m <= L-1; ++m)
1108  {
1109  memcpy(inout,
1110  Fmnb + 0 + bext_stride*(
1111  m + m_offset + m_stride*(
1112  n + n_offset)),
1113  bext_stride*sizeof(*Fmnb));
1114  fftw_execute(plan);
1115 
1116  // Apply spatial shift
1117  for (mm = -L+1; mm <= L-1; ++mm)
1118  {
1119  int mm_shift = mm < 0 ? 2*L-1 : 0;
1120  Fmnm[m + m_offset + m_stride*(
1121  mm + mm_offset + mm_stride*(
1122  n + n_offset))] =
1123  inout[mm + mm_shift];
1124  }
1125  }
1126  fftw_destroy_plan(plan);
1127  free(inout);
1128 
1129  // Apply phase modulation to account for sampling offset.
1130  for (n = n_start; n <= n_stop; n += n_inc)
1131  for (mm = -L+1; mm <= L-1; ++mm)
1132  for (m = -L+1; m <= L-1; ++m)
1133  Fmnm[m + m_offset + m_stride*(
1134  mm + mm_offset + mm_stride*(
1135  n + n_offset))] *=
1136  expsmm[mm + mm_offset];
1137 
1138  // Compute flmn.
1139  double *dl, *dl8 = NULL;
1140  dl = ssht_dl_calloc(L, SSHT_DL_QUARTER);
1141  SO3_ERROR_MEM_ALLOC_CHECK(dl);
1142  if (dl_method == SSHT_DL_RISBO)
1143  {
1144  dl8 = ssht_dl_calloc(L, SSHT_DL_QUARTER_EXTENDED);
1145  SO3_ERROR_MEM_ALLOC_CHECK(dl8);
1146  }
1147  int dl_offset = ssht_dl_get_offset(L, SSHT_DL_QUARTER);
1148  int dl_stride = ssht_dl_get_stride(L, SSHT_DL_QUARTER);
1149  for (n = 0; n <= N-1; ++n)
1150  for (el = n; el < L; ++el)
1151  for (m = -el; m <= el; ++m)
1152  {
1153  int ind;
1154  so3_sampling_elmn2ind_real(&ind, el, m, n, parameters);
1155  flmn[ind] = 0.0;
1156  }
1157 
1158  for (el = L0; el < L; ++el)
1159  {
1160  int eltmp;
1161 
1162  // Compute Wigner plane.
1163  switch (dl_method)
1164  {
1165  case SSHT_DL_RISBO:
1166  if (el != 0 && el == L0)
1167  {
1168  for(eltmp = 0; eltmp <= L0; ++eltmp)
1169  ssht_dl_beta_risbo_eighth_table(dl8, SO3_PION2, L,
1170  SSHT_DL_QUARTER_EXTENDED,
1171  eltmp, sqrt_tbl, signs);
1172  ssht_dl_beta_risbo_fill_eighth2quarter_table(dl,
1173  dl8, L,
1174  SSHT_DL_QUARTER,
1175  SSHT_DL_QUARTER_EXTENDED,
1176  el,
1177  signs);
1178  }
1179  else
1180  {
1181  ssht_dl_beta_risbo_eighth_table(dl8, SO3_PION2, L,
1182  SSHT_DL_QUARTER_EXTENDED,
1183  el, sqrt_tbl, signs);
1184  ssht_dl_beta_risbo_fill_eighth2quarter_table(dl,
1185  dl8, L,
1186  SSHT_DL_QUARTER,
1187  SSHT_DL_QUARTER_EXTENDED,
1188  el,
1189  signs);
1190  }
1191  break;
1192 
1193  case SSHT_DL_TRAPANI:
1194  if (el != 0 && el == L0)
1195  {
1196  for(eltmp = 0; eltmp <= L0; ++eltmp)
1197  ssht_dl_halfpi_trapani_eighth_table(dl, L,
1198  SSHT_DL_QUARTER,
1199  eltmp, sqrt_tbl);
1200  ssht_dl_halfpi_trapani_fill_eighth2quarter_table(dl, L,
1201  SSHT_DL_QUARTER,
1202  el, signs);
1203  }
1204  else
1205  {
1206  ssht_dl_halfpi_trapani_eighth_table(dl, L,
1207  SSHT_DL_QUARTER,
1208  el, sqrt_tbl);
1209  ssht_dl_halfpi_trapani_fill_eighth2quarter_table(dl, L,
1210  SSHT_DL_QUARTER,
1211  el, signs);
1212  }
1213  break;
1214 
1215  default:
1216  SO3_ERROR_GENERIC("Invalid dl method");
1217  }
1218 
1219  // Compute flmn for current el.
1220 
1221 
1222  switch (n_mode)
1223  {
1224  case SO3_N_MODE_ALL:
1225  n_start = 0;
1226  n_stop = MIN( N-1, el);
1227  n_inc = 1;
1228  break;
1229  case SO3_N_MODE_EVEN:
1230  n_start = 0;
1231  n_stop = MIN( N-1, el);
1232  n_stop -= n_stop%2;
1233  n_inc = 2;
1234  break;
1235  case SO3_N_MODE_ODD:
1236  n_start = 1;
1237  n_stop = MIN( N-1, el);
1238  n_stop -= 1-n_stop%2;
1239  n_inc = 2;
1240  break;
1241  case SO3_N_MODE_MAXIMUM:
1242  if (el < N-1)
1243  continue;
1244  n_start = N-1;
1245  n_stop = N-1;
1246  n_inc = 1;
1247  break;
1248  case SO3_N_MODE_L:
1249  if (el >= N)
1250  continue;
1251  n_start = el;
1252  n_stop = el;
1253  n_inc = 1;
1254  break;
1255  default:
1256  SO3_ERROR_GENERIC("Invalid n-mode.");
1257  }
1258 
1259  // Factor which depends only on el.
1260  double elfactor = (2.0*el+1.0)/(8.0*SO3_PI*SO3_PI);
1261 
1262  // TODO: Pull out a few multiplications into precomputations
1263  // or split up loops to avoid conditionals to check signs.
1264  for (mm = -el; mm <= el; ++mm)
1265  {
1266  // These signs are needed for the symmetry relations of
1267  // Wigner symbols.
1268  double elmmsign = signs[el] * signs[abs(mm)];
1269 
1270  for (n = n_start; n <= n_stop; n += n_inc)
1271  {
1272  double mmsign = mm >= 0 ? 1.0 : signs[el] * signs[n];
1273 
1274  // Factor which does not depend on m.
1275  double elnmm_factor = mmsign * elfactor
1276  * dl[n + dl_offset + abs(mm)*dl_stride];
1277 
1278  for (m = -el; m <= el; ++m)
1279  {
1280  mmsign = mm >= 0 ? 1.0 : signs[el] * signs[abs(m)];
1281  double elmsign = m >= 0 ? 1.0 : elmmsign;
1282  int ind;
1283  so3_sampling_elmn2ind_real(&ind, el, m, n, parameters);
1284  int mod = ((m-n)%4 + 4)%4;
1285  flmn[ind] +=
1286  exps[mod]
1287  * elnmm_factor
1288  * mmsign * elmsign
1289  * dl[abs(m) + dl_offset + abs(mm)*dl_stride]
1290  * Fmnm[m + m_offset + m_stride*(
1291  mm + mm_offset + mm_stride*(
1292  n + n_offset))];
1293 
1294  }
1295  }
1296  }
1297  }
1298 
1299  free(dl);
1300  if (dl_method == SSHT_DL_RISBO)
1301  free(dl8);
1302  free(Fmnb);
1303  free(Fmnm);
1304  free(sqrt_tbl);
1305  free(signs);
1306  free(exps);
1307  free(expsmm);
1308 
1309  if (verbosity > 0)
1310  printf("%sAdjoint inverse transform computed!\n", SO3_PROMPT);
1311 }
1312 
1314  double *f, const complex double *flmn,
1315  const so3_parameters_t *parameters
1316 ) {
1317  int L0, L, N;
1318  so3_sampling_t sampling;
1319  so3_storage_t storage;
1320  so3_n_mode_t n_mode;
1321  ssht_dl_method_t dl_method;
1322  int steerable;
1323  int verbosity;
1324 
1325  L0 = parameters->L0;
1326  L = parameters->L;
1327  N = parameters->N;
1328  sampling = parameters->sampling_scheme;
1329  storage = parameters->storage;
1330  // TODO: Add optimisations for all n-modes.
1331  n_mode = parameters->n_mode;
1332  dl_method = parameters->dl_method;
1333  verbosity = parameters->verbosity;
1334  steerable = parameters->steerable;
1335 
1336  // Print messages depending on verbosity level.
1337  if (verbosity > 0)
1338  {
1339  printf("%sComputing forward adjoint transform using MW sampling with\n", SO3_PROMPT);
1340  printf("%sparameters (L, N, reality) = (%d, %d, FALSE)\n", SO3_PROMPT, L, N);
1341  if (verbosity > 1)
1342  printf("%sUsing routine so3_adjoint_forward_direct with storage method %d...\n"
1343  , SO3_PROMPT
1344  , storage);
1345  }
1346 
1347  // Iterators
1348  int el, m, n, mm, a, b, g; // mm for m'
1349 
1350  int m_stride = 2*L-1;
1351  int m_offset = L-1;
1352  int n_offset = 0;
1353  int n_stride = N;
1354  int mm_stride = 2*L-1;
1355  int mm_offset = L-1;
1356  int a_stride = 2*L-1;
1357  int b_stride = L;
1358  int bext_stride = 2*L-1;
1359  int g_stride = 2*N-1;
1360 
1361  complex double* inout; // Used as temporary storage for various FFTs.
1362 
1363  // Allocate memory.
1364  double *sqrt_tbl = calloc(2*(L-1)+2, sizeof(*sqrt_tbl));
1365  SO3_ERROR_MEM_ALLOC_CHECK(sqrt_tbl);
1366  double *signs = calloc(L+1, sizeof(*signs));
1367  SO3_ERROR_MEM_ALLOC_CHECK(signs);
1368  complex double *exps = calloc(4, sizeof(*exps));
1369  SO3_ERROR_MEM_ALLOC_CHECK(exps);
1370  complex double *expsmm = calloc(2*L-1, sizeof(*expsmm));
1371  SO3_ERROR_MEM_ALLOC_CHECK(expsmm);
1372 
1373  // Perform precomputations.
1374  for (el = 0; el <= 2*L-1; ++el)
1375  sqrt_tbl[el] = sqrt((double)el);
1376  for (m = 0; m <= L-1; m += 2)
1377  {
1378  signs[m] = 1.0;
1379  signs[m+1] = -1.0;
1380  }
1381  int i;
1382  for (i = 0; i < 4; ++i)
1383  exps[i] = cexp(I*SO3_PION2*i);
1384  for (mm = -L+1; mm <= L-1; ++mm)
1385  expsmm[mm + mm_offset] = cexp(I*mm*SSHT_PI/(2.0*L-1.0));
1386 
1387  // Compute Gmnm'
1388  // TODO: Currently m is fastest-varying, then n, then m'.
1389  // Should this order be changed to m-m'-n?
1390  complex double *Gmnm = calloc((2*L-1)*(2*L-1)*N, sizeof(*Gmnm));
1391  SO3_ERROR_MEM_ALLOC_CHECK(Gmnm);
1392 
1393  int n_start, n_stop, n_inc;
1394 
1395 
1396  double *dl = ssht_dl_calloc(L, SSHT_DL_QUARTER);
1397  SO3_ERROR_MEM_ALLOC_CHECK(dl);
1398  double *dl8 = NULL;
1399  if (dl_method == SSHT_DL_RISBO)
1400  {
1401  dl8 = ssht_dl_calloc(L, SSHT_DL_QUARTER_EXTENDED);
1402  SO3_ERROR_MEM_ALLOC_CHECK(dl8);
1403  }
1404  int dl_offset = ssht_dl_get_offset(L, SSHT_DL_QUARTER);
1405  int dl_stride = ssht_dl_get_stride(L, SSHT_DL_QUARTER);
1406 
1407  complex double *mn_factors = calloc((2*L-1)*(2*N-1), sizeof *mn_factors);
1408  SO3_ERROR_MEM_ALLOC_CHECK(mn_factors);
1409 
1410  // TODO: SSHT starts this loop from MAX(L0, abs(spin)).
1411  // Can we use a similar optimisation? el can probably
1412  // be limited by n, but then we'd need to switch the
1413  // loop order, which means we'd have to recompute the
1414  // Wigner plane for each n. That seems wrong?
1415  for (el = L0; el <= L-1; ++el)
1416  {
1417  int eltmp;
1418  // Compute Wigner plane.
1419  switch (dl_method)
1420  {
1421  case SSHT_DL_RISBO:
1422  if (el != 0 && el == L0)
1423  {
1424  for(eltmp = 0; eltmp <= L0; ++eltmp)
1425  ssht_dl_beta_risbo_eighth_table(
1426  dl8,
1427  SO3_PION2,
1428  L,
1429  SSHT_DL_QUARTER_EXTENDED,
1430  eltmp,
1431  sqrt_tbl,
1432  signs);
1433  ssht_dl_beta_risbo_fill_eighth2quarter_table(dl,
1434  dl8, L,
1435  SSHT_DL_QUARTER,
1436  SSHT_DL_QUARTER_EXTENDED,
1437  el,
1438  signs);
1439  }
1440  else
1441  {
1442  ssht_dl_beta_risbo_eighth_table(
1443  dl8,
1444  SO3_PION2,
1445  L,
1446  SSHT_DL_QUARTER_EXTENDED,
1447  el,
1448  sqrt_tbl,
1449  signs);
1450  ssht_dl_beta_risbo_fill_eighth2quarter_table(
1451  dl,
1452  dl8, L,
1453  SSHT_DL_QUARTER,
1454  SSHT_DL_QUARTER_EXTENDED,
1455  el,
1456  signs);
1457  }
1458  break;
1459 
1460  case SSHT_DL_TRAPANI:
1461  if (el != 0 && el == L0)
1462  {
1463  for(eltmp = 0; eltmp <= L0; ++eltmp)
1464  ssht_dl_halfpi_trapani_eighth_table(
1465  dl,
1466  L,
1467  SSHT_DL_QUARTER,
1468  eltmp,
1469  sqrt_tbl);
1470  ssht_dl_halfpi_trapani_fill_eighth2quarter_table(
1471  dl,
1472  L,
1473  SSHT_DL_QUARTER,
1474  el,
1475  signs);
1476  }
1477  else
1478  {
1479  ssht_dl_halfpi_trapani_eighth_table(
1480  dl,
1481  L,
1482  SSHT_DL_QUARTER,
1483  el,
1484  sqrt_tbl);
1485  ssht_dl_halfpi_trapani_fill_eighth2quarter_table(
1486  dl,
1487  L,
1488  SSHT_DL_QUARTER,
1489  el,
1490  signs);
1491  }
1492  break;
1493 
1494  default:
1495  SO3_ERROR_GENERIC("Invalid dl method");
1496  }
1497 
1498  // Compute Gmnm' contribution for current el.
1499 
1500  switch (n_mode)
1501  {
1502  case SO3_N_MODE_ALL:
1503  n_start = 0;
1504  n_stop = MIN( N-1, el);
1505  n_inc = 1;
1506  break;
1507  case SO3_N_MODE_EVEN:
1508  n_start = 0;
1509  n_stop = MIN( N-1, el);
1510  n_stop -= n_stop%2;
1511  n_inc = 2;
1512  break;
1513  case SO3_N_MODE_ODD:
1514  n_start = 1;
1515  n_stop = MIN( N-1, el);
1516  n_stop -= 1-n_stop%2;
1517  n_inc = 2;
1518  break;
1519  case SO3_N_MODE_MAXIMUM:
1520  if (el < N-1)
1521  continue;
1522  n_start = N-1;
1523  n_stop = N-1;
1524  n_inc = 1;
1525  break;
1526  case SO3_N_MODE_L:
1527  if (el >= N)
1528  continue;
1529  n_start = el;
1530  n_stop = el;
1531  n_inc = 1;
1532  break;
1533  default:
1534  SO3_ERROR_GENERIC("Invalid n-mode.");
1535  }
1536 
1537  // Factors which do not depend on m'.
1538  for (n = n_start; n <= n_stop; n += n_inc)
1539  for (m = -el; m <= el; ++m)
1540  {
1541  int ind;
1542  so3_sampling_elmn2ind_real(&ind, el, m, n, parameters);
1543  int mod = ((n-m)%4 + 4)%4;
1544  mn_factors[m + m_offset + m_stride*(
1545  n + n_offset)] =
1546  flmn[ind] * exps[mod];
1547  }
1548 
1549  for (mm = 0; mm <= el; ++mm)
1550  {
1551  // These signs are needed for the symmetry relations of
1552  // Wigner symbols.
1553  double elmmsign = signs[el] * signs[mm];
1554 
1555  for (n = n_start; n <= n_stop; n += n_inc)
1556  {
1557  // Factor which does not depend on m.
1558  double elnmm_factor = dl[n + dl_offset + mm*dl_stride];
1559  for (m = -el; m < 0; ++m)
1560  Gmnm[mm + mm_offset + mm_stride*(
1561  m + m_offset + m_stride*(
1562  n + n_offset))] +=
1563  elnmm_factor
1564  * mn_factors[m + m_offset + m_stride*(
1565  n + n_offset)]
1566  * elmmsign * dl[-m + dl_offset + mm*dl_stride];
1567  for (m = 0; m <= el; ++m)
1568  Gmnm[mm + mm_offset + mm_stride*(
1569  m + m_offset + m_stride*(
1570  n + n_offset))] +=
1571  elnmm_factor
1572  * mn_factors[m + m_offset + m_stride*(
1573  n + n_offset)]
1574  * dl[m + dl_offset + mm*dl_stride];
1575  }
1576  }
1577  }
1578 
1579  // Free dl memory.
1580  free(dl);
1581  if (dl_method == SSHT_DL_RISBO)
1582  free(dl8);
1583  free(mn_factors);
1584 
1585  switch (n_mode)
1586  {
1587  case SO3_N_MODE_ALL:
1588  case SO3_N_MODE_L:
1589  n_start = 0;
1590  n_stop = N-1;
1591  n_inc = 1;
1592  break;
1593  case SO3_N_MODE_EVEN:
1594  n_start = 0;
1595  n_stop = ((N-1) % 2 == 0) ? N-1 : N-2;
1596  n_inc = 2;
1597  break;
1598  case SO3_N_MODE_ODD:
1599  n_start = 1;
1600  n_stop = ((N-1) % 2 != 0) ? N-1 : N-2;
1601  n_inc = 2;
1602  break;
1603  case SO3_N_MODE_MAXIMUM:
1604  n_start = N-1;
1605  n_stop = N-1;
1606  n_inc = 1;
1607  break;
1608  default:
1609  SO3_ERROR_GENERIC("Invalid n-mode.");
1610  }
1611 
1612  // Use symmetry to compute Gmnm' for negative m'.
1613  for (n = n_start; n <= n_stop; n += n_inc)
1614  for (m = -L+1; m <= L-1; ++m)
1615  for (mm = -L+1; mm < 0; ++mm)
1616  Gmnm[mm + mm_offset + mm_stride*(
1617  m + m_offset + m_stride*(
1618  n + n_offset))] =
1619  signs[abs(m+n)%2]
1620  * Gmnm[-mm + mm_offset + mm_stride*(
1621  m + m_offset + m_stride*(
1622  n + n_offset))];
1623 
1624 
1625  // Compute weights.
1626  complex double *w = calloc(4*L-3, sizeof(*w));
1627  SO3_ERROR_MEM_ALLOC_CHECK(w);
1628  int w_offset = 2*(L-1);
1629  for (mm = -2*(L-1); mm <= 2*(L-1); ++mm)
1630  w[mm+w_offset] = so3_sampling_weight(parameters, mm);
1631 
1632  // Compute IFFT of w to give wr.
1633  complex double *wr = calloc(4*L-3, sizeof(*w));
1634  SO3_ERROR_MEM_ALLOC_CHECK(wr);
1635  inout = calloc(4*L-3, sizeof(*inout));
1636  SO3_ERROR_MEM_ALLOC_CHECK(inout);
1637  fftw_plan plan_bwd = fftw_plan_dft_1d(
1638  4*L-3,
1639  inout, inout,
1640  FFTW_BACKWARD,
1641  FFTW_MEASURE);
1642  fftw_plan plan_fwd = fftw_plan_dft_1d(
1643  4*L-3,
1644  inout,
1645  inout,
1646  FFTW_FORWARD,
1647  FFTW_MEASURE);
1648 
1649  // Apply spatial shift.
1650  for (mm = 1; mm <= 2*L-2; ++mm)
1651  inout[mm + w_offset] = w[mm - 2*(L-1) - 1 + w_offset];
1652  for (mm = -2*(L-1); mm <= 0; ++mm)
1653  inout[mm + w_offset] = w[mm + 2*(L-1) + w_offset];
1654 
1655  fftw_execute_dft(plan_bwd, inout, inout);
1656 
1657  // Apply spatial shift.
1658  for (mm = 0; mm <= 2*L-2; ++mm)
1659  wr[mm + w_offset] = inout[mm - 2*(L-1) + w_offset];
1660  for (mm = -2*(L-1); mm <= -1; ++mm)
1661  wr[mm + w_offset] = inout[mm + 2*(L-1) + 1 + w_offset];
1662 
1663  // Compute Fmnm'' by convolution implemented as product in real space.
1664  complex double *Gmnm_pad = calloc(4*L-3, sizeof(*Gmnm_pad));
1665  SO3_ERROR_MEM_ALLOC_CHECK(Gmnm_pad);
1666  complex double *Fmnm = calloc((2*L-1)*(2*L-1)*N, sizeof(*Fmnm));
1667  SO3_ERROR_MEM_ALLOC_CHECK(Fmnm);
1668  for (n = n_start; n <= n_stop; n += n_inc)
1669  for (m = -L+1; m <= L-1; ++m)
1670  {
1671 
1672  // Zero-pad Gmnm'.
1673  for (mm = -2*(L-1); mm <= -L; ++mm)
1674  Gmnm_pad[mm+w_offset] = 0.0;
1675  for (mm = L; mm <= 2*(L-1); ++mm)
1676  Gmnm_pad[mm+w_offset] = 0.0;
1677  for (mm = -(L-1); mm <= L-1; ++mm)
1678  Gmnm_pad[mm + w_offset] =
1679  Gmnm[mm + mm_offset + mm_stride*(
1680  m + m_offset + m_stride*(
1681  n + n_offset))];
1682 
1683  // Apply spatial shift.
1684  for (mm = 1; mm <= 2*L-2; ++mm)
1685  inout[mm + w_offset] = Gmnm_pad[mm - 2*(L-1) - 1 + w_offset];
1686  for (mm = -2*(L-1); mm <= 0; ++mm)
1687  inout[mm + w_offset] = Gmnm_pad[mm + 2*(L-1) + w_offset];
1688  // Compute IFFT of Gmnm'.
1689  fftw_execute_dft(plan_bwd, inout, inout);
1690  // Apply spatial shift.
1691  for (mm = 0; mm <= 2*L-2; ++mm)
1692  Gmnm_pad[mm + w_offset] = inout[mm - 2*(L-1) + w_offset];
1693  for (mm = -2*(L-1); mm <= -1; ++mm)
1694  Gmnm_pad[mm + w_offset] = inout[mm + 2*(L-1) + 1 + w_offset];
1695 
1696  // Compute product of Gmnm' and weight in real space.
1697  int r;
1698  for (r = -2*(L-1); r <= 2*(L-1); ++r)
1699  Gmnm_pad[r + w_offset] *= wr[r + w_offset];
1700 
1701  // Apply spatial shift.
1702  for (mm = 1; mm <= 2*L-2; ++mm)
1703  inout[mm + w_offset] = Gmnm_pad[mm - 2*(L-1) - 1 + w_offset];
1704  for (mm = -2*(L-1); mm <= 0; ++mm)
1705  inout[mm + w_offset] = Gmnm_pad[mm + 2*(L-1) + w_offset];
1706  // Compute Fmnm'' by FFT.
1707  fftw_execute_dft(plan_fwd, inout, inout);
1708  // Apply spatial shift.
1709  for (mm = 0; mm <= 2*L-2; ++mm)
1710  Gmnm_pad[mm + w_offset] = inout[mm - 2*(L-1) + w_offset];
1711  for (mm = -2*(L-1); mm <= -1; ++mm)
1712  Gmnm_pad[mm + w_offset] = inout[mm + 2*(L-1) + 1 + w_offset];
1713 
1714  // Extract section of Fmnm'' of interest.
1715  for (mm = -(L-1); mm <= L-1; ++mm)
1716  Fmnm[mm + mm_offset + mm_stride*(
1717  m + m_offset + m_stride*(
1718  n + n_offset))] =
1719  Gmnm_pad[mm + w_offset]
1720  * 4.0 * SSHT_PI * SSHT_PI / (4.0*L-3.0);
1721 
1722  }
1723  fftw_destroy_plan(plan_bwd);
1724  fftw_destroy_plan(plan_fwd);
1725 
1726  // Apply phase modulation to account for sampling offset.
1727  for (n = n_start; n <= n_stop; n += n_inc)
1728  for (m = -L+1; m <= L-1; ++m)
1729  for (mm = -L+1; mm <= L-1; ++mm)
1730  Fmnm[mm + mm_offset + mm_stride*(
1731  m + m_offset + m_stride*(
1732  n + n_offset))] *=
1733  expsmm[mm + mm_offset];
1734 
1735  // Compute Fourier transform over mm, i.e. compute Fmnb.
1736  complex double *Fmnb = calloc((2*L-1)*(2*L-1)*N, sizeof(*Fmnb));
1737  SO3_ERROR_MEM_ALLOC_CHECK(Fmnm);
1738 
1739  fftw_plan plan = fftw_plan_dft_1d(
1740  2*L-1,
1741  inout, inout,
1742  FFTW_BACKWARD,
1743  FFTW_ESTIMATE);
1744  for (n = n_start; n <= n_stop; n += n_inc)
1745  for (m = -L+1; m <= L-1; ++m)
1746  {
1747  // Apply spatial shift and normalisation factor
1748  for (mm = -L+1; mm <= L-1; ++mm)
1749  {
1750  int mm_shift = mm < 0 ? 2*L-1 : 0;
1751  inout[mm + mm_shift] =
1752  Fmnm[mm + mm_offset + mm_stride*(
1753  m + m_offset + m_stride*(
1754  n + n_offset))] / (2.0*L-1.0);
1755  }
1756  fftw_execute_dft(plan, inout, inout);
1757  memcpy(Fmnb + 0 + bext_stride*(
1758  m + m_offset + m_stride*(
1759  n + n_offset)),
1760  inout,
1761  bext_stride*sizeof(*Fmnb));
1762 
1763  }
1764 
1765  fftw_destroy_plan(plan);
1766  free(inout);
1767 
1768 
1769 
1770  // Adjoint of periodic extension of Ftm.
1771  for (n = n_start; n <= n_stop; n += n_inc)
1772  for (m = -L+1; m <= L-1; ++m)
1773  {
1774  int signmn = signs[abs(m+n)%2];
1775  for (b = 0; b <= L-2; ++b)
1776  Fmnb[b + bext_stride*(
1777  m + m_offset + m_stride*(
1778  n + n_offset))] +=
1779  signmn
1780  * Fmnb[(2*L-2-b) + bext_stride*(
1781  m + m_offset + m_stride*(
1782  n + n_offset))];
1783  }
1784 
1785  // Compute Fourier transform over alpha and gamma, i.e. compute f.
1786  complex double *fft_in = calloc((2*L-1)*N, sizeof(*fft_in));
1787  SO3_ERROR_MEM_ALLOC_CHECK(fft_in);
1788  double *fft_out = calloc((2*L-1)*(2*N-1), sizeof(*fft_out));
1789  SO3_ERROR_MEM_ALLOC_CHECK(fft_out);
1790  // Redundant dimension needs to be last
1791  plan = fftw_plan_dft_c2r_2d(
1792  2*L-1, 2*N-1,
1793  fft_in, fft_out,
1794  FFTW_ESTIMATE);
1795 
1796  double norm_factor = 1.0/(2.0*L-1.0)/(2.0*N-1.0);
1797 
1798  for (b = 0; b < L; ++b)
1799  {
1800  // Apply spatial shift and normalisation factor
1801  for (n = n_start; n <= n_stop; n += n_inc)
1802  {
1803  for (m = -L+1; m <= L-1; ++m)
1804  {
1805  int m_shift = m < 0 ? 2*L-1 : 0;
1806  fft_in[n + n_stride*(
1807  m + m_shift)] =
1808  Fmnb[b + bext_stride*(
1809  m + m_offset + m_stride*(
1810  n + n_offset))] * norm_factor;
1811  }
1812  }
1813 
1814  fftw_execute(plan);
1815 
1816  // TODO: This loop could probably be avoided by using
1817  // a more elaborate FFTW plan which performs the FFT directly
1818  // over the 1st and 3rd dimensions of f.
1819  for (a = 0; a < 2*L-1; ++a)
1820  for (g = 0; g < 2*N-1; ++g)
1821  f[a + a_stride*(
1822  b + b_stride*(
1823  g))] = fft_out[g + g_stride*(
1824  a)];
1825 
1826  }
1827  fftw_destroy_plan(plan);
1828 
1829  free(Fmnb);
1830  free(Fmnm);
1831  free(fft_in);
1832  free(fft_out);
1833  free(w);
1834  free(wr);
1835  free(Gmnm_pad);
1836  free(Gmnm);
1837  free(sqrt_tbl);
1838  free(signs);
1839  free(exps);
1840  free(expsmm);
1841 }
so3_adjoint_forward_direct
void so3_adjoint_forward_direct(complex double *f, const complex double *flmn, const so3_parameters_t *parameters)
Definition: so3_adjoint.c:402
MAX
#define MAX(a, b)
Definition: so3_adjoint.c:27
inverse_real_ssht
void(* inverse_real_ssht)(double *, const complex double *, int, int, ssht_dl_method_t, int)
Definition: so3_adjoint.c:30
so3_sampling_elmn2ind_real
void so3_sampling_elmn2ind_real(int *ind, int el, int m, int n, const so3_parameters_t *parameters)
Definition: so3_sampling.c:569
MIN
#define MIN(a, b)
Definition: so3_adjoint.c:26
inverse_complex_ssht
void(* inverse_complex_ssht)(complex double *, const complex double *, int, int, int, ssht_dl_method_t, int)
Definition: so3_adjoint.c:29
forward_real_ssht
void(* forward_real_ssht)(complex double *, const double *, int, int, ssht_dl_method_t, int)
Definition: so3_adjoint.c:32
so3_adjoint_inverse_direct_real
void so3_adjoint_inverse_direct_real(complex double *flmn, const double *f, const so3_parameters_t *parameters)
Definition: so3_adjoint.c:939
forward_complex_ssht
void(* forward_complex_ssht)(complex double *, const complex double *, int, int, int, ssht_dl_method_t, int)
Definition: so3_adjoint.c:31
so3_sampling_weight
complex double so3_sampling_weight(const so3_parameters_t *parameters, int p)
Definition: so3_sampling.c:30
so3_adjoint_forward_direct_real
void so3_adjoint_forward_direct_real(double *f, const complex double *flmn, const so3_parameters_t *parameters)
Definition: so3_adjoint.c:1313
so3_sampling_elmn2ind
void so3_sampling_elmn2ind(int *ind, int el, int m, int n, const so3_parameters_t *parameters)
Definition: so3_sampling.c:367
so3_adjoint_inverse_direct
void so3_adjoint_inverse_direct(complex double *flmn, const complex double *f, const so3_parameters_t *parameters)
Definition: so3_adjoint.c:34