gsl_specfunc__bessel_i.c 8.8 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330
  1. /* specfunc/bessel_i.c
  2. *
  3. * Copyright (C) 1996, 1997, 1998, 1999, 2000 Gerard Jungman
  4. *
  5. * This program is free software; you can redistribute it and/or modify
  6. * it under the terms of the GNU General Public License as published by
  7. * the Free Software Foundation; either version 3 of the License, or (at
  8. * your option) any later version.
  9. *
  10. * This program is distributed in the hope that it will be useful, but
  11. * WITHOUT ANY WARRANTY; without even the implied warranty of
  12. * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
  13. * General Public License for more details.
  14. *
  15. * You should have received a copy of the GNU General Public License
  16. * along with this program; if not, write to the Free Software
  17. * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301, USA.
  18. */
  19. /* Author: G. Jungman */
  20. #include "gsl__config.h"
  21. #include "gsl_math.h"
  22. #include "gsl_errno.h"
  23. #include "gsl_sf_pow_int.h"
  24. #include "gsl_sf_bessel.h"
  25. #include "gsl_specfunc__error.h"
  26. #include "gsl_specfunc__bessel.h"
  27. /* i_{l+1}/i_l
  28. */
  29. static
  30. int
  31. bessel_il_CF1(const int l, const double x, const double threshold, double * ratio)
  32. {
  33. const int kmax = 2000;
  34. double tk = 1.0;
  35. double sum = 1.0;
  36. double rhok = 0.0;
  37. int k;
  38. for(k=1; k<=kmax; k++) {
  39. double ak = (x/(2.0*l+1.0+2.0*k)) * (x/(2.0*l+3.0+2.0*k));
  40. rhok = -ak*(1.0 + rhok)/(1.0 + ak*(1.0 + rhok));
  41. tk *= rhok;
  42. sum += tk;
  43. if(fabs(tk/sum) < threshold) break;
  44. }
  45. *ratio = x/(2.0*l+3.0) * sum;
  46. if(k == kmax)
  47. GSL_ERROR ("error", GSL_EMAXITER);
  48. else
  49. return GSL_SUCCESS;
  50. }
  51. /*-*-*-*-*-*-*-*-*-*-*-* Functions with Error Codes *-*-*-*-*-*-*-*-*-*-*-*/
  52. int gsl_sf_bessel_i0_scaled_e(const double x, gsl_sf_result * result)
  53. {
  54. double ax = fabs(x);
  55. /* CHECK_POINTER(result) */
  56. if(x == 0.0) {
  57. result->val = 1.0;
  58. result->err = 0.0;
  59. return GSL_SUCCESS;
  60. }
  61. else if(ax < 0.2) {
  62. const double eax = exp(-ax);
  63. const double y = ax*ax;
  64. const double c1 = 1.0/6.0;
  65. const double c2 = 1.0/120.0;
  66. const double c3 = 1.0/5040.0;
  67. const double c4 = 1.0/362880.0;
  68. const double c5 = 1.0/39916800.0;
  69. const double sum = 1.0 + y*(c1 + y*(c2 + y*(c3 + y*(c4 + y*c5))));
  70. result->val = eax * sum;
  71. result->err = 2.0 * GSL_DBL_EPSILON * result->val;
  72. }
  73. else if(ax < -0.5*GSL_LOG_DBL_EPSILON) {
  74. result->val = (1.0 - exp(-2.0*ax))/(2.0*ax);
  75. result->err = 2.0 * GSL_DBL_EPSILON * result->val;
  76. }
  77. else {
  78. result->val = 1.0/(2.0*ax);
  79. result->err = 2.0 * GSL_DBL_EPSILON * result->val;
  80. }
  81. return GSL_SUCCESS;
  82. }
  83. int gsl_sf_bessel_i1_scaled_e(const double x, gsl_sf_result * result)
  84. {
  85. double ax = fabs(x);
  86. /* CHECK_POINTER(result) */
  87. if(x == 0.0) {
  88. result->val = 0.0;
  89. result->err = 0.0;
  90. return GSL_SUCCESS;
  91. }
  92. else if(ax < 3.0*GSL_DBL_MIN) {
  93. UNDERFLOW_ERROR(result);
  94. }
  95. else if(ax < 0.25) {
  96. const double eax = exp(-ax);
  97. const double y = x*x;
  98. const double c1 = 1.0/10.0;
  99. const double c2 = 1.0/280.0;
  100. const double c3 = 1.0/15120.0;
  101. const double c4 = 1.0/1330560.0;
  102. const double c5 = 1.0/172972800.0;
  103. const double sum = 1.0 + y*(c1 + y*(c2 + y*(c3 + y*(c4 + y*c5))));
  104. result->val = eax * x/3.0 * sum;
  105. result->err = 2.0 * GSL_DBL_EPSILON * fabs(result->val);
  106. return GSL_SUCCESS;
  107. }
  108. else {
  109. double ex = exp(-2.0*ax);
  110. result->val = 0.5 * (ax*(1.0+ex) - (1.0-ex)) / (ax*ax);
  111. result->err = 2.0 * GSL_DBL_EPSILON * fabs(result->val);
  112. if(x < 0.0) result->val = -result->val;
  113. return GSL_SUCCESS;
  114. }
  115. }
  116. int gsl_sf_bessel_i2_scaled_e(const double x, gsl_sf_result * result)
  117. {
  118. double ax = fabs(x);
  119. /* CHECK_POINTER(result) */
  120. if(x == 0.0) {
  121. result->val = 0.0;
  122. result->err = 0.0;
  123. return GSL_SUCCESS;
  124. }
  125. else if(ax < 4.0*GSL_SQRT_DBL_MIN) {
  126. UNDERFLOW_ERROR(result);
  127. }
  128. else if(ax < 0.25) {
  129. const double y = x*x;
  130. const double c1 = 1.0/14.0;
  131. const double c2 = 1.0/504.0;
  132. const double c3 = 1.0/33264.0;
  133. const double c4 = 1.0/3459456.0;
  134. const double c5 = 1.0/518918400.0;
  135. const double sum = 1.0 + y*(c1 + y*(c2 + y*(c3 + y*(c4 + y*c5))));
  136. const double pre = exp(-ax) * x*x/15.0;
  137. result->val = pre * sum;
  138. result->err = 2.0 * GSL_DBL_EPSILON * fabs(result->val);
  139. return GSL_SUCCESS;
  140. }
  141. else {
  142. double ex = exp(-2.0*ax);
  143. double x2 = x*x;
  144. result->val = 0.5 * ((3.0+x2)*(1.0-ex) - 3.0*ax*(1.0+ex))/(ax*ax*ax);
  145. result->err = 2.0 * GSL_DBL_EPSILON * fabs(result->val);
  146. return GSL_SUCCESS;
  147. }
  148. }
  149. int gsl_sf_bessel_il_scaled_e(const int l, double x, gsl_sf_result * result)
  150. {
  151. double sgn = 1.0;
  152. double ax = fabs(x);
  153. if(x < 0.0) {
  154. /* i_l(-x) = (-1)^l i_l(x) */
  155. sgn = ( GSL_IS_ODD(l) ? -1.0 : 1.0 );
  156. x = -x;
  157. }
  158. if(l < 0) {
  159. DOMAIN_ERROR(result);
  160. }
  161. else if(x == 0.0) {
  162. result->val = ( l == 0 ? 1.0 : 0.0 );
  163. result->err = 0.0;
  164. return GSL_SUCCESS;
  165. }
  166. else if(l == 0) {
  167. gsl_sf_result il;
  168. int stat_il = gsl_sf_bessel_i0_scaled_e(x, &il);
  169. result->val = sgn * il.val;
  170. result->err = il.err;
  171. return stat_il;
  172. }
  173. else if(l == 1) {
  174. gsl_sf_result il;
  175. int stat_il = gsl_sf_bessel_i1_scaled_e(x, &il);
  176. result->val = sgn * il.val;
  177. result->err = il.err;
  178. return stat_il;
  179. }
  180. else if(l == 2) {
  181. gsl_sf_result il;
  182. int stat_il = gsl_sf_bessel_i2_scaled_e(x, &il);
  183. result->val = sgn * il.val;
  184. result->err = il.err;
  185. return stat_il;
  186. }
  187. else if(x*x < 10.0*(l+1.5)/M_E) {
  188. gsl_sf_result b;
  189. int stat = gsl_sf_bessel_IJ_taylor_e(l+0.5, x, 1, 50, GSL_DBL_EPSILON, &b);
  190. double pre = exp(-ax) * sqrt((0.5*M_PI)/x);
  191. result->val = sgn * pre * b.val;
  192. result->err = pre * b.err;
  193. result->err += 2.0 * GSL_DBL_EPSILON * fabs(result->val);
  194. return stat;
  195. }
  196. else if(l < 150) {
  197. gsl_sf_result i0_scaled;
  198. int stat_i0 = gsl_sf_bessel_i0_scaled_e(ax, &i0_scaled);
  199. double rat;
  200. int stat_CF1 = bessel_il_CF1(l, ax, GSL_DBL_EPSILON, &rat);
  201. double iellp1 = rat * GSL_SQRT_DBL_MIN;
  202. double iell = GSL_SQRT_DBL_MIN;
  203. double iellm1;
  204. int ell;
  205. for(ell = l; ell >= 1; ell--) {
  206. iellm1 = iellp1 + (2*ell + 1)/x * iell;
  207. iellp1 = iell;
  208. iell = iellm1;
  209. }
  210. result->val = sgn * i0_scaled.val * (GSL_SQRT_DBL_MIN / iell);
  211. result->err = i0_scaled.err * (GSL_SQRT_DBL_MIN / iell);
  212. result->err += 2.0 * GSL_DBL_EPSILON * fabs(result->val);
  213. return GSL_ERROR_SELECT_2(stat_i0, stat_CF1);
  214. }
  215. else if(GSL_MIN(0.29/(l*l+1.0), 0.5/(l*l+1.0+x*x)) < 0.5*GSL_ROOT3_DBL_EPSILON) {
  216. int status = gsl_sf_bessel_Inu_scaled_asymp_unif_e(l + 0.5, x, result);
  217. double pre = sqrt((0.5*M_PI)/x);
  218. result->val *= sgn * pre;
  219. result->err *= pre;
  220. return status;
  221. }
  222. else {
  223. /* recurse down from safe values */
  224. double rt_term = sqrt((0.5*M_PI)/x);
  225. const int LMAX = 2 + (int) (1.2 / GSL_ROOT6_DBL_EPSILON);
  226. gsl_sf_result r_iellp1;
  227. gsl_sf_result r_iell;
  228. int stat_a1 = gsl_sf_bessel_Inu_scaled_asymp_unif_e(LMAX + 1 + 0.5, x, &r_iellp1);
  229. int stat_a2 = gsl_sf_bessel_Inu_scaled_asymp_unif_e(LMAX + 0.5, x, &r_iell);
  230. double iellp1 = r_iellp1.val;
  231. double iell = r_iell.val;
  232. double iellm1 = 0.0;
  233. int ell;
  234. iellp1 *= rt_term;
  235. iell *= rt_term;
  236. for(ell = LMAX; ell >= l+1; ell--) {
  237. iellm1 = iellp1 + (2*ell + 1)/x * iell;
  238. iellp1 = iell;
  239. iell = iellm1;
  240. }
  241. result->val = sgn * iellm1;
  242. result->err = fabs(result->val)*(GSL_DBL_EPSILON + fabs(r_iellp1.err/r_iellp1.val) + fabs(r_iell.err/r_iell.val));
  243. result->err += 2.0 * GSL_DBL_EPSILON * fabs(result->val);
  244. return GSL_ERROR_SELECT_2(stat_a1, stat_a2);
  245. }
  246. }
  247. int gsl_sf_bessel_il_scaled_array(const int lmax, const double x, double * result_array)
  248. {
  249. if(x == 0.0) {
  250. int ell;
  251. result_array[0] = 1.0;
  252. for (ell = lmax; ell >= 1; ell--) {
  253. result_array[ell] = 0.0;
  254. };
  255. return GSL_SUCCESS;
  256. } else {
  257. int ell;
  258. gsl_sf_result r_iellp1;
  259. gsl_sf_result r_iell;
  260. int stat_0 = gsl_sf_bessel_il_scaled_e(lmax+1, x, &r_iellp1);
  261. int stat_1 = gsl_sf_bessel_il_scaled_e(lmax, x, &r_iell);
  262. double iellp1 = r_iellp1.val;
  263. double iell = r_iell.val;
  264. double iellm1;
  265. result_array[lmax] = iell;
  266. for(ell = lmax; ell >= 1; ell--) {
  267. iellm1 = iellp1 + (2*ell + 1)/x * iell;
  268. iellp1 = iell;
  269. iell = iellm1;
  270. result_array[ell-1] = iellm1;
  271. }
  272. return GSL_ERROR_SELECT_2(stat_0, stat_1);
  273. }
  274. }
  275. /*-*-*-*-*-*-*-*-*-* Functions w/ Natural Prototypes *-*-*-*-*-*-*-*-*-*-*/
  276. #include "gsl_specfunc__eval.h"
  277. double gsl_sf_bessel_i0_scaled(const double x)
  278. {
  279. EVAL_RESULT(gsl_sf_bessel_i0_scaled_e(x, &result));
  280. }
  281. double gsl_sf_bessel_i1_scaled(const double x)
  282. {
  283. EVAL_RESULT(gsl_sf_bessel_i1_scaled_e(x, &result));
  284. }
  285. double gsl_sf_bessel_i2_scaled(const double x)
  286. {
  287. EVAL_RESULT(gsl_sf_bessel_i2_scaled_e(x, &result));
  288. }
  289. double gsl_sf_bessel_il_scaled(const int l, const double x)
  290. {
  291. EVAL_RESULT(gsl_sf_bessel_il_scaled_e(l, x, &result));
  292. }