28 #define MAX(a,b) ((a) > (b) ? (a) : (b))
78 int L0,
int L,
int spin,
85 double *sqrt_tbl, *signs;
87 double ssign, elfactor;
91 int dl_offset, dl_stride;
94 double elmmsign, elssign;
97 int Fmm_offset, Fmm_stride, fext_stride;
101 sqrt_tbl = (
double*)calloc(2*(L-1)+2,
sizeof(double));
103 signs = (
double*)calloc(L+1,
sizeof(
double));
107 m_factors = calloc(2*L-1,
sizeof *m_factors);
111 for (el=0; el<=2*(L-1)+1; el++)
112 sqrt_tbl[el] = sqrt((
double)el);
113 for (m=0; m<=L-1; m=m+2) {
117 ssign = signs[abs(spin)];
118 spinneg = spin <= 0 ? spin : -spin;
120 for (m=-(L-1); m<=L-1; m++)
121 exps[m + exps_offset] = cexp(-I*
SSHT_PION2*(m+spin));
126 "Computing inverse transform using MW sampling with ");
127 printf(
"%s%s%d%s%d%s\n",
SSHT_PROMPT,
"parameters (L,spin,reality) = (",
128 L,
",", spin,
", FALSE)");
131 "Using routine ssht_core_mw_inverse_sov_sym...");
148 for (el=
MAX(L0, abs(spin)); el<=L-1; el++) {
154 if (el!=0 && el==
MAX(L0, abs(spin))) {
155 for(eltmp=0; eltmp<=
MAX(L0, abs(spin)); eltmp++)
158 eltmp, sqrt_tbl, signs);
169 el, sqrt_tbl, signs);
180 if (el!=0 && el==
MAX(L0, abs(spin))) {
181 for(eltmp=0; eltmp<=
MAX(L0, abs(spin)); eltmp++)
204 elfactor = ssign * sqrt((
double)(2.0*el+1.0)/(4.0*
SSHT_PI));
206 for (m=-el; m<=el; m++)
207 m_factors[m + m_offset] = flm[el2pel + m] * exps[m + exps_offset];
209 for (mm=0; mm<=el; mm++) {
211 int mm_offset = mm*dl_stride;
212 elmmsign = signs[el] * signs[mm];
213 elssign = spin <= 0 ? 1.0 : elmmsign;
214 mm_factor = elfactor * elssign * dl[mm_offset - spinneg + dl_offset];
215 for (m=-el; m<=-1; m++) {
216 Fmm[(mm + Fmm_offset)*Fmm_stride + m + Fmm_offset] +=
218 * m_factors[m + m_offset]
219 * elmmsign * dl[mm_offset - m + dl_offset];
221 for (m=0; m<=el; m++) {
222 Fmm[(mm + Fmm_offset)*Fmm_stride + m + Fmm_offset] +=
224 * m_factors[m + m_offset]
225 * dl[mm_offset + m + dl_offset];
238 for (mm=-(L-1); mm<=-1; mm++)
239 for (m=-(L-1); m<=L-1; m++)
240 Fmm[(mm + Fmm_offset)*Fmm_stride + m + Fmm_offset] =
241 signs[abs(m)] * ssign
242 * Fmm[(-mm + Fmm_offset)*Fmm_stride + m + Fmm_offset];
245 for (mm=-(L-1); mm<=L-1; mm++) {
246 mmfactor = cexp(I*mm*
SSHT_PI/(2.0*L-1.0));
247 for (m=-(L-1); m<=L-1; m++)
248 Fmm[(mm + Fmm_offset)*Fmm_stride + m + Fmm_offset] *=
258 for (mm=0; mm<=L-1; mm++)
259 for (m=0; m<=L-1; m++)
260 fext[mm*fext_stride + m] =
261 Fmm[(mm + Fmm_offset)*Fmm_stride + m + Fmm_offset];
262 for (mm=0; mm<=L-1; mm++)
263 for (m=-(L-1); m<=-1; m++)
264 fext[mm*fext_stride + (m+2*L-1)] =
265 Fmm[(mm + Fmm_offset)*Fmm_stride + m + Fmm_offset];
266 for (mm=-(L-1); mm<=-1; mm++)
267 for (m=0; m<=L-1; m++)
268 fext[(mm + 2*L-1)*fext_stride + m] =
269 Fmm[(mm + Fmm_offset)*Fmm_stride + m + Fmm_offset];
270 for (mm=-(L-1); mm<=-1; mm++)
271 for (m=-(L-1); m<=-1; m++)
272 fext[(mm+2*L-1)*fext_stride + m + 2*L-1] =
273 Fmm[(mm + Fmm_offset)*Fmm_stride + m + Fmm_offset];
276 plan = fftw_plan_dft_2d(2*L-1, 2*L-1, Fmm, Fmm,
277 FFTW_BACKWARD, FFTW_ESTIMATE);
278 fftw_execute_dft(plan, fext, fext);
279 fftw_destroy_plan(plan);
285 memcpy(f, fext, L*(2*L-1)*
sizeof(
SSHT_COMPLEX(
double)));
297 printf(
"%s%s",
SSHT_PROMPT,
"Inverse transform computed!");
356 double *sqrt_tbl, *signs;
357 int el2pel, m_offset;
358 double ssign, elfactor;
363 int dl_offset, dl_stride;
366 double elmmsign, elssign;
370 int Fmm_offset, Fmm_stride;
376 sqrt_tbl = (
double*)calloc(2*(L-1)+2,
sizeof(double));
378 signs = (
double*)calloc(L+1,
sizeof(
double));
382 m_factors = calloc(2*L-1,
sizeof *m_factors);
386 for (el=0; el<=2*(L-1)+1; el++)
387 sqrt_tbl[el] = sqrt((
double)el);
388 for (m=0; m<=L-1; m=m+2) {
392 ssign = signs[abs(spin)];
393 spinneg = spin <= 0 ? spin : -spin;
395 for (m=-(L-1); m<=L-1; m++)
396 exps[m + exps_offset] = cexp(-I*
SSHT_PION2*(m+spin));
401 "Computing inverse transform using MW sampling with ");
402 printf(
"%s%s%d%s%d%s\n",
SSHT_PROMPT,
"parameters (L,spin,reality) = (",
403 L,
",", spin,
", TRUE)");
406 "Using routine ssht_core_mw_inverse_sov_sym_real...");
423 for (el=
MAX(L0, abs(spin)); el<=L-1; el++) {
429 if (el!=0 && el==
MAX(L0, abs(spin))) {
430 for(eltmp=0; eltmp<=
MAX(L0, abs(spin)); eltmp++)
433 eltmp, sqrt_tbl, signs);
444 el, sqrt_tbl, signs);
455 if (el!=0 && el==
MAX(L0, abs(spin))) {
456 for(eltmp=0; eltmp<=
MAX(L0, abs(spin)); eltmp++)
480 elfactor = ssign * sqrt((
double)(2.0*el+1.0)/(4.0*
SSHT_PI));
482 for (m=-el; m<=el; m++)
483 m_factors[m + m_offset] = flm[el2pel + m] * exps[m + exps_offset];
485 for (mm=0; mm<=el; mm++) {
487 int mm_offset = mm*dl_stride;
488 elmmsign = signs[el] * signs[mm];
489 elssign = spin <= 0 ? 1.0 : elmmsign;
491 mm_factor = elfactor * elssign * dl[mm_offset - spinneg + dl_offset];
493 for (m=0; m<=el; m++) {
494 Fmm[(mm + Fmm_offset)*Fmm_stride + m] +=
496 * m_factors[m + m_offset]
497 * dl[mm_offset + m + dl_offset];
510 for (mm=-(L-1); mm<=-1; mm++)
511 for (m=0; m<=L-1; m++)
512 Fmm[(mm + Fmm_offset)*Fmm_stride + m] =
513 signs[abs(m)] * ssign
514 * Fmm[(-mm + Fmm_offset)*Fmm_stride + m];
517 for (mm=-(L-1); mm<=L-1; mm++) {
518 mmfactor = cexp(I*mm*
SSHT_PI/(2.0*L-1.0));
519 for (m=0; m<=L-1; m++)
520 Fmm[(mm + Fmm_offset)*Fmm_stride + m] *=
527 for (mm=0; mm<=L-1; mm++)
528 for (m=0; m<=L-1; m++)
529 Fmm_shift[mm*Fmm_stride + m] =
530 Fmm[(mm + Fmm_offset)*Fmm_stride + m];
531 for (mm=-(L-1); mm<=-1; mm++)
532 for (m=0; m<=L-1; m++)
533 Fmm_shift[(mm + 2*L-1)*Fmm_stride + m] =
534 Fmm[(mm + Fmm_offset)*Fmm_stride + m];
537 fext_real = (
double*)calloc((2*L-1)*(2*L-1),
sizeof(
double));
542 plan = fftw_plan_dft_c2r_2d(2*L-1, 2*L-1, Fmm_shift, fext_real,
544 fftw_execute_dft_c2r(plan, Fmm_shift, fext_real);
545 fftw_destroy_plan(plan);
552 memcpy(f, fext_real, L*(2*L-1)*
sizeof(
double));
564 printf(
"%s%s",
SSHT_PROMPT,
"Inverse transform computed!");
591 int L,
int spin,
int verbosity) {
593 int t, p, m, el, ind, eltmp;
596 double theta, phi, elfactor;
598 int dl_offset, dl_stride, f_stride;
601 sqrt_tbl = (
double*)calloc(2*(L-1)+2,
sizeof(double));
605 for (el=0; el<=2*(L-1)+1; el++)
606 sqrt_tbl[el] = sqrt((
double)el);
608 ssign = 1 - ssign - ssign;
613 "Computing inverse transform using MW sampling with ");
614 printf(
"%s%s%d%s%d%s\n",
SSHT_PROMPT,
"parameters (L,spin,reality) = (",
615 L,
",", spin,
", FALSE)");
618 "Using routine ssht_core_mwdirect_inverse...");
623 for (t=0; t<=L-1; t++)
624 for (p=0; p<=2*L-2; p++)
625 f[t*f_stride + p] = 0.0;
632 for (t=0; t<=L-1; t++) {
634 for (el=abs(spin); el<=L-1; el++) {
635 elfactor = sqrt((
double)(2.0*el+1.0)/(4.0*
SSHT_PI));
636 if (el!=0 && el==abs(spin)) {
637 for(eltmp=0; eltmp<=abs(spin); eltmp++)
648 for (m=-el; m<=el; m++) {
649 ssht_sampling_elm2ind(&ind, el, m);
650 for (p=0; p<=2*L-2; p++) {
656 * dl[(m+dl_offset)*dl_stride - spin + dl_offset]
669 printf(
"%s%s",
SSHT_PROMPT,
"Inverse transform computed!");
690 int L,
int spin,
int verbosity) {
692 int t, p, m, el, ind;
693 int ftm_stride, ftm_offset, f_stride;
694 double *dlm1p1_line, *dl_line;
696 double *sqrt_tbl, *signs;
698 double theta, ssign, elfactor;
702 sqrt_tbl = (
double*)calloc(2*(L-1)+2,
sizeof(double));
704 signs = (
double*)calloc(L+1,
sizeof(
double));
708 for (el=0; el<=2*(L-1)+1; el++)
709 sqrt_tbl[el] = sqrt((
double)el);
710 for (m=0; m<=L-1; m=m+2) {
714 ssign = signs[abs(spin)];
719 "Computing inverse transform using MW sampling with ");
720 printf(
"%s%s%d%s%d%s\n",
SSHT_PROMPT,
"parameters (L,spin,reality) = (",
721 L,
",", spin,
", FALSE)");
724 "Using routine ssht_core_mwdirect_inverse_sov...");
732 dlm1p1_line = (
double*)calloc(2*L-1,
sizeof(
double));
734 dl_line = (
double*)calloc(2*L-1,
sizeof(
double));
736 for (t=0; t<=L-1; t++) {
738 for (el=abs(spin); el<=L-1; el++) {
739 elfactor = sqrt((
double)(2.0*el+1.0)/(4.0*
SSHT_PI));
747 dl_line = dlm1p1_line;
748 dlm1p1_line = dl_ptr;
750 for (m=-el; m<=el; m++) {
751 ssht_sampling_elm2ind(&ind, el, m);
752 ftm[t*ftm_stride + m + ftm_offset] +=
769 plan = fftw_plan_dft_1d(2*L-1, inout, inout, FFTW_BACKWARD, FFTW_MEASURE);
770 for (t=0; t<=L-1; t++) {
771 for (m=0; m<=L-1; m++)
772 inout[m] = ftm[t*ftm_stride + m + ftm_offset];
773 for (m=-(L-1); m<=-1; m++)
774 inout[m+2*L-1] = ftm[t*ftm_stride + m + ftm_offset];
775 fftw_execute_dft(plan, inout, inout);
776 for (p=0; p<=2*L-2; p++)
777 f[t*f_stride + p] = inout[p];
779 fftw_destroy_plan(plan);
789 printf(
"%s%s",
SSHT_PROMPT,
"Inverse transform computed!");
838 int L0,
int L,
int spin,
842 int el, m, mm, ind, t, r;
844 double *sqrt_tbl, *signs;
846 double ssign, elfactor;
847 fftw_plan plan, plan_bwd, plan_fwd;
852 int f_stride, Fmt_stride, Fmt_offset, Fmm_stride, Fmm_offset;
855 int dl_offset, dl_stride;
859 int elmmsign, elssign;
863 sqrt_tbl = (
double*)calloc(2*(L-1)+2,
sizeof(double));
865 signs = (
double*)calloc(L+1,
sizeof(
double));
873 for (el=0; el<=2*(L-1)+1; el++)
874 sqrt_tbl[el] = sqrt((
double)el);
875 for (m=0; m<=L-1; m=m+2) {
879 ssign = signs[abs(spin)];
880 spinneg = spin <= 0 ? spin : -spin;
882 for (m=-(L-1); m<=L-1; m++)
883 expsm[m + exps_offset] = cexp(I*
SSHT_PION2*(m+spin));
884 for (mm=-(L-1); mm<=L-1; mm++)
885 expsmm[mm + exps_offset] = cexp(-I*mm*
SSHT_PI/(2.0*L-1.0));
890 "Computing forward transform using MW sampling with ");
891 printf(
"%s%s%d%s%d%s\n",
SSHT_PROMPT,
"parameters (L,spin,reality) = (",
892 L,
",", spin,
", FALSE)");
895 "Using routine ssht_core_mw_forward_sov_conv_sym...");
906 plan = fftw_plan_dft_1d(2*L-1, inout, inout, FFTW_FORWARD, FFTW_ESTIMATE);
907 for (t=0; t<=L-1; t++) {
908 memcpy(inout, &f[t*f_stride], f_stride*
sizeof(
SSHT_COMPLEX(
double)));
909 fftw_execute_dft(plan, inout, inout);
910 for(m=0; m<=L-1; m++)
911 Fmt[(m+Fmt_offset)*Fmt_stride + t] = inout[m] / (2.0*L-1.0);
912 for(m=-(L-1); m<=-1; m++)
913 Fmt[(m+Fmt_offset)*Fmt_stride + t] = inout[m+2*L-1] / (2.0*L-1.0);
917 for (m=-(L-1); m<=L-1; m++)
918 for (t=L; t<=2*L-2; t++)
919 Fmt[(m+Fmt_offset)*Fmt_stride + t] =
920 signs[abs(m)] * ssign * Fmt[(m+Fmt_offset)*Fmt_stride + (2*L-2-t)];
927 for (m=-(L-1); m<=L-1; m++) {
928 memcpy(inout, &Fmt[(m+Fmt_offset)*Fmt_stride], Fmt_stride*
sizeof(
SSHT_COMPLEX(
double)));
929 fftw_execute_dft(plan, inout, inout);
930 for(mm=0; mm<=L-1; mm++)
931 Fmm[(m+Fmm_offset)*Fmm_stride + mm + Fmm_offset] =
932 inout[mm] / (2.0*L-1.0);
933 for(mm=-(L-1); mm<=-1; mm++)
934 Fmm[(m+Fmm_offset)*Fmm_stride + mm + Fmm_offset] =
935 inout[mm+2*L-1] / (2.0*L-1.0);
937 fftw_destroy_plan(plan);
941 for (m=-(L-1); m<=L-1; m++)
942 for (mm=-(L-1); mm<=L-1; mm++)
943 Fmm[(m+Fmm_offset)*Fmm_stride + mm + Fmm_offset] *=
944 expsmm[mm + exps_offset];
950 for (mm=-2*(L-1); mm<=2*(L-1); mm++)
951 w[mm+w_offset] = ssht_sampling_weight_mw(mm);
958 plan_bwd = fftw_plan_dft_1d(4*L-3, inout, inout, FFTW_BACKWARD, FFTW_MEASURE);
959 plan_fwd = fftw_plan_dft_1d(4*L-3, inout, inout, FFTW_FORWARD, FFTW_MEASURE);
960 for (mm=1; mm<=2*L-2; mm++)
961 inout[mm + w_offset] = w[mm - 2*(L-1) - 1 + w_offset];
962 for (mm=-2*(L-1); mm<=0; mm++)
963 inout[mm + w_offset] = w[mm + 2*(L-1) + w_offset];
964 fftw_execute_dft(plan_bwd, inout, inout);
965 for (mm=0; mm<=2*L-2; mm++)
966 wr[mm + w_offset] = inout[mm - 2*(L-1) + w_offset];
967 for (mm=-2*(L-1); mm<=-1; mm++)
968 wr[mm + w_offset] = inout[mm + 2*(L-1) + 1 + w_offset];
977 for (m=-(L-1); m<=L-1; m++) {
980 for (mm=-2*(L-1); mm<=-L; mm++)
981 Fmm_pad[mm+w_offset] = 0.0;
982 for (mm=L; mm<=2*(L-1); mm++)
983 Fmm_pad[mm+w_offset] = 0.0;
984 for (mm=-(L-1); mm<=L-1; mm++)
985 Fmm_pad[mm + w_offset] =
986 Fmm[(m+Fmm_offset)*Fmm_stride + mm + Fmm_offset];
989 for (mm=1; mm<=2*L-2; mm++)
990 inout[mm + w_offset] = Fmm_pad[mm - 2*(L-1) - 1 + w_offset];
991 for (mm=-2*(L-1); mm<=0; mm++)
992 inout[mm + w_offset] = Fmm_pad[mm + 2*(L-1) + w_offset];
993 fftw_execute_dft(plan_bwd, inout, inout);
994 for (mm=0; mm<=2*L-2; mm++)
995 Fmm_pad[mm + w_offset] = inout[mm - 2*(L-1) + w_offset];
996 for (mm=-2*(L-1); mm<=-1; mm++)
997 Fmm_pad[mm + w_offset] = inout[mm + 2*(L-1) + 1 + w_offset];
1000 for (r=-2*(L-1); r<=2*(L-1); r++)
1001 Fmm_pad[r + w_offset] *= wr[-r + w_offset];
1004 for (mm=1; mm<=2*L-2; mm++)
1005 inout[mm + w_offset] = Fmm_pad[mm - 2*(L-1) - 1 + w_offset];
1006 for (mm=-2*(L-1); mm<=0; mm++)
1007 inout[mm + w_offset] = Fmm_pad[mm + 2*(L-1) + w_offset];
1008 fftw_execute_dft(plan_fwd, inout, inout);
1009 for (mm=0; mm<=2*L-2; mm++)
1010 Fmm_pad[mm + w_offset] = inout[mm - 2*(L-1) + w_offset];
1011 for (mm=-2*(L-1); mm<=-1; mm++)
1012 Fmm_pad[mm + w_offset] = inout[mm + 2*(L-1) + 1 + w_offset];
1015 for (mm=-(L-1); mm<=L-1; mm++)
1016 Gmm[(mm+Fmm_offset)*Fmm_stride + m + Fmm_offset] =
1017 Fmm_pad[mm + w_offset] * 2.0 *
SSHT_PI / (4.0*L-3.0);
1020 fftw_destroy_plan(plan_bwd);
1021 fftw_destroy_plan(plan_fwd);
1025 m_mm_factor = calloc(L*(2*L-1),
sizeof *m_mm_factor);
1028 for (mm = 1; mm < L; mm++) {
1029 for (m = -L+1; m < L; m++) {
1030 m_mm_factor[mm*Fmm_stride + m + Fmm_offset] =
1031 expsm[m + exps_offset]
1032 * ( Gmm[(mm+Fmm_offset)*Fmm_stride + m + Fmm_offset]
1033 + signs[abs(m)] * ssign
1034 * Gmm[(-mm+Fmm_offset)*Fmm_stride + m + Fmm_offset] );
1047 for (el=0; el<=L-1; el++) {
1048 for (m=-el; m<=el; m++) {
1049 ssht_sampling_elm2ind(&ind, el, m);
1053 for (el=
MAX(L0, abs(spin)); el<=L-1; el++) {
1057 switch (dl_method) {
1060 if (el!=0 && el==
MAX(L0, abs(spin))) {
1061 for(eltmp=0; eltmp<=
MAX(L0, abs(spin)); eltmp++)
1064 eltmp, sqrt_tbl, signs);
1075 el, sqrt_tbl, signs);
1086 if (el!=0 && el==
MAX(L0, abs(spin))) {
1087 for(eltmp=0; eltmp<=
MAX(L0, abs(spin)); eltmp++)
1110 elfactor = sqrt((
double)(2.0*el+1.0)/(4.0*
SSHT_PI));
1111 el2pel = el *el + el;
1112 elssign = spin <= 0 ? 1.0 : signs[el];
1114 mm_factor = ssign * elfactor * elssign * dl[- spinneg + dl_offset];
1116 for (m=-el; m<=-1; m++) {
1121 * expsm[m + exps_offset]
1122 * dl[0*dl_stride - m + dl_offset]
1123 * Gmm[(0+Fmm_offset)*Fmm_stride + m + Fmm_offset];
1125 for (m=0; m<=el; m++) {
1129 * expsm[m + exps_offset]
1130 * dl[0*dl_stride + m + dl_offset]
1131 * Gmm[(0+Fmm_offset)*Fmm_stride + m + Fmm_offset];
1134 for (mm=1; mm<=el; mm++) {
1135 int mm_offset = mm * dl_stride;
1136 elmmsign = signs[el] * signs[mm];
1137 elssign = spin <= 0 ? 1.0 : elmmsign;
1143 * dl[mm_offset - spinneg + dl_offset];
1145 for (m=-el; m<=-1; m++) {
1148 * elmmsign * dl[mm_offset - m + dl_offset]
1149 * m_mm_factor[mm*Fmm_stride + m + Fmm_offset];
1151 for (m=0; m<=el; m++) {
1154 * dl[mm_offset + m + dl_offset]
1155 * m_mm_factor[mm*Fmm_stride + m + Fmm_offset];
1182 printf(
"%s%s",
SSHT_PROMPT,
"Forward transform computed!");
1235 int el, m, mm, ind, ind_nm, t, r;
1237 double *sqrt_tbl, *signs;
1239 double ssign, elfactor;
1240 fftw_plan plan, plan_bwd, plan_fwd;
1246 int f_stride, Fmt_stride, Fmt_offset, Fmm_stride, Fmm_offset, Gmm_stride;
1249 int dl_offset, dl_stride;
1253 int elmmsign, elssign;
1258 sqrt_tbl = (
double*)calloc(2*(L-1)+2,
sizeof(double));
1260 signs = (
double*)calloc(L+1,
sizeof(
double));
1268 for (el=0; el<=2*(L-1)+1; el++)
1269 sqrt_tbl[el] = sqrt((
double)el);
1270 for (m=0; m<=L-1; m=m+2) {
1274 ssign = signs[abs(spin)];
1275 spinneg = spin <= 0 ? spin : -spin;
1277 for (m=0; m<=L-1; m++)
1279 for (mm=-(L-1); mm<=L-1; mm++)
1280 expsmm[mm + exps_offset] = cexp(-I*mm*
SSHT_PI/(2.0*L-1.0));
1283 if (verbosity > 0) {
1285 "Computing forward transform using MW sampling with ");
1286 printf(
"%s%s%d%s%d%s\n",
SSHT_PROMPT,
"parameters (L,spin,reality) = (",
1287 L,
",", spin,
", TRUE)");
1290 "Using routine ssht_core_mw_forward_sov_conv_sym_real...");
1298 in_real = (
double*)calloc(2*L-1,
sizeof(
double));
1302 plan = fftw_plan_dft_r2c_1d(2*L-1, in_real, out, FFTW_MEASURE);
1303 for (t=0; t<=L-1; t++) {
1304 memcpy(in_real, &f[t*f_stride], f_stride*
sizeof(
double));
1305 fftw_execute_dft_r2c(plan, in_real, out);
1306 for(m=0; m<=L-1; m++)
1307 Fmt[m*Fmt_stride + t] = out[m] / (2.0*L-1.0);
1311 fftw_destroy_plan(plan);
1314 for (m=0; m<=L-1; m++)
1315 for (t=L; t<=2*L-2; t++)
1316 Fmt[m*Fmt_stride + t] =
1317 signs[abs(m)] * ssign * Fmt[m*Fmt_stride + (2*L-2-t)];
1326 plan = fftw_plan_dft_1d(2*L-1, inout, inout, FFTW_FORWARD, FFTW_MEASURE);
1327 for (m=0; m<=L-1; m++) {
1328 memcpy(inout, &Fmt[m*Fmt_stride], Fmt_stride*
sizeof(
SSHT_COMPLEX(
double)));
1329 fftw_execute_dft(plan, inout, inout);
1330 for(mm=0; mm<=L-1; mm++)
1331 Fmm[m*Fmm_stride + mm + Fmm_offset] =
1332 inout[mm] / (2.0*L-1.0);
1333 for(mm=-(L-1); mm<=-1; mm++)
1334 Fmm[m*Fmm_stride + mm + Fmm_offset] =
1335 inout[mm+2*L-1] / (2.0*L-1.0);
1337 fftw_destroy_plan(plan);
1341 for (m=0; m<=L-1; m++)
1342 for (mm=-(L-1); mm<=L-1; mm++)
1343 Fmm[m*Fmm_stride + mm + Fmm_offset] *=
1344 expsmm[mm + exps_offset];
1350 for (mm=-2*(L-1); mm<=2*(L-1); mm++)
1351 w[mm+w_offset] = ssht_sampling_weight_mw(mm);
1358 plan_bwd = fftw_plan_dft_1d(4*L-3, inout, inout, FFTW_BACKWARD, FFTW_MEASURE);
1359 plan_fwd = fftw_plan_dft_1d(4*L-3, inout, inout, FFTW_FORWARD, FFTW_MEASURE);
1360 for (mm=1; mm<=2*L-2; mm++)
1361 inout[mm + w_offset] = w[mm - 2*(L-1) - 1 + w_offset];
1362 for (mm=-2*(L-1); mm<=0; mm++)
1363 inout[mm + w_offset] = w[mm + 2*(L-1) + w_offset];
1364 fftw_execute_dft(plan_bwd, inout, inout);
1365 for (mm=0; mm<=2*L-2; mm++)
1366 wr[mm + w_offset] = inout[mm - 2*(L-1) + w_offset];
1367 for (mm=-2*(L-1); mm<=-1; mm++)
1368 wr[mm + w_offset] = inout[mm + 2*(L-1) + 1 + w_offset];
1378 for (m=0; m<=L-1; m++) {
1381 for (mm=-2*(L-1); mm<=-L; mm++)
1382 Fmm_pad[mm+w_offset] = 0.0;
1383 for (mm=L; mm<=2*(L-1); mm++)
1384 Fmm_pad[mm+w_offset] = 0.0;
1385 for (mm=-(L-1); mm<=L-1; mm++)
1386 Fmm_pad[mm + w_offset] =
1387 Fmm[m*Fmm_stride + mm + Fmm_offset];
1390 for (mm=1; mm<=2*L-2; mm++)
1391 inout[mm + w_offset] = Fmm_pad[mm - 2*(L-1) - 1 + w_offset];
1392 for (mm=-2*(L-1); mm<=0; mm++)
1393 inout[mm + w_offset] = Fmm_pad[mm + 2*(L-1) + w_offset];
1394 fftw_execute_dft(plan_bwd, inout, inout);
1395 for (mm=0; mm<=2*L-2; mm++)
1396 Fmm_pad[mm + w_offset] = inout[mm - 2*(L-1) + w_offset];
1397 for (mm=-2*(L-1); mm<=-1; mm++)
1398 Fmm_pad[mm + w_offset] = inout[mm + 2*(L-1) + 1 + w_offset];
1401 for (r=-2*(L-1); r<=2*(L-1); r++)
1402 Fmm_pad[r + w_offset] *= wr[-r + w_offset];
1405 for (mm=1; mm<=2*L-2; mm++)
1406 inout[mm + w_offset] = Fmm_pad[mm - 2*(L-1) - 1 + w_offset];
1407 for (mm=-2*(L-1); mm<=0; mm++)
1408 inout[mm + w_offset] = Fmm_pad[mm + 2*(L-1) + w_offset];
1409 fftw_execute_dft(plan_fwd, inout, inout);
1410 for (mm=0; mm<=2*L-2; mm++)
1411 Fmm_pad[mm + w_offset] = inout[mm - 2*(L-1) + w_offset];
1412 for (mm=-2*(L-1); mm<=-1; mm++)
1413 Fmm_pad[mm + w_offset] = inout[mm + 2*(L-1) + 1 + w_offset];
1416 for (mm=-(L-1); mm<=L-1; mm++)
1417 Gmm[(mm+Fmm_offset)*Gmm_stride + m] =
1418 Fmm_pad[mm + w_offset] * 2.0 *
SSHT_PI / (4.0*L-3.0);
1421 fftw_destroy_plan(plan_bwd);
1422 fftw_destroy_plan(plan_fwd);
1426 m_mm_factor = calloc(L*L,
sizeof *m_mm_factor);
1429 for (mm = 1; mm < L; mm++) {
1430 for (m = 0; m < L; m++) {
1431 m_mm_factor[mm*Gmm_stride + m] =
1433 * ( Gmm[(mm+Fmm_offset)*Gmm_stride + m]
1435 * Gmm[(-mm+Fmm_offset)*Gmm_stride + m] );
1448 for (el=0; el<=L-1; el++) {
1449 for (m=0; m<=el; m++) {
1450 ssht_sampling_elm2ind(&ind, el, m);
1454 for (el=
MAX(L0, abs(spin)); el<=L-1; el++) {
1458 switch (dl_method) {
1461 if (el!=0 && el==
MAX(L0, abs(spin))) {
1462 for(eltmp=0; eltmp<=
MAX(L0, abs(spin)); eltmp++)
1465 eltmp, sqrt_tbl, signs);
1476 el, sqrt_tbl, signs);
1487 if (el!=0 && el==
MAX(L0, abs(spin))) {
1488 for(eltmp=0; eltmp<=
MAX(L0, abs(spin)); eltmp++)
1512 elfactor = sqrt((
double)(2.0*el+1.0)/(4.0*
SSHT_PI));
1513 el2pel = el *el + el;
1514 elssign = spin <= 0 ? 1.0 : signs[el];
1516 mm_factor = ssign * elfactor * elssign * dl[- spinneg + dl_offset];
1518 for (m=0; m<=el; m++) {
1523 * dl[0*dl_stride + m + dl_offset]
1524 * Gmm[(0+Fmm_offset)*Gmm_stride + m];
1527 for (mm=1; mm<=el; mm++) {
1528 int mm_offset = mm * dl_stride;
1529 elmmsign = signs[el] * signs[mm];
1530 elssign = spin <= 0 ? 1.0 : elmmsign;
1536 * dl[mm_offset - spinneg + dl_offset];
1538 for (m=0; m<=el; m++) {
1541 * dl[mm_offset + m + dl_offset]
1542 * m_mm_factor[mm*Gmm_stride + m];
1550 for (el=abs(spin); el<=L-1; el++) {
1551 for (m=1; m<=el; m++) {
1552 ssht_sampling_elm2ind(&ind, el, m);
1553 ssht_sampling_elm2ind(&ind_nm, el, -m);
1554 flm[ind_nm] = signs[m] * conj(flm[ind]);
1578 printf(
"%s%s",
SSHT_PROMPT,
"Forward transform computed!");
1616 int f_stride = 2*L-1;
1624 dl_method, verbosity);
1627 memcpy(f, f_full, (L-1)*(2*L-1)*
sizeof(
SSHT_COMPLEX(
double)));
1628 *f_sp = f_full[(L-1)*f_stride + 0];
1662 int f_stride = 2*L-1;
1665 f_full = (
double*)calloc(L*(2*L-1),
sizeof(double));
1670 dl_method, verbosity);
1673 memcpy(f, f_full, (L-1)*(2*L-1)*
sizeof(
double));
1674 *f_sp = f_full[(L-1)*f_stride + 0];
1709 int p, f_stride = 2*L-1;
1715 memcpy(f_full, f, (L-1)*(2*L-1)*
sizeof(
SSHT_COMPLEX(
double)));
1718 for (p=0; p<=2*L-2; p++) {
1720 f_full[(L-1)*f_stride + p] = f_sp * cexp(I*spin*(phi-phi_sp));
1725 dl_method, verbosity);
1758 int p, f_stride = 2*L-1;
1761 f_full = (
double*)calloc(L*(2*L-1),
sizeof(double));
1763 memcpy(f_full, f, (L-1)*(2*L-1)*
sizeof(
double));
1766 for (p=0; p<=2*L-2; p++)
1767 f_full[(L-1)*f_stride + p] = f_sp;
1771 dl_method, verbosity);
1825 int L0,
int L,
int spin,
1832 double *sqrt_tbl, *signs;
1833 int el2pel, inds_offset;
1835 double ssign, elfactor;
1838 int dl_offset, dl_stride;
1841 double elmmsign, elssign;
1844 int Fmm_offset, Fmm_stride, fext_stride;
1848 sqrt_tbl = (
double*)calloc(2*(L-1)+2,
sizeof(double));
1850 signs = (
double*)calloc(L+1,
sizeof(
double));
1854 inds = (
int*)calloc(2*L-1,
sizeof(
int));
1858 for (el=0; el<=2*(L-1)+1; el++)
1859 sqrt_tbl[el] = sqrt((
double)el);
1860 for (m=0; m<=L-1; m=m+2) {
1864 ssign = signs[abs(spin)];
1865 spinneg = spin <= 0 ? spin : -spin;
1867 for (m=-(L-1); m<=L-1; m++)
1868 exps[m + exps_offset] = cexp(-I*
SSHT_PION2*(m+spin));
1871 if (verbosity > 0) {
1873 "Computing inverse transform using MW symmetric sampling with ");
1874 printf(
"%s%s%d%s%d%s\n",
SSHT_PROMPT,
"parameters (L,spin,reality) = (",
1875 L,
",", spin,
", FALSE)");
1878 "Using routine ssht_core_mw_inverse_sov_sym_ss...");
1897 for (el=
MAX(L0, abs(spin)); el<=L-1; el++) {
1900 switch (dl_method) {
1903 if (el!=0 && el==
MAX(L0, abs(spin))) {
1904 for(eltmp=0; eltmp<=
MAX(L0, abs(spin)); eltmp++)
1907 eltmp, sqrt_tbl, signs);
1918 el, sqrt_tbl, signs);
1929 if (el!=0 && el==
MAX(L0, abs(spin))) {
1930 for(eltmp=0; eltmp<=
MAX(L0, abs(spin)); eltmp++)
1953 elfactor = sqrt((
double)(2.0*el+1.0)/(4.0*
SSHT_PI));
1954 el2pel = el *el + el;
1955 for (m=-el; m<=el; m++)
1956 inds[m + inds_offset] = el2pel + m;
1957 for (mm=0; mm<=el; mm++) {
1958 elmmsign = signs[el] * signs[mm];
1959 elssign = spin <= 0 ? 1.0 : elmmsign;
1961 for (m=-el; m<=-1; m++) {
1962 ind = inds[m + inds_offset];
1963 Fmm[(mm + Fmm_offset)*Fmm_stride + m + Fmm_offset] +=
1966 * exps[m + exps_offset]
1967 * elmmsign * dl[mm*dl_stride - m + dl_offset]
1968 * elssign * dl[mm*dl_stride - spinneg + dl_offset]
1971 for (m=0; m<=el; m++) {
1972 ind = inds[m + inds_offset];
1973 Fmm[(mm + Fmm_offset)*Fmm_stride + m + Fmm_offset] +=
1976 * exps[m + exps_offset]
1977 * dl[mm*dl_stride + m + dl_offset]
1978 * elssign * dl[mm*dl_stride - spinneg + dl_offset]
1992 for (mm=-(L-1); mm<=-1; mm++)
1993 for (m=-(L-1); m<=L; m++)
1994 Fmm[(mm + Fmm_offset)*Fmm_stride + m + Fmm_offset] =
1995 signs[abs(m)] * ssign
1996 * Fmm[(-mm + Fmm_offset)*Fmm_stride + m + Fmm_offset];
2006 for (mm=0; mm<=L; mm++)
2007 for (m=0; m<=L; m++)
2008 fext[mm*fext_stride + m] =
2009 Fmm[(mm + Fmm_offset)*Fmm_stride + m + Fmm_offset];
2010 for (mm=0; mm<=L; mm++)
2011 for (m=-(L-1); m<=-1; m++)
2012 fext[mm*fext_stride + (m+2*L-1+1)] =
2013 Fmm[(mm + Fmm_offset)*Fmm_stride + m + Fmm_offset];
2014 for (mm=-(L-1); mm<=-1; mm++)
2015 for (m=0; m<=L; m++)
2016 fext[(mm + 2*L-1+1)*fext_stride + m] =
2017 Fmm[(mm + Fmm_offset)*Fmm_stride + m + Fmm_offset];
2018 for (mm=-(L-1); mm<=-1; mm++)
2019 for (m=-(L-1); m<=-1; m++)
2020 fext[(mm+2*L-1+1)*fext_stride + m + 2*L-1+1] =
2021 Fmm[(mm + Fmm_offset)*Fmm_stride + m + Fmm_offset];
2024 plan = fftw_plan_dft_2d(2*L, 2*L, Fmm, Fmm,
2025 FFTW_BACKWARD, FFTW_ESTIMATE);
2026 fftw_execute_dft(plan, fext, fext);
2027 fftw_destroy_plan(plan);
2035 memcpy(f, fext, (L+1)*(2*L)*
sizeof(
SSHT_COMPLEX(
double)));
2042 printf(
"%s%s",
SSHT_PROMPT,
"Inverse transform computed!");
2103 double *sqrt_tbl, *signs;
2104 int el2pel, inds_offset;
2106 double ssign, elfactor;
2109 int dl_offset, dl_stride;
2112 double elmmsign, elssign;
2116 int Fmm_offset, Fmm_stride;
2121 sqrt_tbl = (
double*)calloc(2*(L-1)+2,
sizeof(double));
2123 signs = (
double*)calloc(L+1,
sizeof(
double));
2127 inds = (
int*)calloc(2*L-1,
sizeof(
int));
2131 for (el=0; el<=2*(L-1)+1; el++)
2132 sqrt_tbl[el] = sqrt((
double)el);
2133 for (m=0; m<=L-1; m=m+2) {
2137 ssign = signs[abs(spin)];
2138 spinneg = spin <= 0 ? spin : -spin;
2140 for (m=-(L-1); m<=L-1; m++)
2141 exps[m + exps_offset] = cexp(-I*
SSHT_PION2*(m+spin));
2144 if (verbosity > 0) {
2146 "Computing inverse transform using MW symmetric sampling with ");
2147 printf(
"%s%s%d%s%d%s\n",
SSHT_PROMPT,
"parameters (L,spin,reality) = (",
2148 L,
",", spin,
", TRUE)");
2151 "Using routine ssht_core_mw_inverse_sov_sym_ss_real...");
2170 for (el=
MAX(L0, abs(spin)); el<=L-1; el++) {
2173 switch (dl_method) {
2176 if (el!=0 && el==
MAX(L0, abs(spin))) {
2177 for(eltmp=0; eltmp<=
MAX(L0, abs(spin)); eltmp++)
2180 eltmp, sqrt_tbl, signs);
2191 el, sqrt_tbl, signs);
2202 if (el!=0 && el==
MAX(L0, abs(spin))) {
2203 for(eltmp=0; eltmp<=
MAX(L0, abs(spin)); eltmp++)
2226 elfactor = sqrt((
double)(2.0*el+1.0)/(4.0*
SSHT_PI));
2227 el2pel = el *el + el;
2228 for (m=0; m<=el; m++)
2229 inds[m + inds_offset] = el2pel + m;
2230 for (mm=0; mm<=el; mm++) {
2231 elmmsign = signs[el] * signs[mm];
2232 elssign = spin <= 0 ? 1.0 : elmmsign;
2234 for (m=0; m<=el; m++) {
2235 ind = inds[m + inds_offset];
2236 Fmm[(mm + Fmm_offset)*Fmm_stride + m] +=
2239 * exps[m + exps_offset]
2240 * dl[mm*dl_stride + m + dl_offset]
2241 * elssign * dl[mm*dl_stride - spinneg + dl_offset]
2255 for (mm=-(L-1); mm<=-1; mm++)
2256 for (m=0; m<=L; m++)
2257 Fmm[(mm + Fmm_offset)*Fmm_stride + m] =
2258 signs[abs(m)] * ssign
2259 * Fmm[(-mm + Fmm_offset)*Fmm_stride + m];
2264 for (mm=0; mm<=L-1; mm++)
2265 for (m=0; m<=L-1; m++)
2266 Fmm_shift[mm*Fmm_stride + m] =
2267 Fmm[(mm + Fmm_offset)*Fmm_stride + m];
2268 for (mm=-(L-1); mm<=-1; mm++)
2269 for (m=0; m<=L-1; m++)
2270 Fmm_shift[(mm + 2*L-1+1)*Fmm_stride + m] =
2271 Fmm[(mm + Fmm_offset)*Fmm_stride + m];
2276 fext_real = (
double*)calloc((2*L)*(2*L),
sizeof(
double));
2280 plan = fftw_plan_dft_c2r_2d(2*L, 2*L, Fmm_shift, fext_real,
2282 fftw_execute_dft_c2r(plan, Fmm_shift, fext_real);
2283 fftw_destroy_plan(plan);
2292 memcpy(f, fext_real, (L+1)*(2*L)*
sizeof(
double));
2299 printf(
"%s%s",
SSHT_PROMPT,
"Inverse transform computed!");
2327 int L,
int spin,
int verbosity) {
2329 int t, p, m, el, ind, eltmp;
2332 double theta, phi, elfactor;
2334 int dl_offset, dl_stride, f_stride;
2337 sqrt_tbl = (
double*)calloc(2*(L-1)+2,
sizeof(double));
2341 for (el=0; el<=2*(L-1)+1; el++)
2342 sqrt_tbl[el] = sqrt((
double)el);
2344 ssign = 1 - ssign - ssign;
2347 if (verbosity > 0) {
2349 "Computing inverse transform using MW sampling with ");
2350 printf(
"%s%s%d%s%d%s\n",
SSHT_PROMPT,
"parameters (L,spin,reality) = (",
2351 L,
",", spin,
", FALSE)");
2354 "Using routine ssht_core_mwdirect_inverse_ss...");
2359 for (t=0; t<=L; t++)
2360 for (p=0; p<=2*L-1; p++)
2361 f[t*f_stride + p] = 0.0;
2368 for (t=0; t<=L; t++) {
2370 for (el=abs(spin); el<=L-1; el++) {
2371 elfactor = sqrt((
double)(2.0*el+1.0)/(4.0*
SSHT_PI));
2372 if (el!=0 && el==abs(spin)) {
2373 for(eltmp=0; eltmp<=abs(spin); eltmp++)
2384 for (m=-el; m<=el; m++) {
2385 ssht_sampling_elm2ind(&ind, el, m);
2386 for (p=0; p<=2*L-1; p++) {
2388 f[t*f_stride + p] +=
2392 * dl[(m+dl_offset)*dl_stride - spin + dl_offset]
2405 printf(
"%s%s",
SSHT_PROMPT,
"Inverse transform computed!");
2453 int L0,
int L,
int spin,
2457 int el, m, mm, ind, t, r;
2459 double *sqrt_tbl, *signs;
2460 int el2pel, inds_offset;
2462 double ssign, elfactor;
2463 fftw_plan plan, plan_bwd, plan_fwd;
2468 int f_stride, Fmt_stride, Fmt_offset, Fmm_stride, Fmm_offset;
2471 int dl_offset, dl_stride;
2475 int elmmsign, elssign;
2477 int Gmm_stride, Gmm_offset;
2480 sqrt_tbl = (
double*)calloc(2*(L-1)+2,
sizeof(double));
2482 signs = (
double*)calloc(L+1,
sizeof(
double));
2488 inds = (
int*)calloc(2*L-1,
sizeof(
int));
2492 for (el=0; el<=2*(L-1)+1; el++)
2493 sqrt_tbl[el] = sqrt((
double)el);
2494 for (m=0; m<=L-1; m=m+2) {
2498 ssign = signs[abs(spin)];
2499 spinneg = spin <= 0 ? spin : -spin;
2501 for (m=-(L-1); m<=L-1; m++)
2502 expsm[m + exps_offset] = cexp(I*
SSHT_PION2*(m+spin));
2503 for (mm=-(L-1); mm<=L-1; mm++)
2504 expsmm[mm + exps_offset] = cexp(-I*mm*
SSHT_PI/(2.0*L-1.0));
2507 if (verbosity > 0) {
2509 "Computing forward transform using MW symmetric sampling with ");
2510 printf(
"%s%s%d%s%d%s\n",
SSHT_PROMPT,
"parameters (L,spin,reality) = (",
2511 L,
",", spin,
", FALSE)");
2514 "Using routine ssht_core_mw_forward_sov_conv_sym_ss...");
2527 plan = fftw_plan_dft_1d(2*L, inout, inout, FFTW_FORWARD, FFTW_MEASURE);
2528 for (t=0; t<=L; t++) {
2529 memcpy(inout, &f[t*f_stride], f_stride*
sizeof(
SSHT_COMPLEX(
double)));
2530 fftw_execute_dft(plan, inout, inout);
2532 Fmt[(m+Fmt_offset)*Fmt_stride + t] = inout[m] / (2.0*L);
2533 for(m=-(L-1); m<=-1; m++)
2534 Fmt[(m+Fmt_offset)*Fmt_stride + t] = inout[m+2*L-1+1] / (2.0*L);
2536 fftw_destroy_plan(plan);
2540 for (m=-(L-1); m<=L; m++)
2541 for (t=L+1; t<=2*L-1; t++)
2542 Fmt[(m+Fmt_offset)*Fmt_stride + t] =
2543 signs[abs(m)] * ssign * Fmt[(m+Fmt_offset)*Fmt_stride + (2*L-t)];
2553 plan = fftw_plan_dft_1d(2*L, inout, inout, FFTW_FORWARD, FFTW_MEASURE);
2554 for (m=-(L-1); m<=L; m++) {
2555 memcpy(inout, &Fmt[(m+Fmt_offset)*Fmt_stride], Fmt_stride*
sizeof(
SSHT_COMPLEX(
double)));
2556 fftw_execute_dft(plan, inout, inout);
2557 for(mm=0; mm<=L; mm++)
2558 Fmm[(m+Fmm_offset)*Fmm_stride + mm + Fmm_offset] =
2559 inout[mm] / (2.0*L);
2560 for(mm=-(L-1); mm<=-1; mm++)
2561 Fmm[(m+Fmm_offset)*Fmm_stride + mm + Fmm_offset] =
2562 inout[mm+2*L-1+1] / (2.0*L);
2564 fftw_destroy_plan(plan);
2571 for (mm=-2*(L-1); mm<=2*(L-1); mm++)
2572 w[mm+w_offset] = ssht_sampling_weight_mw(mm);
2579 plan_bwd = fftw_plan_dft_1d(4*L-3, inout, inout, FFTW_BACKWARD, FFTW_MEASURE);
2580 plan_fwd = fftw_plan_dft_1d(4*L-3, inout, inout, FFTW_FORWARD, FFTW_MEASURE);
2581 for (mm=1; mm<=2*L-2; mm++)
2582 inout[mm + w_offset] = w[mm - 2*(L-1) - 1 + w_offset];
2583 for (mm=-2*(L-1); mm<=0; mm++)
2584 inout[mm + w_offset] = w[mm + 2*(L-1) + w_offset];
2585 fftw_execute_dft(plan_bwd, inout, inout);
2586 for (mm=0; mm<=2*L-2; mm++)
2587 wr[mm + w_offset] = inout[mm - 2*(L-1) + w_offset];
2588 for (mm=-2*(L-1); mm<=-1; mm++)
2589 wr[mm + w_offset] = inout[mm + 2*(L-1) + 1 + w_offset];
2600 for (m=-(L-1); m<=L-1; m++) {
2603 for (mm=-2*(L-1); mm<=-L; mm++)
2604 Fmm_pad[mm+w_offset] = 0.0;
2605 for (mm=L; mm<=2*(L-1); mm++)
2606 Fmm_pad[mm+w_offset] = 0.0;
2607 for (mm=-(L-1); mm<=L-1; mm++)
2608 Fmm_pad[mm + w_offset] =
2609 Fmm[(m+Fmm_offset)*Fmm_stride + mm + Fmm_offset];
2612 for (mm=1; mm<=2*L-2; mm++)
2613 inout[mm + w_offset] = Fmm_pad[mm - 2*(L-1) - 1 + w_offset];
2614 for (mm=-2*(L-1); mm<=0; mm++)
2615 inout[mm + w_offset] = Fmm_pad[mm + 2*(L-1) + w_offset];
2616 fftw_execute_dft(plan_bwd, inout, inout);
2617 for (mm=0; mm<=2*L-2; mm++)
2618 Fmm_pad[mm + w_offset] = inout[mm - 2*(L-1) + w_offset];
2619 for (mm=-2*(L-1); mm<=-1; mm++)
2620 Fmm_pad[mm + w_offset] = inout[mm + 2*(L-1) + 1 + w_offset];
2623 for (r=-2*(L-1); r<=2*(L-1); r++)
2624 Fmm_pad[r + w_offset] *= wr[-r + w_offset];
2627 for (mm=1; mm<=2*L-2; mm++)
2628 inout[mm + w_offset] = Fmm_pad[mm - 2*(L-1) - 1 + w_offset];
2629 for (mm=-2*(L-1); mm<=0; mm++)
2630 inout[mm + w_offset] = Fmm_pad[mm + 2*(L-1) + w_offset];
2631 fftw_execute_dft(plan_fwd, inout, inout);
2632 for (mm=0; mm<=2*L-2; mm++)
2633 Fmm_pad[mm + w_offset] = inout[mm - 2*(L-1) + w_offset];
2634 for (mm=-2*(L-1); mm<=-1; mm++)
2635 Fmm_pad[mm + w_offset] = inout[mm + 2*(L-1) + 1 + w_offset];
2638 for (mm=-(L-1); mm<=L-1; mm++)
2639 Gmm[(mm+Gmm_offset)*Gmm_stride + m + Gmm_offset] =
2640 Fmm_pad[mm + w_offset] * 2.0 *
SSHT_PI / (4.0*L-3.0);
2643 fftw_destroy_plan(plan_bwd);
2644 fftw_destroy_plan(plan_fwd);
2656 for (el=0; el<=L-1; el++) {
2657 for (m=-el; m<=el; m++) {
2658 ssht_sampling_elm2ind(&ind, el, m);
2662 for (el=
MAX(L0, abs(spin)); el<=L-1; el++) {
2665 switch (dl_method) {
2668 if (el!=0 && el==
MAX(L0, abs(spin))) {
2669 for(eltmp=0; eltmp<=
MAX(L0, abs(spin)); eltmp++)
2672 eltmp, sqrt_tbl, signs);
2683 el, sqrt_tbl, signs);
2694 if (el!=0 && el==
MAX(L0, abs(spin))) {
2695 for(eltmp=0; eltmp<=
MAX(L0, abs(spin)); eltmp++)
2718 elfactor = sqrt((
double)(2.0*el+1.0)/(4.0*
SSHT_PI));
2719 el2pel = el *el + el;
2720 for (m=-el; m<=el; m++)
2721 inds[m + inds_offset] = el2pel + m;
2722 elssign = spin <= 0 ? 1.0 : signs[el];
2724 for (m=-el; m<=-1; m++) {
2726 ind = inds[m + inds_offset];
2730 * expsm[m + exps_offset]
2731 * signs[el] * dl[0*dl_stride - m + dl_offset]
2732 * elssign * dl[0*dl_stride - spinneg + dl_offset]
2733 * Gmm[(0+Gmm_offset)*Gmm_stride + m + Gmm_offset];
2735 for (m=0; m<=el; m++) {
2737 ind = inds[m + inds_offset];
2741 * expsm[m + exps_offset]
2742 * dl[0*dl_stride + m + dl_offset]
2743 * elssign * dl[0*dl_stride - spinneg + dl_offset]
2744 * Gmm[(0+Gmm_offset)*Gmm_stride + m + Gmm_offset];
2747 for (mm=1; mm<=el; mm++) {
2748 elmmsign = signs[el] * signs[mm];
2749 elssign = spin <= 0 ? 1.0 : elmmsign;
2751 for (m=-el; m<=-1; m++) {
2752 ind = inds[m + inds_offset];
2756 * expsm[m + exps_offset]
2757 * elmmsign * dl[mm*dl_stride - m + dl_offset]
2758 * elssign * dl[mm*dl_stride - spinneg + dl_offset]
2759 * ( Gmm[(mm+Gmm_offset)*Gmm_stride + m + Gmm_offset]
2761 * Gmm[(-mm+Gmm_offset)*Gmm_stride + m + Gmm_offset]);
2763 for (m=0; m<=el; m++) {
2764 ind = inds[m + inds_offset];
2768 * expsm[m + exps_offset]
2769 * dl[mm*dl_stride + m + dl_offset]
2770 * elssign * dl[mm*dl_stride - spinneg + dl_offset]
2771 * ( Gmm[(mm+Gmm_offset)*Gmm_stride + m + Gmm_offset]
2773 * Gmm[(-mm+Gmm_offset)*Gmm_stride + m + Gmm_offset]);
2800 printf(
"%s%s",
SSHT_PROMPT,
"Forward transform computed!");
2852 int el, m, mm, ind, ind_nm, t, r;
2854 double *sqrt_tbl, *signs;
2855 int el2pel, inds_offset;
2857 double ssign, elfactor;
2858 fftw_plan plan, plan_bwd, plan_fwd;
2864 int f_stride, Fmt_stride, Fmt_offset, Fmm_stride, Fmm_offset, Gmm_stride;
2867 int dl_offset, dl_stride;
2871 int elmmsign, elssign;
2876 sqrt_tbl = (
double*)calloc(2*(L-1)+2,
sizeof(double));
2878 signs = (
double*)calloc(L+1,
sizeof(
double));
2882 inds = (
int*)calloc(L,
sizeof(
int));
2886 for (el=0; el<=2*(L-1)+1; el++)
2887 sqrt_tbl[el] = sqrt((
double)el);
2888 for (m=0; m<=L-1; m=m+2) {
2892 ssign = signs[abs(spin)];
2893 spinneg = spin <= 0 ? spin : -spin;
2894 for (m=0; m<=L-1; m++)
2898 if (verbosity > 0) {
2900 "Computing forward transform using MW symmetric sampling with ");
2901 printf(
"%s%s%d%s%d%s\n",
SSHT_PROMPT,
"parameters (L,spin,reality) = (",
2902 L,
",", spin,
", TRUE)");
2905 "Using routine ssht_core_mw_forward_sov_conv_sym_ss_real...");
2915 in_real = (
double*)calloc(2*L,
sizeof(
double));
2919 plan = fftw_plan_dft_r2c_1d(2*L, in_real, out, FFTW_MEASURE);
2920 for (t=0; t<=L; t++) {
2921 memcpy(in_real, &f[t*f_stride], f_stride*
sizeof(
double));
2922 fftw_execute_dft_r2c(plan, in_real, out);
2924 Fmt[m*Fmt_stride + t] = out[m] / (2.0*L);
2929 fftw_destroy_plan(plan);
2932 for (m=0; m<=L; m++)
2933 for (t=L+1; t<=2*L-1; t++)
2934 Fmt[m*Fmt_stride + t] =
2935 signs[abs(m)] * ssign * Fmt[m*Fmt_stride + (2*L-t)];
2945 plan = fftw_plan_dft_1d(2*L, inout, inout, FFTW_FORWARD, FFTW_MEASURE);
2946 for (m=0; m<=L; m++) {
2947 memcpy(inout, &Fmt[m*Fmt_stride], Fmt_stride*
sizeof(
SSHT_COMPLEX(
double)));
2948 fftw_execute_dft(plan, inout, inout);
2949 for(mm=0; mm<=L; mm++)
2950 Fmm[m*Fmm_stride + mm + Fmm_offset] =
2951 inout[mm] / (2.0*L);
2952 for(mm=-(L-1); mm<=-1; mm++)
2953 Fmm[m*Fmm_stride + mm + Fmm_offset] =
2954 inout[mm+2*L-1+1] / (2.0*L);
2956 fftw_destroy_plan(plan);
2963 for (mm=-2*(L-1); mm<=2*(L-1); mm++)
2964 w[mm+w_offset] = ssht_sampling_weight_mw(mm);
2971 plan_bwd = fftw_plan_dft_1d(4*L-3, inout, inout, FFTW_BACKWARD, FFTW_MEASURE);
2972 plan_fwd = fftw_plan_dft_1d(4*L-3, inout, inout, FFTW_FORWARD, FFTW_MEASURE);
2973 for (mm=1; mm<=2*L-2; mm++)
2974 inout[mm + w_offset] = w[mm - 2*(L-1) - 1 + w_offset];
2975 for (mm=-2*(L-1); mm<=0; mm++)
2976 inout[mm + w_offset] = w[mm + 2*(L-1) + w_offset];
2977 fftw_execute_dft(plan_bwd, inout, inout);
2978 for (mm=0; mm<=2*L-2; mm++)
2979 wr[mm + w_offset] = inout[mm - 2*(L-1) + w_offset];
2980 for (mm=-2*(L-1); mm<=-1; mm++)
2981 wr[mm + w_offset] = inout[mm + 2*(L-1) + 1 + w_offset];
2991 for (m=0; m<=L-1; m++) {
2994 for (mm=-2*(L-1); mm<=-L; mm++)
2995 Fmm_pad[mm+w_offset] = 0.0;
2996 for (mm=L; mm<=2*(L-1); mm++)
2997 Fmm_pad[mm+w_offset] = 0.0;
2998 for (mm=-(L-1); mm<=L-1; mm++)
2999 Fmm_pad[mm + w_offset] =
3000 Fmm[m*Fmm_stride + mm + Fmm_offset];
3003 for (mm=1; mm<=2*L-2; mm++)
3004 inout[mm + w_offset] = Fmm_pad[mm - 2*(L-1) - 1 + w_offset];
3005 for (mm=-2*(L-1); mm<=0; mm++)
3006 inout[mm + w_offset] = Fmm_pad[mm + 2*(L-1) + w_offset];
3007 fftw_execute_dft(plan_bwd, inout, inout);
3008 for (mm=0; mm<=2*L-2; mm++)
3009 Fmm_pad[mm + w_offset] = inout[mm - 2*(L-1) + w_offset];
3010 for (mm=-2*(L-1); mm<=-1; mm++)
3011 Fmm_pad[mm + w_offset] = inout[mm + 2*(L-1) + 1 + w_offset];
3014 for (r=-2*(L-1); r<=2*(L-1); r++)
3015 Fmm_pad[r + w_offset] *= wr[-r + w_offset];
3018 for (mm=1; mm<=2*L-2; mm++)
3019 inout[mm + w_offset] = Fmm_pad[mm - 2*(L-1) - 1 + w_offset];
3020 for (mm=-2*(L-1); mm<=0; mm++)
3021 inout[mm + w_offset] = Fmm_pad[mm + 2*(L-1) + w_offset];
3022 fftw_execute_dft(plan_fwd, inout, inout);
3023 for (mm=0; mm<=2*L-2; mm++)
3024 Fmm_pad[mm + w_offset] = inout[mm - 2*(L-1) + w_offset];
3025 for (mm=-2*(L-1); mm<=-1; mm++)
3026 Fmm_pad[mm + w_offset] = inout[mm + 2*(L-1) + 1 + w_offset];
3029 for (mm=-(L-1); mm<=L-1; mm++)
3030 Gmm[(mm+Fmm_offset)*Gmm_stride + m] =
3031 Fmm_pad[mm + w_offset] * 2.0 *
SSHT_PI / (4.0*L-3.0);
3034 fftw_destroy_plan(plan_bwd);
3035 fftw_destroy_plan(plan_fwd);
3047 for (el=0; el<=L-1; el++) {
3048 for (m=0; m<=el; m++) {
3049 ssht_sampling_elm2ind(&ind, el, m);
3053 for (el=
MAX(L0, abs(spin)); el<=L-1; el++) {
3056 switch (dl_method) {
3059 if (el!=0 && el==
MAX(L0, abs(spin))) {
3060 for(eltmp=0; eltmp<=
MAX(L0, abs(spin)); eltmp++)
3063 eltmp, sqrt_tbl, signs);
3074 el, sqrt_tbl, signs);
3085 if (el!=0 && el==
MAX(L0, abs(spin))) {
3086 for(eltmp=0; eltmp<=
MAX(L0, abs(spin)); eltmp++)
3109 elfactor = sqrt((
double)(2.0*el+1.0)/(4.0*
SSHT_PI));
3110 el2pel = el *el + el;
3111 for (m=0; m<=el; m++)
3112 inds[m + inds_offset] = el2pel + m;
3113 elssign = spin <= 0 ? 1.0 : signs[el];
3115 for (m=0; m<=el; m++) {
3117 ind = inds[m + inds_offset];
3122 * dl[0*dl_stride + m + dl_offset]
3123 * elssign * dl[0*dl_stride - spinneg + dl_offset]
3124 * Gmm[(0+Fmm_offset)*Gmm_stride + m];
3127 for (mm=1; mm<=el; mm++) {
3128 elmmsign = signs[el] * signs[mm];
3129 elssign = spin <= 0 ? 1.0 : elmmsign;
3131 for (m=0; m<=el; m++) {
3132 ind = inds[m + inds_offset];
3137 * dl[mm*dl_stride + m + dl_offset]
3138 * elssign * dl[mm*dl_stride - spinneg + dl_offset]
3139 * ( Gmm[(mm+Fmm_offset)*Gmm_stride + m]
3141 * Gmm[(-mm+Fmm_offset)*Gmm_stride + m]);
3149 for (el=abs(spin); el<=L-1; el++) {
3150 for (m=1; m<=el; m++) {
3151 ssht_sampling_elm2ind(&ind, el, m);
3152 ssht_sampling_elm2ind(&ind_nm, el, -m);
3153 flm[ind_nm] = signs[m] * conj(flm[ind]);
3176 printf(
"%s%s",
SSHT_PROMPT,
"Forward transform computed!");
3218 int t, f_stride = 2*L;
3226 dl_method, verbosity);
3229 for (t=1; t<=L-1; t++)
3230 memcpy(&f[(t-1)*f_stride], &f_full[t*f_stride],
3234 *f_sp = f_full[L*f_stride + 0];
3270 int t, f_stride = 2*L;
3273 f_full = (
double*)calloc((L+1)*(2*L),
sizeof(
double));
3278 dl_method, verbosity);
3281 for (t=1; t<=L-1; t++)
3282 memcpy(&f[(t-1)*f_stride], &f_full[t*f_stride],
3283 (2*L)*
sizeof(
double));
3285 *f_sp = f_full[L*f_stride + 0];
3324 int t, p, f_stride = 2*L;
3330 for (t=1; t<=L-1; t++)
3331 memcpy(&f_full[t*f_stride], &f[(t-1)*f_stride],
3335 for (p=0; p<=2*L-1; p++) {
3337 f_full[0*f_stride + p] = f_np * cexp(-I*spin*(phi-phi_np));
3338 f_full[L*f_stride + p] = f_sp * cexp(I*spin*(phi-phi_sp));
3343 dl_method, verbosity);
3378 int t, p, f_stride = 2*L;
3381 f_full = (
double*)calloc((L+1)*(2*L),
sizeof(
double));
3383 for (t=1; t<=L-1; t++)
3384 memcpy(&f_full[t*f_stride], &f[(t-1)*f_stride],
3385 (2*L)*
sizeof(
double));
3388 for (p=0; p<=2*L-1; p++) {
3389 f_full[0*f_stride + p] = f_np;
3390 f_full[L*f_stride + p] = f_sp;
3395 dl_method, verbosity);
3422 int L,
int spin,
int verbosity) {
3424 int t, p, m, el, ind;
3425 int ftm_stride, ftm_offset, f_stride;
3426 double *dlm1p1_line, *dl_line;
3428 double *sqrt_tbl, *signs;
3430 double theta, ssign, elfactor;
3432 double *thetas, *weights;
3435 sqrt_tbl = (
double*)calloc(2*(L-1)+2,
sizeof(double));
3437 signs = (
double*)calloc(L+1,
sizeof(
double));
3441 for (el=0; el<=2*(L-1)+1; el++)
3442 sqrt_tbl[el] = sqrt((
double)el);
3443 for (m=0; m<=L-1; m=m+2) {
3447 ssign = signs[abs(spin)];
3450 if (verbosity > 0) {
3452 "Computing inverse transform using GL sampling with ");
3453 printf(
"%s%s%d%s%d%s\n",
SSHT_PROMPT,
"parameters (L,spin,reality) = (",
3454 L,
",", spin,
", FALSE)");
3457 "Using routine ssht_core_gl_inverse_sov...");
3461 thetas = (
double*)calloc(L,
sizeof(
double));
3463 weights = (
double*)calloc(L,
sizeof(
double));
3472 dlm1p1_line = (
double*)calloc(2*L-1,
sizeof(
double));
3474 dl_line = (
double*)calloc(2*L-1,
sizeof(
double));
3476 for (t=0; t<=L-1; t++) {
3478 for (el=abs(spin); el<=L-1; el++) {
3479 elfactor = sqrt((
double)(2.0*el+1.0)/(4.0*
SSHT_PI));
3483 theta, L, -spin, el,
3487 dl_line = dlm1p1_line;
3488 dlm1p1_line = dl_ptr;
3490 for (m=-el; m<=el; m++) {
3491 ssht_sampling_elm2ind(&ind, el, m);
3492 ftm[t*ftm_stride + m + ftm_offset] +=
3513 plan = fftw_plan_dft_1d(2*L-1, inout, inout, FFTW_BACKWARD, FFTW_MEASURE);
3514 for (t=0; t<=L-1; t++) {
3515 for (m=0; m<=L-1; m++)
3516 inout[m] = ftm[t*ftm_stride + m + ftm_offset];
3517 for (m=-(L-1); m<=-1; m++)
3518 inout[m+2*L-1] = ftm[t*ftm_stride + m + ftm_offset];
3519 fftw_execute_dft(plan, inout, inout);
3520 for (p=0; p<=2*L-2; p++)
3521 f[t*f_stride + p] = inout[p];
3523 fftw_destroy_plan(plan);
3533 printf(
"%s%s",
SSHT_PROMPT,
"Inverse transform computed!");
3553 int L,
int verbosity) {
3555 int t, p, m, el, ind;
3556 int ftm_stride, ftm_offset, f_stride;
3557 double *dlm1p1_line, *dl_line;
3559 double *sqrt_tbl, *signs;
3563 double theta, ssign, elfactor;
3565 double *thetas, *weights;
3569 sqrt_tbl = (
double*)calloc(2*(L-1)+2,
sizeof(double));
3571 signs = (
double*)calloc(L+1,
sizeof(
double));
3575 for (el=0; el<=2*(L-1)+1; el++)
3576 sqrt_tbl[el] = sqrt((
double)el);
3577 for (m=0; m<=L-1; m=m+2) {
3581 ssign = signs[abs(spin)];
3584 if (verbosity > 0) {
3586 "Computing inverse transform using GL sampling with ");
3587 printf(
"%s%s%d%s%d%s\n",
SSHT_PROMPT,
"parameters (L,spin,reality) = (",
3588 L,
",", spin,
", TRUE)");
3591 "Using routine ssht_core_gl_inverse_sov_real...");
3595 thetas = (
double*)calloc(L,
sizeof(
double));
3597 weights = (
double*)calloc(L,
sizeof(
double));
3606 dlm1p1_line = (
double*)calloc(L,
sizeof(
double));
3608 dl_line = (
double*)calloc(L,
sizeof(
double));
3610 for (t=0; t<=L-1; t++) {
3612 for (el=abs(spin); el<=L-1; el++) {
3613 elfactor = sqrt((
double)(2.0*el+1.0)/(4.0*
SSHT_PI));
3617 theta, L, -spin, el,
3621 dl_line = dlm1p1_line;
3622 dlm1p1_line = dl_ptr;
3624 for (m=0; m<=el; m++) {
3625 ssht_sampling_elm2ind(&ind, el, m);
3626 ftm[t*ftm_stride + m + ftm_offset] +=
3646 out = (
double*)calloc(2*L-1,
sizeof(
double));
3648 plan = fftw_plan_dft_c2r_1d(2*L-1, in, out, FFTW_MEASURE);
3650 for (t=0; t<=L-1; t++) {
3651 memcpy(in, &ftm[t*ftm_stride], L*
sizeof(
SSHT_COMPLEX(
double)));
3652 fftw_execute_dft_c2r(plan, in, out);
3653 for (p=0; p<=2*L-2; p++)
3654 f[t*f_stride + p] = out[p];
3656 fftw_destroy_plan(plan);
3667 printf(
"%s%s",
SSHT_PROMPT,
"Inverse transform computed!");
3686 int L,
int spin,
int verbosity) {
3690 double *dlm1p1_line, *dl_line;
3692 int el2pel, inds_offset;
3694 double *sqrt_tbl, *signs;
3695 int Ftm_stride, Ftm_offset;
3697 double theta, ssign, elfactor;
3699 double *thetas, *weights;
3703 sqrt_tbl = (
double*)calloc(2*(L-1)+2,
sizeof(double));
3705 signs = (
double*)calloc(L+1,
sizeof(
double));
3707 inds = (
int*)calloc(2*L-1,
sizeof(
int));
3711 for (el=0; el<=2*(L-1)+1; el++)
3712 sqrt_tbl[el] = sqrt((
double)el);
3713 for (m=0; m<=L-1; m=m+2) {
3717 ssign = signs[abs(spin)];
3720 if (verbosity > 0) {
3722 "Computing forward transform using GL sampling with ");
3723 printf(
"%s%s%d%s%d%s\n",
SSHT_PROMPT,
"parameters (L,spin,reality) = (",
3724 L,
",", spin,
", FALSE)");
3727 "Using routine ssht_core_gl_forward_sov...");
3731 thetas = (
double*)calloc(L,
sizeof(
double));
3733 weights = (
double*)calloc(L,
sizeof(
double));
3745 plan = fftw_plan_dft_1d(2*L-1, inout, inout, FFTW_FORWARD, FFTW_MEASURE);
3746 for (t=0; t<=L-1; t++) {
3747 memcpy(inout, &f[t*f_stride], f_stride*
sizeof(
SSHT_COMPLEX(
double)));
3748 fftw_execute_dft(plan, inout, inout);
3749 for(m=0; m<=L-1; m++)
3750 Ftm[t*Ftm_stride + m + Ftm_offset] =
3751 inout[m] * 2.0 *
SSHT_PI / (2.0*L-1.0) ;
3752 for(m=-(L-1); m<=-1; m++)
3753 Ftm[t*Ftm_stride + m + Ftm_offset] =
3754 inout[m+2*L-1] * 2.0 *
SSHT_PI / (2.0*L-1.0);
3756 fftw_destroy_plan(plan);
3759 dlm1p1_line = (
double*)calloc(2*L-1,
sizeof(
double));
3761 dl_line = (
double*)calloc(2*L-1,
sizeof(
double));
3764 for (el=0; el<=L-1; el++) {
3765 for (m=-el; m<=el; m++) {
3766 ssht_sampling_elm2ind(&ind, el, m);
3770 for (t=0; t<=L-1; t++) {
3773 for (el=abs(spin); el<=L-1; el++) {
3774 elfactor = sqrt((
double)(2.0*el+1.0)/(4.0*
SSHT_PI));
3775 el2pel = el *el + el;
3776 for (m=-el; m<=el; m++)
3777 inds[m + inds_offset] = el2pel + m;
3781 theta, L, -spin, el,
3785 dl_line = dlm1p1_line;
3786 dlm1p1_line = dl_ptr;
3788 for (m=-el; m<=el; m++) {
3789 ind = inds[m + inds_offset];
3795 * Ftm[t*Ftm_stride + m + Ftm_offset];
3813 printf(
"%s%s",
SSHT_PROMPT,
"Forward transform computed!");
3833 int L,
int verbosity) {
3835 int t, m, el, ind, ind_nm;
3837 double *dlm1p1_line, *dl_line;
3839 int el2pel, inds_offset;
3841 double *sqrt_tbl, *signs;
3842 int Ftm_stride, Ftm_offset;
3846 double theta, ssign, elfactor;
3848 double *thetas, *weights;
3853 sqrt_tbl = (
double*)calloc(2*(L-1)+2,
sizeof(double));
3855 signs = (
double*)calloc(L+1,
sizeof(
double));
3857 inds = (
int*)calloc(L,
sizeof(
int));
3861 for (el=0; el<=2*(L-1)+1; el++)
3862 sqrt_tbl[el] = sqrt((
double)el);
3863 for (m=0; m<=L-1; m=m+2) {
3867 ssign = signs[abs(spin)];
3870 if (verbosity > 0) {
3872 "Computing forward transform using GL sampling with ");
3873 printf(
"%s%s%d%s%d%s\n",
SSHT_PROMPT,
"parameters (L,spin,reality) = (",
3874 L,
",", spin,
", TRUE)");
3877 "Using routine ssht_core_gl_forward_sov_real...");
3881 thetas = (
double*)calloc(L,
sizeof(
double));
3883 weights = (
double*)calloc(L,
sizeof(
double));
3893 in_real = (
double*)calloc(2*L-1,
sizeof(
double));
3897 plan = fftw_plan_dft_r2c_1d(2*L-1, in_real, out, FFTW_MEASURE);
3898 for (t=0; t<=L-1; t++) {
3899 memcpy(in_real, &f[t*f_stride], f_stride*
sizeof(
double));
3900 fftw_execute_dft_r2c(plan, in_real, out);
3901 for(m=0; m<=L-1; m++)
3902 Ftm[t*Ftm_stride + m + Ftm_offset] =
3903 out[m] * 2.0 *
SSHT_PI / (2.0*L-1.0);
3907 fftw_destroy_plan(plan);
3910 dlm1p1_line = (
double*)calloc(L,
sizeof(
double));
3912 dl_line = (
double*)calloc(L,
sizeof(
double));
3915 for (el=0; el<=L-1; el++) {
3916 for (m=0; m<=el; m++) {
3917 ssht_sampling_elm2ind(&ind, el, m);
3921 for (t=0; t<=L-1; t++) {
3924 for (el=abs(spin); el<=L-1; el++) {
3925 elfactor = sqrt((
double)(2.0*el+1.0)/(4.0*
SSHT_PI));
3926 el2pel = el *el + el;
3927 for (m=0; m<=el; m++)
3928 inds[m + inds_offset] = el2pel + m;
3932 theta, L, -spin, el,
3937 dl_line = dlm1p1_line;
3938 dlm1p1_line = dl_ptr;
3940 for (m=0; m<=el; m++) {
3941 ind = inds[m + inds_offset];
3947 * Ftm[t*Ftm_stride + m + Ftm_offset];
3953 for (el=abs(spin); el<=L-1; el++) {
3954 for (m=1; m<=el; m++) {
3955 ssht_sampling_elm2ind(&ind, el, m);
3956 ssht_sampling_elm2ind(&ind_nm, el, -m);
3957 flm[ind_nm] = signs[m] * conj(flm[ind]);
3973 printf(
"%s%s",
SSHT_PROMPT,
"Forward transform computed!");
3997 int L,
int spin,
int verbosity) {
3999 int t, p, m, el, ind;
4000 int ftm_stride, ftm_offset, f_stride;
4001 double *dlm1p1_line, *dl_line;
4003 double *sqrt_tbl, *signs;
4005 double theta, ssign, elfactor;
4009 sqrt_tbl = (
double*)calloc(2*(L-1)+2,
sizeof(double));
4011 signs = (
double*)calloc(L+1,
sizeof(
double));
4015 for (el=0; el<=2*(L-1)+1; el++)
4016 sqrt_tbl[el] = sqrt((
double)el);
4017 for (m=0; m<=L-1; m=m+2) {
4021 ssign = signs[abs(spin)];
4024 if (verbosity > 0) {
4026 "Computing inverse transform using DH sampling with ");
4027 printf(
"%s%s%d%s%d%s\n",
SSHT_PROMPT,
"parameters (L,spin,reality) = (",
4028 L,
",", spin,
", FALSE)");
4031 "Using routine ssht_core_dh_inverse_sov...");
4039 dlm1p1_line = (
double*)calloc(2*L-1,
sizeof(
double));
4041 dl_line = (
double*)calloc(2*L-1,
sizeof(
double));
4043 for (t=0; t<=2*L-1; t++) {
4045 for (el=abs(spin); el<=L-1; el++) {
4046 elfactor = sqrt((
double)(2.0*el+1.0)/(4.0*
SSHT_PI));
4050 theta, L, -spin, el,
4054 dl_line = dlm1p1_line;
4055 dlm1p1_line = dl_ptr;
4057 for (m=-el; m<=el; m++) {
4058 ssht_sampling_elm2ind(&ind, el, m);
4059 ftm[t*ftm_stride + m + ftm_offset] +=
4076 plan = fftw_plan_dft_1d(2*L-1, inout, inout, FFTW_BACKWARD, FFTW_MEASURE);
4077 for (t=0; t<=2*L-1; t++) {
4078 for (m=0; m<=L-1; m++)
4079 inout[m] = ftm[t*ftm_stride + m + ftm_offset];
4080 for (m=-(L-1); m<=-1; m++)
4081 inout[m+2*L-1] = ftm[t*ftm_stride + m + ftm_offset];
4082 fftw_execute_dft(plan, inout, inout);
4083 for (p=0; p<=2*L-2; p++)
4084 f[t*f_stride + p] = inout[p];
4086 fftw_destroy_plan(plan);
4096 printf(
"%s%s",
SSHT_PROMPT,
"Inverse transform computed!");
4116 int L,
int verbosity) {
4118 int t, p, m, el, ind;
4119 int ftm_stride, ftm_offset, f_stride;
4120 double *dlm1p1_line, *dl_line;
4122 double *sqrt_tbl, *signs;
4126 double theta, ssign, elfactor;
4131 sqrt_tbl = (
double*)calloc(2*(L-1)+2,
sizeof(double));
4133 signs = (
double*)calloc(L+1,
sizeof(
double));
4137 for (el=0; el<=2*(L-1)+1; el++)
4138 sqrt_tbl[el] = sqrt((
double)el);
4139 for (m=0; m<=L-1; m=m+2) {
4143 ssign = signs[abs(spin)];
4146 if (verbosity > 0) {
4148 "Computing inverse transform using DH sampling with ");
4149 printf(
"%s%s%d%s%d%s\n",
SSHT_PROMPT,
"parameters (L,spin,reality) = (",
4150 L,
",", spin,
", TRUE)");
4153 "Using routine ssht_core_dh_inverse_sov_real...");
4161 dlm1p1_line = (
double*)calloc(L,
sizeof(
double));
4163 dl_line = (
double*)calloc(L,
sizeof(
double));
4165 for (t=0; t<=2*L-1; t++) {
4167 for (el=abs(spin); el<=L-1; el++) {
4168 elfactor = sqrt((
double)(2.0*el+1.0)/(4.0*
SSHT_PI));
4172 theta, L, -spin, el,
4176 dl_line = dlm1p1_line;
4177 dlm1p1_line = dl_ptr;
4179 for (m=0; m<=el; m++) {
4180 ssht_sampling_elm2ind(&ind, el, m);
4181 ftm[t*ftm_stride + m + ftm_offset] +=
4197 out = (
double*)calloc(2*L-1,
sizeof(
double));
4199 plan = fftw_plan_dft_c2r_1d(2*L-1, in, out, FFTW_MEASURE);
4201 for (t=0; t<=2*L-1; t++) {
4202 memcpy(in, &ftm[t*ftm_stride], L*
sizeof(
SSHT_COMPLEX(
double)));
4203 fftw_execute_dft_c2r(plan, in, out);
4204 for (p=0; p<=2*L-2; p++)
4205 f[t*f_stride + p] = out[p];
4207 fftw_destroy_plan(plan);
4218 printf(
"%s%s",
SSHT_PROMPT,
"Inverse transform computed!");
4237 int L,
int spin,
int verbosity) {
4241 double *dlm1p1_line, *dl_line;
4243 int el2pel, inds_offset;
4245 double *sqrt_tbl, *signs;
4246 int Ftm_stride, Ftm_offset;
4248 double theta, ssign, elfactor;
4253 sqrt_tbl = (
double*)calloc(2*(L-1)+2,
sizeof(double));
4255 signs = (
double*)calloc(L+1,
sizeof(
double));
4257 inds = (
int*)calloc(2*L-1,
sizeof(
int));
4261 for (el=0; el<=2*(L-1)+1; el++)
4262 sqrt_tbl[el] = sqrt((
double)el);
4263 for (m=0; m<=L-1; m=m+2) {
4267 ssign = signs[abs(spin)];
4270 if (verbosity > 0) {
4272 "Computing forward transform using DH sampling with ");
4273 printf(
"%s%s%d%s%d%s\n",
SSHT_PROMPT,
"parameters (L,spin,reality) = (",
4274 L,
",", spin,
", FALSE)");
4277 "Using routine ssht_core_dh_forward_sov...");
4288 plan = fftw_plan_dft_1d(2*L-1, inout, inout, FFTW_FORWARD, FFTW_MEASURE);
4289 for (t=0; t<=2*L-1; t++) {
4290 memcpy(inout, &f[t*f_stride], f_stride*
sizeof(
SSHT_COMPLEX(
double)));
4291 fftw_execute_dft(plan, inout, inout);
4292 for(m=0; m<=L-1; m++)
4293 Ftm[t*Ftm_stride + m + Ftm_offset] =
4294 inout[m] * 2.0 *
SSHT_PI / (2.0*L-1.0) ;
4295 for(m=-(L-1); m<=-1; m++)
4296 Ftm[t*Ftm_stride + m + Ftm_offset] =
4297 inout[m+2*L-1] * 2.0 *
SSHT_PI / (2.0*L-1.0);
4299 fftw_destroy_plan(plan);
4302 dlm1p1_line = (
double*)calloc(2*L-1,
sizeof(
double));
4304 dl_line = (
double*)calloc(2*L-1,
sizeof(
double));
4307 for (el=0; el<=L-1; el++) {
4308 for (m=-el; m<=el; m++) {
4309 ssht_sampling_elm2ind(&ind, el, m);
4313 for (t=0; t<=2*L-1; t++) {
4316 for (el=abs(spin); el<=L-1; el++) {
4317 elfactor = sqrt((
double)(2.0*el+1.0)/(4.0*
SSHT_PI));
4318 el2pel = el *el + el;
4319 for (m=-el; m<=el; m++)
4320 inds[m + inds_offset] = el2pel + m;
4324 theta, L, -spin, el,
4328 dl_line = dlm1p1_line;
4329 dlm1p1_line = dl_ptr;
4331 for (m=-el; m<=el; m++) {
4332 ind = inds[m + inds_offset];
4338 * Ftm[t*Ftm_stride + m + Ftm_offset];
4354 printf(
"%s%s",
SSHT_PROMPT,
"Forward transform computed!");
4374 int L,
int verbosity) {
4376 int t, m, el, ind, ind_nm;
4378 double *dlm1p1_line, *dl_line;
4380 int el2pel, inds_offset;
4382 double *sqrt_tbl, *signs;
4383 int Ftm_stride, Ftm_offset;
4387 double theta, ssign, elfactor;
4393 sqrt_tbl = (
double*)calloc(2*(L-1)+2,
sizeof(double));
4395 signs = (
double*)calloc(L+1,
sizeof(
double));
4397 inds = (
int*)calloc(L,
sizeof(
int));
4401 for (el=0; el<=2*(L-1)+1; el++)
4402 sqrt_tbl[el] = sqrt((
double)el);
4403 for (m=0; m<=L-1; m=m+2) {
4407 ssign = signs[abs(spin)];
4410 if (verbosity > 0) {
4412 "Computing forward transform using GL sampling with ");
4413 printf(
"%s%s%d%s%d%s\n",
SSHT_PROMPT,
"parameters (L,spin,reality) = (",
4414 L,
",", spin,
", TRUE)");
4417 "Using routine ssht_core_gl_forward_sov_real...");
4426 in_real = (
double*)calloc(2*L-1,
sizeof(
double));
4430 plan = fftw_plan_dft_r2c_1d(2*L-1, in_real, out, FFTW_MEASURE);
4431 for (t=0; t<=2*L-1; t++) {
4432 memcpy(in_real, &f[t*f_stride], f_stride*
sizeof(
double));
4433 fftw_execute_dft_r2c(plan, in_real, out);
4434 for(m=0; m<=L-1; m++)
4435 Ftm[t*Ftm_stride + m + Ftm_offset] =
4436 out[m] * 2.0 *
SSHT_PI / (2.0*L-1.0);
4440 fftw_destroy_plan(plan);
4443 dlm1p1_line = (
double*)calloc(L,
sizeof(
double));
4445 dl_line = (
double*)calloc(L,
sizeof(
double));
4448 for (el=0; el<=L-1; el++) {
4449 for (m=0; m<=el; m++) {
4450 ssht_sampling_elm2ind(&ind, el, m);
4454 for (t=0; t<=2*L-1; t++) {
4457 for (el=abs(spin); el<=L-1; el++) {
4458 elfactor = sqrt((
double)(2.0*el+1.0)/(4.0*
SSHT_PI));
4459 el2pel = el *el + el;
4460 for (m=0; m<=el; m++)
4461 inds[m + inds_offset] = el2pel + m;
4465 theta, L, -spin, el,
4470 dl_line = dlm1p1_line;
4471 dlm1p1_line = dl_ptr;
4473 for (m=0; m<=el; m++) {
4474 ind = inds[m + inds_offset];
4480 * Ftm[t*Ftm_stride + m + Ftm_offset];
4486 for (el=abs(spin); el<=L-1; el++) {
4487 for (m=1; m<=el; m++) {
4488 ssht_sampling_elm2ind(&ind, el, m);
4489 ssht_sampling_elm2ind(&ind_nm, el, -m);
4490 flm[ind_nm] = signs[m] * conj(flm[ind]);
4504 printf(
"%s%s",
SSHT_PROMPT,
"Forward transform computed!");