WIP hexarray modes example
[qpms.git] / qpms / bessel.c
blob000ee49781c665948c275880ebf3f577be501b15
1 #include <assert.h>
2 #include "qpms_specfunc.h"
3 #include <stdlib.h>
4 #include <stddef.h>
5 #include <string.h>
6 #include "kahansum.h"
7 #include <gsl/gsl_sf_bessel.h>
8 #include <complex.h>
9 #include "qpms_error.h"
10 #include <camos.h>
11 #include <math.h>
13 #ifndef M_LN2
14 #define M_LN2 0.69314718055994530942 /* log_e 2 */
15 #endif
17 static inline complex double ipow(int x) {
18 return cpow(I,x);
21 #if 0
22 // Inspired by scipy/special/_spherical_bessel.pxd
23 static inline complex double spherical_jn(qpms_l_t l, complex double z) {
24 if (isnan(creal(z)) || isnan(cimag(z))) return NAN+I*NAN;
25 if (l < 0) QPMS_WTF;
26 if (fpclassify(creal(z)) == FP_INFINITE)
27 if(0 == cimag(z))
28 return 0;
29 else
30 return INFINITY + I * INFINITY;
31 if (z == 0)
32 if (l == 0) return 1;
33 else return 0;
34 return csqrt(M_PI_2/z) * zbesj(l + .5, z);
37 static inline complex double spherical_yn(qpms_l_t l, complex double z) {
38 if (isnan(creal(z)) || isnan(cimag(z))) return NAN+I*NAN;
39 if (l < 0) QPMS_WTF;
40 if (fpclassify(creal(z)) == FP_INFINITE)
41 if(0 == cimag(z))
42 return 0;
43 else
44 return INFINITY + I * INFINITY;
45 if (z == 0)
46 return NAN;
47 return csqrt(M_PI_2/z) * zbesy(l + .5, z);
49 #endif
51 // There is a big issue with gsl's precision of spherical bessel function; these have to be implemented differently
52 qpms_errno_t qpms_sph_bessel_realx_fill(qpms_bessel_t typ, qpms_l_t lmax, double x, complex double *result_array) {
53 int retval;
54 double tmparr[lmax+1];
55 switch(typ) {
56 case QPMS_BESSEL_REGULAR:
57 retval = gsl_sf_bessel_jl_steed_array(lmax, x, tmparr);
58 for (int l = 0; l <= lmax; ++l) result_array[l] = tmparr[l];
59 return retval;
60 break;
61 case QPMS_BESSEL_SINGULAR: //FIXME: is this precise enough? Would it be better to do it one-by-one?
62 retval = gsl_sf_bessel_yl_array(lmax,x,tmparr);
63 for (int l = 0; l <= lmax; ++l) result_array[l] = tmparr[l];
64 return retval;
65 break;
66 case QPMS_HANKEL_PLUS:
67 case QPMS_HANKEL_MINUS:
68 retval = gsl_sf_bessel_jl_steed_array(lmax, x, tmparr);
69 for (int l = 0; l <= lmax; ++l) result_array[l] = tmparr[l];
70 if(retval) return retval;
71 retval = gsl_sf_bessel_yl_array(lmax, x, tmparr);
72 if (typ==QPMS_HANKEL_PLUS)
73 for (int l = 0; l <= lmax; ++l) result_array[l] += I * tmparr[l];
74 else
75 for (int l = 0; l <= lmax; ++l) result_array[l] +=-I * tmparr[l];
76 return retval;
77 break;
78 default:
79 abort();
80 //return GSL_EDOM;
82 assert(0);
85 // TODO DOC
86 qpms_errno_t qpms_sph_bessel_fill(qpms_bessel_t typ, qpms_l_t lmax, complex double x, complex double *res) {
87 if(!cimag(x))
88 return qpms_sph_bessel_realx_fill(typ, lmax, creal(x), res);
89 else if (isnan(creal(x)) || isnan(cimag(x)))
90 for(qpms_l_t l = 0; l <= lmax; ++l) res[l] = NAN + I*NAN;
91 else if (lmax < 0) QPMS_WTF;
92 else if (fpclassify(creal(x)) == FP_INFINITE)
93 for(qpms_l_t l = 0; l <= lmax; ++l) res[l] = INFINITY + I * INFINITY;
94 else {
95 try_again: ;
96 int retry_counter = 0;
97 const double zr = creal(x), zi = cimag(x), fnu = 0.5;
98 const int n = lmax + 1, kode = 1 /* No exponential scaling */;
99 double cyr[n], cyi[n];
100 int ierr, nz;
101 unsigned int kindchar; // Only for error output
102 const complex double prefac = csqrt(M_PI_2/x);
103 switch(typ) {
104 case QPMS_BESSEL_REGULAR:
105 kindchar = 'j';
106 ierr = camos_zbesj(zr, zi, fnu, kode, n, cyr, cyi, &nz);
107 break;
108 case QPMS_BESSEL_SINGULAR:
109 kindchar = 'y';
111 double cwrkr[lmax + 1], cwrki[lmax + 1];
112 ierr = camos_zbesy(zr, zi, fnu, kode, n, cyr, cyi, &nz, cwrkr, cwrki);
114 break;
115 case QPMS_HANKEL_PLUS:
116 case QPMS_HANKEL_MINUS:
117 kindchar = 'h';
119 const int m = (typ == QPMS_HANKEL_PLUS) ? 1 : 2;
120 ierr = camos_zbesh(zr, zi, fnu, kode, m, n, cyr, cyi, &nz);
122 break;
123 default:
124 QPMS_WTF;
126 // TODO check for underflows? (nz != 0)
127 if (ierr == 0 || ierr == 3) {
128 for (qpms_l_t l = 0; l <= lmax; ++l)
129 res[l] = prefac * (cyr[l] + I * cyi[l]);
130 if (ierr == 3)
131 QPMS_WARN("Amos's zbes%c computation done but losses of significance "
132 "by argument reduction produce less than half of machine accuracy.",
133 kindchar);
134 return QPMS_SUCCESS; //TODO maybe something else if ierr == 3
136 else {
137 if (retry_counter < 5) {
138 QPMS_WARN("Amos's zbes%c failed with ierr == %d (lMax = %d, x = %+.16g%+.16gi). Retrying.\n",
139 kindchar, (int) ierr, lmax, creal(x), cimag(x));
140 ++retry_counter;
141 goto try_again;
142 } else QPMS_PR_ERROR("Amos's zbes%c failed with ierr == %d (lMax = %d, x = %+.16g%+.16gi).",
143 kindchar, (int) ierr, lmax, creal(x), cimag(x));
146 return QPMS_SUCCESS;
149 static inline ptrdiff_t akn_index(qpms_l_t n, qpms_l_t k) {
150 assert(k <= n);
151 return ((ptrdiff_t) n + 1) * n / 2 + k;
153 static inline ptrdiff_t bkn_index(qpms_l_t n, qpms_l_t k) {
154 assert(k <= n+1);
155 return ((ptrdiff_t) n + 2) * (n + 1) / 2 - 1 + k;
158 static inline qpms_errno_t qpms_sbessel_calculator_ensure_lMax(qpms_sbessel_calculator_t *c, qpms_l_t lMax) {
159 if (lMax <= c->lMax)
160 return QPMS_SUCCESS;
161 else {
162 if ( NULL == (c->akn = realloc(c->akn, sizeof(double) * akn_index(lMax + 2, 0))))
163 abort();
164 //if ( NULL == (c->bkn = realloc(c->bkn, sizeof(complex double) * bkn_index(lMax + 1, 0))))
165 // abort();
166 for(qpms_l_t n = c->lMax+1; n <= lMax + 1; ++n)
167 for(qpms_l_t k = 0; k <= n; ++k)
168 c->akn[akn_index(n,k)] = exp(lgamma(n + k + 1) - k*M_LN2 - lgamma(k + 1) - lgamma(n - k + 1));
169 // ... TODO derivace
170 c->lMax = lMax;
171 return QPMS_SUCCESS;
175 complex double qpms_sbessel_calc_h1(qpms_sbessel_calculator_t *c, qpms_l_t n, complex double x) {
176 if(QPMS_SUCCESS != qpms_sbessel_calculator_ensure_lMax(c, n))
177 abort();
178 complex double z = I/x;
179 complex double result = 0;
180 for (qpms_l_t k = n; k >= 0; --k)
181 // can we use fma for complex?
182 //result = fma(result, z, c->akn(n, k));
183 result = result * z + c->akn[akn_index(n,k)];
184 result *= z * ipow(-n-2) * cexp(I * x);
185 return result;
188 qpms_errno_t qpms_sbessel_calc_h1_fill(qpms_sbessel_calculator_t * const c,
189 const qpms_l_t lMax, const complex double x, complex double * const target) {
190 if(QPMS_SUCCESS != qpms_sbessel_calculator_ensure_lMax(c, lMax))
191 abort();
192 memset(target, 0, sizeof(complex double) * lMax);
193 complex double kahancomp[lMax];
194 memset(kahancomp, 0, sizeof(complex double) * lMax);
195 for(qpms_l_t k = 0; k <= lMax; ++k){
196 double xp = cpow(x, -k-1);
197 for(qpms_l_t l = k; l <= lMax; ++l)
198 ckahanadd(target + l, kahancomp + l, c->akn[akn_index(l,k)] * xp * ipow(k-l-1));
200 complex double eix = cexp(I * x);
201 for (qpms_l_t l = 0; l <= lMax; ++l)
202 target[l] *= eix;
203 return QPMS_SUCCESS;
206 qpms_sbessel_calculator_t *qpms_sbessel_calculator_init() {
207 qpms_sbessel_calculator_t *c = malloc(sizeof(qpms_sbessel_calculator_t));
208 c->akn = NULL;
209 //c->bkn = NULL;
210 c->lMax = -1;
211 return c;
214 void qpms_sbessel_calculator_pfree(qpms_sbessel_calculator_t *c) {
215 if(c->akn) free(c->akn);
216 //if(c->bkn) free(c->bkn);
217 free(c);