My Project
so3_core.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 
49  complex double *f, const complex double *flmn,
50  const so3_parameters_t *parameters
51 ) {
52 
53  int L0, L, N;
54  so3_sampling_t sampling;
55  so3_storage_t storage;
56  so3_n_mode_t n_mode;
57  ssht_dl_method_t dl_method;
58  int steerable;
59  int verbosity;
60 
61  // Iterator
62  int n;
63  // Intermediate results
64  complex double *fn, *ftemp;
65  // Stride for several arrays
66  int fn_n_stride;
67  // FFTW-related variables
68  int fftw_rank, fftw_howmany;
69  int fftw_idist, fftw_odist;
70  int fftw_istride, fftw_ostride;
71  int fftw_n;
72  complex double *fftw_target;
73  fftw_plan plan;
74 
76 
77  L0 = parameters->L0;
78  L = parameters->L;
79  N = parameters->N;
80  sampling = parameters->sampling_scheme;
81  storage = parameters->storage;
82  // TODO: Support N_MODE_L
83  n_mode = parameters->n_mode;
84  dl_method = parameters->dl_method;
85  verbosity = parameters->verbosity;
86  steerable = parameters->steerable;
87 
88  // Print messages depending on verbosity level.
89  if (verbosity > 0) {
90  printf("%sComputing inverse transform using MW sampling with\n", SO3_PROMPT);
91  printf("%sparameters (L, N, reality) = (%d, %d, FALSE)\n", SO3_PROMPT, L, N);
92  if (verbosity > 1)
93  printf("%sUsing routine so3_core_mw_inverse_via_ssht with storage method %d...\n"
94  , SO3_PROMPT
95  , storage);
96  }
97 
98  switch (sampling)
99  {
100  case SO3_SAMPLING_MW:
101  fn_n_stride = L * (2*L-1);
102  ssht = ssht_core_mw_lb_inverse_sov_sym;
103  break;
104  case SO3_SAMPLING_MW_SS:
105  fn_n_stride = (L+1) * 2*L;
106  ssht = ssht_core_mw_lb_inverse_sov_sym_ss;
107  break;
108  default:
109  SO3_ERROR_GENERIC("Invalid sampling scheme.");
110  }
111 
112  // Compute fn(a,b)
113 
114  if (steerable)
115  {
116  // For steerable signals, we need to supersample in n/gamma,
117  // in order to create a symmetric sampling.
118  fftw_n = 2*N; // Each transform is over 2*N
119 
120  // We need to perform the FFT into a temporary buffer, because
121  // the result will be twice as large as the output we need.
122  ftemp = malloc(2*N*fn_n_stride * sizeof *ftemp);
123  SO3_ERROR_MEM_ALLOC_CHECK(ftemp);
124 
125  fftw_target = ftemp;
126  }
127  else
128  {
129  fftw_n = 2*N-1; // Each transform is over 2*N-1
130 
131  fftw_target = f;
132  }
133 
134  fn = calloc(fftw_n*fn_n_stride, sizeof *fn);
135  SO3_ERROR_MEM_ALLOC_CHECK(fn);
136 
137 
138 
139  // Initialize fftw_plan first. With FFTW_ESTIMATE this is technically not
140  // necessary but still good practice.
141  fftw_rank = 1; // We compute 1d transforms
142  fftw_howmany = fn_n_stride; // We need L*(2*L-1) of these transforms
143 
144  // We want to transform columns
145  fftw_idist = fftw_odist = 1; // The starts of the columns are contiguous in memory
146  fftw_istride = fftw_ostride = fn_n_stride; // Distance between two elements of the same column
147 
148  plan = fftw_plan_many_dft(
149  fftw_rank, &fftw_n, fftw_howmany,
150  fn, NULL, fftw_istride, fftw_idist,
151  fftw_target, NULL, fftw_ostride, fftw_odist,
152  FFTW_BACKWARD, FFTW_ESTIMATE
153  );
154 
155  for(n = -N+1; n <= N-1; ++n)
156  {
157  int ind, offset, i, el;
158  int L0e = MAX(L0, abs(n)); // 'e' for 'effective'
159  double factor;
160  complex double *flm;
161 
162  if ((n_mode == SO3_N_MODE_EVEN && n % 2)
163  || (n_mode == SO3_N_MODE_ODD && !(n % 2))
164  || (n_mode == SO3_N_MODE_MAXIMUM && abs(n) < N-1)
165  ) {
166  continue;
167  }
168 
169  flm = malloc(L*L * sizeof *flm);
170  SO3_ERROR_MEM_ALLOC_CHECK(flm);
171 
172  switch (storage)
173  {
174  case SO3_STORAGE_PADDED:
175  so3_sampling_elmn2ind(&ind, 0, 0, n, parameters);
176  memcpy(flm, flmn + ind, L*L * sizeof(complex double));
177  break;
178  case SO3_STORAGE_COMPACT:
179  so3_sampling_elmn2ind(&ind, abs(n), -abs(n), n, parameters);
180  memcpy(flm + n*n, flmn + ind, (L*L - n*n) * sizeof(complex double));
181  for(i = 0; i < n*n; ++i)
182  flm[i] = 0.0;
183  break;
184  default:
185  SO3_ERROR_GENERIC("Invalid storage method.");
186  }
187 
188  el = L0e;
189  i = offset = el*el;
190  for(; el < L; ++el)
191  {
192  factor = sqrt((double)(2*el+1)/(16.*pow(SO3_PI, 3.)));
193  for (; i < offset + 2*el+1; ++i)
194  flm[i] *= factor;
195 
196  offset = i;
197  }
198 
199  // The conditional applies the spatial transform, so that we store
200  // the results in n-order 0, 1, 2, -2, -1
201  offset = (n < 0 ? n + fftw_n : n);
202 
203  (*ssht)(
204  fn + offset*fn_n_stride, flm,
205  L0e, L, -n,
206  dl_method,
207  verbosity
208  );
209 
210  if(n % 2)
211  for(i = 0; i < fn_n_stride; ++i)
212  fn[offset*fn_n_stride + i] = -fn[offset*fn_n_stride + i];
213 
214  if (verbosity > 0)
215  printf("\n");
216 
217  free(flm);
218  }
219 
220 
221  fftw_execute(plan);
222  fftw_destroy_plan(plan);
223 
224  if (steerable)
225  {
226  memcpy(f, ftemp, N*fn_n_stride * sizeof(complex double));
227  free(ftemp);
228  }
229 
230  free(fn);
231 
232  if (verbosity > 0)
233  printf("%sInverse transform computed!\n", SO3_PROMPT);
234 }
235 
236 
254  complex double *flmn, const complex double *f,
255  const so3_parameters_t *parameters
256 ) {
257  int L0, L, N;
258  so3_sampling_t sampling;
259  so3_storage_t storage;
260  so3_n_mode_t n_mode;
261  ssht_dl_method_t dl_method;
262  int steerable;
263  int verbosity;
264 
265  // Iterator
266  int i, n;
267  // Intermediate results
268  complex double *ftemp, *fn;
269  // Stride for several arrays
270  int fn_n_stride;
271  // FFTW-related variables
272  int fftw_rank, fftw_howmany;
273  int fftw_idist, fftw_odist;
274  int fftw_istride, fftw_ostride;
275  int fftw_n;
276  fftw_plan plan;
277 
279 
280  // for precomputation
281  double factor;
282 
283  L0 = parameters->L0;
284  L = parameters->L;
285  N = parameters->N;
286  sampling = parameters->sampling_scheme;
287  storage = parameters->storage;
288  // TODO: Support N_MODE_L
289  n_mode = parameters->n_mode;
290  dl_method = parameters->dl_method;
291  verbosity = parameters->verbosity;
292  steerable = parameters->steerable;
293 
294  // Print messages depending on verbosity level.
295  if (verbosity > 0) {
296  printf("%sComputing forward transform using MW sampling with\n", SO3_PROMPT);
297  printf("%sparameters (L, N, reality) = (%d, %d, FALSE)\n", SO3_PROMPT, L, N);
298  if (verbosity > 1)
299  printf("%sUsing routine so3_core_mw_forward_via_ssht with storage method %d...\n"
300  , SO3_PROMPT
301  , storage);
302  }
303 
304  switch (sampling)
305  {
306  case SO3_SAMPLING_MW:
307  fn_n_stride = L * (2*L-1);
308  ssht = ssht_core_mw_lb_forward_sov_conv_sym;
309  break;
310  case SO3_SAMPLING_MW_SS:
311  fn_n_stride = (L+1) * 2*L;
312  ssht = ssht_core_mw_lb_forward_sov_conv_sym_ss;
313  break;
314  default:
315  SO3_ERROR_GENERIC("Invalid sampling scheme.");
316  }
317 
318  if (steerable)
319  {
320  int g, offset;
321 
322  fn = calloc((2*N-1)*fn_n_stride, sizeof *fn);
323  SO3_ERROR_MEM_ALLOC_CHECK(fn);
324 
325  for (n = -N+1; n < N; n+=2)
326  {
327  // The conditional applies the spatial transform, because the fn
328  // are to be stored in n-order 0, 1, 2, -2, -1
329  offset = (n < 0 ? n + 2*N-1 : n);
330 
331  for (g = 0; g < N; ++g)
332  {
333  double gamma = g * SO3_PI / N;
334  for (i = 0; i < fn_n_stride; ++i)
335  {
336  double weight = 2*SO3_PI/N;
337  fn[offset * fn_n_stride + i] += weight*f[g * fn_n_stride + i]*cexp(-I*n*gamma);
338  }
339  }
340  }
341  }
342  else
343  {
344  // Make a copy of the input, because input is const
345  // This could potentially be avoided by copying the input into fn and using an
346  // in-place FFTW. The performance impact has to be profiled, though.
347  ftemp = malloc((2*N-1)*fn_n_stride * sizeof *ftemp);
348  SO3_ERROR_MEM_ALLOC_CHECK(ftemp);
349  memcpy(ftemp, f, (2*N-1)*fn_n_stride * sizeof(complex double));
350 
351  fn = malloc((2*N-1)*fn_n_stride * sizeof *fn);
352  SO3_ERROR_MEM_ALLOC_CHECK(fn);
353 
354  // Initialize fftw_plan first. With FFTW_ESTIMATE this is technically not
355  // necessary but still good practice.
356  fftw_rank = 1;
357  fftw_n = 2*N-1;
358  fftw_howmany = fn_n_stride;
359  fftw_idist = fftw_odist = 1;
360  fftw_istride = fftw_ostride = fn_n_stride;
361 
362  plan = fftw_plan_many_dft(
363  fftw_rank, &fftw_n, fftw_howmany,
364  ftemp, NULL, fftw_istride, fftw_idist,
365  fn, NULL, fftw_ostride, fftw_odist,
366  FFTW_FORWARD, FFTW_ESTIMATE
367  );
368 
369  fftw_execute(plan);
370  fftw_destroy_plan(plan);
371 
372  free(ftemp);
373 
374  factor = 2*SO3_PI/(double)(2*N-1);
375  for(i = 0; i < (2*N-1)*fn_n_stride; ++i)
376  fn[i] *= factor;
377  }
378 
379  for(n = -N+1; n <= N-1; ++n)
380  {
381  int ind, offset, el, sign;
382  int L0e = MAX(L0, abs(n)); // 'e' for 'effective'
383 
384  complex double *flm = NULL;
385 
386  if ((n_mode == SO3_N_MODE_EVEN && n % 2)
387  || (n_mode == SO3_N_MODE_ODD && !(n % 2))
388  || (n_mode == SO3_N_MODE_MAXIMUM && abs(n) < N-1)
389  ) {
390  continue;
391  }
392 
393  if (storage == SO3_STORAGE_COMPACT)
394  flm = malloc(L*L * sizeof *flm);
395 
396  // The conditional applies the spatial transform, because the fn
397  // are stored in n-order 0, 1, 2, -2, -1
398  offset = (n < 0 ? n + 2*N-1 : n);
399 
400  complex double *flm_block;
401  complex double *fn_block = fn + offset*fn_n_stride;
402 
403  el = L0e;
404  switch (storage)
405  {
406  case SO3_STORAGE_PADDED:
407  so3_sampling_elmn2ind(&ind, 0, 0, n, parameters);
408  flm_block = flmn + ind;
409  i = offset = el*el;
410  break;
411  case SO3_STORAGE_COMPACT:
412  flm_block = flm;
413  i = offset = el*el-n*n;
414  break;
415  default:
416  SO3_ERROR_GENERIC("Invalid storage method.");
417  }
418 
419  (*ssht)(
420  flm_block, fn_block,
421  L0e, L, -n,
422  dl_method,
423  verbosity
424  );
425 
426  if (storage == SO3_STORAGE_COMPACT)
427  {
428  so3_sampling_elmn2ind(&ind, abs(n), -abs(n), n, parameters);
429  memcpy(flmn + ind, flm + n*n, (L*L - n*n) * sizeof(complex double));
430  }
431 
432  if (n % 2)
433  sign = -1;
434  else
435  sign = 1;
436 
437  for(; el < L; ++el)
438  {
439  factor = sign*sqrt(4.0*SO3_PI/(double)(2*el+1));
440  for (; i < offset + 2*el+1; ++i)
441  flmn[ind + i] *= factor;
442 
443  offset = i;
444  }
445 
446  if (storage == SO3_STORAGE_COMPACT)
447  free(flm);
448 
449  if (verbosity > 0)
450  printf("\n");
451  }
452 
453  free(fn);
454 
455  if (verbosity > 0)
456  printf("%sForward transform computed!\n", SO3_PROMPT);
457 
458 }
459 
476  double *f, const complex double *flmn,
477  const so3_parameters_t *parameters
478 ) {
479  int L0, L, N;
480  so3_sampling_t sampling;
481  so3_storage_t storage;
482  so3_n_mode_t n_mode;
483  ssht_dl_method_t dl_method;
484  int steerable;
485  int verbosity;
486 
487  // Iterator
488  int n;
489  // Intermediate results
490  complex double *fn, *flm;
491  double *ftemp;
492  // Stride for several arrays
493  int fn_n_stride;
494  // FFTW-related variables
495  int fftw_rank, fftw_howmany;
496  int fftw_idist, fftw_odist;
497  int fftw_istride, fftw_ostride;
498  int fftw_n;
499  double *fftw_target;
500  fftw_plan plan;
501 
502  inverse_complex_ssht complex_ssht;
503  inverse_real_ssht real_ssht;
504 
505  L0 = parameters->L0;
506  L = parameters->L;
507  N = parameters->N;
508  sampling = parameters->sampling_scheme;
509  storage = parameters->storage;
510  // TODO: Support N_MODE_L
511  n_mode = parameters->n_mode;
512  dl_method = parameters->dl_method;
513  verbosity = parameters->verbosity;
514  steerable = parameters->steerable;
515 
516  // Print messages depending on verbosity level.
517  if (verbosity > 0) {
518  printf("%sComputing inverse transform using MW sampling with\n", SO3_PROMPT);
519  printf("%sparameters (L, N, reality) = (%d, %d, FALSE)\n", SO3_PROMPT, L, N);
520  if (verbosity > 1)
521  printf("%sUsing routine so3_core_mw_inverse_via_ssht_real with storage method %d...\n"
522  , SO3_PROMPT
523  , storage);
524  }
525 
526  switch (sampling)
527  {
528  case SO3_SAMPLING_MW:
529  fn_n_stride = L * (2*L-1);
530  complex_ssht = ssht_core_mw_lb_inverse_sov_sym;
531  real_ssht = ssht_core_mw_lb_inverse_sov_sym_real;
532  break;
533  case SO3_SAMPLING_MW_SS:
534  fn_n_stride = (L+1) * 2*L;
535  complex_ssht = ssht_core_mw_lb_inverse_sov_sym_ss;
536  real_ssht = ssht_core_mw_lb_inverse_sov_sym_ss_real;
537  break;
538  default:
539  SO3_ERROR_GENERIC("Invalid sampling scheme.");
540  }
541 
542  // Compute fn(a,b)
543 
544  // For steerable signals, we need to supersample in n/gamma,
545  // in order to create a symmetric sampling.
546  // Each transform is over fftw_n samples (logically; physically, fn for negative n will be omitted)
547  //if (steerable)
548  //{
549  // For steerable signals, we need to supersample in n/gamma,
550  // in order to create a symmetric sampling.
551  // fftw_n = 2*N;
552 
553  // We need to perform the FFT into a temporary buffer, because
554  // the result will be twice as large as the output we need.
555  // ftemp = malloc(2*N*fn_n_stride * sizeof *ftemp);
556  // SO3_ERROR_MEM_ALLOC_CHECK(ftemp);
557 
558  // fftw_target = ftemp;
559  //}
560  //else
561  //{
562  fftw_n = 2*N-1;
563 
564  fftw_target = f;
565  //}
566 
567  // Only need to store for non-negative n
568  fn = calloc((fftw_n/2+1)*fn_n_stride, sizeof *fn);
569  SO3_ERROR_MEM_ALLOC_CHECK(fn);
570 
571  // Initialize fftw_plan first. With FFTW_ESTIMATE this is technically not
572  // necessary but still good practice.
573  fftw_rank = 1; // We compute 1d transforms
574  fftw_howmany = fn_n_stride; // We need L*(2*L-1) of these transforms
575 
576  // We want to transform columns
577  fftw_idist = fftw_odist = 1; // The starts of the columns are contiguous in memory
578  fftw_istride = fftw_ostride = fn_n_stride; // Distance between two elements of the same column
579 
580  plan = fftw_plan_many_dft_c2r(
581  fftw_rank, &fftw_n, fftw_howmany,
582  fn, NULL, fftw_istride, fftw_idist,
583  fftw_target, NULL, fftw_ostride, fftw_odist,
584  FFTW_ESTIMATE
585  );
586 
587  flm = malloc(L*L * sizeof *flm);
588  SO3_ERROR_MEM_ALLOC_CHECK(flm);
589 
590  for(n = 0; n <= N-1; ++n)
591  {
592  int ind, offset, i, el;
593  int L0e = MAX(L0, abs(n)); // 'e' for 'effective'
594  double factor;
595 
596  if ((n_mode == SO3_N_MODE_EVEN && n % 2)
597  || (n_mode == SO3_N_MODE_ODD && !(n % 2))
598  || (n_mode == SO3_N_MODE_MAXIMUM && abs(n) < N-1)
599  ) {
600  continue;
601  }
602 
603  switch (storage)
604  {
605  case SO3_STORAGE_PADDED:
606  so3_sampling_elmn2ind_real(&ind, 0, 0, n, parameters); //L, N, SO3_STORE_NEG_FIRST_PAD);
607  memcpy(flm, flmn + ind, L*L * sizeof(complex double));
608  break;
609  case SO3_STORAGE_COMPACT:
610  so3_sampling_elmn2ind_real(&ind, n, -n, n, parameters); //L, N, SO3_STORE_NEG_FIRST_COMPACT);
611  memcpy(flm + n*n, flmn + ind, (L*L - n*n) * sizeof(complex double));
612  for(i = 0; i < n*n; ++i)
613  flm[i] = 0.0;
614  break;
615  default:
616  SO3_ERROR_GENERIC("Invalid storage method.");
617  }
618 
619  el = L0e;
620  i = offset = el*el;
621  for(; el < L; ++el)
622  {
623  factor = sqrt((double)(2*el+1)/(16.*pow(SO3_PI, 3.)));
624  for (; i < offset + 2*el+1; ++i)
625  flm[i] *= factor;
626 
627  offset = i;
628  }
629 
630  if (N > 1 || n)
631  {
632  (*complex_ssht)(
633  fn + n*fn_n_stride, flm,
634  L0e, L, -n,
635  dl_method,
636  verbosity
637  );
638  }
639  else
640  {
641  double *fn_r;
642 
643  // Create an array of real doubles for n = 0
644  fn_r = calloc(fn_n_stride, sizeof *fn_r);
645  SO3_ERROR_MEM_ALLOC_CHECK(fn_r);
646 
647  (*real_ssht)(
648  fn_r, flm,
649  L0e, L,
650  dl_method,
651  verbosity
652  );
653 
654  for (i = 0; i < fn_n_stride; ++i)
655  fn[i] = fn_r[i];
656  }
657 
658 
659  if(n % 2)
660  for(i = 0; i < fn_n_stride; ++i)
661  fn[n*fn_n_stride + i] = -fn[n*fn_n_stride + i];
662 
663  if (verbosity > 0)
664  printf("\n");
665  }
666 
667  free(flm);
668 
669  fftw_execute(plan);
670  fftw_destroy_plan(plan);
671 
672  //if (steerable)
673  //{
674  // memcpy(f, ftemp, N*fn_n_stride * sizeof *f);
675  // free(ftemp);
676  //}
677 
678  free(fn);
679 
680  if (verbosity > 0)
681  printf("%sInverse transform computed!\n", SO3_PROMPT);
682 }
683 
701  complex double *flmn, const double *f,
702  const so3_parameters_t *parameters
703 ) {
704  int L0, L, N;
705  so3_sampling_t sampling;
706  so3_storage_t storage;
707  so3_n_mode_t n_mode;
708  ssht_dl_method_t dl_method;
709  int steerable;
710  int verbosity;
711 
712  // Iterator
713  int i, n;
714  // Intermediate results
715  double *ftemp;
716  complex double *flm = NULL, *fn;
717  // Stride for several arrays
718  int fn_n_stride;
719  // FFTW-related variables
720  int fftw_rank, fftw_howmany;
721  int fftw_idist, fftw_odist;
722  int fftw_istride, fftw_ostride;
723  int fftw_n;
724  fftw_plan plan;
725 
726  forward_complex_ssht complex_ssht;
727  forward_real_ssht real_ssht;
728 
729  // for precomputation
730  double factor;
731 
732  L0 = parameters->L0;
733  L = parameters->L;
734  N = parameters->N;
735  sampling = parameters->sampling_scheme;
736  storage = parameters->storage;
737  // TODO: Support N_MODE_L
738  n_mode = parameters->n_mode;
739  dl_method = parameters->dl_method;
740  steerable = parameters->steerable;
741  verbosity = parameters->verbosity;
742 
743  // Print messages depending on verbosity level.
744  if (verbosity > 0) {
745  printf("%sComputing forward transform using MW sampling with\n", SO3_PROMPT);
746  printf("%sparameters (L, N, reality) = (%d, %d, FALSE)\n", SO3_PROMPT, L, N);
747  if (verbosity > 1)
748  printf("%sUsing routine so3_core_mw_forward_via_ssht_real with storage method %d...\n"
749  , SO3_PROMPT
750  , storage);
751  }
752 
753  switch (sampling)
754  {
755  case SO3_SAMPLING_MW:
756  fn_n_stride = L * (2*L-1);
757  complex_ssht = ssht_core_mw_lb_forward_sov_conv_sym;
758  real_ssht = ssht_core_mw_lb_forward_sov_conv_sym_real;
759  break;
760  case SO3_SAMPLING_MW_SS:
761  fn_n_stride = (L+1) * 2*L;
762  complex_ssht = ssht_core_mw_lb_forward_sov_conv_sym_ss;
763  real_ssht = ssht_core_mw_lb_forward_sov_conv_sym_ss_real;
764  break;
765  default:
766  SO3_ERROR_GENERIC("Invalid sampling scheme.");
767  }
768 
769 
770  //if (steerable)
771  //{
772  // int g, offset;
773 
774  // fn = calloc((2*N-1)*fn_n_stride, sizeof *fn);
775  // SO3_ERROR_MEM_ALLOC_CHECK(fn);
776 
777  // for (n = -N+1; n < N; n+=2)
778  // {
779  // The conditional applies the spatial transform, because the fn
780  // are to be stored in n-order 0, 1, 2, -2, -1
781  // offset = (n < 0 ? n + 2*N-1 : n);
782 
783  // for (g = 0; g < N; ++g)
784  // {
785  // double gamma = g * SO3_PI / N;
786  // for (i = 0; i < fn_n_stride; ++i)
787  // {
788  // double weight = 2*SO3_PI/N;
789  // fn[offset * fn_n_stride + i] += weight*f[g * fn_n_stride + i]*cexp(-I*n*gamma);
790  // }
791  // }
792  // }
793  //}
794  //else
795  //{
796  // Make a copy of the input, because input is const
797  // This could potentially be avoided by copying the input into fn and using an
798  // in-place FFTW. The performance impact has to be profiled, though.
799  ftemp = malloc((2*N-1)*fn_n_stride * sizeof *ftemp);
800  SO3_ERROR_MEM_ALLOC_CHECK(ftemp);
801  memcpy(ftemp, f, (2*N-1)*fn_n_stride * sizeof(double));
802 
803  fn = malloc(N*fn_n_stride * sizeof *fn);
804  SO3_ERROR_MEM_ALLOC_CHECK(fn);
805  // Initialize fftw_plan first. With FFTW_ESTIMATE this is technically not
806  // necessary but still good practice.
807  fftw_rank = 1; // We compute 1d transforms
808  fftw_n = 2*N-1; // Each transform is over 2*N-1 (logically; physically, fn for negative n will be omitted)
809  fftw_howmany = fn_n_stride; // We need L*(2*L-1) of these transforms
810 
811  // We want to transform columns
812  fftw_idist = fftw_odist = 1; // The starts of the columns are contiguous in memory
813  fftw_istride = fftw_ostride = fn_n_stride; // Distance between two elements of the same column
814 
815  plan = fftw_plan_many_dft_r2c(
816  fftw_rank, &fftw_n, fftw_howmany,
817  ftemp, NULL, fftw_istride, fftw_idist,
818  fn, NULL, fftw_ostride, fftw_odist,
819  FFTW_ESTIMATE
820  );
821 
822  fftw_execute(plan);
823  fftw_destroy_plan(plan);
824 
825  free(ftemp);
826 
827  factor = 2*SO3_PI/(double)(2*N-1);
828  for(i = 0; i < N*fn_n_stride; ++i)
829  fn[i] *= factor;
830  //}
831 
832  if (storage == SO3_STORAGE_COMPACT)
833  flm = malloc(L*L * sizeof *flm);
834 
835  for(n = 0; n <= N-1; ++n)
836  {
837  int ind, offset, el, sign;
838  int L0e = MAX(L0, abs(n)); // 'e' for 'effective'
839 
840  complex double* flm_block;
841 
842  if ((n_mode == SO3_N_MODE_EVEN && n % 2)
843  || (n_mode == SO3_N_MODE_ODD && !(n % 2))
844  || (n_mode == SO3_N_MODE_MAXIMUM && abs(n) < N-1)
845  ) {
846  continue;
847  }
848 
849  el = L0e;
850  i = offset = el*el;
851  switch (storage)
852  {
853  case SO3_STORAGE_PADDED:
854  so3_sampling_elmn2ind_real(&ind, 0, 0, n, parameters);
855  flm_block = flmn + ind;
856  break;
857  case SO3_STORAGE_COMPACT:
858  flm_block = flm;
859 
860  offset -= n*n;
861  i = offset;
862  break;
863  default:
864  SO3_ERROR_GENERIC("Invalid storage method.");
865  }
866 
867 
868  if (N > 1)
869  {
870  (*complex_ssht)(
871  flm_block, fn + n*fn_n_stride,
872  L0e, L, -n,
873  dl_method,
874  verbosity
875  );
876  }
877  else
878  {
879  // Now we know n = 0 in which case the reality conditions
880  // for SO3 and SSHT coincide.
881  int j;
882  double *fn_r;
883 
884  // Create an array of real doubles for n = 0
885  fn_r = malloc(fn_n_stride * sizeof *fn_r);
886  SO3_ERROR_MEM_ALLOC_CHECK(fn_r);
887  for (j = 0; j < fn_n_stride; ++j)
888  fn_r[j] = creal(fn[j]);
889 
890  // Now use real SSHT transforms
891  (*real_ssht)(
892  flm_block, fn_r,
893  L0e, L,
894  dl_method,
895  verbosity
896  );
897 
898  free(fn_r);
899  }
900 
901  if (storage == SO3_STORAGE_COMPACT)
902  {
903  so3_sampling_elmn2ind_real(&ind, n, -n, n, parameters);
904  memcpy(flmn + ind, flm + n*n, (L*L - n*n) * sizeof(complex double));
905  }
906 
907  if (n % 2)
908  sign = -1;
909  else
910  sign = 1;
911 
912  for(; el < L; ++el)
913  {
914  factor = sign*sqrt(4.0*SO3_PI/(double)(2*el+1));
915  for (; i < offset + 2*el+1; ++i)
916  flmn[ind + i] *= factor;
917 
918  offset = i;
919  }
920 
921  if (verbosity > 0)
922  printf("\n");
923  }
924 
925  if (storage == SO3_STORAGE_COMPACT)
926  free(flm);
927 
928  free(fn);
929 
930  if (verbosity > 0)
931  printf("%sForward transform computed!\n", SO3_PROMPT);
932 
933 }
934 
935 
952  complex double *f, const complex double *flmn,
953  const so3_parameters_t *parameters
954 ) {
955 
956  int L0, L, N;
957  so3_sampling_t sampling;
958  so3_storage_t storage;
959  so3_n_mode_t n_mode;
960  ssht_dl_method_t dl_method;
961  int steerable;
962  int verbosity;
963 
964  L0 = parameters->L0;
965  L = parameters->L;
966  N = parameters->N;
967  sampling = parameters->sampling_scheme;
968  storage = parameters->storage;
969  // TODO: Add optimisations for all n-modes.
970  n_mode = parameters->n_mode;
971  dl_method = parameters->dl_method;
972  verbosity = parameters->verbosity;
973  steerable = parameters->steerable;
974 
975  // Print messages depending on verbosity level.
976  if (verbosity > 0)
977  {
978  printf("%sComputing inverse transform using MW sampling with\n", SO3_PROMPT);
979  printf("%sparameters (L, N, reality) = (%d, %d, FALSE)\n", SO3_PROMPT, L, N);
980  if (verbosity > 1)
981  printf("%sUsing routine so3_core_mw_inverse_direct with storage method %d...\n"
982  , SO3_PROMPT
983  , storage);
984  }
985 
986  // Iterators
987  int el, m, n, mm; // mm for m'
988 
989  // Allocate memory.
990  double *sqrt_tbl = calloc(2*(L-1)+2, sizeof(*sqrt_tbl));
991  SO3_ERROR_MEM_ALLOC_CHECK(sqrt_tbl);
992  double *signs = calloc(L+1, sizeof(*signs));
993  SO3_ERROR_MEM_ALLOC_CHECK(signs);
994  complex double *exps = calloc(4, sizeof(*exps));
995  SO3_ERROR_MEM_ALLOC_CHECK(exps);
996 
997  // Perform precomputations.
998  for (el = 0; el <= 2*L-1; ++el)
999  sqrt_tbl[el] = sqrt((double)el);
1000  for (m = 0; m <= L-1; m += 2)
1001  {
1002  signs[m] = 1.0;
1003  signs[m+1] = -1.0;
1004  }
1005  int i;
1006  for (i = 0; i < 4; ++i)
1007  exps[i] = cexp(I*SO3_PION2*i);
1008 
1009  // Compute Fmnm'
1010  // TODO: Currently m is fastest-varying, then n, then m'.
1011  // Should this order be changed to m-m'-n?
1012  complex double *Fmnm = calloc((2*L-1)*(2*L-1)*(2*N-1), sizeof(*Fmnm));
1013  SO3_ERROR_MEM_ALLOC_CHECK(Fmnm);
1014  int m_offset = L-1;
1015  int m_stride = 2*L-1;
1016  int n_offset = N-1;
1017  int n_stride = 2*N-1;
1018  int mm_offset = L-1;
1019  int mm_stride = 2*L-1;
1020 
1021  int n_start, n_stop, n_inc;
1022 
1023  double *dl = ssht_dl_calloc(L, SSHT_DL_QUARTER);
1024  SO3_ERROR_MEM_ALLOC_CHECK(dl);
1025  double *dl8 = NULL;
1026  if (dl_method == SSHT_DL_RISBO)
1027  {
1028  dl8 = ssht_dl_calloc(L, SSHT_DL_QUARTER_EXTENDED);
1029  SO3_ERROR_MEM_ALLOC_CHECK(dl8);
1030  }
1031  int dl_offset = ssht_dl_get_offset(L, SSHT_DL_QUARTER);
1032  int dl_stride = ssht_dl_get_stride(L, SSHT_DL_QUARTER);
1033 
1034  complex double *mn_factors = calloc((2*L-1)*(2*N-1), sizeof *mn_factors);
1035  SO3_ERROR_MEM_ALLOC_CHECK(mn_factors);
1036 
1037  // TODO: SSHT starts this loop from MAX(L0, abs(spin)).
1038  // Can we use a similar optimisation? el can probably
1039  // be limited by n, but then we'd need to switch the
1040  // loop order, which means we'd have to recompute the
1041  // Wigner plane for each n. That seems wrong?
1042  for (el = L0; el <= L-1; ++el)
1043  {
1044  int eltmp;
1045  // Compute Wigner plane.
1046  switch (dl_method)
1047  {
1048  case SSHT_DL_RISBO:
1049  if (el != 0 && el == L0)
1050  {
1051  for(eltmp = 0; eltmp <= L0; ++eltmp)
1052  ssht_dl_beta_risbo_eighth_table(
1053  dl8,
1054  SO3_PION2,
1055  L,
1056  SSHT_DL_QUARTER_EXTENDED,
1057  eltmp,
1058  sqrt_tbl,
1059  signs);
1060  ssht_dl_beta_risbo_fill_eighth2quarter_table(dl,
1061  dl8, L,
1062  SSHT_DL_QUARTER,
1063  SSHT_DL_QUARTER_EXTENDED,
1064  el,
1065  signs);
1066  }
1067  else
1068  {
1069  ssht_dl_beta_risbo_eighth_table(
1070  dl8,
1071  SO3_PION2,
1072  L,
1073  SSHT_DL_QUARTER_EXTENDED,
1074  el,
1075  sqrt_tbl,
1076  signs);
1077  ssht_dl_beta_risbo_fill_eighth2quarter_table(
1078  dl,
1079  dl8, L,
1080  SSHT_DL_QUARTER,
1081  SSHT_DL_QUARTER_EXTENDED,
1082  el,
1083  signs);
1084  }
1085  break;
1086 
1087  case SSHT_DL_TRAPANI:
1088  if (el != 0 && el == L0)
1089  {
1090  for(eltmp = 0; eltmp <= L0; ++eltmp)
1091  ssht_dl_halfpi_trapani_eighth_table(
1092  dl,
1093  L,
1094  SSHT_DL_QUARTER,
1095  eltmp,
1096  sqrt_tbl);
1097  ssht_dl_halfpi_trapani_fill_eighth2quarter_table(
1098  dl,
1099  L,
1100  SSHT_DL_QUARTER,
1101  el,
1102  signs);
1103  }
1104  else
1105  {
1106  ssht_dl_halfpi_trapani_eighth_table(
1107  dl,
1108  L,
1109  SSHT_DL_QUARTER,
1110  el,
1111  sqrt_tbl);
1112  ssht_dl_halfpi_trapani_fill_eighth2quarter_table(
1113  dl,
1114  L,
1115  SSHT_DL_QUARTER,
1116  el,
1117  signs);
1118  }
1119  break;
1120 
1121  default:
1122  SO3_ERROR_GENERIC("Invalid dl method");
1123  }
1124 
1125  // Compute Fmnm' contribution for current el.
1126 
1127  // Factor which depends only on el.
1128  double elfactor = (2.0*el+1.0)/(8.0*SO3_PI*SO3_PI);
1129 
1130  switch (n_mode)
1131  {
1132  case SO3_N_MODE_ALL:
1133  n_start = MAX(-N+1,-el);
1134  n_stop = MIN( N-1, el);
1135  n_inc = 1;
1136  break;
1137  case SO3_N_MODE_EVEN:
1138  n_start = MAX(-N+1,-el);
1139  n_start += (-n_start)%2;
1140  n_stop = MIN( N-1, el);
1141  n_stop -= n_stop%2;
1142  n_inc = 2;
1143  break;
1144  case SO3_N_MODE_ODD:
1145  n_start = MAX(-N+1,-el);
1146  n_start += 1+n_start%2;
1147  n_stop = MIN( N-1, el);
1148  n_stop -= 1-n_stop%2;
1149  n_inc = 2;
1150  break;
1151  case SO3_N_MODE_MAXIMUM:
1152  if (el < N-1)
1153  continue;
1154  n_start = -N+1;
1155  n_stop = N-1;
1156  n_inc = MAX(1,2*N-2);
1157  break;
1158  case SO3_N_MODE_L:
1159  if (el >= N)
1160  continue;
1161  n_start = -el;
1162  n_stop = el;
1163  n_inc = MAX(1,2*el);
1164  break;
1165  default:
1166  SO3_ERROR_GENERIC("Invalid n-mode.");
1167  }
1168 
1169  // Factors which do not depend on m'.
1170  for (n = n_start; n <= n_stop; n += n_inc)
1171  for (m = -el; m <= el; ++m)
1172  {
1173  int ind;
1174  so3_sampling_elmn2ind(&ind, el, m, n, parameters);
1175  int mod = ((n-m)%4 + 4)%4;
1176  mn_factors[m + m_offset + m_stride*(
1177  n + n_offset)] =
1178  flmn[ind] * exps[mod];
1179  }
1180 
1181  for (mm = 0; mm <= el; ++mm)
1182  {
1183  // These signs are needed for the symmetry relations of
1184  // Wigner symbols.
1185  double elmmsign = signs[el] * signs[mm];
1186 
1187  // TODO: If the conditional for elnsign is a bottleneck
1188  // this loop can be split up just like the inner loop.
1189  for (n = n_start; n <= n_stop; n += n_inc)
1190  {
1191  double elnsign = n >= 0 ? 1.0 : elmmsign;
1192  // Factor which does not depend on m.
1193  double elnmm_factor = elfactor * elnsign
1194  * dl[abs(n) + dl_offset + mm*dl_stride];
1195  for (m = -el; m < 0; ++m)
1196  Fmnm[m + m_offset + m_stride*(
1197  n + n_offset + n_stride*(
1198  mm + mm_offset))] +=
1199  elnmm_factor
1200  * mn_factors[m + m_offset + m_stride*(
1201  n + n_offset)]
1202  * elmmsign * dl[-m + dl_offset + mm*dl_stride];
1203  for (m = 0; m <= el; ++m)
1204  Fmnm[m + m_offset + m_stride*(
1205  n + n_offset + n_stride*(
1206  mm + mm_offset))] +=
1207  elnmm_factor
1208  * mn_factors[m + m_offset + m_stride*(
1209  n + n_offset)]
1210  * dl[m + dl_offset + mm*dl_stride];
1211  }
1212  }
1213  }
1214 
1215  // Free dl memory.
1216  free(dl);
1217  if (dl_method == SSHT_DL_RISBO)
1218  free(dl8);
1219 
1220  switch (n_mode)
1221  {
1222  case SO3_N_MODE_ALL:
1223  case SO3_N_MODE_L:
1224  n_start = -N+1;
1225  n_stop = N-1;
1226  n_inc = 1;
1227  break;
1228  case SO3_N_MODE_EVEN:
1229  n_start = ((N-1) % 2 == 0) ? -N+1 : -N+2;
1230  n_stop = ((N-1) % 2 == 0) ? N-1 : N-2;
1231  n_inc = 2;
1232  break;
1233  case SO3_N_MODE_ODD:
1234  n_start = ((N-1) % 2 != 0) ? -N+1 : -N+2;
1235  n_stop = ((N-1) % 2 != 0) ? N-1 : N-2;
1236  n_inc = 2;
1237  break;
1238  case SO3_N_MODE_MAXIMUM:
1239  n_start = -N+1;
1240  n_stop = N-1;
1241  n_inc = MAX(1,2*N - 2);
1242  break;
1243  default:
1244  SO3_ERROR_GENERIC("Invalid n-mode.");
1245  }
1246 
1247  // Use symmetry to compute Fmnm' for negative m'.
1248  for (mm = -L+1; mm < 0; ++mm)
1249  for (n = n_start; n <= n_stop; n += n_inc)
1250  for (m = -L+1; m <= L-1; ++m)
1251  Fmnm[m + m_offset + m_stride*(
1252  n + n_offset + n_stride*(
1253  mm + mm_offset))] =
1254  signs[abs(m+n)%2]
1255  * Fmnm[m + m_offset + m_stride*(
1256  n + n_offset + n_stride*(
1257  -mm + mm_offset))];
1258 
1259  // Apply phase modulation to account for sampling offset.
1260  for (mm = -L+1; mm <= L-1; ++mm)
1261  {
1262  complex double mmfactor = cexp(I*mm*SO3_PI/(2.0*L-1.0));
1263  for (n = n_start; n <= n_stop; n += n_inc)
1264  for (m = -L+1; m <= L-1; ++m)
1265  Fmnm[m + m_offset + m_stride*(
1266  n + n_offset + n_stride*(
1267  mm + mm_offset))] *= mmfactor;
1268  }
1269 
1270  // Allocate space for function values.
1271  complex double *fext = calloc((2*L-1)*(2*L-1)*(2*N-1), sizeof(*fext));
1272  SO3_ERROR_MEM_ALLOC_CHECK(fext);
1273 
1274 
1275  // Set up plan before initialising array.
1276  fftw_plan plan = fftw_plan_dft_3d(
1277  2*N-1, 2*L-1, 2*L-1,
1278  fext, fext,
1279  FFTW_BACKWARD,
1280  FFTW_ESTIMATE);
1281 
1282  // Apply spatial shift.
1283  for (mm = -L+1; mm <= L-1; ++mm)
1284  {
1285  int mm_shift = mm < 0 ? 2*L-1 : 0;
1286  for (n = n_start; n <= n_stop; n += n_inc)
1287  {
1288  int n_shift = n < 0 ? 2*N-1 : 0;
1289  for (m = -L+1; m <= L-1; ++m)
1290  {
1291  int m_shift = m < 0 ? 2*L-1 : 0;
1292  fext[m + m_shift + m_stride*(
1293  mm + mm_shift + mm_stride*(
1294  n + n_shift))] =
1295  Fmnm[m + m_offset + m_stride*(
1296  n + n_offset + n_stride*(
1297  mm + mm_offset))];
1298  }
1299  }
1300  }
1301 
1302  // Free Fmnm' memory.
1303  free(Fmnm);
1304 
1305  // Perform 3D FFT.
1306  fftw_execute(plan);
1307  fftw_destroy_plan(plan);
1308 
1309  // Extract f from the extended torus.
1310  int a,b,g;
1311  int a_stride = 2*L-1;
1312  int b_ext_stride = 2*L-1;
1313  int b_stride = L;
1314  for (g = 0; g < 2*N-1; ++g)
1315  for (b = 0; b < L; ++b)
1316  for (a = 0; a < 2*L-1; ++a)
1317  f[a + a_stride*(
1318  b + b_stride*(
1319  g))] = fext[a + a_stride*(
1320  b + b_ext_stride*(
1321  g))];
1322 
1323  // Free fext memory.
1324  free(fext);
1325 
1326  if (verbosity > 0)
1327  printf("%sInverse transform computed!\n", SO3_PROMPT);
1328 
1329  // Free precomputation memory.
1330  free(sqrt_tbl);
1331  free(signs);
1332  free(exps);
1333  free(mn_factors);
1334 }
1335 
1336 
1355  complex double *flmn, const complex double *f,
1356  const so3_parameters_t *parameters
1357 ) {
1358  int L0, L, N;
1359  so3_sampling_t sampling;
1360  so3_storage_t storage;
1361  so3_n_mode_t n_mode;
1362  ssht_dl_method_t dl_method;
1363  int steerable;
1364  int verbosity;
1365 
1366  L0 = parameters->L0;
1367  L = parameters->L;
1368  N = parameters->N;
1369  sampling = parameters->sampling_scheme;
1370  storage = parameters->storage;
1371  // TODO: Add optimisations for all n-modes.
1372  n_mode = parameters->n_mode;
1373  dl_method = parameters->dl_method;
1374  verbosity = parameters->verbosity;
1375  steerable = parameters->steerable;
1376 
1377  // Print messages depending on verbosity level.
1378  if (verbosity > 0) {
1379  printf("%sComputing forward transform using MW sampling with\n", SO3_PROMPT);
1380  printf("%sparameters (L, N, reality) = (%d, %d, FALSE)\n", SO3_PROMPT, L, N);
1381  if (verbosity > 1)
1382  printf("%sUsing routine so3_core_mw_forward_direct with storage method %d...\n"
1383  , SO3_PROMPT
1384  , storage);
1385  }
1386 
1387  int m_stride = 2*L-1;
1388  int m_offset = L-1;
1389  // unused: int n_stride = 2*N-1;
1390  int n_offset = N-1;
1391  int mm_stride = 2*L-1;
1392  int mm_offset = L-1;
1393  int a_stride = 2*L-1;
1394  int b_stride = L;
1395  int bext_stride = 2*L-1;
1396  // unused: int g_stride = 2*N-1;
1397 
1398  int n_start, n_stop, n_inc;
1399 
1400  switch (n_mode)
1401  {
1402  case SO3_N_MODE_ALL:
1403  case SO3_N_MODE_L:
1404  n_start = -N+1;
1405  n_stop = N-1;
1406  n_inc = 1;
1407  break;
1408  case SO3_N_MODE_EVEN:
1409  n_start = ((N-1) % 2 == 0) ? -N+1 : -N+2;
1410  n_stop = ((N-1) % 2 == 0) ? N-1 : N-2;
1411  n_inc = 2;
1412  break;
1413  case SO3_N_MODE_ODD:
1414  n_start = ((N-1) % 2 != 0) ? -N+1 : -N+2;
1415  n_stop = ((N-1) % 2 != 0) ? N-1 : N-2;
1416  n_inc = 2;
1417  break;
1418  case SO3_N_MODE_MAXIMUM:
1419  n_start = -N+1;
1420  n_stop = N-1;
1421  n_inc = MAX(1,2*N - 2);
1422  break;
1423  default:
1424  SO3_ERROR_GENERIC("Invalid n-mode.");
1425  }
1426 
1427  double *sqrt_tbl = calloc(2*(L-1)+2, sizeof(*sqrt_tbl));
1428  SO3_ERROR_MEM_ALLOC_CHECK(sqrt_tbl);
1429  double *signs = calloc(L+1, sizeof(*signs));
1430  SO3_ERROR_MEM_ALLOC_CHECK(signs);
1431  complex double *exps = calloc(4, sizeof(*exps));
1432  SO3_ERROR_MEM_ALLOC_CHECK(exps);
1433  complex double *expsmm = calloc(2*L-1, sizeof(*expsmm));
1434  SO3_ERROR_MEM_ALLOC_CHECK(expsmm);
1435 
1436  int el, m, n, mm; // mm is for m'
1437  // Perform precomputations.
1438  for (el = 0; el <= 2*(L-1)+1; ++el)
1439  sqrt_tbl[el] = sqrt((double)el);
1440  for (m = 0; m <= L-1; m += 2)
1441  {
1442  signs[m] = 1.0;
1443  signs[m+1] = -1.0;
1444  }
1445  int i;
1446  for (i = 0; i < 4; ++i)
1447  exps[i] = cexp(I*SO3_PION2*i);
1448  for (mm = -L+1; mm <= L-1; ++mm)
1449  expsmm[mm + mm_offset] = cexp(-I*mm*SSHT_PI/(2.0*L-1.0));
1450 
1451  double norm_factor = 1.0/(2.0*L-1.0)/(2.0*N-1.0);
1452 
1453  // Compute Fourier transform over alpha and gamma, i.e. compute Fmn(b).
1454  complex double *Fmnb = calloc((2*L-1)*(2*L-1)*(2*N-1), sizeof(*Fmnb));
1455  SO3_ERROR_MEM_ALLOC_CHECK(Fmnb);
1456  complex double *inout = calloc((2*L-1)*(2*N-1), sizeof(*inout));
1457  SO3_ERROR_MEM_ALLOC_CHECK(inout);
1458  fftw_plan plan = fftw_plan_dft_2d(
1459  2*N-1, 2*L-1,
1460  inout, inout,
1461  FFTW_FORWARD,
1462  FFTW_ESTIMATE);
1463 
1464  int b, g;
1465  for (b = 0; b < L; ++b)
1466  {
1467  // TODO: This memcpy loop could probably be avoided by using
1468  // a more elaborate FFTW plan which performs the FFT directly
1469  // over the 1st and 3rd dimensions of f.
1470  for (g = 0; g < 2*N-1; ++g)
1471  memcpy(
1472  inout + g*a_stride,
1473  f + 0 + a_stride*(
1474  b + b_stride*(
1475  g)),
1476  a_stride*sizeof(*f));
1477  fftw_execute_dft(plan, inout, inout);
1478 
1479  // Apply spatial shift and normalisation factor
1480  for (n = n_start; n <= n_stop; n += n_inc)
1481  {
1482  int n_shift = n < 0 ? 2*N-1 : 0;
1483  for (m = -L+1; m <= L-1; ++m)
1484  {
1485  int m_shift = m < 0 ? 2*L-1 : 0;
1486  Fmnb[b + bext_stride*(
1487  m + m_offset + m_stride*(
1488  n + n_offset))] =
1489  inout[m + m_shift + m_stride*(
1490  n + n_shift)] * norm_factor;
1491  }
1492  }
1493  }
1494  fftw_destroy_plan(plan);
1495 
1496  // Extend Fmnb periodically.
1497  for (n = n_start; n <= n_stop; n += n_inc)
1498  for (m = -L+1; m <= L-1; ++m)
1499  {
1500  int signmn = signs[abs(m+n)%2];
1501  for (b = L; b < 2*L-1; ++b)
1502  Fmnb[b + bext_stride*(
1503  m + m_offset + m_stride*(
1504  n + n_offset))] =
1505  signmn
1506  * Fmnb[(2*L-2-b) + bext_stride*(
1507  m + m_offset + m_stride*(
1508  n + n_offset))];
1509  }
1510 
1511 
1512  // Compute Fourier transform over beta, i.e. compute Fmnm'.
1513  complex double *Fmnm = calloc((2*L-1)*(2*L-1)*(2*N-1), sizeof(*Fmnm));
1514  SO3_ERROR_MEM_ALLOC_CHECK(Fmnm);
1515 
1516  plan = fftw_plan_dft_1d(
1517  2*L-1,
1518  inout, inout,
1519  FFTW_FORWARD,
1520  FFTW_ESTIMATE);
1521  for (n = n_start; n <= n_stop; n += n_inc)
1522  for (m = -L+1; m <= L-1; ++m)
1523  {
1524  memcpy(inout,
1525  Fmnb + 0 + bext_stride*(
1526  m + m_offset + m_stride*(
1527  n + n_offset)),
1528  bext_stride*sizeof(*Fmnb));
1529  fftw_execute_dft(plan, inout, inout);
1530 
1531  // Apply spatial shift and normalisation factor
1532  for (mm = -L+1; mm <= L-1; ++mm)
1533  {
1534  int mm_shift = mm < 0 ? 2*L-1 : 0;
1535  Fmnm[mm + mm_offset + mm_stride*(
1536  m + m_offset + m_stride*(
1537  n + n_offset))] =
1538  inout[mm + mm_shift] / (2.0*L-1.0);
1539  }
1540  }
1541  fftw_destroy_plan(plan);
1542  free(inout);
1543 
1544  // Apply phase modulation to account for sampling offset.
1545  for (n = n_start; n <= n_stop; n += n_inc)
1546  for (m = -L+1; m <= L-1; ++m)
1547  for (mm = -L+1; mm <= L-1; ++mm)
1548  Fmnm[mm + mm_offset + mm_stride*(
1549  m + m_offset + m_stride*(
1550  n + n_offset))] *=
1551  expsmm[mm + mm_offset];
1552 
1553  // Compute weights.
1554  complex double *w = calloc(4*L-3, sizeof(*w));
1555  SO3_ERROR_MEM_ALLOC_CHECK(w);
1556  int w_offset = 2*(L-1);
1557  for (mm = -2*(L-1); mm <= 2*(L-1); ++mm)
1558  w[mm+w_offset] = so3_sampling_weight(parameters, mm);
1559 
1560  // Compute IFFT of w to give wr.
1561  complex double *wr = calloc(4*L-3, sizeof(*w));
1562  SO3_ERROR_MEM_ALLOC_CHECK(wr);
1563  inout = calloc(4*L-3, sizeof(*inout));
1564  SO3_ERROR_MEM_ALLOC_CHECK(inout);
1565  fftw_plan plan_bwd = fftw_plan_dft_1d(
1566  4*L-3,
1567  inout, inout,
1568  FFTW_BACKWARD,
1569  FFTW_MEASURE);
1570  fftw_plan plan_fwd = fftw_plan_dft_1d(
1571  4*L-3,
1572  inout,
1573  inout,
1574  FFTW_FORWARD,
1575  FFTW_MEASURE);
1576 
1577  // Apply spatial shift.
1578  for (mm = 1; mm <= 2*L-2; ++mm)
1579  inout[mm + w_offset] = w[mm - 2*(L-1) - 1 + w_offset];
1580  for (mm = -2*(L-1); mm <= 0; ++mm)
1581  inout[mm + w_offset] = w[mm + 2*(L-1) + w_offset];
1582 
1583  fftw_execute_dft(plan_bwd, inout, inout);
1584 
1585  // Apply spatial shift.
1586  for (mm = 0; mm <= 2*L-2; ++mm)
1587  wr[mm + w_offset] = inout[mm - 2*(L-1) + w_offset];
1588  for (mm = -2*(L-1); mm <= -1; ++mm)
1589  wr[mm + w_offset] = inout[mm + 2*(L-1) + 1 + w_offset];
1590 
1591  // Compute Gmnm' by convolution implemented as product in real space.
1592  complex double *Fmnm_pad = calloc(4*L-3, sizeof(*Fmnm_pad));
1593  SO3_ERROR_MEM_ALLOC_CHECK(Fmnm_pad);
1594  complex double *Gmnm = calloc((2*L-1)*(2*L-1)*(2*N-1), sizeof(*Gmnm));
1595  SO3_ERROR_MEM_ALLOC_CHECK(Gmnm);
1596  for (n = n_start; n <= n_stop; n += n_inc)
1597  for (m = -L+1; m <= L-1; ++m)
1598  {
1599 
1600  // Zero-pad Fmnm'.
1601  for (mm = -2*(L-1); mm <= -L; ++mm)
1602  Fmnm_pad[mm+w_offset] = 0.0;
1603  for (mm = L; mm <= 2*(L-1); ++mm)
1604  Fmnm_pad[mm+w_offset] = 0.0;
1605  for (mm = -(L-1); mm <= L-1; ++mm)
1606  Fmnm_pad[mm + w_offset] =
1607  Fmnm[mm + mm_offset + mm_stride*(
1608  m + m_offset + m_stride*(
1609  n + n_offset))];
1610 
1611  // Apply spatial shift.
1612  for (mm = 1; mm <= 2*L-2; ++mm)
1613  inout[mm + w_offset] = Fmnm_pad[mm - 2*(L-1) - 1 + w_offset];
1614  for (mm = -2*(L-1); mm <= 0; ++mm)
1615  inout[mm + w_offset] = Fmnm_pad[mm + 2*(L-1) + w_offset];
1616  // Compute IFFT of Fmnm'.
1617  fftw_execute_dft(plan_bwd, inout, inout);
1618  // Apply spatial shift.
1619  for (mm = 0; mm <= 2*L-2; ++mm)
1620  Fmnm_pad[mm + w_offset] = inout[mm - 2*(L-1) + w_offset];
1621  for (mm = -2*(L-1); mm <= -1; ++mm)
1622  Fmnm_pad[mm + w_offset] = inout[mm + 2*(L-1) + 1 + w_offset];
1623 
1624  // Compute product of Fmnm' and weight in real space.
1625  int r;
1626  for (r = -2*(L-1); r <= 2*(L-1); ++r)
1627  Fmnm_pad[r + w_offset] *= wr[r + w_offset];
1628 
1629  // Apply spatial shift.
1630  for (mm = 1; mm <= 2*L-2; ++mm)
1631  inout[mm + w_offset] = Fmnm_pad[mm - 2*(L-1) - 1 + w_offset];
1632  for (mm = -2*(L-1); mm <= 0; ++mm)
1633  inout[mm + w_offset] = Fmnm_pad[mm + 2*(L-1) + w_offset];
1634  // Compute Gmnm' by FFT.
1635  fftw_execute_dft(plan_fwd, inout, inout);
1636  // Apply spatial shift.
1637  for (mm = 0; mm <= 2*L-2; ++mm)
1638  Fmnm_pad[mm + w_offset] = inout[mm - 2*(L-1) + w_offset];
1639  for (mm = -2*(L-1); mm <= -1; ++mm)
1640  Fmnm_pad[mm + w_offset] = inout[mm + 2*(L-1) + 1 + w_offset];
1641 
1642  // Extract section of Gmnm' of interest.
1643  for (mm = -(L-1); mm <= L-1; ++mm)
1644  Gmnm[m + m_offset + m_stride*(
1645  mm + mm_offset + mm_stride*(
1646  n + n_offset))] =
1647  Fmnm_pad[mm + w_offset]
1648  * 4.0 * SSHT_PI * SSHT_PI / (4.0*L-3.0);
1649 
1650  }
1651  fftw_destroy_plan(plan_bwd);
1652  fftw_destroy_plan(plan_fwd);
1653 
1654  // Compute flmn.
1655  double *dl, *dl8 = NULL;
1656  dl = ssht_dl_calloc(L, SSHT_DL_QUARTER);
1657  SO3_ERROR_MEM_ALLOC_CHECK(dl);
1658  if (dl_method == SSHT_DL_RISBO)
1659  {
1660  dl8 = ssht_dl_calloc(L, SSHT_DL_QUARTER_EXTENDED);
1661  SO3_ERROR_MEM_ALLOC_CHECK(dl8);
1662  }
1663  int dl_offset = ssht_dl_get_offset(L, SSHT_DL_QUARTER);
1664  int dl_stride = ssht_dl_get_stride(L, SSHT_DL_QUARTER);
1665  for (n = -N+1; n <= N-1; ++n)
1666  for (el = abs(n); el < L; ++el)
1667  for (m = -el; m <= el; ++m)
1668  {
1669  int ind;
1670  so3_sampling_elmn2ind(&ind, el, m, n, parameters);
1671  flmn[ind] = 0.0;
1672  }
1673 
1674  for (el = L0; el < L; ++el)
1675  {
1676  int eltmp;
1677 
1678  // Compute Wigner plane.
1679  switch (dl_method)
1680  {
1681  case SSHT_DL_RISBO:
1682  if (el != 0 && el == L0)
1683  {
1684  for(eltmp = 0; eltmp <= L0; ++eltmp)
1685  ssht_dl_beta_risbo_eighth_table(dl8, SO3_PION2, L,
1686  SSHT_DL_QUARTER_EXTENDED,
1687  eltmp, sqrt_tbl, signs);
1688  ssht_dl_beta_risbo_fill_eighth2quarter_table(dl,
1689  dl8, L,
1690  SSHT_DL_QUARTER,
1691  SSHT_DL_QUARTER_EXTENDED,
1692  el,
1693  signs);
1694  }
1695  else
1696  {
1697  ssht_dl_beta_risbo_eighth_table(dl8, SO3_PION2, L,
1698  SSHT_DL_QUARTER_EXTENDED,
1699  el, sqrt_tbl, signs);
1700  ssht_dl_beta_risbo_fill_eighth2quarter_table(dl,
1701  dl8, L,
1702  SSHT_DL_QUARTER,
1703  SSHT_DL_QUARTER_EXTENDED,
1704  el,
1705  signs);
1706  }
1707  break;
1708 
1709  case SSHT_DL_TRAPANI:
1710  if (el != 0 && el == L0)
1711  {
1712  for(eltmp = 0; eltmp <= L0; ++eltmp)
1713  ssht_dl_halfpi_trapani_eighth_table(dl, L,
1714  SSHT_DL_QUARTER,
1715  eltmp, sqrt_tbl);
1716  ssht_dl_halfpi_trapani_fill_eighth2quarter_table(dl, L,
1717  SSHT_DL_QUARTER,
1718  el, signs);
1719  }
1720  else
1721  {
1722  ssht_dl_halfpi_trapani_eighth_table(dl, L,
1723  SSHT_DL_QUARTER,
1724  el, sqrt_tbl);
1725  ssht_dl_halfpi_trapani_fill_eighth2quarter_table(dl, L,
1726  SSHT_DL_QUARTER,
1727  el, signs);
1728  }
1729  break;
1730 
1731  default:
1732  SO3_ERROR_GENERIC("Invalid dl method");
1733  }
1734 
1735  // Compute flmn for current el.
1736 
1737  switch (n_mode)
1738  {
1739  case SO3_N_MODE_ALL:
1740  n_start = MAX(-N+1,-el);
1741  n_stop = MIN( N-1, el);
1742  n_inc = 1;
1743  break;
1744  case SO3_N_MODE_EVEN:
1745  n_start = MAX(-N+1,-el);
1746  n_start += (-n_start)%2;
1747  n_stop = MIN( N-1, el);
1748  n_stop -= n_stop%2;
1749  n_inc = 2;
1750  break;
1751  case SO3_N_MODE_ODD:
1752  n_start = MAX(-N+1,-el);
1753  n_start += 1+n_start%2;
1754  n_stop = MIN( N-1, el);
1755  n_stop -= 1-n_stop%2;
1756  n_inc = 2;
1757  break;
1758  case SO3_N_MODE_MAXIMUM:
1759  if (el < N-1)
1760  continue;
1761  n_start = -N+1;
1762  n_stop = N-1;
1763  n_inc = MAX(1,2*N-2);
1764  break;
1765  case SO3_N_MODE_L:
1766  if (el >= N)
1767  continue;
1768  n_start = -el;
1769  n_stop = el;
1770  n_inc = MAX(1,2*el);
1771  break;
1772  default:
1773  SO3_ERROR_GENERIC("Invalid n-mode.");
1774  }
1775 
1776  // TODO: Pull out a few multiplications into precomputations
1777  // or split up loops to avoid conditionals to check signs.
1778  for (mm = -el; mm <= el; ++mm)
1779  {
1780  // These signs are needed for the symmetry relations of
1781  // Wigner symbols.
1782  double elmmsign = signs[el] * signs[abs(mm)];
1783 
1784  for (n = n_start; n <= n_stop; n += n_inc)
1785  {
1786  double mmsign = mm >= 0 ? 1.0 : signs[el] * signs[abs(n)];
1787  double elnsign = n >= 0 ? 1.0 : elmmsign;
1788 
1789  // Factor which does not depend on m.
1790  double elnmm_factor = mmsign * elnsign
1791  * dl[abs(n) + dl_offset + abs(mm)*dl_stride];
1792 
1793  for (m = -el; m <= el; ++m)
1794  {
1795  mmsign = mm >= 0 ? 1.0 : signs[el] * signs[abs(m)];
1796  double elmsign = m >= 0 ? 1.0 : elmmsign;
1797  int ind;
1798  so3_sampling_elmn2ind(&ind, el, m, n, parameters);
1799  int mod = ((m-n)%4 + 4)%4;
1800  flmn[ind] +=
1801  exps[mod]
1802  * elnmm_factor
1803  * mmsign * elmsign
1804  * dl[abs(m) + dl_offset + abs(mm)*dl_stride]
1805  * Gmnm[m + m_offset + m_stride*(
1806  mm + mm_offset + mm_stride*(
1807  n + n_offset))];
1808 
1809  }
1810  }
1811  }
1812  }
1813 
1814  free(dl);
1815  if (dl_method == SSHT_DL_RISBO)
1816  free(dl8);
1817  free(Fmnb);
1818  free(Fmnm);
1819  free(inout);
1820  free(w);
1821  free(wr);
1822  free(Fmnm_pad);
1823  free(Gmnm);
1824  free(sqrt_tbl);
1825  free(signs);
1826  free(exps);
1827  free(expsmm);
1828 
1829  if (verbosity > 0)
1830  printf("%sForward transform computed!\n", SO3_PROMPT);
1831 
1832 }
1833 
1851  double *f, const complex double *flmn,
1852  const so3_parameters_t *parameters
1853 ) {
1854 
1855  int L0, L, N;
1856  so3_sampling_t sampling;
1857  so3_storage_t storage;
1858  so3_n_mode_t n_mode;
1859  ssht_dl_method_t dl_method;
1860  int steerable;
1861  int verbosity;
1862 
1863  L0 = parameters->L0;
1864  L = parameters->L;
1865  N = parameters->N;
1866  sampling = parameters->sampling_scheme;
1867  storage = parameters->storage;
1868  // TODO: Add optimisations for all n-modes.
1869  n_mode = parameters->n_mode;
1870  dl_method = parameters->dl_method;
1871  verbosity = parameters->verbosity;
1872  steerable = parameters->steerable;
1873 
1874  // Print messages depending on verbosity level.
1875  if (verbosity > 0)
1876  {
1877  printf("%sComputing inverse transform using MW sampling with\n", SO3_PROMPT);
1878  printf("%sparameters (L, N, reality) = (%d, %d, FALSE)\n", SO3_PROMPT, L, N);
1879  if (verbosity > 1)
1880  printf("%sUsing routine so3_core_mw_inverse_direct with storage method %d...\n"
1881  , SO3_PROMPT
1882  , storage);
1883  }
1884 
1885  // Iterators
1886  int el, m, n, mm; // mm for m'
1887 
1888  // Allocate memory.
1889  double *sqrt_tbl = calloc(2*(L-1)+2, sizeof(*sqrt_tbl));
1890  SO3_ERROR_MEM_ALLOC_CHECK(sqrt_tbl);
1891  double *signs = calloc(L+1, sizeof(*signs));
1892  SO3_ERROR_MEM_ALLOC_CHECK(signs);
1893  complex double *exps = calloc(4, sizeof(*exps));
1894  SO3_ERROR_MEM_ALLOC_CHECK(exps);
1895 
1896  // Perform precomputations.
1897  for (el = 0; el <= 2*L-1; ++el)
1898  sqrt_tbl[el] = sqrt((double)el);
1899  for (m = 0; m <= L-1; m += 2)
1900  {
1901  signs[m] = 1.0;
1902  signs[m+1] = -1.0;
1903  }
1904  int i;
1905  for (i = 0; i < 4; ++i)
1906  exps[i] = cexp(I*SO3_PION2*i);
1907 
1908  // Compute Fmnm'
1909  // TODO: Currently m is fastest-varying, then n, then m'.
1910  // Should this order be changed to m-m'-n?
1911  complex double *Fmnm = calloc((2*L-1)*(2*L-1)*N, sizeof(*Fmnm));
1912  SO3_ERROR_MEM_ALLOC_CHECK(Fmnm);
1913  int m_offset = L-1;
1914  int m_stride = 2*L-1;
1915  int n_offset = 0;
1916  int n_stride = N;
1917  int mm_offset = L-1;
1918  // unused: int mm_stride = 2*L-1;
1919 
1920  int n_start, n_stop, n_inc;
1921 
1922  double *dl = ssht_dl_calloc(L, SSHT_DL_QUARTER);
1923  SO3_ERROR_MEM_ALLOC_CHECK(dl);
1924  double *dl8 = NULL;
1925  if (dl_method == SSHT_DL_RISBO)
1926  {
1927  dl8 = ssht_dl_calloc(L, SSHT_DL_QUARTER_EXTENDED);
1928  SO3_ERROR_MEM_ALLOC_CHECK(dl8);
1929  }
1930  int dl_offset = ssht_dl_get_offset(L, SSHT_DL_QUARTER);
1931  int dl_stride = ssht_dl_get_stride(L, SSHT_DL_QUARTER);
1932 
1933  complex double *mn_factors = calloc((2*L-1)*N, sizeof *mn_factors);
1934  SO3_ERROR_MEM_ALLOC_CHECK(mn_factors);
1935 
1936  // TODO: SSHT starts this loop from MAX(L0, abs(spin)).
1937  // Can we use a similar optimisation? el can probably
1938  // be limited by n, but then we'd need to switch the
1939  // loop order, which means we'd have to recompute the
1940  // Wigner plane for each n. That seems wrong?
1941  for (el = L0; el <= L-1; ++el)
1942  {
1943  int eltmp;
1944  // Compute Wigner plane.
1945  switch (dl_method)
1946  {
1947  case SSHT_DL_RISBO:
1948  if (el != 0 && el == L0)
1949  {
1950  for(eltmp = 0; eltmp <= L0; ++eltmp)
1951  ssht_dl_beta_risbo_eighth_table(
1952  dl8,
1953  SO3_PION2,
1954  L,
1955  SSHT_DL_QUARTER_EXTENDED,
1956  eltmp,
1957  sqrt_tbl,
1958  signs);
1959  ssht_dl_beta_risbo_fill_eighth2quarter_table(dl,
1960  dl8, L,
1961  SSHT_DL_QUARTER,
1962  SSHT_DL_QUARTER_EXTENDED,
1963  el,
1964  signs);
1965  }
1966  else
1967  {
1968  ssht_dl_beta_risbo_eighth_table(
1969  dl8,
1970  SO3_PION2,
1971  L,
1972  SSHT_DL_QUARTER_EXTENDED,
1973  el,
1974  sqrt_tbl,
1975  signs);
1976  ssht_dl_beta_risbo_fill_eighth2quarter_table(
1977  dl,
1978  dl8, L,
1979  SSHT_DL_QUARTER,
1980  SSHT_DL_QUARTER_EXTENDED,
1981  el,
1982  signs);
1983  }
1984  break;
1985 
1986  case SSHT_DL_TRAPANI:
1987  if (el != 0 && el == L0)
1988  {
1989  for(eltmp = 0; eltmp <= L0; ++eltmp)
1990  ssht_dl_halfpi_trapani_eighth_table(
1991  dl,
1992  L,
1993  SSHT_DL_QUARTER,
1994  eltmp,
1995  sqrt_tbl);
1996  ssht_dl_halfpi_trapani_fill_eighth2quarter_table(
1997  dl,
1998  L,
1999  SSHT_DL_QUARTER,
2000  el,
2001  signs);
2002  }
2003  else
2004  {
2005  ssht_dl_halfpi_trapani_eighth_table(
2006  dl,
2007  L,
2008  SSHT_DL_QUARTER,
2009  el,
2010  sqrt_tbl);
2011  ssht_dl_halfpi_trapani_fill_eighth2quarter_table(
2012  dl,
2013  L,
2014  SSHT_DL_QUARTER,
2015  el,
2016  signs);
2017  }
2018  break;
2019 
2020  default:
2021  SO3_ERROR_GENERIC("Invalid dl method");
2022  }
2023 
2024  // Compute Fmnm' contribution for current el.
2025 
2026  // Factor which depends only on el.
2027  double elfactor = (2.0*el+1.0)/(8.0*SO3_PI*SO3_PI);
2028 
2029 
2030  switch (n_mode)
2031  {
2032  case SO3_N_MODE_ALL:
2033  n_start = 0;
2034  n_stop = MIN( N-1, el);
2035  n_inc = 1;
2036  break;
2037  case SO3_N_MODE_EVEN:
2038  n_start = 0;
2039  n_stop = MIN( N-1, el);
2040  n_stop -= n_stop%2;
2041  n_inc = 2;
2042  break;
2043  case SO3_N_MODE_ODD:
2044  n_start = 1;
2045  n_stop = MIN( N-1, el);
2046  n_stop -= 1-n_stop%2;
2047  n_inc = 2;
2048  break;
2049  case SO3_N_MODE_MAXIMUM:
2050  if (el < N-1)
2051  continue;
2052  n_start = N-1;
2053  n_stop = N-1;
2054  n_inc = 1;
2055  break;
2056  case SO3_N_MODE_L:
2057  if (el >= N)
2058  continue;
2059  n_start = el;
2060  n_stop = el;
2061  n_inc = 1;
2062  break;
2063  default:
2064  SO3_ERROR_GENERIC("Invalid n-mode.");
2065  }
2066 
2067  // Factors which do not depend on m'.
2068  for (n = n_start; n <= n_stop; n += n_inc)
2069  for (m = -el; m <= el; ++m)
2070  {
2071  int ind;
2072  so3_sampling_elmn2ind_real(&ind, el, m, n, parameters);
2073  int mod = ((n-m)%4 + 4)%4;
2074  mn_factors[m + m_offset + m_stride*(
2075  n + n_offset)] =
2076  flmn[ind] * exps[mod];
2077  }
2078 
2079  for (mm = 0; mm <= el; ++mm)
2080  {
2081  // These signs are needed for the symmetry relations of
2082  // Wigner symbols.
2083  double elmmsign = signs[el] * signs[mm];
2084 
2085  for (n = n_start; n <= n_stop; n += n_inc)
2086  {
2087  // Factor which does not depend on m.
2088  double elnmm_factor = elfactor
2089  * dl[n + dl_offset + mm*dl_stride];
2090  for (m = -el; m < 0; ++m)
2091  Fmnm[m + m_offset + m_stride*(
2092  n + n_offset + n_stride*(
2093  mm + mm_offset))] +=
2094  elnmm_factor
2095  * mn_factors[m + m_offset + m_stride*(
2096  n + n_offset)]
2097  * elmmsign * dl[-m + dl_offset + mm*dl_stride];
2098  for (m = 0; m <= el; ++m)
2099  Fmnm[m + m_offset + m_stride*(
2100  n + n_offset + n_stride*(
2101  mm + mm_offset))] +=
2102  elnmm_factor
2103  * mn_factors[m + m_offset + m_stride*(
2104  n + n_offset)]
2105  * dl[m + dl_offset + mm*dl_stride];
2106  }
2107  }
2108  }
2109 
2110  // Free dl memory.
2111  free(dl);
2112  if (dl_method == SSHT_DL_RISBO)
2113  free(dl8);
2114 
2115  switch (n_mode)
2116  {
2117  case SO3_N_MODE_ALL:
2118  case SO3_N_MODE_L:
2119  n_start = 0;
2120  n_stop = N-1;
2121  n_inc = 1;
2122  break;
2123  case SO3_N_MODE_EVEN:
2124  n_start = 0;
2125  n_stop = ((N-1) % 2 == 0) ? N-1 : N-2;
2126  n_inc = 2;
2127  break;
2128  case SO3_N_MODE_ODD:
2129  n_start = 1;
2130  n_stop = ((N-1) % 2 != 0) ? N-1 : N-2;
2131  n_inc = 2;
2132  break;
2133  case SO3_N_MODE_MAXIMUM:
2134  n_start = N-1;
2135  n_stop = N-1;
2136  n_inc = 1;
2137  break;
2138  default:
2139  SO3_ERROR_GENERIC("Invalid n-mode.");
2140  }
2141 
2142  // Use symmetry to compute Fmnm' for negative m'.
2143  for (mm = -L+1; mm < 0; ++mm)
2144  for (n = n_start; n <= n_stop; n += n_inc)
2145  for (m = -L+1; m <= L-1; ++m)
2146  Fmnm[m + m_offset + m_stride*(
2147  n + n_offset + n_stride*(
2148  mm + mm_offset))] =
2149  signs[abs(m+n)%2]
2150  * Fmnm[m + m_offset + m_stride*(
2151  n + n_offset + n_stride*(
2152  -mm + mm_offset))];
2153 
2154  // Apply phase modulation to account for sampling offset.
2155  for (mm = -L+1; mm <= L-1; ++mm)
2156  {
2157  complex double mmfactor = cexp(I*mm*SO3_PI/(2.0*L-1.0));
2158  for (n = n_start; n <= n_stop; n += n_inc)
2159  for (m = -L+1; m <= L-1; ++m)
2160  Fmnm[m + m_offset + m_stride*(
2161  n + n_offset + n_stride*(
2162  mm + mm_offset))] *= mmfactor;
2163  }
2164 
2165  // Allocate space for shifted Fmnm'.
2166  complex double *Fmnm_shift = calloc((2*L-1)*(2*L-1)*N, sizeof(*Fmnm_shift));
2167  SO3_ERROR_MEM_ALLOC_CHECK(Fmnm_shift);
2168 
2169  // Allocate space for function values.
2170  double *fext = calloc((2*L-1)*(2*L-1)*(2*N-1), sizeof(*fext));
2171  SO3_ERROR_MEM_ALLOC_CHECK(fext);
2172 
2173  // Set up plan before initialising array.
2174  // The redundant dimension needs to be the last one.
2175  fftw_plan plan = fftw_plan_dft_c2r_3d(
2176  2*L-1, 2*L-1, 2*N-1,
2177  Fmnm_shift, fext,
2178  FFTW_ESTIMATE);
2179 
2180  // Apply spatial shift.
2181  // This also reshapes the array to make n the inner dimension.
2182  for (mm = -L+1; mm <= L-1; ++mm)
2183  {
2184  int mm_shift = mm < 0 ? 2*L-1 : 0;
2185  for (n = n_start; n <= n_stop; n += n_inc)
2186  {
2187  for (m = -L+1; m <= L-1; ++m)
2188  {
2189  int m_shift = m < 0 ? 2*L-1 : 0;
2190  Fmnm_shift[n + n_stride*(
2191  m + m_shift + m_stride*(
2192  mm + mm_shift))] =
2193  Fmnm[m + m_offset + m_stride*(
2194  n + n_offset + n_stride*(
2195  mm + mm_offset))];
2196  }
2197  }
2198  }
2199 
2200  // Free Fmnm' memory.
2201  free(Fmnm);
2202 
2203  // Perform 3D FFT.
2204  fftw_execute(plan);
2205  fftw_destroy_plan(plan);
2206 
2207  free(Fmnm_shift);
2208 
2209  // Extract f from the extended torus.
2210  // Again, we reshape the array in the process.
2211  int a,b,g;
2212  int a_stride = 2*L-1;
2213  // unused: int b_ext_stride = 2*L-1;
2214  int b_stride = L;
2215  int g_stride = 2*N-1;
2216  for (g = 0; g < 2*N-1; ++g)
2217  for (b = 0; b < L; ++b)
2218  for (a = 0; a < 2*L-1; ++a)
2219  f[a + a_stride*(
2220  b + b_stride*(
2221  g))] = fext[g + g_stride*(
2222  a + a_stride*(
2223  b))];
2224 
2225  // Free fext memory.
2226  free(fext);
2227 
2228  if (verbosity > 0)
2229  printf("%sInverse transform computed!\n", SO3_PROMPT);
2230 
2231  // Free precomputation memory.
2232  free(sqrt_tbl);
2233  free(signs);
2234  free(exps);
2235  free(mn_factors);
2236 }
2237 
2256  complex double *flmn, const double *f,
2257  const so3_parameters_t *parameters
2258 ) {
2259  int L0, L, N;
2260  so3_sampling_t sampling;
2261  so3_storage_t storage;
2262  so3_n_mode_t n_mode;
2263  ssht_dl_method_t dl_method;
2264  int steerable;
2265  int verbosity;
2266 
2267  L0 = parameters->L0;
2268  L = parameters->L;
2269  N = parameters->N;
2270  sampling = parameters->sampling_scheme;
2271  storage = parameters->storage;
2272  // TODO: Add optimisations for all n-modes.
2273  n_mode = parameters->n_mode;
2274  dl_method = parameters->dl_method;
2275  verbosity = parameters->verbosity;
2276  steerable = parameters->steerable;
2277 
2278  // Print messages depending on verbosity level.
2279  if (verbosity > 0) {
2280  printf("%sComputing forward transform using MW sampling with\n", SO3_PROMPT);
2281  printf("%sparameters (L, N, reality) = (%d, %d, FALSE)\n", SO3_PROMPT, L, N);
2282  if (verbosity > 1)
2283  printf("%sUsing routine so3_core_mw_forward_direct with storage method %d...\n"
2284  , SO3_PROMPT
2285  , storage);
2286  }
2287 
2288  int m_stride = 2*L-1;
2289  int m_offset = L-1;
2290  int n_offset = 0;
2291  int n_stride = N;
2292  int mm_stride = 2*L-1;
2293  int mm_offset = L-1;
2294  int a_stride = 2*L-1;
2295  int b_stride = L;
2296  int bext_stride = 2*L-1;
2297  int g_stride = 2*N-1;
2298 
2299  int n_start, n_stop, n_inc;
2300 
2301  switch (n_mode)
2302  {
2303  case SO3_N_MODE_ALL:
2304  case SO3_N_MODE_L:
2305  n_start = 0;
2306  n_stop = N-1;
2307  n_inc = 1;
2308  break;
2309  case SO3_N_MODE_EVEN:
2310  n_start = 0;
2311  n_stop = ((N-1) % 2 == 0) ? N-1 : N-2;
2312  n_inc = 2;
2313  break;
2314  case SO3_N_MODE_ODD:
2315  n_start = 1;
2316  n_stop = ((N-1) % 2 != 0) ? N-1 : N-2;
2317  n_inc = 2;
2318  break;
2319  case SO3_N_MODE_MAXIMUM:
2320  n_start = N-1;
2321  n_stop = N-1;
2322  n_inc = 1;
2323  break;
2324  default:
2325  SO3_ERROR_GENERIC("Invalid n-mode.");
2326  }
2327 
2328  double *sqrt_tbl = calloc(2*(L-1)+2, sizeof(*sqrt_tbl));
2329  SO3_ERROR_MEM_ALLOC_CHECK(sqrt_tbl);
2330  double *signs = calloc(L+1, sizeof(*signs));
2331  SO3_ERROR_MEM_ALLOC_CHECK(signs);
2332  complex double *exps = calloc(4, sizeof(*exps));
2333  SO3_ERROR_MEM_ALLOC_CHECK(exps);
2334  complex double *expsmm = calloc(2*L-1, sizeof(*expsmm));
2335  SO3_ERROR_MEM_ALLOC_CHECK(expsmm);
2336 
2337  int el, m, n, mm; // mm is for m'
2338  // Perform precomputations.
2339  for (el = 0; el <= 2*(L-1)+1; ++el)
2340  sqrt_tbl[el] = sqrt((double)el);
2341  for (m = 0; m <= L-1; m += 2)
2342  {
2343  signs[m] = 1.0;
2344  signs[m+1] = -1.0;
2345  }
2346  int i;
2347  for (i = 0; i < 4; ++i)
2348  exps[i] = cexp(I*SO3_PION2*i);
2349  for (mm = -L+1; mm <= L-1; ++mm)
2350  expsmm[mm + mm_offset] = cexp(-I*mm*SSHT_PI/(2.0*L-1.0));
2351 
2352  double norm_factor = 1.0/(2.0*L-1.0)/(2.0*N-1.0);
2353 
2354  // Compute Fourier transform over alpha and gamma, i.e. compute Fmn(b).
2355  complex double *Fmnb = calloc((2*L-1)*(2*L-1)*N, sizeof(*Fmnb));
2356  SO3_ERROR_MEM_ALLOC_CHECK(Fmnb);
2357  double *fft_in = calloc((2*L-1)*(2*N-1), sizeof(*fft_in));
2358  SO3_ERROR_MEM_ALLOC_CHECK(fft_in);
2359  complex double *fft_out = calloc((2*L-1)*N, sizeof(*fft_out));
2360  SO3_ERROR_MEM_ALLOC_CHECK(fft_out);
2361  // Redundant dimension needs to be last
2362  fftw_plan plan = fftw_plan_dft_r2c_2d(
2363  2*L-1, 2*N-1,
2364  fft_in, fft_out,
2365  FFTW_ESTIMATE);
2366 
2367  int a, b, g;
2368  for (b = 0; b < L; ++b)
2369  {
2370  // TODO: This loop could probably be avoided by using
2371  // a more elaborate FFTW plan which performs the FFT directly
2372  // over the 1st and 3rd dimensions of f.
2373  // Instead, for each index in the 2nd dimension, we copy the
2374  // corresponding values in the 1st and 3rd dimension into a
2375  // new 2D array, to perform a standard 2D FFT there. While
2376  // we're at it, we also reshape that array such that gamma
2377  // is the inner dimension, as required by FFTW.
2378  for (a = 0; a < 2*L-1; ++a)
2379  for (g = 0; g < 2*N-1; ++g)
2380  fft_in[g + g_stride*(
2381  a)] =
2382  f[a + a_stride*(
2383  b + b_stride*(
2384  g))];
2385 
2386  fftw_execute(plan);
2387 
2388  // Apply spatial shift and normalisation factor, while
2389  // reshaping the dimensions once more.
2390  for (n = n_start; n <= n_stop; n += n_inc)
2391  {
2392  for (m = -L+1; m <= L-1; ++m)
2393  {
2394  int m_shift = m < 0 ? 2*L-1 : 0;
2395  Fmnb[b + bext_stride*(
2396  m + m_offset + m_stride*(
2397  n + n_offset))] =
2398  fft_out[n + n_stride*(
2399  m + m_shift)] * norm_factor;
2400  }
2401  }
2402  }
2403  fftw_destroy_plan(plan);
2404 
2405  // Extend Fmnb periodically.
2406  for (n = n_start; n <= n_stop; n += n_inc)
2407  for (m = -L+1; m <= L-1; ++m)
2408  {
2409  int signmn = signs[abs(m+n)%2];
2410  for (b = L; b < 2*L-1; ++b)
2411  Fmnb[b + bext_stride*(
2412  m + m_offset + m_stride*(
2413  n + n_offset))] =
2414  signmn
2415  * Fmnb[(2*L-2-b) + bext_stride*(
2416  m + m_offset + m_stride*(
2417  n + n_offset))];
2418  }
2419 
2420 
2421  // Compute Fourier transform over beta, i.e. compute Fmnm'.
2422  complex double *Fmnm = calloc((2*L-1)*(2*L-1)*(2*N-1), sizeof(*Fmnm));
2423  SO3_ERROR_MEM_ALLOC_CHECK(Fmnm);
2424  complex double *inout = calloc(2*L-1, sizeof(*inout));
2425  SO3_ERROR_MEM_ALLOC_CHECK(inout);
2426 
2427  plan = fftw_plan_dft_1d(
2428  2*L-1,
2429  inout, inout,
2430  FFTW_FORWARD,
2431  FFTW_ESTIMATE);
2432  for (n = n_start; n <= n_stop; n += n_inc)
2433  for (m = -L+1; m <= L-1; ++m)
2434  {
2435  memcpy(inout,
2436  Fmnb + 0 + bext_stride*(
2437  m + m_offset + m_stride*(
2438  n + n_offset)),
2439  bext_stride*sizeof(*Fmnb));
2440  fftw_execute(plan);
2441 
2442  // Apply spatial shift and normalisation factor
2443  for (mm = -L+1; mm <= L-1; ++mm)
2444  {
2445  int mm_shift = mm < 0 ? 2*L-1 : 0;
2446  Fmnm[mm + mm_offset + mm_stride*(
2447  m + m_offset + m_stride*(
2448  n + n_offset))] =
2449  inout[mm + mm_shift] / (2.0*L-1.0);
2450  }
2451  }
2452  fftw_destroy_plan(plan);
2453  free(inout);
2454 
2455  // Apply phase modulation to account for sampling offset.
2456  for (n = n_start; n <= n_stop; n += n_inc)
2457  for (m = -L+1; m <= L-1; ++m)
2458  for (mm = -L+1; mm <= L-1; ++mm)
2459  Fmnm[mm + mm_offset + mm_stride*(
2460  m + m_offset + m_stride*(
2461  n + n_offset))] *=
2462  expsmm[mm + mm_offset];
2463 
2464  // Compute weights.
2465  complex double *w = calloc(4*L-3, sizeof(*w));
2466  SO3_ERROR_MEM_ALLOC_CHECK(w);
2467  int w_offset = 2*(L-1);
2468  for (mm = -2*(L-1); mm <= 2*(L-1); ++mm)
2469  w[mm+w_offset] = so3_sampling_weight(parameters, mm);
2470 
2471  // Compute IFFT of w to give wr.
2472  complex double *wr = calloc(4*L-3, sizeof(*w));
2473  SO3_ERROR_MEM_ALLOC_CHECK(wr);
2474  inout = calloc(4*L-3, sizeof(*inout));
2475  SO3_ERROR_MEM_ALLOC_CHECK(inout);
2476  fftw_plan plan_bwd = fftw_plan_dft_1d(
2477  4*L-3,
2478  inout, inout,
2479  FFTW_BACKWARD,
2480  FFTW_MEASURE);
2481  fftw_plan plan_fwd = fftw_plan_dft_1d(
2482  4*L-3,
2483  inout,
2484  inout,
2485  FFTW_FORWARD,
2486  FFTW_MEASURE);
2487 
2488  // Apply spatial shift.
2489  for (mm = 1; mm <= 2*L-2; ++mm)
2490  inout[mm + w_offset] = w[mm - 2*(L-1) - 1 + w_offset];
2491  for (mm = -2*(L-1); mm <= 0; ++mm)
2492  inout[mm + w_offset] = w[mm + 2*(L-1) + w_offset];
2493 
2494  fftw_execute_dft(plan_bwd, inout, inout);
2495 
2496  // Apply spatial shift.
2497  for (mm = 0; mm <= 2*L-2; ++mm)
2498  wr[mm + w_offset] = inout[mm - 2*(L-1) + w_offset];
2499  for (mm = -2*(L-1); mm <= -1; ++mm)
2500  wr[mm + w_offset] = inout[mm + 2*(L-1) + 1 + w_offset];
2501 
2502  // Compute Gmnm' by convolution implemented as product in real space.
2503  complex double *Fmnm_pad = calloc(4*L-3, sizeof(*Fmnm_pad));
2504  SO3_ERROR_MEM_ALLOC_CHECK(Fmnm_pad);
2505  complex double *Gmnm = calloc((2*L-1)*(2*L-1)*N, sizeof(*Gmnm));
2506  SO3_ERROR_MEM_ALLOC_CHECK(Gmnm);
2507  for (n = n_start; n <= n_stop; n += n_inc)
2508  for (m = -L+1; m <= L-1; ++m)
2509  {
2510 
2511  // Zero-pad Fmnm'.
2512  for (mm = -2*(L-1); mm <= -L; ++mm)
2513  Fmnm_pad[mm+w_offset] = 0.0;
2514  for (mm = L; mm <= 2*(L-1); ++mm)
2515  Fmnm_pad[mm+w_offset] = 0.0;
2516  for (mm = -(L-1); mm <= L-1; ++mm)
2517  Fmnm_pad[mm + w_offset] =
2518  Fmnm[mm + mm_offset + mm_stride*(
2519  m + m_offset + m_stride*(
2520  n + n_offset))];
2521 
2522  // Apply spatial shift.
2523  for (mm = 1; mm <= 2*L-2; ++mm)
2524  inout[mm + w_offset] = Fmnm_pad[mm - 2*(L-1) - 1 + w_offset];
2525  for (mm = -2*(L-1); mm <= 0; ++mm)
2526  inout[mm + w_offset] = Fmnm_pad[mm + 2*(L-1) + w_offset];
2527  // Compute IFFT of Fmnm'.
2528  fftw_execute_dft(plan_bwd, inout, inout);
2529  // Apply spatial shift.
2530  for (mm = 0; mm <= 2*L-2; ++mm)
2531  Fmnm_pad[mm + w_offset] = inout[mm - 2*(L-1) + w_offset];
2532  for (mm = -2*(L-1); mm <= -1; ++mm)
2533  Fmnm_pad[mm + w_offset] = inout[mm + 2*(L-1) + 1 + w_offset];
2534 
2535  // Compute product of Fmnm' and weight in real space.
2536  int r;
2537  for (r = -2*(L-1); r <= 2*(L-1); ++r)
2538  Fmnm_pad[r + w_offset] *= wr[r + w_offset];
2539 
2540  // Apply spatial shift.
2541  for (mm = 1; mm <= 2*L-2; ++mm)
2542  inout[mm + w_offset] = Fmnm_pad[mm - 2*(L-1) - 1 + w_offset];
2543  for (mm = -2*(L-1); mm <= 0; ++mm)
2544  inout[mm + w_offset] = Fmnm_pad[mm + 2*(L-1) + w_offset];
2545  // Compute Gmnm' by FFT.
2546  fftw_execute_dft(plan_fwd, inout, inout);
2547  // Apply spatial shift.
2548  for (mm = 0; mm <= 2*L-2; ++mm)
2549  Fmnm_pad[mm + w_offset] = inout[mm - 2*(L-1) + w_offset];
2550  for (mm = -2*(L-1); mm <= -1; ++mm)
2551  Fmnm_pad[mm + w_offset] = inout[mm + 2*(L-1) + 1 + w_offset];
2552 
2553  // Extract section of Gmnm' of interest.
2554  for (mm = -(L-1); mm <= L-1; ++mm)
2555  Gmnm[m + m_offset + m_stride*(
2556  mm + mm_offset + mm_stride*(
2557  n + n_offset))] =
2558  Fmnm_pad[mm + w_offset]
2559  * 4.0 * SSHT_PI * SSHT_PI / (4.0*L-3.0);
2560 
2561  }
2562  fftw_destroy_plan(plan_bwd);
2563  fftw_destroy_plan(plan_fwd);
2564 
2565  // Compute flmn.
2566  double *dl, *dl8 = NULL;
2567  dl = ssht_dl_calloc(L, SSHT_DL_QUARTER);
2568  SO3_ERROR_MEM_ALLOC_CHECK(dl);
2569  if (dl_method == SSHT_DL_RISBO)
2570  {
2571  dl8 = ssht_dl_calloc(L, SSHT_DL_QUARTER_EXTENDED);
2572  SO3_ERROR_MEM_ALLOC_CHECK(dl8);
2573  }
2574  int dl_offset = ssht_dl_get_offset(L, SSHT_DL_QUARTER);
2575  int dl_stride = ssht_dl_get_stride(L, SSHT_DL_QUARTER);
2576  for (n = 0; n <= N-1; ++n)
2577  for (el = n; el < L; ++el)
2578  for (m = -el; m <= el; ++m)
2579  {
2580  int ind;
2581  so3_sampling_elmn2ind_real(&ind, el, m, n, parameters);
2582  flmn[ind] = 0.0;
2583  }
2584 
2585  for (el = L0; el < L; ++el)
2586  {
2587  int eltmp;
2588 
2589  // Compute Wigner plane.
2590  switch (dl_method)
2591  {
2592  case SSHT_DL_RISBO:
2593  if (el != 0 && el == L0)
2594  {
2595  for(eltmp = 0; eltmp <= L0; ++eltmp)
2596  ssht_dl_beta_risbo_eighth_table(dl8, SO3_PION2, L,
2597  SSHT_DL_QUARTER_EXTENDED,
2598  eltmp, sqrt_tbl, signs);
2599  ssht_dl_beta_risbo_fill_eighth2quarter_table(dl,
2600  dl8, L,
2601  SSHT_DL_QUARTER,
2602  SSHT_DL_QUARTER_EXTENDED,
2603  el,
2604  signs);
2605  }
2606  else
2607  {
2608  ssht_dl_beta_risbo_eighth_table(dl8, SO3_PION2, L,
2609  SSHT_DL_QUARTER_EXTENDED,
2610  el, sqrt_tbl, signs);
2611  ssht_dl_beta_risbo_fill_eighth2quarter_table(dl,
2612  dl8, L,
2613  SSHT_DL_QUARTER,
2614  SSHT_DL_QUARTER_EXTENDED,
2615  el,
2616  signs);
2617  }
2618  break;
2619 
2620  case SSHT_DL_TRAPANI:
2621  if (el != 0 && el == L0)
2622  {
2623  for(eltmp = 0; eltmp <= L0; ++eltmp)
2624  ssht_dl_halfpi_trapani_eighth_table(dl, L,
2625  SSHT_DL_QUARTER,
2626  eltmp, sqrt_tbl);
2627  ssht_dl_halfpi_trapani_fill_eighth2quarter_table(dl, L,
2628  SSHT_DL_QUARTER,
2629  el, signs);
2630  }
2631  else
2632  {
2633  ssht_dl_halfpi_trapani_eighth_table(dl, L,
2634  SSHT_DL_QUARTER,
2635  el, sqrt_tbl);
2636  ssht_dl_halfpi_trapani_fill_eighth2quarter_table(dl, L,
2637  SSHT_DL_QUARTER,
2638  el, signs);
2639  }
2640  break;
2641 
2642  default:
2643  SO3_ERROR_GENERIC("Invalid dl method");
2644  }
2645 
2646  // Compute flmn for current el.
2647 
2648 
2649  switch (n_mode)
2650  {
2651  case SO3_N_MODE_ALL:
2652  n_start = 0;
2653  n_stop = MIN( N-1, el);
2654  n_inc = 1;
2655  break;
2656  case SO3_N_MODE_EVEN:
2657  n_start = 0;
2658  n_stop = MIN( N-1, el);
2659  n_stop -= n_stop%2;
2660  n_inc = 2;
2661  break;
2662  case SO3_N_MODE_ODD:
2663  n_start = 1;
2664  n_stop = MIN( N-1, el);
2665  n_stop -= 1-n_stop%2;
2666  n_inc = 2;
2667  break;
2668  case SO3_N_MODE_MAXIMUM:
2669  if (el < N-1)
2670  continue;
2671  n_start = N-1;
2672  n_stop = N-1;
2673  n_inc = 1;
2674  break;
2675  case SO3_N_MODE_L:
2676  if (el >= N)
2677  continue;
2678  n_start = el;
2679  n_stop = el;
2680  n_inc = 1;
2681  break;
2682  default:
2683  SO3_ERROR_GENERIC("Invalid n-mode.");
2684  }
2685 
2686  // TODO: Pull out a few multiplications into precomputations
2687  // or split up loops to avoid conditionals to check signs.
2688  for (mm = -el; mm <= el; ++mm)
2689  {
2690  // These signs are needed for the symmetry relations of
2691  // Wigner symbols.
2692  double elmmsign = signs[el] * signs[abs(mm)];
2693 
2694  for (n = n_start; n <= n_stop; n += n_inc)
2695  {
2696  double mmsign = mm >= 0 ? 1.0 : signs[el] * signs[n];
2697 
2698  // Factor which does not depend on m.
2699  double elnmm_factor = mmsign
2700  * dl[n + dl_offset + abs(mm)*dl_stride];
2701 
2702  for (m = -el; m <= el; ++m)
2703  {
2704  mmsign = mm >= 0 ? 1.0 : signs[el] * signs[abs(m)];
2705  double elmsign = m >= 0 ? 1.0 : elmmsign;
2706  int ind;
2707  so3_sampling_elmn2ind_real(&ind, el, m, n, parameters);
2708  int mod = ((m-n)%4 + 4)%4;
2709  flmn[ind] +=
2710  exps[mod]
2711  * elnmm_factor
2712  * mmsign * elmsign
2713  * dl[abs(m) + dl_offset + abs(mm)*dl_stride]
2714  * Gmnm[m + m_offset + m_stride*(
2715  mm + mm_offset + mm_stride*(
2716  n + n_offset))];
2717 
2718  }
2719  }
2720  }
2721  }
2722 
2723  free(dl);
2724  if (dl_method == SSHT_DL_RISBO)
2725  free(dl8);
2726  free(Fmnb);
2727  free(Fmnm);
2728  free(inout);
2729  free(w);
2730  free(wr);
2731  free(Fmnm_pad);
2732  free(Gmnm);
2733  free(sqrt_tbl);
2734  free(signs);
2735  free(exps);
2736  free(expsmm);
2737 
2738  if (verbosity > 0)
2739  printf("%sForward transform computed!\n", SO3_PROMPT);
2740 }
so3_core_inverse_direct
void so3_core_inverse_direct(complex double *f, const complex double *flmn, const so3_parameters_t *parameters)
Definition: so3_core.c:951
MAX
#define MAX(a, b)
Definition: so3_core.c:27
so3_sampling_elmn2ind_real
void so3_sampling_elmn2ind_real(int *ind, int el, int m, int n, const so3_parameters_t *parameters)
Definition: so3_sampling.c:569
so3_core_forward_direct_real
void so3_core_forward_direct_real(complex double *flmn, const double *f, const so3_parameters_t *parameters)
Definition: so3_core.c:2255
so3_core_inverse_via_ssht
void so3_core_inverse_via_ssht(complex double *f, const complex double *flmn, const so3_parameters_t *parameters)
Definition: so3_core.c:48
so3_core_inverse_direct_real
void so3_core_inverse_direct_real(double *f, const complex double *flmn, const so3_parameters_t *parameters)
Definition: so3_core.c:1850
so3_core_forward_via_ssht_real
void so3_core_forward_via_ssht_real(complex double *flmn, const double *f, const so3_parameters_t *parameters)
Definition: so3_core.c:700
so3_core_forward_via_ssht
void so3_core_forward_via_ssht(complex double *flmn, const complex double *f, const so3_parameters_t *parameters)
Definition: so3_core.c:253
forward_real_ssht
void(* forward_real_ssht)(complex double *, const double *, int, int, ssht_dl_method_t, int)
Definition: so3_core.c:32
inverse_complex_ssht
void(* inverse_complex_ssht)(complex double *, const complex double *, int, int, int, ssht_dl_method_t, int)
Definition: so3_core.c:29
so3_core_inverse_via_ssht_real
void so3_core_inverse_via_ssht_real(double *f, const complex double *flmn, const so3_parameters_t *parameters)
Definition: so3_core.c:475
MIN
#define MIN(a, b)
Definition: so3_core.c:26
so3_core_forward_direct
void so3_core_forward_direct(complex double *flmn, const complex double *f, const so3_parameters_t *parameters)
Definition: so3_core.c:1354
so3_sampling_weight
complex double so3_sampling_weight(const so3_parameters_t *parameters, int p)
Definition: so3_sampling.c:30
forward_complex_ssht
void(* forward_complex_ssht)(complex double *, const complex double *, int, int, int, ssht_dl_method_t, int)
Definition: so3_core.c:31
inverse_real_ssht
void(* inverse_real_ssht)(double *, const complex double *, int, int, ssht_dl_method_t, int)
Definition: so3_core.c:30
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