My Project
ssht_dl.c
Go to the documentation of this file.
1 // SSHT package to perform spin spherical harmonic transforms
2 // Copyright (C) 2011 Jason McEwen
3 // See LICENSE.txt for license details
4 
5 
42 #include <stdlib.h>
43 #include <math.h>
44 #include "ssht_types.h"
45 #include "ssht_error.h"
46 #include "ssht_dl.h"
47 
48 double logfact(int n);
49 
50 
64 double* ssht_dl_calloc(int L, ssht_dl_size_t dl_size) {
65 
66  double *dl;
67 
68  // Allocate memory.
69  switch (dl_size) {
70 
71  case SSHT_DL_QUARTER:
72  dl = (double*)calloc(L*L, sizeof(double));
73  break;
74 
76  dl = (double*)calloc((L+2)*(L+2), sizeof(double));
77  break;
78 
79  case SSHT_DL_HALF:
80  dl = (double*)calloc(L*(2*L-1), sizeof(double));
81  break;
82 
83  case SSHT_DL_FULL:
84  dl = (double*)calloc((2*L-1)*(2*L-1), sizeof(double));
85  break;
86 
87  default:
88  SSHT_ERROR_GENERIC("Invalid dl size")
89 
90  }
91 
93  return dl;
94 
95 }
96 
97 
111 int ssht_dl_get_offset(int L, ssht_dl_size_t dl_size) {
112 
113  switch (dl_size) {
114  case SSHT_DL_QUARTER:
115  return 0;
116 
118  return L - 1;
119 
120  case SSHT_DL_HALF:
121  return L - 1;
122 
123  case SSHT_DL_FULL:
124  return L - 1;
125 
126  default:
127  SSHT_ERROR_GENERIC("Invalid dl size")
128  }
129 
130 }
131 
132 
146 int ssht_dl_get_stride(int L, ssht_dl_size_t dl_size) {
147 
148  switch (dl_size) {
149  case SSHT_DL_QUARTER:
150  return L;
151 
153  return L + 2;
154 
155  case SSHT_DL_HALF:
156  return 2*L - 1;
157 
158  case SSHT_DL_FULL:
159  return 2*L - 1;
160 
161  default:
162  SSHT_ERROR_GENERIC("Invalid dl size")
163  }
164 
165 }
166 
167 
190 void ssht_dl_beta_risbo_full_table(double *dl, double beta, int L,
191  ssht_dl_size_t dl_size,
192  int el, double *sqrt_tbl) {
193 
194  int offset, stride;
195  double cosb, sinb, coshb, sinhb;
196  int i, j, k;
197  double rj, dlj, ddj;
198  double *dd;
199 
200  // Get mm offset and stride for accessing dl data.
201  offset = ssht_dl_get_offset(L, dl_size);
202  stride = ssht_dl_get_stride(L, dl_size);
203 
204  // Compute Wigner plane.
205  if (el == 0) {
206 
207  dl[(0+offset)*stride + 0 + offset] = 1.0;
208 
209  }
210  else if (el == 1) {
211 
212  cosb = cos(beta);
213  sinb = sin(beta);
214  coshb = cos(beta / 2.0);
215  sinhb = sin(beta / 2.0);
216 
217  dl[(-1 + offset)*stride - 1 + offset] = coshb * coshb;
218  dl[(-1 + offset)*stride + 0 + offset] = sinb / SSHT_SQRT2;
219  dl[(-1 + offset)*stride + 1 + offset] = sinhb * sinhb;
220 
221  dl[(0 + offset)*stride - 1 + offset] = -sinb / SSHT_SQRT2;
222  dl[(0 + offset)*stride + 0 + offset] = cosb;
223  dl[(0 + offset)*stride + 1 + offset] = sinb / SSHT_SQRT2;
224 
225  dl[(1 + offset)*stride - 1 + offset] = sinhb * sinhb;
226  dl[(1 + offset)*stride + 0 + offset] = -sinb / SSHT_SQRT2;
227  dl[(1 + offset)*stride + 1 + offset] = coshb * coshb;
228 
229  }
230  else {
231 
232  coshb = -cos(beta / 2.0);
233  sinhb = sin(beta / 2.0);
234 
235  // Initialise the plane of the dl-matrix to 0.0 for the recursion
236  // from l - 1 to l - 1/2.
237  dd = (double*)calloc((2*el+2)*(2*el+2), sizeof(double));
239 
240  j = 2*el - 1;
241  rj = (double) j;
242  for (k=0; k<=j-1; k++) {
243  for (i=0; i<=j-1; i++) {
244  dlj = dl[(k-(el-1)+offset)*stride + i-(el-1) + offset] / rj;
245  dd[k*(2*el+2) + i] +=
246  sqrt_tbl[j-i] * sqrt_tbl[j-k] * dlj * coshb;
247  dd[k*(2*el+2) + i+1] -=
248  sqrt_tbl[i+1] * sqrt_tbl[j-k] * dlj * sinhb;
249  dd[(k+1)*(2*el+2) + i] +=
250  sqrt_tbl[j-i] * sqrt_tbl[k+1] * dlj * sinhb;
251  dd[(k+1)*(2*el+2) + i+1] +=
252  sqrt_tbl[i+1] * sqrt_tbl[k+1] * dlj * coshb;
253  }
254  }
255 
256  // Having constructed the d^(l+1/2) matrix in dd, do the second
257  // half-step recursion from dd to dl. Start by initilalising
258  // the plane of the dl-matrix to 0.0.
259  for (k=-el; k<=el; k++)
260  for (i=-el; i<=el; i++)
261  dl[(k+offset)*stride + i + offset] = 0.0;
262 
263  j = 2*el;
264  rj = (double) j;
265  for (k=0; k<=j-1; k++) {
266  for (i=0; i<=j-1; i++) {
267  ddj = dd[k*(2*el+2) + i] / rj;
268  dl[(k-el+offset)*stride + i-el + offset] +=
269  sqrt_tbl[j-i] * sqrt_tbl[j-k] * ddj * coshb;
270  dl[(k-el+offset)*stride + i+1-el + offset] -=
271  sqrt_tbl[i+1] * sqrt_tbl[j-k] * ddj * sinhb;
272  dl[(k+1-el+offset)*stride + i-el + offset] +=
273  sqrt_tbl[j-i] * sqrt_tbl[k+1] * ddj * sinhb;
274  dl[(k+1-el+offset)*stride + i+1-el + offset] +=
275  sqrt_tbl[i+1] * sqrt_tbl[k+1] * ddj * coshb;
276  }
277  }
278 
279  // Free temporary memory.
280  free(dd);
281  }
282 
283 }
284 
285 
313 void ssht_dl_beta_risbo_half_table(double *dl, double beta, int L,
314  ssht_dl_size_t dl_size,
315  int el, double *sqrt_tbl,
316  double *signs) {
317 
318  int offset, stride;
319  double cosb, sinb, coshb, sinhb;
320  int i, j, k;
321  double rj, dlj, ddj;
322  double *dd;
323  int m, mm;
324 
325  // Get mm offset and stride for accessing dl data.
326  offset = ssht_dl_get_offset(L, dl_size);
327  stride = ssht_dl_get_stride(L, dl_size);
328 
329  // Compute Wigner plane.
330  if (el == 0) {
331 
332  dl[(0+offset)*stride + 0 + offset] = 1.0;
333 
334  }
335  else if (el == 1) {
336 
337  cosb = cos(beta);
338  sinb = sin(beta);
339  coshb = cos(beta / 2.0);
340  sinhb = sin(beta / 2.0);
341 
342  dl[(-1 + offset)*stride - 1 + offset] = coshb * coshb;
343  dl[(-1 + offset)*stride + 0 + offset] = sinb / SSHT_SQRT2;
344  dl[(-1 + offset)*stride + 1 + offset] = sinhb * sinhb;
345 
346  dl[(0 + offset)*stride - 1 + offset] = -sinb / SSHT_SQRT2;
347  dl[(0 + offset)*stride + 0 + offset] = cosb;
348  dl[(0 + offset)*stride + 1 + offset] = sinb / SSHT_SQRT2;
349 
350  dl[(1 + offset)*stride - 1 + offset] = sinhb * sinhb;
351  dl[(1 + offset)*stride + 0 + offset] = -sinb / SSHT_SQRT2;
352  dl[(1 + offset)*stride + 1 + offset] = coshb * coshb;
353 
354  }
355  else {
356 
357  coshb = -cos(beta / 2.0);
358  sinhb = sin(beta / 2.0);
359 
360  // Initialise the plane of the dl-matrix to 0.0 for the recursion
361  // from l - 1 to l - 1/2.
362  dd = (double*)calloc((2*el+2)*(2*el+2), sizeof(double));
364 
365  j = 2*el - 1;
366  rj = (double) j;
367  for (k=0; k<=j-1; k++) {
368  for (i=0; i<=el; i++) {
369  dlj = dl[(k-(el-1)+offset)*stride + i-(el-1) + offset] / rj;
370  dd[k*(2*el+2) + i] +=
371  sqrt_tbl[j-i] * sqrt_tbl[j-k] * dlj * coshb;
372  dd[k*(2*el+2) + i+1] -=
373  sqrt_tbl[i+1] * sqrt_tbl[j-k] * dlj * sinhb;
374  dd[(k+1)*(2*el+2) + i] +=
375  sqrt_tbl[j-i] * sqrt_tbl[k+1] * dlj * sinhb;
376  dd[(k+1)*(2*el+2) + i+1] +=
377  sqrt_tbl[i+1] * sqrt_tbl[k+1] * dlj * coshb;
378  }
379  }
380 
381  // Having constructed the d^(l+1/2) matrix in dd, do the second
382  // half-step recursion from dd to dl. Start by initilalising
383  // the plane of the dl-matrix to 0.0.
384  for (k=-el; k<=el; k++)
385  for (i=-el; i<=el; i++)
386  dl[(k+offset)*stride + i + offset] = 0.0;
387 
388  j = 2*el;
389  rj = (double) j;
390  for (k=0; k<=j-1; k++) {
391  for (i=0; i<=el; i++) {
392  ddj = dd[k*(2*el+2) + i] / rj;
393  dl[(k-el+offset)*stride + i-el + offset] +=
394  sqrt_tbl[j-i] * sqrt_tbl[j-k] * ddj * coshb;
395  dl[(k-el+offset)*stride + i+1-el + offset] -=
396  sqrt_tbl[i+1] * sqrt_tbl[j-k] * ddj * sinhb;
397  dl[(k+1-el+offset)*stride + i-el + offset] +=
398  sqrt_tbl[j-i] * sqrt_tbl[k+1] * ddj * sinhb;
399  dl[(k+1-el+offset)*stride + i+1-el + offset] +=
400  sqrt_tbl[i+1] * sqrt_tbl[k+1] * ddj * coshb;
401  }
402  }
403 
404  // Fill top half of plane using symmetry.
405  for (m=-el; m<=el; m++)
406  for (mm=1; mm<=el; mm++)
407  dl[(m+offset)*stride + mm + offset] =
408  signs[abs(m)] * signs[abs(mm)]
409  * dl[(-m+offset)*stride - mm + offset];
410 
411  // Free temporary memory.
412  free(dd);
413  }
414 
415 }
416 
417 
449 void ssht_dl_beta_risbo_eighth_table(double *dl, double beta, int L,
450  ssht_dl_size_t dl_size,
451  int el, double *sqrt_tbl,
452  double *signs) {
453 
454  int offset, stride, imax;
455  double cosb, sinb, coshb, sinhb;
456  int i, j, k;
457  double rj, dlj, ddj;
458  double *dd;
459  int m, mm;
460 
461  // Get mm offset and stride for accessing dl data.
462  offset = ssht_dl_get_offset(L, dl_size);
463  stride = ssht_dl_get_stride(L, dl_size);
464 
465  // Compute Wigner plane.
466  if (el == 0) {
467 
468  dl[(0+offset)*stride + 0 + offset] = 1.0;
469 
470  }
471  else if (el == 1) {
472 
473  cosb = cos(beta);
474  sinb = sin(beta);
475  coshb = cos(beta / 2.0);
476  sinhb = sin(beta / 2.0);
477 
478  dl[(-1 + offset)*stride - 1 + offset] = coshb * coshb;
479  dl[(-1 + offset)*stride + 0 + offset] = sinb / SSHT_SQRT2;
480  dl[(-1 + offset)*stride + 1 + offset] = sinhb * sinhb;
481 
482  dl[(0 + offset)*stride - 1 + offset] = -sinb / SSHT_SQRT2;
483  dl[(0 + offset)*stride + 0 + offset] = cosb;
484  dl[(0 + offset)*stride + 1 + offset] = sinb / SSHT_SQRT2;
485 
486  dl[(1 + offset)*stride - 1 + offset] = sinhb * sinhb;
487  dl[(1 + offset)*stride + 0 + offset] = -sinb / SSHT_SQRT2;
488  dl[(1 + offset)*stride + 1 + offset] = coshb * coshb;
489 
490  }
491  else {
492 
493  coshb = -cos(beta / 2.0);
494  sinhb = sin(beta / 2.0);
495 
496  // Initialise the plane of the dl-matrix to 0.0 for the recursion
497  // from l - 1 to l - 1/2.
498  dd = (double*)calloc((el+3)*(el+3), sizeof(double));
500 
501  j = 2*el - 1;
502  rj = (double) j;
503  for (k=0; k<=el; k++) {
504  if (k==el)
505  imax = k+1;
506  else
507  imax = k+2;
508  for (i=0; i<=imax; i++) {
509  dlj = dl[(k-(el-1)+offset)*stride + i-(el-1) + offset] / rj;
510  dd[k*(el+3) + i] +=
511  sqrt_tbl[j-i] * sqrt_tbl[j-k] * dlj * coshb;
512  dd[k*(el+3) + i+1] -=
513  sqrt_tbl[i+1] * sqrt_tbl[j-k] * dlj * sinhb;
514  dd[(k+1)*(el+3) + i] +=
515  sqrt_tbl[j-i] * sqrt_tbl[k+1] * dlj * sinhb;
516  dd[(k+1)*(el+3) + i+1] +=
517  sqrt_tbl[i+1] * sqrt_tbl[k+1] * dlj * coshb;
518  }
519  }
520 
521  // Having constructed the d^(l+1/2) matrix in dd, do the second
522  // half-step recursion from dd to dl. Start by initilalising
523  // the plane of the dl-matrix to 0.0.
524  for (k=-el; k<=1; k++)
525  for (i=-el; i<=3; i++)
526  dl[(k+offset)*stride + i + offset] = 0.0;
527 
528  j = 2*el;
529  rj = (double) j;
530  for (k=0; k<=el; k++) {
531  for (i=0; i<=k+1; i++) {
532  ddj = dd[k*(el+3) + i] / rj;
533  dl[(k-el+offset)*stride + i-el + offset] +=
534  sqrt_tbl[j-i] * sqrt_tbl[j-k] * ddj * coshb;
535  dl[(k-el+offset)*stride + i+1-el + offset] -=
536  sqrt_tbl[i+1] * sqrt_tbl[j-k] * ddj * sinhb;
537  dl[(k+1-el+offset)*stride + i-el + offset] +=
538  sqrt_tbl[j-i] * sqrt_tbl[k+1] * ddj * sinhb;
539  dl[(k+1-el+offset)*stride + i+1-el + offset] +=
540  sqrt_tbl[i+1] * sqrt_tbl[k+1] * ddj * coshb;
541  }
542  }
543 
544  // Extend dl plane about boundaries, using symmetries, to the
545  // extents required for recursion.
546 
547  // Extend above diagonal by two in mm.
548  for (m=-el; m<=0; m++)
549  for (mm=m+1; mm<=m+2; mm++)
550  dl[(m+offset)*stride + mm + offset] =
551  signs[abs(m)] * signs[abs(mm)]
552  * dl[(mm+offset)*stride + m + offset];
553 
554  // Extend right by one in m.
555  for (m=1; m<=1; m++)
556  for (mm=-el; mm<=0; mm++)
557  dl[(m+offset)*stride + mm + offset] =
558  signs[abs(el)] * signs[abs(mm)]
559  * dl[(-m+offset)*stride + mm + offset];
560 
561  // Extend up by one in mm.
562  for (m=-el; m<=1; m++)
563  for (mm=1; mm<=1; mm++)
564  dl[(m+offset)*stride + mm + offset] =
565  signs[abs(el)] * signs[abs(m)]
566  * dl[(m+offset)*stride - mm + offset];
567 
568  // Free temporary memory.
569  free(dd);
570  }
571 
572 }
573 
574 
600  double *dl8,
601  int L,
602  ssht_dl_size_t dl_size,
603  ssht_dl_size_t dl8_size,
604  int el,
605  double *signs) {
606 
607  int offset, stride, offset8, stride8;
608  int m, mm;
609 
610  // Get mm offset and stride for accessing dl data.
611  offset = ssht_dl_get_offset(L, dl_size);
612  stride = ssht_dl_get_stride(L, dl_size);
613  offset8 = ssht_dl_get_offset(L, dl8_size);
614  stride8 = ssht_dl_get_stride(L, dl8_size);
615 
616  // Symmetry through origin to get eighth of the required quarter
617  // (first quadrant).
618  for (m=0; m<=el; m++)
619  for (mm=m; mm<=el; mm++)
620  dl[(m+offset)*stride + mm + offset] =
621  signs[m] * signs[mm]
622  * dl8[(-m+offset8)*stride8 - mm + offset8];
623 
624  // Diagonal symmetry to fill remaining quarter.
625  for (m=0; m<=el; m++)
626  for (mm=0; mm<=m-1; mm++)
627  dl[(m+offset)*stride + mm + offset] =
628  signs[m] * signs[mm]
629  * dl[(mm+offset)*stride + m + offset];
630 
631 }
632 
633 
662 void ssht_dl_beta_kostelec_full_table(double *dlm1p1, double *dl,
663  double beta, int L,
664  ssht_dl_size_t dl_size,
665  int el,
666  double *sqrt_tbl, double *signs) {
667 
668  int offset, stride;
669  double cosb, sinb, coshb, sinhb;
670  double lnAlm, lnfact2el;
671  int m, mm, elm1;
672  double elr, elm1r;
673 
674  // Get mm offset and stride for accessing dl data.
675  offset = ssht_dl_get_offset(L, dl_size);
676  stride = ssht_dl_get_stride(L, dl_size);
677 
678  // Compute Wigner plane.
679  if (el == 0) {
680 
681  dlm1p1[(0+offset)*stride + 0 + offset] = 1.0;
682 
683  }
684  else if (el == 1) {
685 
686  cosb = cos(beta);
687  sinb = sin(beta);
688  coshb = cos(beta / 2.0);
689  sinhb = sin(beta / 2.0);
690 
691  // These terms follow directly from the boundary conditions and
692  // one recursion to get the central term.
693 
694  dlm1p1[(-1 + offset)*stride - 1 + offset] = coshb * coshb;
695  dlm1p1[(-1 + offset)*stride + 0 + offset] = sinb / SSHT_SQRT2;
696  dlm1p1[(-1 + offset)*stride + 1 + offset] = sinhb * sinhb;
697 
698  dlm1p1[(0 + offset)*stride - 1 + offset] = -sinb / SSHT_SQRT2;
699  dlm1p1[(0 + offset)*stride + 0 + offset] = cosb;
700  dlm1p1[(0 + offset)*stride + 1 + offset] = sinb / SSHT_SQRT2;
701 
702  dlm1p1[(1 + offset)*stride - 1 + offset] = sinhb * sinhb;
703  dlm1p1[(1 + offset)*stride + 0 + offset] = -sinb / SSHT_SQRT2;
704  dlm1p1[(1 + offset)*stride + 1 + offset] = coshb * coshb;
705 
706  }
707  else {
708 
709  cosb = cos(beta);
710  coshb = cos(beta / 2.0);
711  sinhb = sin(beta / 2.0);
712 
713  elr = (double) el;
714  elm1 = el - 1;
715  elm1r = (double) elm1;
716 
717  // Recurse over all of plane except boundaries.
718  for (m=-(el-1); m<=el-1; m++) {
719  for (mm=-(el-1); mm<=el-1; mm++) {
720 
721  // Compute 3-term recursion.
722  dlm1p1[(m + offset)*stride + mm + offset] =
723  (cosb - m*mm/(elm1r*elr)) * dl[(m + offset)*stride + mm + offset]
724  -
725  sqrt_tbl[elm1+m] * sqrt_tbl[elm1-m] * sqrt_tbl[elm1+mm] * sqrt_tbl[elm1-mm]
726  / (elm1r * (2.0*elm1r + 1.0))
727  * dlm1p1[(m + offset)*stride + mm + offset];
728 
729  // Perform scaling.
730  dlm1p1[(m + offset)*stride + mm + offset] *=
731  el * (2*elm1 + 1.0)
732  / (sqrt_tbl[el-m] * sqrt_tbl[el+m] * sqrt_tbl[el-mm] * sqrt_tbl[el+mm]);
733 
734  }
735  }
736 
737  // Compute boundaries.
738  lnfact2el = logfact(2*el);
739  for (m=-el; m<=el; m++) {
740 
741  lnAlm = (lnfact2el - logfact(el+m) - logfact(el-m)) / 2.0;
742 
743  // Right line.
744  dlm1p1[(el + offset)*stride + m + offset] =
745  signs[el] * signs[abs(m)]
746  * exp(lnAlm + (el+m)*log(coshb) + (el-m)*log(sinhb));
747 
748  // Left line.
749  dlm1p1[(-el + offset)*stride + m + offset] =
750  exp(lnAlm + (el-m)*log(coshb) + (el+m)*log(sinhb));
751 
752  // Top line.
753  dlm1p1[(m + offset)*stride + el + offset] =
754  exp(lnAlm + (el+m)*log(coshb) + (el-m)*log(sinhb));
755 
756  // Bottom line.
757  dlm1p1[(m + offset)*stride - el + offset] =
758  signs[el] * signs[abs(m)]
759  * exp(lnAlm + (el-m)*log(coshb) + (el+m)*log(sinhb));
760  }
761 
762  }
763 
764 }
765 
766 
795 void ssht_dl_beta_kostelec_line_table(double *dlm1p1_line, double *dl_line,
796  double beta, int L, int mm, int el,
797  double *sqrt_tbl, double *signs) {
798 
799  int offset;
800  double cosb, sinb, coshb, sinhb;
801  double lnAlm, lnAlmm, lnfact2el;
802  int m, elm1;
803  double elr, elm1r;
804 
805  // Compute m offset for accessing dl line.
806  offset = L-1;
807 
808  // Compute Wigner plane.
809  if (el < abs(mm)) {
810  // Do nothing (dl line should remain zero).
811  return;
812  }
813  else if (el == 1) {
814 
815  cosb = cos(beta);
816  sinb = sin(beta);
817  coshb = cos(beta / 2.0);
818  sinhb = sin(beta / 2.0);
819 
820  if (mm == -1) {
821  dlm1p1_line[-1 + offset] = coshb * coshb;
822  dlm1p1_line[ 0 + offset] = -sinb / SSHT_SQRT2;
823  dlm1p1_line[ 1 + offset] = sinhb * sinhb;
824  }
825  else if (mm == 0) {
826  dlm1p1_line[-1 + offset] = sinb / SSHT_SQRT2;
827  dlm1p1_line[ 0 + offset] = cosb;
828  dlm1p1_line[ 1 + offset] = -sinb / SSHT_SQRT2;
829  }
830  else {
831  // mm == +1
832  dlm1p1_line[-1 + offset] = sinhb * sinhb;
833  dlm1p1_line[ 0 + offset] = sinb / SSHT_SQRT2;
834  dlm1p1_line[ 1 + offset] = coshb * coshb;
835  }
836 
837  }
838  else if (el == abs(mm)) {
839 
840  coshb = cos(beta / 2.0);
841  sinhb = sin(beta / 2.0);
842 
843  // Initalise line using equation (4.8) or (4.9) from K&R (2010).
844  if (mm >= 0) {
845 
846  // Initialise using equation (4.8), i.e. top line.
847  lnfact2el = logfact(2*el);
848  for (m=-el; m<=el; m++) {
849  lnAlm = (lnfact2el - logfact(el+m) - logfact(el-m)) / 2.0;
850  dlm1p1_line[m + offset] =
851  exp(lnAlm + (el+m)*log(coshb) + (el-m)*log(sinhb));
852  }
853 
854  }
855  else {
856 
857  // Initialise using equation (4.9), i.e. bottom line.
858  lnfact2el = logfact(2*el);
859  for (m=-el; m<=el; m++) {
860  lnAlm = (lnfact2el - logfact(el+m) - logfact(el-m)) / 2.0;
861  dlm1p1_line[m + offset] =
862  signs[el] * signs[abs(m)]
863  * exp(lnAlm + (el-m)*log(coshb) + (el+m)*log(sinhb));
864  }
865 
866  }
867 
868  }
869  else {
870 
871  cosb = cos(beta);
872  coshb = cos(beta / 2.0);
873  sinhb = sin(beta / 2.0);
874 
875  elr = (double) el;
876  elm1 = el - 1;
877  elm1r = (double) elm1;
878 
879  // Recuse over line.
880  for (m=-(el-1); m<=el-1; m++) {
881 
882  // Compute 3-term recursion.
883  dlm1p1_line[m + offset] =
884  (cosb - m*mm/(elm1r*elr)) * dl_line[m + offset]
885  -
886  sqrt_tbl[elm1+m] * sqrt_tbl[elm1-m] * sqrt_tbl[elm1+mm] * sqrt_tbl[elm1-mm]
887  / (elm1r * (2.0*elm1r + 1.0))
888  * dlm1p1_line[m + offset];
889 
890  // Perform scaling.
891  dlm1p1_line[m + offset] *=
892  el * (2*elm1 + 1.0)
893  / (sqrt_tbl[el-m] * sqrt_tbl[el+m] * sqrt_tbl[el-mm] * sqrt_tbl[el+mm]);
894 
895  }
896 
897  // Compute edges...
898  lnfact2el = logfact(2*el);
899  lnAlmm = (lnfact2el - logfact(el+mm) - logfact(el-mm)) / 2.0;
900 
901  // Left edge.
902  dlm1p1_line[-el + offset] =
903  exp(lnAlmm + (el-mm)*log(coshb) + (el+mm)*log(sinhb));
904 
905  // Right edge.
906  dlm1p1_line[el + offset] =
907  signs[el] * signs[abs(mm)]
908  * exp(lnAlmm + (el+mm)*log(coshb) + (el-mm)*log(sinhb));
909 
910  }
911 
912 }
913 
914 
943 void ssht_dl_beta_kostelec_halfline_table(double *dlm1p1_line, double *dl_line,
944  double beta, int L, int mm, int el,
945  double *sqrt_tbl, double *signs) {
946 
947  int offset;
948  double cosb, sinb, coshb, sinhb;
949  double lnAlm, lnAlmm, lnfact2el;
950  int m, elm1;
951  double elr, elm1r;
952 
953  // Compute m offset for accessing dl line.
954  offset = 0;
955 
956  // Compute Wigner plane.
957  if (el < abs(mm)) {
958  // Do nothing (dl line should remain zero).
959  return;
960  }
961  else if (el == 1) {
962 
963  cosb = cos(beta);
964  sinb = sin(beta);
965  coshb = cos(beta / 2.0);
966  sinhb = sin(beta / 2.0);
967 
968  if (mm == -1) {
969  dlm1p1_line[ 0 + offset] = -sinb / SSHT_SQRT2;
970  dlm1p1_line[ 1 + offset] = sinhb * sinhb;
971  }
972  else if (mm == 0) {
973  dlm1p1_line[ 0 + offset] = cosb;
974  dlm1p1_line[ 1 + offset] = -sinb / SSHT_SQRT2;
975  }
976  else {
977  // mm == +1
978  dlm1p1_line[ 0 + offset] = sinb / SSHT_SQRT2;
979  dlm1p1_line[ 1 + offset] = coshb * coshb;
980  }
981 
982  }
983  else if (el == abs(mm)) {
984 
985  coshb = cos(beta / 2.0);
986  sinhb = sin(beta / 2.0);
987 
988  // Initalise line using equation (4.8) or (4.9) from K&R (2010).
989  if (mm >= 0) {
990 
991  // Initialise using equation (4.8), i.e. top line.
992  lnfact2el = logfact(2*el);
993  for (m=0; m<=el; m++) {
994  lnAlm = (lnfact2el - logfact(el+m) - logfact(el-m)) / 2.0;
995  dlm1p1_line[m + offset] =
996  exp(lnAlm + (el+m)*log(coshb) + (el-m)*log(sinhb));
997  }
998 
999  }
1000  else {
1001 
1002  // Initialise using equation (4.9), i.e. bottom line.
1003  lnfact2el = logfact(2*el);
1004  for (m=0; m<=el; m++) {
1005  lnAlm = (lnfact2el - logfact(el+m) - logfact(el-m)) / 2.0;
1006  dlm1p1_line[m + offset] =
1007  signs[el] * signs[abs(m)]
1008  * exp(lnAlm + (el-m)*log(coshb) + (el+m)*log(sinhb));
1009  }
1010 
1011  }
1012 
1013  }
1014  else {
1015 
1016  cosb = cos(beta);
1017  coshb = cos(beta / 2.0);
1018  sinhb = sin(beta / 2.0);
1019 
1020  elr = (double) el;
1021  elm1 = el - 1;
1022  elm1r = (double) elm1;
1023 
1024  // Recuse over line.
1025  for (m=0; m<=el-1; m++) {
1026 
1027  // Compute 3-term recursion.
1028  dlm1p1_line[m + offset] =
1029  (cosb - m*mm/(elm1r*elr)) * dl_line[m + offset]
1030  -
1031  sqrt_tbl[elm1+m] * sqrt_tbl[elm1-m] * sqrt_tbl[elm1+mm] * sqrt_tbl[elm1-mm]
1032  / (elm1r * (2.0*elm1r + 1.0))
1033  * dlm1p1_line[m + offset];
1034 
1035  // Perform scaling.
1036  dlm1p1_line[m + offset] *=
1037  el * (2*elm1 + 1.0)
1038  / (sqrt_tbl[el-m] * sqrt_tbl[el+m] * sqrt_tbl[el-mm] * sqrt_tbl[el+mm]);
1039 
1040  }
1041 
1042  // Compute edges...
1043  lnfact2el = logfact(2*el);
1044  lnAlmm = (lnfact2el - logfact(el+mm) - logfact(el-mm)) / 2.0;
1045 
1046  // Right edge.
1047  dlm1p1_line[el + offset] =
1048  signs[el] * signs[abs(mm)]
1049  * exp(lnAlmm + (el+mm)*log(coshb) + (el-mm)*log(sinhb));
1050 
1051  }
1052 
1053 }
1054 
1055 
1077 void ssht_dl_halfpi_trapani_eighth_table(double *dl, int L,
1078  ssht_dl_size_t dl_size,
1079  int el, double *sqrt_tbl) {
1080 
1081  int m, mm, offset, stride;
1082  double *dmm;
1083  double t1, t2, s1, s2;
1084 
1085  // Allocate temporary memory.
1086  dmm = (double*)calloc(el+1, sizeof(double));
1088 
1089  // Get mm offset and stride for accessing dl data.
1090  offset = ssht_dl_get_offset(L, dl_size);
1091  stride = ssht_dl_get_stride(L, dl_size);
1092 
1093  // Compute Wigner plane.
1094  if (el == 0) {
1095 
1096  dl[0*stride + 0 + offset] = 1.0;
1097 
1098  }
1099  else {
1100 
1101  // Eqn (9) of T&N (2006).
1102  dmm[0] = - sqrt_tbl[2*el-1] / sqrt_tbl[2*el]
1103  * dl[(el-1)*stride + 0 + offset];
1104 
1105  // Eqn (10) of T&N (2006).
1106  for (mm=1; mm<=el; mm++) {
1107  dmm[mm] = sqrt_tbl[el] / SSHT_SQRT2
1108  * sqrt_tbl[2*el-1] / sqrt_tbl[el+mm] / sqrt_tbl[el+mm-1]
1109  * dl[(el-1)*stride + (mm-1) + offset];
1110  }
1111 
1112  // Initialise dl for next el.
1113  for (mm=0; mm<=el; mm++) {
1114  dl[el*stride + mm + offset] = dmm[mm];
1115  }
1116 
1117 /* LOGICAL BUT *NOT* MOST EFFICIENT ALGORITHM
1118  // Eqn (11) of T&N (2006).
1119  for (mm=0; mm<=el; mm++) {
1120 
1121  // m = el-1 case (t2 = 0).
1122  m = el-1;
1123  dl[m*stride + mm + offset] = 2e0 * mm / sqrt_tbl[el-m] / sqrt_tbl[el+m+1]
1124  * dl[(m+1)*stride + mm + offset];
1125 
1126  // Remaining m cases.
1127  for (m=el-2; m>=mm; m--) {
1128  t1 = 2e0 * mm / sqrt_tbl[el-m] / sqrt_tbl[el+m+1]
1129  * dl[(m+1)*stride + mm + offset];
1130  t2 = sqrt_tbl[el-m-1] * sqrt_tbl[el+m+2] / sqrt_tbl[el-m] / sqrt_tbl[el+m+1]
1131  * dl[(m+2)*stride + mm + offset];
1132  dl[m*stride + mm + offset] = t1 - t2;
1133  }
1134 
1135  }
1136 */
1137 
1138  // Eqn (11) of T&N (2006).
1139  // OPTIMISED FOR MEMORY ACCESS.
1140  m = el-1;
1141  s1 = sqrt_tbl[el-m] * sqrt_tbl[el+m+1];
1142  for (mm=0; mm<=el; mm++) {
1143  // m = el-1 case (t2 = 0).
1144  dl[m*stride + mm + offset] = 2e0 * mm / s1
1145  * dl[(m+1)*stride + mm + offset];
1146  }
1147 
1148  for (m=el-2; m>=0; m--) {
1149  s1 = sqrt_tbl[el-m] * sqrt_tbl[el+m+1];
1150  s2 = sqrt_tbl[el-m-1] * sqrt_tbl[el+m+2] / sqrt_tbl[el-m] / sqrt_tbl[el+m+1];
1151  for (mm=0; mm<=m; mm++) {
1152  t1 = 2e0 * mm / s1
1153  * dl[(m+1)*stride + mm + offset];
1154  t2 = s2
1155  * dl[(m+2)*stride + mm + offset];
1156  dl[m*stride + mm + offset] = t1 - t2;
1157  }
1158 
1159  }
1160 
1161  }
1162 
1163  // Free memory.
1164  free(dmm);
1165 
1166 }
1167 
1168 
1196  ssht_dl_size_t dl_size,
1197  int el, double *sqrt_tbl) {
1198 
1199  int m, mm, offset, stride;
1200  double *dmm;
1201  double t1, t2, s1, s2;
1202 
1203  // Allocate temporary memory.
1204  dmm = (double*)calloc(el+1, sizeof(double));
1206 
1207  // Get mm offset and stride for accessing dl data.
1208  offset = ssht_dl_get_offset(L, dl_size);
1209  stride = ssht_dl_get_stride(L, dl_size);
1210 
1211  // Compute Wigner plane.
1212  if (el == 0) {
1213 
1214  dl[0*stride + 0 + offset] = 1.0;
1215 
1216  }
1217  else {
1218 
1219  // Eqn (9) of T&N (2006).
1220  dmm[0] = - sqrt_tbl[2*el-1] / sqrt_tbl[2*el]
1221  * dl[(el-1)*stride + 0 + offset];
1222 
1223  // Eqn (10) of T&N (2006).
1224  for (mm=1; mm<=el; mm++) {
1225  dmm[mm] = sqrt_tbl[el] / SSHT_SQRT2
1226  * sqrt_tbl[2*el-1] / sqrt_tbl[el+mm] / sqrt_tbl[el+mm-1]
1227  * dl[(el-1)*stride + (mm-1) + offset];
1228  }
1229 
1230  // Initialise dl for next el.
1231  for (mm=0; mm<=el; mm++) {
1232  dl[el*stride + mm + offset] = dmm[mm];
1233  }
1234 
1235 /* LOGICAL BUT *NOT* MOST EFFICIENT ALGORITHM
1236  // Eqn (11) of T&N (2006).
1237  for (mm=0; mm<=el; mm++) {
1238 
1239  // m = el-1 case (t2 = 0).
1240  m = el-1;
1241  dl[m*stride + mm + offset] = 2e0 * mm / sqrt_tbl[el-m] / sqrt_tbl[el+m+1]
1242  * dl[(m+1)*stride + mm + offset];
1243 
1244  // Remaining m cases.
1245  for (m=el-2; m>=mm; m--) {
1246  t1 = 2e0 * mm / sqrt_tbl[el-m] / sqrt_tbl[el+m+1]
1247  * dl[(m+1)*stride + mm + offset];
1248  t2 = sqrt_tbl[el-m-1] * sqrt_tbl[el+m+2] / sqrt_tbl[el-m] / sqrt_tbl[el+m+1]
1249  * dl[(m+2)*stride + mm + offset];
1250  dl[m*stride + mm + offset] = t1 - t2;
1251  }
1252 
1253  }
1254 */
1255 
1256  // Eqn (11) of T&N (2006).
1257  // OPTIMISED FOR MEMORY ACCESS.
1258  m = el-1;
1259  s1 = sqrt_tbl[el-m] * sqrt_tbl[el+m+1];
1260  for (mm=0; mm<=el; mm++) {
1261  // m = el-1 case (t2 = 0).
1262  dl[m*stride + mm + offset] = 2e0 * mm / s1
1263  * dl[(m+1)*stride + mm + offset];
1264  }
1265 
1266  for (m=el-2; m>=0; m--) {
1267  s1 = sqrt_tbl[el-m] * sqrt_tbl[el+m+1];
1268  s2 = sqrt_tbl[el-m-1] * sqrt_tbl[el+m+2] / sqrt_tbl[el-m] / sqrt_tbl[el+m+1];
1269  for (mm=0; mm<=el; mm++) {
1270  t1 = 2e0 * mm / s1
1271  * dl[(m+1)*stride + mm + offset];
1272  t2 = s2
1273  * dl[(m+2)*stride + mm + offset];
1274  dl[m*stride + mm + offset] = t1 - t2;
1275  }
1276 
1277  }
1278 
1279  }
1280 
1281  // Free memory.
1282  free(dmm);
1283 
1284 }
1285 
1286 
1308  ssht_dl_size_t dl_size,
1309  int el, double *signs) {
1310 
1311  int m, mm, offset, stride;
1312 
1313  // Get mm offset and stride for accessing dl data.
1314  offset = ssht_dl_get_offset(L, dl_size);
1315  stride = ssht_dl_get_stride(L, dl_size);
1316 
1317  // Diagonal symmetry to fill in quarter.
1318  for (m=0; m<=el; m++)
1319  for (mm=m+1; mm<=el; mm++)
1320  dl[m*stride + mm + offset] =
1321  signs[m] * signs[mm] * dl[mm*stride + m + offset];
1322 
1323  // Symmetry in mm to fill in half.
1324  for (m=0; m<=el; m++)
1325  for (mm=-el; mm<=-1; mm++)
1326  dl[m*stride + mm + offset] =
1327  signs[el] * signs[m] * dl[m*stride - mm + offset];
1328 
1329 }
1330 
1331 
1353  ssht_dl_size_t dl_size,
1354  int el, double *signs) {
1355 
1356  int m, mm, offset, stride;
1357 
1358  // Get mm offset and stride for accessing dl data.
1359  offset = ssht_dl_get_offset(L, dl_size);
1360  stride = ssht_dl_get_stride(L, dl_size);
1361 
1362  // Diagonal symmetry to fill in quarter.
1363  for (m=0; m<=el; m++)
1364  for (mm=m+1; mm<=el; mm++)
1365  dl[m*stride + mm + offset] =
1366  signs[m] * signs[mm] * dl[mm*stride + m + offset];
1367 
1368 }
1369 
1370 
1379 double logfact(int n) {
1380 
1381  double y, temp, sum, c[6], loggamma, x;
1382  int nn;
1383 
1384  if (n < 0) {
1385 
1386  SSHT_ERROR_GENERIC("Factorial argument negative")
1387 
1388  }
1389  else {
1390 
1391  // The engine of this function actually calculates the gamma function,
1392  // for which the real argument is x = n + 1.
1393 
1394  x = (double) (n) + 1.0;
1395 
1396  // Table of fitting constants.
1397 
1398  c[0] = 76.18009172947146;
1399  c[1] = - 86.50532032941677;
1400  c[2] = 24.01409824083091;
1401  c[3] = - 1.231739572450155;
1402  c[4] = 0.1208650973866179e-2;
1403  c[5] = - 0.5395239384953e-5;
1404 
1405  // Add up fit.
1406 
1407  temp = x + 5.5 - (x + 0.5) * log(x + 5.5);
1408  sum = 1.000000000190015;
1409  y = x;
1410 
1411  for (nn=0; nn<=5; nn++) {
1412  y = y + 1.0;
1413  sum = sum + c[nn] / y;
1414  }
1415 
1416  loggamma = - temp + log(2.5066282746310005 * sum / x);
1417 
1418  }
1419 
1420  // Finally make explicit the conversion back to log of the factorial.
1421  return loggamma;
1422 
1423 }
ssht_dl_beta_risbo_full_table
void ssht_dl_beta_risbo_full_table(double *dl, double beta, int L, ssht_dl_size_t dl_size, int el, double *sqrt_tbl)
Definition: ssht_dl.c:190
ssht_dl_beta_risbo_half_table
void ssht_dl_beta_risbo_half_table(double *dl, double beta, int L, ssht_dl_size_t dl_size, int el, double *sqrt_tbl, double *signs)
Definition: ssht_dl.c:313
ssht_error.h
SSHT_DL_QUARTER_EXTENDED
@ SSHT_DL_QUARTER_EXTENDED
Definition: ssht_dl.h:16
SSHT_ERROR_GENERIC
#define SSHT_ERROR_GENERIC(comment)
Definition: ssht_error.h:19
ssht_dl.h
ssht_dl_get_stride
int ssht_dl_get_stride(int L, ssht_dl_size_t dl_size)
Definition: ssht_dl.c:146
SSHT_DL_FULL
@ SSHT_DL_FULL
Definition: ssht_dl.h:18
ssht_dl_calloc
double * ssht_dl_calloc(int L, ssht_dl_size_t dl_size)
Definition: ssht_dl.c:64
ssht_dl_beta_kostelec_halfline_table
void ssht_dl_beta_kostelec_halfline_table(double *dlm1p1_line, double *dl_line, double beta, int L, int mm, int el, double *sqrt_tbl, double *signs)
Definition: ssht_dl.c:943
ssht_dl_halfpi_trapani_eighth_table
void ssht_dl_halfpi_trapani_eighth_table(double *dl, int L, ssht_dl_size_t dl_size, int el, double *sqrt_tbl)
Definition: ssht_dl.c:1077
SSHT_SQRT2
#define SSHT_SQRT2
Definition: ssht_types.h:41
ssht_dl_size_t
ssht_dl_size_t
Definition: ssht_dl.h:15
ssht_dl_beta_risbo_eighth_table
void ssht_dl_beta_risbo_eighth_table(double *dl, double beta, int L, ssht_dl_size_t dl_size, int el, double *sqrt_tbl, double *signs)
Definition: ssht_dl.c:449
logfact
double logfact(int n)
Definition: ssht_dl.c:1379
ssht_dl_halfpi_trapani_fill_eighth2quarter_table
void ssht_dl_halfpi_trapani_fill_eighth2quarter_table(double *dl, int L, ssht_dl_size_t dl_size, int el, double *signs)
Definition: ssht_dl.c:1352
ssht_dl_halfpi_trapani_fill_eighth2righthalf_table
void ssht_dl_halfpi_trapani_fill_eighth2righthalf_table(double *dl, int L, ssht_dl_size_t dl_size, int el, double *signs)
Definition: ssht_dl.c:1307
ssht_dl_beta_kostelec_full_table
void ssht_dl_beta_kostelec_full_table(double *dlm1p1, double *dl, double beta, int L, ssht_dl_size_t dl_size, int el, double *sqrt_tbl, double *signs)
Definition: ssht_dl.c:662
SSHT_ERROR_MEM_ALLOC_CHECK
#define SSHT_ERROR_MEM_ALLOC_CHECK(pointer)
Definition: ssht_error.h:28
ssht_dl_beta_risbo_fill_eighth2quarter_table
void ssht_dl_beta_risbo_fill_eighth2quarter_table(double *dl, double *dl8, int L, ssht_dl_size_t dl_size, ssht_dl_size_t dl8_size, int el, double *signs)
Definition: ssht_dl.c:599
SSHT_DL_HALF
@ SSHT_DL_HALF
Definition: ssht_dl.h:17
ssht_dl_get_offset
int ssht_dl_get_offset(int L, ssht_dl_size_t dl_size)
Definition: ssht_dl.c:111
ssht_dl_beta_kostelec_line_table
void ssht_dl_beta_kostelec_line_table(double *dlm1p1_line, double *dl_line, double beta, int L, int mm, int el, double *sqrt_tbl, double *signs)
Definition: ssht_dl.c:795
ssht_dl_halfpi_trapani_quarter_table
void ssht_dl_halfpi_trapani_quarter_table(double *dl, int L, ssht_dl_size_t dl_size, int el, double *sqrt_tbl)
Definition: ssht_dl.c:1195
ssht_types.h
SSHT_DL_QUARTER
@ SSHT_DL_QUARTER
Definition: ssht_dl.h:15