gsl_specfunc__bessel_Yn.c 5.4 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218
  1. /* specfunc/bessel_Yn.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_gamma.h"
  24. #include "gsl_sf_psi.h"
  25. #include "gsl_sf_bessel.h"
  26. #include "gsl_specfunc__error.h"
  27. #include "gsl_specfunc__bessel.h"
  28. #include "gsl_specfunc__bessel_amp_phase.h"
  29. #include "gsl_specfunc__bessel_olver.h"
  30. /*-*-*-*-*-*-*-*-*-*-*-* Private Section *-*-*-*-*-*-*-*-*-*-*-*/
  31. /* assumes n >= 1 */
  32. static int bessel_Yn_small_x(const int n, const double x, gsl_sf_result * result)
  33. {
  34. int k;
  35. double y = 0.25 * x * x;
  36. double ln_x_2 = log(0.5*x);
  37. gsl_sf_result ln_nm1_fact;
  38. double k_term;
  39. double term1, sum1, ln_pre1;
  40. double term2, sum2, pre2;
  41. gsl_sf_lnfact_e((unsigned int)(n-1), &ln_nm1_fact);
  42. ln_pre1 = -n*ln_x_2 + ln_nm1_fact.val;
  43. if(ln_pre1 > GSL_LOG_DBL_MAX - 3.0) GSL_ERROR ("error", GSL_EOVRFLW);
  44. sum1 = 1.0;
  45. k_term = 1.0;
  46. for(k=1; k<=n-1; k++) {
  47. k_term *= y/(k * (n-k));
  48. sum1 += k_term;
  49. }
  50. term1 = -exp(ln_pre1) * sum1 / M_PI;
  51. pre2 = -exp(n*ln_x_2) / M_PI;
  52. if(fabs(pre2) > 0.0) {
  53. const int KMAX = 20;
  54. gsl_sf_result psi_n;
  55. gsl_sf_result npk_fact;
  56. double yk = 1.0;
  57. double k_fact = 1.0;
  58. double psi_kp1 = -M_EULER;
  59. double psi_npkp1;
  60. gsl_sf_psi_int_e(n, &psi_n);
  61. gsl_sf_fact_e((unsigned int)n, &npk_fact);
  62. psi_npkp1 = psi_n.val + 1.0/n;
  63. sum2 = (psi_kp1 + psi_npkp1 - 2.0*ln_x_2)/npk_fact.val;
  64. for(k=1; k<KMAX; k++) {
  65. psi_kp1 += 1./k;
  66. psi_npkp1 += 1./(n+k);
  67. k_fact *= k;
  68. npk_fact.val *= n+k;
  69. yk *= -y;
  70. k_term = yk*(psi_kp1 + psi_npkp1 - 2.0*ln_x_2)/(k_fact*npk_fact.val);
  71. sum2 += k_term;
  72. }
  73. term2 = pre2 * sum2;
  74. }
  75. else {
  76. term2 = 0.0;
  77. }
  78. result->val = term1 + term2;
  79. result->err = GSL_DBL_EPSILON * (fabs(ln_pre1)*fabs(term1) + fabs(term2));
  80. result->err += 2.0 * GSL_DBL_EPSILON * fabs(result->val);
  81. return GSL_SUCCESS;
  82. }
  83. /*-*-*-*-*-*-*-*-*-*-*-* Functions with Error Codes *-*-*-*-*-*-*-*-*-*-*-*/
  84. int
  85. gsl_sf_bessel_Yn_e(int n, const double x, gsl_sf_result * result)
  86. {
  87. int sign = 1;
  88. if(n < 0) {
  89. /* reduce to case n >= 0 */
  90. n = -n;
  91. if(GSL_IS_ODD(n)) sign = -1;
  92. }
  93. /* CHECK_POINTER(result) */
  94. if(n == 0) {
  95. int status = gsl_sf_bessel_Y0_e(x, result);
  96. result->val *= sign;
  97. return status;
  98. }
  99. else if(n == 1) {
  100. int status = gsl_sf_bessel_Y1_e(x, result);
  101. result->val *= sign;
  102. return status;
  103. }
  104. else {
  105. if(x <= 0.0) {
  106. DOMAIN_ERROR(result);
  107. }
  108. if(x < 5.0) {
  109. int status = bessel_Yn_small_x(n, x, result);
  110. result->val *= sign;
  111. return status;
  112. }
  113. else if(GSL_ROOT3_DBL_EPSILON * x > (n*n + 1.0)) {
  114. int status = gsl_sf_bessel_Ynu_asympx_e((double)n, x, result);
  115. result->val *= sign;
  116. return status;
  117. }
  118. else if(n > 50) {
  119. int status = gsl_sf_bessel_Ynu_asymp_Olver_e((double)n, x, result);
  120. result->val *= sign;
  121. return status;
  122. }
  123. else {
  124. double two_over_x = 2.0/x;
  125. gsl_sf_result r_by;
  126. gsl_sf_result r_bym;
  127. int stat_1 = gsl_sf_bessel_Y1_e(x, &r_by);
  128. int stat_0 = gsl_sf_bessel_Y0_e(x, &r_bym);
  129. double bym = r_bym.val;
  130. double by = r_by.val;
  131. double byp;
  132. int j;
  133. for(j=1; j<n; j++) {
  134. byp = j*two_over_x*by - bym;
  135. bym = by;
  136. by = byp;
  137. }
  138. result->val = sign * by;
  139. result->err = fabs(result->val) * (fabs(r_by.err/r_by.val) + fabs(r_bym.err/r_bym.val));
  140. result->err += 2.0 * GSL_DBL_EPSILON * fabs(result->val);
  141. return GSL_ERROR_SELECT_2(stat_1, stat_0);
  142. }
  143. }
  144. }
  145. int
  146. gsl_sf_bessel_Yn_array(const int nmin, const int nmax, const double x, double * result_array)
  147. {
  148. /* CHECK_POINTER(result_array) */
  149. if(nmin < 0 || nmax < nmin || x <= 0.0) {
  150. int j;
  151. for(j=0; j<=nmax-nmin; j++) result_array[j] = 0.0;
  152. GSL_ERROR ("error", GSL_EDOM);
  153. }
  154. else {
  155. gsl_sf_result r_Ynm1;
  156. gsl_sf_result r_Yn;
  157. int stat_nm1 = gsl_sf_bessel_Yn_e(nmin, x, &r_Ynm1);
  158. int stat_n = gsl_sf_bessel_Yn_e(nmin+1, x, &r_Yn);
  159. double Ynp1;
  160. double Yn = r_Yn.val;
  161. double Ynm1 = r_Ynm1.val;
  162. int n;
  163. int stat = GSL_ERROR_SELECT_2(stat_nm1, stat_n);
  164. if(stat == GSL_SUCCESS) {
  165. for(n=nmin+1; n<=nmax+1; n++) {
  166. result_array[n-nmin-1] = Ynm1;
  167. Ynp1 = -Ynm1 + 2.0*n/x * Yn;
  168. Ynm1 = Yn;
  169. Yn = Ynp1;
  170. }
  171. }
  172. else {
  173. for(n=nmin; n<=nmax; n++) {
  174. result_array[n-nmin] = 0.0;
  175. }
  176. }
  177. return stat;
  178. }
  179. }
  180. /*-*-*-*-*-*-*-*-*-* Functions w/ Natural Prototypes *-*-*-*-*-*-*-*-*-*-*/
  181. #include "gsl_specfunc__eval.h"
  182. double gsl_sf_bessel_Yn(const int n, const double x)
  183. {
  184. EVAL_RESULT(gsl_sf_bessel_Yn_e(n, x, &result));
  185. }