20 #include "ssht/ssht.h"
22 #include "so3_types.h"
23 #include "so3_error.h"
24 #include "so3_sampling.h"
26 #define MIN(a,b) ((a < b) ? (a) : (b))
27 #define MAX(a,b) ((a > b) ? (a) : (b))
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);
49 complex
double *f,
const complex
double *flmn,
50 const so3_parameters_t *parameters
54 so3_sampling_t sampling;
55 so3_storage_t storage;
57 ssht_dl_method_t dl_method;
64 complex
double *fn, *ftemp;
68 int fftw_rank, fftw_howmany;
69 int fftw_idist, fftw_odist;
70 int fftw_istride, fftw_ostride;
72 complex
double *fftw_target;
80 sampling = parameters->sampling_scheme;
81 storage = parameters->storage;
83 n_mode = parameters->n_mode;
84 dl_method = parameters->dl_method;
85 verbosity = parameters->verbosity;
86 steerable = parameters->steerable;
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);
93 printf(
"%sUsing routine so3_core_mw_inverse_via_ssht with storage method %d...\n"
100 case SO3_SAMPLING_MW:
101 fn_n_stride = L * (2*L-1);
102 ssht = ssht_core_mw_lb_inverse_sov_sym;
104 case SO3_SAMPLING_MW_SS:
105 fn_n_stride = (L+1) * 2*L;
106 ssht = ssht_core_mw_lb_inverse_sov_sym_ss;
109 SO3_ERROR_GENERIC(
"Invalid sampling scheme.");
122 ftemp = malloc(2*N*fn_n_stride *
sizeof *ftemp);
123 SO3_ERROR_MEM_ALLOC_CHECK(ftemp);
134 fn = calloc(fftw_n*fn_n_stride,
sizeof *fn);
135 SO3_ERROR_MEM_ALLOC_CHECK(fn);
142 fftw_howmany = fn_n_stride;
145 fftw_idist = fftw_odist = 1;
146 fftw_istride = fftw_ostride = fn_n_stride;
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
155 for(n = -N+1; n <= N-1; ++n)
157 int ind, offset, i, el;
158 int L0e =
MAX(L0, abs(n));
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)
169 flm = malloc(L*L *
sizeof *flm);
170 SO3_ERROR_MEM_ALLOC_CHECK(flm);
174 case SO3_STORAGE_PADDED:
176 memcpy(flm, flmn + ind, L*L *
sizeof(complex
double));
178 case SO3_STORAGE_COMPACT:
180 memcpy(flm + n*n, flmn + ind, (L*L - n*n) *
sizeof(complex
double));
181 for(i = 0; i < n*n; ++i)
185 SO3_ERROR_GENERIC(
"Invalid storage method.");
192 factor = sqrt((
double)(2*el+1)/(16.*pow(SO3_PI, 3.)));
193 for (; i < offset + 2*el+1; ++i)
201 offset = (n < 0 ? n + fftw_n : n);
204 fn + offset*fn_n_stride, flm,
211 for(i = 0; i < fn_n_stride; ++i)
212 fn[offset*fn_n_stride + i] = -fn[offset*fn_n_stride + i];
222 fftw_destroy_plan(plan);
226 memcpy(f, ftemp, N*fn_n_stride *
sizeof(complex
double));
233 printf(
"%sInverse transform computed!\n", SO3_PROMPT);
254 complex
double *flmn,
const complex
double *f,
255 const so3_parameters_t *parameters
258 so3_sampling_t sampling;
259 so3_storage_t storage;
261 ssht_dl_method_t dl_method;
268 complex
double *ftemp, *fn;
272 int fftw_rank, fftw_howmany;
273 int fftw_idist, fftw_odist;
274 int fftw_istride, fftw_ostride;
286 sampling = parameters->sampling_scheme;
287 storage = parameters->storage;
289 n_mode = parameters->n_mode;
290 dl_method = parameters->dl_method;
291 verbosity = parameters->verbosity;
292 steerable = parameters->steerable;
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);
299 printf(
"%sUsing routine so3_core_mw_forward_via_ssht with storage method %d...\n"
306 case SO3_SAMPLING_MW:
307 fn_n_stride = L * (2*L-1);
308 ssht = ssht_core_mw_lb_forward_sov_conv_sym;
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;
315 SO3_ERROR_GENERIC(
"Invalid sampling scheme.");
322 fn = calloc((2*N-1)*fn_n_stride,
sizeof *fn);
323 SO3_ERROR_MEM_ALLOC_CHECK(fn);
325 for (n = -N+1; n < N; n+=2)
329 offset = (n < 0 ? n + 2*N-1 : n);
331 for (g = 0; g < N; ++g)
333 double gamma = g * SO3_PI / N;
334 for (i = 0; i < fn_n_stride; ++i)
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);
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));
351 fn = malloc((2*N-1)*fn_n_stride *
sizeof *fn);
352 SO3_ERROR_MEM_ALLOC_CHECK(fn);
358 fftw_howmany = fn_n_stride;
359 fftw_idist = fftw_odist = 1;
360 fftw_istride = fftw_ostride = fn_n_stride;
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
370 fftw_destroy_plan(plan);
374 factor = 2*SO3_PI/(double)(2*N-1);
375 for(i = 0; i < (2*N-1)*fn_n_stride; ++i)
379 for(n = -N+1; n <= N-1; ++n)
381 int ind, offset, el, sign;
382 int L0e =
MAX(L0, abs(n));
384 complex
double *flm = NULL;
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)
393 if (storage == SO3_STORAGE_COMPACT)
394 flm = malloc(L*L *
sizeof *flm);
398 offset = (n < 0 ? n + 2*N-1 : n);
400 complex
double *flm_block;
401 complex
double *fn_block = fn + offset*fn_n_stride;
406 case SO3_STORAGE_PADDED:
408 flm_block = flmn + ind;
411 case SO3_STORAGE_COMPACT:
413 i = offset = el*el-n*n;
416 SO3_ERROR_GENERIC(
"Invalid storage method.");
426 if (storage == SO3_STORAGE_COMPACT)
429 memcpy(flmn + ind, flm + n*n, (L*L - n*n) *
sizeof(complex
double));
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;
446 if (storage == SO3_STORAGE_COMPACT)
456 printf(
"%sForward transform computed!\n", SO3_PROMPT);
476 double *f,
const complex
double *flmn,
477 const so3_parameters_t *parameters
480 so3_sampling_t sampling;
481 so3_storage_t storage;
483 ssht_dl_method_t dl_method;
490 complex
double *fn, *flm;
495 int fftw_rank, fftw_howmany;
496 int fftw_idist, fftw_odist;
497 int fftw_istride, fftw_ostride;
508 sampling = parameters->sampling_scheme;
509 storage = parameters->storage;
511 n_mode = parameters->n_mode;
512 dl_method = parameters->dl_method;
513 verbosity = parameters->verbosity;
514 steerable = parameters->steerable;
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);
521 printf(
"%sUsing routine so3_core_mw_inverse_via_ssht_real with storage method %d...\n"
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;
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;
539 SO3_ERROR_GENERIC(
"Invalid sampling scheme.");
568 fn = calloc((fftw_n/2+1)*fn_n_stride,
sizeof *fn);
569 SO3_ERROR_MEM_ALLOC_CHECK(fn);
574 fftw_howmany = fn_n_stride;
577 fftw_idist = fftw_odist = 1;
578 fftw_istride = fftw_ostride = fn_n_stride;
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,
587 flm = malloc(L*L *
sizeof *flm);
588 SO3_ERROR_MEM_ALLOC_CHECK(flm);
590 for(n = 0; n <= N-1; ++n)
592 int ind, offset, i, el;
593 int L0e =
MAX(L0, abs(n));
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)
605 case SO3_STORAGE_PADDED:
607 memcpy(flm, flmn + ind, L*L *
sizeof(complex
double));
609 case SO3_STORAGE_COMPACT:
611 memcpy(flm + n*n, flmn + ind, (L*L - n*n) *
sizeof(complex
double));
612 for(i = 0; i < n*n; ++i)
616 SO3_ERROR_GENERIC(
"Invalid storage method.");
623 factor = sqrt((
double)(2*el+1)/(16.*pow(SO3_PI, 3.)));
624 for (; i < offset + 2*el+1; ++i)
633 fn + n*fn_n_stride, flm,
644 fn_r = calloc(fn_n_stride,
sizeof *fn_r);
645 SO3_ERROR_MEM_ALLOC_CHECK(fn_r);
654 for (i = 0; i < fn_n_stride; ++i)
660 for(i = 0; i < fn_n_stride; ++i)
661 fn[n*fn_n_stride + i] = -fn[n*fn_n_stride + i];
670 fftw_destroy_plan(plan);
681 printf(
"%sInverse transform computed!\n", SO3_PROMPT);
701 complex
double *flmn,
const double *f,
702 const so3_parameters_t *parameters
705 so3_sampling_t sampling;
706 so3_storage_t storage;
708 ssht_dl_method_t dl_method;
716 complex
double *flm = NULL, *fn;
720 int fftw_rank, fftw_howmany;
721 int fftw_idist, fftw_odist;
722 int fftw_istride, fftw_ostride;
735 sampling = parameters->sampling_scheme;
736 storage = parameters->storage;
738 n_mode = parameters->n_mode;
739 dl_method = parameters->dl_method;
740 steerable = parameters->steerable;
741 verbosity = parameters->verbosity;
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);
748 printf(
"%sUsing routine so3_core_mw_forward_via_ssht_real with storage method %d...\n"
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;
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;
766 SO3_ERROR_GENERIC(
"Invalid sampling scheme.");
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));
803 fn = malloc(N*fn_n_stride *
sizeof *fn);
804 SO3_ERROR_MEM_ALLOC_CHECK(fn);
809 fftw_howmany = fn_n_stride;
812 fftw_idist = fftw_odist = 1;
813 fftw_istride = fftw_ostride = fn_n_stride;
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,
823 fftw_destroy_plan(plan);
827 factor = 2*SO3_PI/(double)(2*N-1);
828 for(i = 0; i < N*fn_n_stride; ++i)
832 if (storage == SO3_STORAGE_COMPACT)
833 flm = malloc(L*L *
sizeof *flm);
835 for(n = 0; n <= N-1; ++n)
837 int ind, offset, el, sign;
838 int L0e =
MAX(L0, abs(n));
840 complex
double* flm_block;
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)
853 case SO3_STORAGE_PADDED:
855 flm_block = flmn + ind;
857 case SO3_STORAGE_COMPACT:
864 SO3_ERROR_GENERIC(
"Invalid storage method.");
871 flm_block, fn + n*fn_n_stride,
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]);
901 if (storage == SO3_STORAGE_COMPACT)
904 memcpy(flmn + ind, flm + n*n, (L*L - n*n) *
sizeof(complex
double));
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;
925 if (storage == SO3_STORAGE_COMPACT)
931 printf(
"%sForward transform computed!\n", SO3_PROMPT);
952 complex
double *f,
const complex
double *flmn,
953 const so3_parameters_t *parameters
957 so3_sampling_t sampling;
958 so3_storage_t storage;
960 ssht_dl_method_t dl_method;
967 sampling = parameters->sampling_scheme;
968 storage = parameters->storage;
970 n_mode = parameters->n_mode;
971 dl_method = parameters->dl_method;
972 verbosity = parameters->verbosity;
973 steerable = parameters->steerable;
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);
981 printf(
"%sUsing routine so3_core_mw_inverse_direct with storage method %d...\n"
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);
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)
1006 for (i = 0; i < 4; ++i)
1007 exps[i] = cexp(I*SO3_PION2*i);
1012 complex
double *Fmnm = calloc((2*L-1)*(2*L-1)*(2*N-1),
sizeof(*Fmnm));
1013 SO3_ERROR_MEM_ALLOC_CHECK(Fmnm);
1015 int m_stride = 2*L-1;
1017 int n_stride = 2*N-1;
1018 int mm_offset = L-1;
1019 int mm_stride = 2*L-1;
1021 int n_start, n_stop, n_inc;
1023 double *dl = ssht_dl_calloc(L, SSHT_DL_QUARTER);
1024 SO3_ERROR_MEM_ALLOC_CHECK(dl);
1026 if (dl_method == SSHT_DL_RISBO)
1028 dl8 = ssht_dl_calloc(L, SSHT_DL_QUARTER_EXTENDED);
1029 SO3_ERROR_MEM_ALLOC_CHECK(dl8);
1031 int dl_offset = ssht_dl_get_offset(L, SSHT_DL_QUARTER);
1032 int dl_stride = ssht_dl_get_stride(L, SSHT_DL_QUARTER);
1034 complex
double *mn_factors = calloc((2*L-1)*(2*N-1),
sizeof *mn_factors);
1035 SO3_ERROR_MEM_ALLOC_CHECK(mn_factors);
1042 for (el = L0; el <= L-1; ++el)
1049 if (el != 0 && el == L0)
1051 for(eltmp = 0; eltmp <= L0; ++eltmp)
1052 ssht_dl_beta_risbo_eighth_table(
1056 SSHT_DL_QUARTER_EXTENDED,
1060 ssht_dl_beta_risbo_fill_eighth2quarter_table(dl,
1063 SSHT_DL_QUARTER_EXTENDED,
1069 ssht_dl_beta_risbo_eighth_table(
1073 SSHT_DL_QUARTER_EXTENDED,
1077 ssht_dl_beta_risbo_fill_eighth2quarter_table(
1081 SSHT_DL_QUARTER_EXTENDED,
1087 case SSHT_DL_TRAPANI:
1088 if (el != 0 && el == L0)
1090 for(eltmp = 0; eltmp <= L0; ++eltmp)
1091 ssht_dl_halfpi_trapani_eighth_table(
1097 ssht_dl_halfpi_trapani_fill_eighth2quarter_table(
1106 ssht_dl_halfpi_trapani_eighth_table(
1112 ssht_dl_halfpi_trapani_fill_eighth2quarter_table(
1122 SO3_ERROR_GENERIC(
"Invalid dl method");
1128 double elfactor = (2.0*el+1.0)/(8.0*SO3_PI*SO3_PI);
1132 case SO3_N_MODE_ALL:
1133 n_start =
MAX(-N+1,-el);
1134 n_stop =
MIN( N-1, el);
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);
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;
1151 case SO3_N_MODE_MAXIMUM:
1156 n_inc =
MAX(1,2*N-2);
1163 n_inc =
MAX(1,2*el);
1166 SO3_ERROR_GENERIC(
"Invalid n-mode.");
1170 for (n = n_start; n <= n_stop; n += n_inc)
1171 for (m = -el; m <= el; ++m)
1175 int mod = ((n-m)%4 + 4)%4;
1176 mn_factors[m + m_offset + m_stride*(
1178 flmn[ind] * exps[mod];
1181 for (mm = 0; mm <= el; ++mm)
1185 double elmmsign = signs[el] * signs[mm];
1189 for (n = n_start; n <= n_stop; n += n_inc)
1191 double elnsign = n >= 0 ? 1.0 : elmmsign;
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))] +=
1200 * mn_factors[m + m_offset + m_stride*(
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))] +=
1208 * mn_factors[m + m_offset + m_stride*(
1210 * dl[m + dl_offset + mm*dl_stride];
1217 if (dl_method == SSHT_DL_RISBO)
1222 case SO3_N_MODE_ALL:
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;
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;
1238 case SO3_N_MODE_MAXIMUM:
1241 n_inc =
MAX(1,2*N - 2);
1244 SO3_ERROR_GENERIC(
"Invalid n-mode.");
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*(
1255 * Fmnm[m + m_offset + m_stride*(
1256 n + n_offset + n_stride*(
1260 for (mm = -L+1; mm <= L-1; ++mm)
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;
1271 complex
double *fext = calloc((2*L-1)*(2*L-1)*(2*N-1),
sizeof(*fext));
1272 SO3_ERROR_MEM_ALLOC_CHECK(fext);
1276 fftw_plan plan = fftw_plan_dft_3d(
1277 2*N-1, 2*L-1, 2*L-1,
1283 for (mm = -L+1; mm <= L-1; ++mm)
1285 int mm_shift = mm < 0 ? 2*L-1 : 0;
1286 for (n = n_start; n <= n_stop; n += n_inc)
1288 int n_shift = n < 0 ? 2*N-1 : 0;
1289 for (m = -L+1; m <= L-1; ++m)
1291 int m_shift = m < 0 ? 2*L-1 : 0;
1292 fext[m + m_shift + m_stride*(
1293 mm + mm_shift + mm_stride*(
1295 Fmnm[m + m_offset + m_stride*(
1296 n + n_offset + n_stride*(
1307 fftw_destroy_plan(plan);
1311 int a_stride = 2*L-1;
1312 int b_ext_stride = 2*L-1;
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)
1319 g))] = fext[a + a_stride*(
1327 printf(
"%sInverse transform computed!\n", SO3_PROMPT);
1355 complex
double *flmn,
const complex
double *f,
1356 const so3_parameters_t *parameters
1359 so3_sampling_t sampling;
1360 so3_storage_t storage;
1361 so3_n_mode_t n_mode;
1362 ssht_dl_method_t dl_method;
1366 L0 = parameters->L0;
1369 sampling = parameters->sampling_scheme;
1370 storage = parameters->storage;
1372 n_mode = parameters->n_mode;
1373 dl_method = parameters->dl_method;
1374 verbosity = parameters->verbosity;
1375 steerable = parameters->steerable;
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);
1382 printf(
"%sUsing routine so3_core_mw_forward_direct with storage method %d...\n"
1387 int m_stride = 2*L-1;
1391 int mm_stride = 2*L-1;
1392 int mm_offset = L-1;
1393 int a_stride = 2*L-1;
1395 int bext_stride = 2*L-1;
1398 int n_start, n_stop, n_inc;
1402 case SO3_N_MODE_ALL:
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;
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;
1418 case SO3_N_MODE_MAXIMUM:
1421 n_inc =
MAX(1,2*N - 2);
1424 SO3_ERROR_GENERIC(
"Invalid n-mode.");
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);
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)
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));
1451 double norm_factor = 1.0/(2.0*L-1.0)/(2.0*N-1.0);
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(
1465 for (b = 0; b < L; ++b)
1470 for (g = 0; g < 2*N-1; ++g)
1476 a_stride*
sizeof(*f));
1477 fftw_execute_dft(plan, inout, inout);
1480 for (n = n_start; n <= n_stop; n += n_inc)
1482 int n_shift = n < 0 ? 2*N-1 : 0;
1483 for (m = -L+1; m <= L-1; ++m)
1485 int m_shift = m < 0 ? 2*L-1 : 0;
1486 Fmnb[b + bext_stride*(
1487 m + m_offset + m_stride*(
1489 inout[m + m_shift + m_stride*(
1490 n + n_shift)] * norm_factor;
1494 fftw_destroy_plan(plan);
1497 for (n = n_start; n <= n_stop; n += n_inc)
1498 for (m = -L+1; m <= L-1; ++m)
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*(
1506 * Fmnb[(2*L-2-b) + bext_stride*(
1507 m + m_offset + m_stride*(
1513 complex
double *Fmnm = calloc((2*L-1)*(2*L-1)*(2*N-1),
sizeof(*Fmnm));
1514 SO3_ERROR_MEM_ALLOC_CHECK(Fmnm);
1516 plan = fftw_plan_dft_1d(
1521 for (n = n_start; n <= n_stop; n += n_inc)
1522 for (m = -L+1; m <= L-1; ++m)
1525 Fmnb + 0 + bext_stride*(
1526 m + m_offset + m_stride*(
1528 bext_stride*
sizeof(*Fmnb));
1529 fftw_execute_dft(plan, inout, inout);
1532 for (mm = -L+1; mm <= L-1; ++mm)
1534 int mm_shift = mm < 0 ? 2*L-1 : 0;
1535 Fmnm[mm + mm_offset + mm_stride*(
1536 m + m_offset + m_stride*(
1538 inout[mm + mm_shift] / (2.0*L-1.0);
1541 fftw_destroy_plan(plan);
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*(
1551 expsmm[mm + mm_offset];
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)
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(
1570 fftw_plan plan_fwd = fftw_plan_dft_1d(
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];
1583 fftw_execute_dft(plan_bwd, inout, inout);
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];
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)
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*(
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];
1617 fftw_execute_dft(plan_bwd, inout, inout);
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];
1626 for (r = -2*(L-1); r <= 2*(L-1); ++r)
1627 Fmnm_pad[r + w_offset] *= wr[r + w_offset];
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];
1635 fftw_execute_dft(plan_fwd, inout, inout);
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];
1643 for (mm = -(L-1); mm <= L-1; ++mm)
1644 Gmnm[m + m_offset + m_stride*(
1645 mm + mm_offset + mm_stride*(
1647 Fmnm_pad[mm + w_offset]
1648 * 4.0 * SSHT_PI * SSHT_PI / (4.0*L-3.0);
1651 fftw_destroy_plan(plan_bwd);
1652 fftw_destroy_plan(plan_fwd);
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)
1660 dl8 = ssht_dl_calloc(L, SSHT_DL_QUARTER_EXTENDED);
1661 SO3_ERROR_MEM_ALLOC_CHECK(dl8);
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)
1674 for (el = L0; el < L; ++el)
1682 if (el != 0 && el == L0)
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,
1691 SSHT_DL_QUARTER_EXTENDED,
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,
1703 SSHT_DL_QUARTER_EXTENDED,
1709 case SSHT_DL_TRAPANI:
1710 if (el != 0 && el == L0)
1712 for(eltmp = 0; eltmp <= L0; ++eltmp)
1713 ssht_dl_halfpi_trapani_eighth_table(dl, L,
1716 ssht_dl_halfpi_trapani_fill_eighth2quarter_table(dl, L,
1722 ssht_dl_halfpi_trapani_eighth_table(dl, L,
1725 ssht_dl_halfpi_trapani_fill_eighth2quarter_table(dl, L,
1732 SO3_ERROR_GENERIC(
"Invalid dl method");
1739 case SO3_N_MODE_ALL:
1740 n_start =
MAX(-N+1,-el);
1741 n_stop =
MIN( N-1, el);
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);
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;
1758 case SO3_N_MODE_MAXIMUM:
1763 n_inc =
MAX(1,2*N-2);
1770 n_inc =
MAX(1,2*el);
1773 SO3_ERROR_GENERIC(
"Invalid n-mode.");
1778 for (mm = -el; mm <= el; ++mm)
1782 double elmmsign = signs[el] * signs[abs(mm)];
1784 for (n = n_start; n <= n_stop; n += n_inc)
1786 double mmsign = mm >= 0 ? 1.0 : signs[el] * signs[abs(n)];
1787 double elnsign = n >= 0 ? 1.0 : elmmsign;
1790 double elnmm_factor = mmsign * elnsign
1791 * dl[abs(n) + dl_offset + abs(mm)*dl_stride];
1793 for (m = -el; m <= el; ++m)
1795 mmsign = mm >= 0 ? 1.0 : signs[el] * signs[abs(m)];
1796 double elmsign = m >= 0 ? 1.0 : elmmsign;
1799 int mod = ((m-n)%4 + 4)%4;
1804 * dl[abs(m) + dl_offset + abs(mm)*dl_stride]
1805 * Gmnm[m + m_offset + m_stride*(
1806 mm + mm_offset + mm_stride*(
1815 if (dl_method == SSHT_DL_RISBO)
1830 printf(
"%sForward transform computed!\n", SO3_PROMPT);
1851 double *f,
const complex
double *flmn,
1852 const so3_parameters_t *parameters
1856 so3_sampling_t sampling;
1857 so3_storage_t storage;
1858 so3_n_mode_t n_mode;
1859 ssht_dl_method_t dl_method;
1863 L0 = parameters->L0;
1866 sampling = parameters->sampling_scheme;
1867 storage = parameters->storage;
1869 n_mode = parameters->n_mode;
1870 dl_method = parameters->dl_method;
1871 verbosity = parameters->verbosity;
1872 steerable = parameters->steerable;
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);
1880 printf(
"%sUsing routine so3_core_mw_inverse_direct with storage method %d...\n"
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);
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)
1905 for (i = 0; i < 4; ++i)
1906 exps[i] = cexp(I*SO3_PION2*i);
1911 complex
double *Fmnm = calloc((2*L-1)*(2*L-1)*N,
sizeof(*Fmnm));
1912 SO3_ERROR_MEM_ALLOC_CHECK(Fmnm);
1914 int m_stride = 2*L-1;
1917 int mm_offset = L-1;
1920 int n_start, n_stop, n_inc;
1922 double *dl = ssht_dl_calloc(L, SSHT_DL_QUARTER);
1923 SO3_ERROR_MEM_ALLOC_CHECK(dl);
1925 if (dl_method == SSHT_DL_RISBO)
1927 dl8 = ssht_dl_calloc(L, SSHT_DL_QUARTER_EXTENDED);
1928 SO3_ERROR_MEM_ALLOC_CHECK(dl8);
1930 int dl_offset = ssht_dl_get_offset(L, SSHT_DL_QUARTER);
1931 int dl_stride = ssht_dl_get_stride(L, SSHT_DL_QUARTER);
1933 complex
double *mn_factors = calloc((2*L-1)*N,
sizeof *mn_factors);
1934 SO3_ERROR_MEM_ALLOC_CHECK(mn_factors);
1941 for (el = L0; el <= L-1; ++el)
1948 if (el != 0 && el == L0)
1950 for(eltmp = 0; eltmp <= L0; ++eltmp)
1951 ssht_dl_beta_risbo_eighth_table(
1955 SSHT_DL_QUARTER_EXTENDED,
1959 ssht_dl_beta_risbo_fill_eighth2quarter_table(dl,
1962 SSHT_DL_QUARTER_EXTENDED,
1968 ssht_dl_beta_risbo_eighth_table(
1972 SSHT_DL_QUARTER_EXTENDED,
1976 ssht_dl_beta_risbo_fill_eighth2quarter_table(
1980 SSHT_DL_QUARTER_EXTENDED,
1986 case SSHT_DL_TRAPANI:
1987 if (el != 0 && el == L0)
1989 for(eltmp = 0; eltmp <= L0; ++eltmp)
1990 ssht_dl_halfpi_trapani_eighth_table(
1996 ssht_dl_halfpi_trapani_fill_eighth2quarter_table(
2005 ssht_dl_halfpi_trapani_eighth_table(
2011 ssht_dl_halfpi_trapani_fill_eighth2quarter_table(
2021 SO3_ERROR_GENERIC(
"Invalid dl method");
2027 double elfactor = (2.0*el+1.0)/(8.0*SO3_PI*SO3_PI);
2032 case SO3_N_MODE_ALL:
2034 n_stop =
MIN( N-1, el);
2037 case SO3_N_MODE_EVEN:
2039 n_stop =
MIN( N-1, el);
2043 case SO3_N_MODE_ODD:
2045 n_stop =
MIN( N-1, el);
2046 n_stop -= 1-n_stop%2;
2049 case SO3_N_MODE_MAXIMUM:
2064 SO3_ERROR_GENERIC(
"Invalid n-mode.");
2068 for (n = n_start; n <= n_stop; n += n_inc)
2069 for (m = -el; m <= el; ++m)
2073 int mod = ((n-m)%4 + 4)%4;
2074 mn_factors[m + m_offset + m_stride*(
2076 flmn[ind] * exps[mod];
2079 for (mm = 0; mm <= el; ++mm)
2083 double elmmsign = signs[el] * signs[mm];
2085 for (n = n_start; n <= n_stop; n += n_inc)
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))] +=
2095 * mn_factors[m + m_offset + m_stride*(
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))] +=
2103 * mn_factors[m + m_offset + m_stride*(
2105 * dl[m + dl_offset + mm*dl_stride];
2112 if (dl_method == SSHT_DL_RISBO)
2117 case SO3_N_MODE_ALL:
2123 case SO3_N_MODE_EVEN:
2125 n_stop = ((N-1) % 2 == 0) ? N-1 : N-2;
2128 case SO3_N_MODE_ODD:
2130 n_stop = ((N-1) % 2 != 0) ? N-1 : N-2;
2133 case SO3_N_MODE_MAXIMUM:
2139 SO3_ERROR_GENERIC(
"Invalid n-mode.");
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*(
2150 * Fmnm[m + m_offset + m_stride*(
2151 n + n_offset + n_stride*(
2155 for (mm = -L+1; mm <= L-1; ++mm)
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;
2166 complex
double *Fmnm_shift = calloc((2*L-1)*(2*L-1)*N,
sizeof(*Fmnm_shift));
2167 SO3_ERROR_MEM_ALLOC_CHECK(Fmnm_shift);
2170 double *fext = calloc((2*L-1)*(2*L-1)*(2*N-1),
sizeof(*fext));
2171 SO3_ERROR_MEM_ALLOC_CHECK(fext);
2175 fftw_plan plan = fftw_plan_dft_c2r_3d(
2176 2*L-1, 2*L-1, 2*N-1,
2182 for (mm = -L+1; mm <= L-1; ++mm)
2184 int mm_shift = mm < 0 ? 2*L-1 : 0;
2185 for (n = n_start; n <= n_stop; n += n_inc)
2187 for (m = -L+1; m <= L-1; ++m)
2189 int m_shift = m < 0 ? 2*L-1 : 0;
2190 Fmnm_shift[n + n_stride*(
2191 m + m_shift + m_stride*(
2193 Fmnm[m + m_offset + m_stride*(
2194 n + n_offset + n_stride*(
2205 fftw_destroy_plan(plan);
2212 int a_stride = 2*L-1;
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)
2221 g))] = fext[g + g_stride*(
2229 printf(
"%sInverse transform computed!\n", SO3_PROMPT);
2256 complex
double *flmn,
const double *f,
2257 const so3_parameters_t *parameters
2260 so3_sampling_t sampling;
2261 so3_storage_t storage;
2262 so3_n_mode_t n_mode;
2263 ssht_dl_method_t dl_method;
2267 L0 = parameters->L0;
2270 sampling = parameters->sampling_scheme;
2271 storage = parameters->storage;
2273 n_mode = parameters->n_mode;
2274 dl_method = parameters->dl_method;
2275 verbosity = parameters->verbosity;
2276 steerable = parameters->steerable;
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);
2283 printf(
"%sUsing routine so3_core_mw_forward_direct with storage method %d...\n"
2288 int m_stride = 2*L-1;
2292 int mm_stride = 2*L-1;
2293 int mm_offset = L-1;
2294 int a_stride = 2*L-1;
2296 int bext_stride = 2*L-1;
2297 int g_stride = 2*N-1;
2299 int n_start, n_stop, n_inc;
2303 case SO3_N_MODE_ALL:
2309 case SO3_N_MODE_EVEN:
2311 n_stop = ((N-1) % 2 == 0) ? N-1 : N-2;
2314 case SO3_N_MODE_ODD:
2316 n_stop = ((N-1) % 2 != 0) ? N-1 : N-2;
2319 case SO3_N_MODE_MAXIMUM:
2325 SO3_ERROR_GENERIC(
"Invalid n-mode.");
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);
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)
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));
2352 double norm_factor = 1.0/(2.0*L-1.0)/(2.0*N-1.0);
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);
2362 fftw_plan plan = fftw_plan_dft_r2c_2d(
2368 for (b = 0; b < L; ++b)
2378 for (a = 0; a < 2*L-1; ++a)
2379 for (g = 0; g < 2*N-1; ++g)
2380 fft_in[g + g_stride*(
2390 for (n = n_start; n <= n_stop; n += n_inc)
2392 for (m = -L+1; m <= L-1; ++m)
2394 int m_shift = m < 0 ? 2*L-1 : 0;
2395 Fmnb[b + bext_stride*(
2396 m + m_offset + m_stride*(
2398 fft_out[n + n_stride*(
2399 m + m_shift)] * norm_factor;
2403 fftw_destroy_plan(plan);
2406 for (n = n_start; n <= n_stop; n += n_inc)
2407 for (m = -L+1; m <= L-1; ++m)
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*(
2415 * Fmnb[(2*L-2-b) + bext_stride*(
2416 m + m_offset + m_stride*(
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);
2427 plan = fftw_plan_dft_1d(
2432 for (n = n_start; n <= n_stop; n += n_inc)
2433 for (m = -L+1; m <= L-1; ++m)
2436 Fmnb + 0 + bext_stride*(
2437 m + m_offset + m_stride*(
2439 bext_stride*
sizeof(*Fmnb));
2443 for (mm = -L+1; mm <= L-1; ++mm)
2445 int mm_shift = mm < 0 ? 2*L-1 : 0;
2446 Fmnm[mm + mm_offset + mm_stride*(
2447 m + m_offset + m_stride*(
2449 inout[mm + mm_shift] / (2.0*L-1.0);
2452 fftw_destroy_plan(plan);
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*(
2462 expsmm[mm + mm_offset];
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)
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(
2481 fftw_plan plan_fwd = fftw_plan_dft_1d(
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];
2494 fftw_execute_dft(plan_bwd, inout, inout);
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];
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)
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*(
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];
2528 fftw_execute_dft(plan_bwd, inout, inout);
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];
2537 for (r = -2*(L-1); r <= 2*(L-1); ++r)
2538 Fmnm_pad[r + w_offset] *= wr[r + w_offset];
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];
2546 fftw_execute_dft(plan_fwd, inout, inout);
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];
2554 for (mm = -(L-1); mm <= L-1; ++mm)
2555 Gmnm[m + m_offset + m_stride*(
2556 mm + mm_offset + mm_stride*(
2558 Fmnm_pad[mm + w_offset]
2559 * 4.0 * SSHT_PI * SSHT_PI / (4.0*L-3.0);
2562 fftw_destroy_plan(plan_bwd);
2563 fftw_destroy_plan(plan_fwd);
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)
2571 dl8 = ssht_dl_calloc(L, SSHT_DL_QUARTER_EXTENDED);
2572 SO3_ERROR_MEM_ALLOC_CHECK(dl8);
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)
2585 for (el = L0; el < L; ++el)
2593 if (el != 0 && el == L0)
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,
2602 SSHT_DL_QUARTER_EXTENDED,
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,
2614 SSHT_DL_QUARTER_EXTENDED,
2620 case SSHT_DL_TRAPANI:
2621 if (el != 0 && el == L0)
2623 for(eltmp = 0; eltmp <= L0; ++eltmp)
2624 ssht_dl_halfpi_trapani_eighth_table(dl, L,
2627 ssht_dl_halfpi_trapani_fill_eighth2quarter_table(dl, L,
2633 ssht_dl_halfpi_trapani_eighth_table(dl, L,
2636 ssht_dl_halfpi_trapani_fill_eighth2quarter_table(dl, L,
2643 SO3_ERROR_GENERIC(
"Invalid dl method");
2651 case SO3_N_MODE_ALL:
2653 n_stop =
MIN( N-1, el);
2656 case SO3_N_MODE_EVEN:
2658 n_stop =
MIN( N-1, el);
2662 case SO3_N_MODE_ODD:
2664 n_stop =
MIN( N-1, el);
2665 n_stop -= 1-n_stop%2;
2668 case SO3_N_MODE_MAXIMUM:
2683 SO3_ERROR_GENERIC(
"Invalid n-mode.");
2688 for (mm = -el; mm <= el; ++mm)
2692 double elmmsign = signs[el] * signs[abs(mm)];
2694 for (n = n_start; n <= n_stop; n += n_inc)
2696 double mmsign = mm >= 0 ? 1.0 : signs[el] * signs[n];
2699 double elnmm_factor = mmsign
2700 * dl[n + dl_offset + abs(mm)*dl_stride];
2702 for (m = -el; m <= el; ++m)
2704 mmsign = mm >= 0 ? 1.0 : signs[el] * signs[abs(m)];
2705 double elmsign = m >= 0 ? 1.0 : elmmsign;
2708 int mod = ((m-n)%4 + 4)%4;
2713 * dl[abs(m) + dl_offset + abs(mm)*dl_stride]
2714 * Gmnm[m + m_offset + m_stride*(
2715 mm + mm_offset + mm_stride*(
2724 if (dl_method == SSHT_DL_RISBO)
2739 printf(
"%sForward transform computed!\n", SO3_PROMPT);