gsl_specfunc__legendre_con.c 42 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813814815816817818819820821822823824825826827828829830831832833834835836837838839840841842843844845846847848849850851852853854855856857858859860861862863864865866867868869870871872873874875876877878879880881882883884885886887888889890891892893894895896897898899900901902903904905906907908909910911912913914915916917918919920921922923924925926927928929930931932933934935936937938939940941942943944945946947948949950951952953954955956957958959960961962963964965966967968969970971972973974975976977978979980981982983984985986987988989990991992993994995996997998999100010011002100310041005100610071008100910101011101210131014101510161017101810191020102110221023102410251026102710281029103010311032103310341035103610371038103910401041104210431044104510461047104810491050105110521053105410551056105710581059106010611062106310641065106610671068106910701071107210731074107510761077107810791080108110821083108410851086108710881089109010911092109310941095109610971098109911001101110211031104110511061107110811091110111111121113111411151116111711181119112011211122112311241125112611271128112911301131113211331134113511361137113811391140114111421143114411451146114711481149115011511152115311541155115611571158115911601161116211631164116511661167116811691170117111721173117411751176117711781179118011811182118311841185118611871188118911901191119211931194119511961197119811991200120112021203120412051206120712081209121012111212121312141215121612171218121912201221122212231224122512261227122812291230123112321233123412351236123712381239124012411242124312441245124612471248124912501251125212531254125512561257125812591260126112621263126412651266126712681269127012711272127312741275127612771278127912801281128212831284128512861287128812891290129112921293129412951296129712981299130013011302130313041305130613071308130913101311131213131314131513161317131813191320132113221323132413251326132713281329133013311332133313341335133613371338133913401341134213431344134513461347134813491350135113521353135413551356135713581359136013611362136313641365136613671368136913701371137213731374
  1. /* specfunc/legendre_con.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_poly.h"
  24. #include "gsl_sf_exp.h"
  25. #include "gsl_sf_trig.h"
  26. #include "gsl_sf_gamma.h"
  27. #include "gsl_sf_ellint.h"
  28. #include "gsl_sf_pow_int.h"
  29. #include "gsl_sf_bessel.h"
  30. #include "gsl_sf_hyperg.h"
  31. #include "gsl_sf_legendre.h"
  32. #include "gsl_specfunc__error.h"
  33. #include "gsl_specfunc__legendre.h"
  34. #define Root_2OverPi_ 0.797884560802865355879892
  35. #define locEPS (1000.0*GSL_DBL_EPSILON)
  36. /*-*-*-*-*-*-*-*-*-*-*-* Private Section *-*-*-*-*-*-*-*-*-*-*-*/
  37. #define RECURSE_LARGE (1.0e-5*GSL_DBL_MAX)
  38. #define RECURSE_SMALL (1.0e+5*GSL_DBL_MIN)
  39. /* Continued fraction for f_{ell+1}/f_ell
  40. * f_ell := P^{-mu-ell}_{-1/2 + I tau}(x), x < 1.0
  41. *
  42. * Uses standard CF method from Temme's book.
  43. */
  44. static
  45. int
  46. conicalP_negmu_xlt1_CF1(const double mu, const int ell, const double tau,
  47. const double x, gsl_sf_result * result)
  48. {
  49. const double RECUR_BIG = GSL_SQRT_DBL_MAX;
  50. const int maxiter = 5000;
  51. int n = 1;
  52. double xi = x/(sqrt(1.0-x)*sqrt(1.0+x));
  53. double Anm2 = 1.0;
  54. double Bnm2 = 0.0;
  55. double Anm1 = 0.0;
  56. double Bnm1 = 1.0;
  57. double a1 = 1.0;
  58. double b1 = 2.0*(mu + ell + 1.0) * xi;
  59. double An = b1*Anm1 + a1*Anm2;
  60. double Bn = b1*Bnm1 + a1*Bnm2;
  61. double an, bn;
  62. double fn = An/Bn;
  63. while(n < maxiter) {
  64. double old_fn;
  65. double del;
  66. n++;
  67. Anm2 = Anm1;
  68. Bnm2 = Bnm1;
  69. Anm1 = An;
  70. Bnm1 = Bn;
  71. an = tau*tau + (mu - 0.5 + ell + n)*(mu - 0.5 + ell + n);
  72. bn = 2.0*(ell + mu + n) * xi;
  73. An = bn*Anm1 + an*Anm2;
  74. Bn = bn*Bnm1 + an*Bnm2;
  75. if(fabs(An) > RECUR_BIG || fabs(Bn) > RECUR_BIG) {
  76. An /= RECUR_BIG;
  77. Bn /= RECUR_BIG;
  78. Anm1 /= RECUR_BIG;
  79. Bnm1 /= RECUR_BIG;
  80. Anm2 /= RECUR_BIG;
  81. Bnm2 /= RECUR_BIG;
  82. }
  83. old_fn = fn;
  84. fn = An/Bn;
  85. del = old_fn/fn;
  86. if(fabs(del - 1.0) < 2.0*GSL_DBL_EPSILON) break;
  87. }
  88. result->val = fn;
  89. result->err = 4.0 * GSL_DBL_EPSILON * (sqrt(n) + 1.0) * fabs(fn);
  90. if(n >= maxiter)
  91. GSL_ERROR ("error", GSL_EMAXITER);
  92. else
  93. return GSL_SUCCESS;
  94. }
  95. /* Continued fraction for f_{ell+1}/f_ell
  96. * f_ell := P^{-mu-ell}_{-1/2 + I tau}(x), x >= 1.0
  97. *
  98. * Uses Gautschi (Euler) equivalent series.
  99. */
  100. static
  101. int
  102. conicalP_negmu_xgt1_CF1(const double mu, const int ell, const double tau,
  103. const double x, gsl_sf_result * result)
  104. {
  105. const int maxk = 20000;
  106. const double gamma = 1.0-1.0/(x*x);
  107. const double pre = sqrt(x-1.0)*sqrt(x+1.0) / (x*(2.0*(ell+mu+1.0)));
  108. double tk = 1.0;
  109. double sum = 1.0;
  110. double rhok = 0.0;
  111. int k;
  112. for(k=1; k<maxk; k++) {
  113. double tlk = 2.0*(ell + mu + k);
  114. double l1k = (ell + mu - 0.5 + 1.0 + k);
  115. double ak = -(tau*tau + l1k*l1k)/(tlk*(tlk+2.0)) * gamma;
  116. rhok = -ak*(1.0 + rhok)/(1.0 + ak*(1.0 + rhok));
  117. tk *= rhok;
  118. sum += tk;
  119. if(fabs(tk/sum) < GSL_DBL_EPSILON) break;
  120. }
  121. result->val = pre * sum;
  122. result->err = fabs(pre * tk);
  123. result->err += 2.0 * GSL_DBL_EPSILON * (sqrt(k) + 1.0) * fabs(pre*sum);
  124. if(k >= maxk)
  125. GSL_ERROR ("error", GSL_EMAXITER);
  126. else
  127. return GSL_SUCCESS;
  128. }
  129. /* Implementation of large negative mu asymptotic
  130. * [Dunster, Proc. Roy. Soc. Edinburgh 119A, 311 (1991), p. 326]
  131. */
  132. inline
  133. static double olver_U1(double beta2, double p)
  134. {
  135. return (p-1.0)/(24.0*(1.0+beta2)) * (3.0 + beta2*(2.0 + 5.0*p*(1.0+p)));
  136. }
  137. inline
  138. static double olver_U2(double beta2, double p)
  139. {
  140. double beta4 = beta2*beta2;
  141. double p2 = p*p;
  142. double poly1 = 4.0*beta4 + 84.0*beta2 - 63.0;
  143. double poly2 = 16.0*beta4 + 90.0*beta2 - 81.0;
  144. double poly3 = beta2*p2*(97.0*beta2 - 432.0 + 77.0*p*(beta2-6.0) - 385.0*beta2*p2*(1.0 + p));
  145. return (1.0-p)/(1152.0*(1.0+beta2)) * (poly1 + poly2 + poly3);
  146. }
  147. static const double U3c1[] = { -1307.0, -1647.0, 3375.0, 3675.0 };
  148. static const double U3c2[] = { 29366.0, 35835.0, -252360.0, -272630.0,
  149. 276810.0, 290499.0 };
  150. static const double U3c3[] = { -29748.0, -8840.0, 1725295.0, 1767025.0,
  151. -7313470.0, -754778.0, 6309875.0, 6480045.0 };
  152. static const double U3c4[] = { 2696.0, -16740.0, -524250.0, -183975.0,
  153. 14670540.0, 14172939.0, -48206730.0, -48461985.0,
  154. 36756720.0, 37182145.0 };
  155. static const double U3c5[] = { 9136.0, 22480.0, 12760.0,
  156. -252480.0, -4662165.0, -1705341.0,
  157. 92370135.0, 86244015.0, -263678415.0,
  158. -260275015.0, 185910725.0, 185910725.0 };
  159. #if 0
  160. static double olver_U3(double beta2, double p)
  161. {
  162. double beta4 = beta2*beta2;
  163. double beta6 = beta4*beta2;
  164. double opb2s = (1.0+beta2)*(1.0+beta2);
  165. double den = 39813120.0 * opb2s*opb2s;
  166. double poly1 = gsl_poly_eval(U3c1, 4, p);
  167. double poly2 = gsl_poly_eval(U3c2, 6, p);
  168. double poly3 = gsl_poly_eval(U3c3, 8, p);
  169. double poly4 = gsl_poly_eval(U3c4, 10, p);
  170. double poly5 = gsl_poly_eval(U3c5, 12, p);
  171. return (p-1.0)*( 1215.0*poly1 + 324.0*beta2*poly2
  172. + 54.0*beta4*poly3 + 12.0*beta6*poly4
  173. + beta4*beta4*poly5
  174. ) / den;
  175. }
  176. #endif /* 0 */
  177. /* Large negative mu asymptotic
  178. * P^{-mu}_{-1/2 + I tau}, mu -> Inf
  179. * |x| < 1
  180. *
  181. * [Dunster, Proc. Roy. Soc. Edinburgh 119A, 311 (1991), p. 326]
  182. */
  183. int
  184. gsl_sf_conicalP_xlt1_large_neg_mu_e(double mu, double tau, double x,
  185. gsl_sf_result * result, double * ln_multiplier)
  186. {
  187. double beta = tau/mu;
  188. double beta2 = beta*beta;
  189. double S = beta * acos((1.0-beta2)/(1.0+beta2));
  190. double p = x/sqrt(beta2*(1.0-x*x) + 1.0);
  191. gsl_sf_result lg_mup1;
  192. int lg_stat = gsl_sf_lngamma_e(mu+1.0, &lg_mup1);
  193. double ln_pre_1 = 0.5*mu*(S - log(1.0+beta2) + log((1.0-p)/(1.0+p))) - lg_mup1.val;
  194. double ln_pre_2 = -0.25 * log(1.0 + beta2*(1.0-x));
  195. double ln_pre_3 = -tau * atan(p*beta);
  196. double ln_pre = ln_pre_1 + ln_pre_2 + ln_pre_3;
  197. double sum = 1.0 - olver_U1(beta2, p)/mu + olver_U2(beta2, p)/(mu*mu);
  198. if(sum == 0.0) {
  199. result->val = 0.0;
  200. result->err = 0.0;
  201. *ln_multiplier = 0.0;
  202. return GSL_SUCCESS;
  203. }
  204. else {
  205. int stat_e = gsl_sf_exp_mult_e(ln_pre, sum, result);
  206. if(stat_e != GSL_SUCCESS) {
  207. result->val = sum;
  208. result->err = 2.0 * GSL_DBL_EPSILON * fabs(sum);
  209. *ln_multiplier = ln_pre;
  210. }
  211. else {
  212. *ln_multiplier = 0.0;
  213. }
  214. return lg_stat;
  215. }
  216. }
  217. /* Implementation of large tau asymptotic
  218. *
  219. * A_n^{-mu}, B_n^{-mu} [Olver, p.465, 469]
  220. */
  221. inline
  222. static double olver_B0_xi(double mu, double xi)
  223. {
  224. return (1.0 - 4.0*mu*mu)/(8.0*xi) * (1.0/tanh(xi) - 1.0/xi);
  225. }
  226. static double olver_A1_xi(double mu, double xi, double x)
  227. {
  228. double B = olver_B0_xi(mu, xi);
  229. double psi;
  230. if(fabs(x - 1.0) < GSL_ROOT4_DBL_EPSILON) {
  231. double y = x - 1.0;
  232. double s = -1.0/3.0 + y*(2.0/15.0 - y *(61.0/945.0 - 452.0/14175.0*y));
  233. psi = (4.0*mu*mu - 1.0)/16.0 * s;
  234. }
  235. else {
  236. psi = (4.0*mu*mu - 1.0)/16.0 * (1.0/(x*x-1.0) - 1.0/(xi*xi));
  237. }
  238. return 0.5*xi*xi*B*B + (mu+0.5)*B - psi + mu/6.0*(0.25 - mu*mu);
  239. }
  240. inline
  241. static double olver_B0_th(double mu, double theta)
  242. {
  243. return -(1.0 - 4.0*mu*mu)/(8.0*theta) * (1.0/tan(theta) - 1.0/theta);
  244. }
  245. static double olver_A1_th(double mu, double theta, double x)
  246. {
  247. double B = olver_B0_th(mu, theta);
  248. double psi;
  249. if(fabs(x - 1.0) < GSL_ROOT4_DBL_EPSILON) {
  250. double y = 1.0 - x;
  251. double s = -1.0/3.0 + y*(2.0/15.0 - y *(61.0/945.0 - 452.0/14175.0*y));
  252. psi = (4.0*mu*mu - 1.0)/16.0 * s;
  253. }
  254. else {
  255. psi = (4.0*mu*mu - 1.0)/16.0 * (1.0/(x*x-1.0) + 1.0/(theta*theta));
  256. }
  257. return -0.5*theta*theta*B*B + (mu+0.5)*B - psi + mu/6.0*(0.25 - mu*mu);
  258. }
  259. /* Large tau uniform asymptotics
  260. * P^{-mu}_{-1/2 + I tau}
  261. * 1 < x
  262. * tau -> Inf
  263. * [Olver, p. 469]
  264. */
  265. int
  266. gsl_sf_conicalP_xgt1_neg_mu_largetau_e(const double mu, const double tau,
  267. const double x, double acosh_x,
  268. gsl_sf_result * result, double * ln_multiplier)
  269. {
  270. double xi = acosh_x;
  271. double ln_xi_pre;
  272. double ln_pre;
  273. double sumA, sumB, sum;
  274. double arg;
  275. gsl_sf_result J_mup1;
  276. gsl_sf_result J_mu;
  277. double J_mum1;
  278. if(xi < GSL_ROOT4_DBL_EPSILON) {
  279. ln_xi_pre = -xi*xi/6.0; /* log(1.0 - xi*xi/6.0) */
  280. }
  281. else {
  282. gsl_sf_result lnshxi;
  283. gsl_sf_lnsinh_e(xi, &lnshxi);
  284. ln_xi_pre = log(xi) - lnshxi.val; /* log(xi/sinh(xi) */
  285. }
  286. ln_pre = 0.5*ln_xi_pre - mu*log(tau);
  287. arg = tau*xi;
  288. gsl_sf_bessel_Jnu_e(mu + 1.0, arg, &J_mup1);
  289. gsl_sf_bessel_Jnu_e(mu, arg, &J_mu);
  290. J_mum1 = -J_mup1.val + 2.0*mu/arg*J_mu.val; /* careful of mu < 1 */
  291. sumA = 1.0 - olver_A1_xi(-mu, xi, x)/(tau*tau);
  292. sumB = olver_B0_xi(-mu, xi);
  293. sum = J_mu.val * sumA - xi/tau * J_mum1 * sumB;
  294. if(sum == 0.0) {
  295. result->val = 0.0;
  296. result->err = 0.0;
  297. *ln_multiplier = 0.0;
  298. return GSL_SUCCESS;
  299. }
  300. else {
  301. int stat_e = gsl_sf_exp_mult_e(ln_pre, sum, result);
  302. if(stat_e != GSL_SUCCESS) {
  303. result->val = sum;
  304. result->err = 2.0 * GSL_DBL_EPSILON * fabs(sum);
  305. *ln_multiplier = ln_pre;
  306. }
  307. else {
  308. *ln_multiplier = 0.0;
  309. }
  310. return GSL_SUCCESS;
  311. }
  312. }
  313. /* Large tau uniform asymptotics
  314. * P^{-mu}_{-1/2 + I tau}
  315. * -1 < x < 1
  316. * tau -> Inf
  317. * [Olver, p. 473]
  318. */
  319. int
  320. gsl_sf_conicalP_xlt1_neg_mu_largetau_e(const double mu, const double tau,
  321. const double x, const double acos_x,
  322. gsl_sf_result * result, double * ln_multiplier)
  323. {
  324. double theta = acos_x;
  325. double ln_th_pre;
  326. double ln_pre;
  327. double sumA, sumB, sum, sumerr;
  328. double arg;
  329. gsl_sf_result I_mup1, I_mu;
  330. double I_mum1;
  331. if(theta < GSL_ROOT4_DBL_EPSILON) {
  332. ln_th_pre = theta*theta/6.0; /* log(1.0 + theta*theta/6.0) */
  333. }
  334. else {
  335. ln_th_pre = log(theta/sin(theta));
  336. }
  337. ln_pre = 0.5 * ln_th_pre - mu * log(tau);
  338. arg = tau*theta;
  339. gsl_sf_bessel_Inu_e(mu + 1.0, arg, &I_mup1);
  340. gsl_sf_bessel_Inu_e(mu, arg, &I_mu);
  341. I_mum1 = I_mup1.val + 2.0*mu/arg * I_mu.val; /* careful of mu < 1 */
  342. sumA = 1.0 - olver_A1_th(-mu, theta, x)/(tau*tau);
  343. sumB = olver_B0_th(-mu, theta);
  344. sum = I_mu.val * sumA - theta/tau * I_mum1 * sumB;
  345. sumerr = fabs(I_mu.err * sumA);
  346. sumerr += fabs(I_mup1.err * theta/tau * sumB);
  347. sumerr += fabs(I_mu.err * theta/tau * sumB * 2.0 * mu/arg);
  348. if(sum == 0.0) {
  349. result->val = 0.0;
  350. result->err = 0.0;
  351. *ln_multiplier = 0.0;
  352. return GSL_SUCCESS;
  353. }
  354. else {
  355. int stat_e = gsl_sf_exp_mult_e(ln_pre, sum, result);
  356. if(stat_e != GSL_SUCCESS) {
  357. result->val = sum;
  358. result->err = sumerr;
  359. result->err += GSL_DBL_EPSILON * fabs(sum);
  360. *ln_multiplier = ln_pre;
  361. }
  362. else {
  363. *ln_multiplier = 0.0;
  364. }
  365. return GSL_SUCCESS;
  366. }
  367. }
  368. /* Hypergeometric function which appears in the
  369. * large x expansion below:
  370. *
  371. * 2F1(1/4 - mu/2 - I tau/2, 3/4 - mu/2 - I tau/2, 1 - I tau, y)
  372. *
  373. * Note that for the usage below y = 1/x^2;
  374. */
  375. static
  376. int
  377. conicalP_hyperg_large_x(const double mu, const double tau, const double y,
  378. double * reF, double * imF)
  379. {
  380. const int kmax = 1000;
  381. const double re_a = 0.25 - 0.5*mu;
  382. const double re_b = 0.75 - 0.5*mu;
  383. const double re_c = 1.0;
  384. const double im_a = -0.5*tau;
  385. const double im_b = -0.5*tau;
  386. const double im_c = -tau;
  387. double re_sum = 1.0;
  388. double im_sum = 0.0;
  389. double re_term = 1.0;
  390. double im_term = 0.0;
  391. int k;
  392. for(k=1; k<=kmax; k++) {
  393. double re_ak = re_a + k - 1.0;
  394. double re_bk = re_b + k - 1.0;
  395. double re_ck = re_c + k - 1.0;
  396. double im_ak = im_a;
  397. double im_bk = im_b;
  398. double im_ck = im_c;
  399. double den = re_ck*re_ck + im_ck*im_ck;
  400. double re_multiplier = ((re_ak*re_bk - im_ak*im_bk)*re_ck + im_ck*(im_ak*re_bk + re_ak*im_bk)) / den;
  401. double im_multiplier = ((im_ak*re_bk + re_ak*im_bk)*re_ck - im_ck*(re_ak*re_bk - im_ak*im_bk)) / den;
  402. double re_tmp = re_multiplier*re_term - im_multiplier*im_term;
  403. double im_tmp = im_multiplier*re_term + re_multiplier*im_term;
  404. double asum = fabs(re_sum) + fabs(im_sum);
  405. re_term = y/k * re_tmp;
  406. im_term = y/k * im_tmp;
  407. if(fabs(re_term/asum) < GSL_DBL_EPSILON && fabs(im_term/asum) < GSL_DBL_EPSILON) break;
  408. re_sum += re_term;
  409. im_sum += im_term;
  410. }
  411. *reF = re_sum;
  412. *imF = im_sum;
  413. if(k == kmax)
  414. GSL_ERROR ("error", GSL_EMAXITER);
  415. else
  416. return GSL_SUCCESS;
  417. }
  418. /* P^{mu}_{-1/2 + I tau}
  419. * x->Inf
  420. */
  421. int
  422. gsl_sf_conicalP_large_x_e(const double mu, const double tau, const double x,
  423. gsl_sf_result * result, double * ln_multiplier)
  424. {
  425. /* 2F1 term
  426. */
  427. double y = ( x < 0.5*GSL_SQRT_DBL_MAX ? 1.0/(x*x) : 0.0 );
  428. double reF, imF;
  429. int stat_F = conicalP_hyperg_large_x(mu, tau, y, &reF, &imF);
  430. /* f = Gamma(+i tau)/Gamma(1/2 - mu + i tau)
  431. * FIXME: shift so it's better for tau-> 0
  432. */
  433. gsl_sf_result lgr_num, lgth_num;
  434. gsl_sf_result lgr_den, lgth_den;
  435. int stat_gn = gsl_sf_lngamma_complex_e(0.0,tau,&lgr_num,&lgth_num);
  436. int stat_gd = gsl_sf_lngamma_complex_e(0.5-mu,tau,&lgr_den,&lgth_den);
  437. double angle = lgth_num.val - lgth_den.val + atan2(imF,reF);
  438. double lnx = log(x);
  439. double lnxp1 = log(x+1.0);
  440. double lnxm1 = log(x-1.0);
  441. double lnpre_const = 0.5*M_LN2 - 0.5*M_LNPI;
  442. double lnpre_comm = (mu-0.5)*lnx - 0.5*mu*(lnxp1 + lnxm1);
  443. double lnpre_err = GSL_DBL_EPSILON * (0.5*M_LN2 + 0.5*M_LNPI)
  444. + GSL_DBL_EPSILON * fabs((mu-0.5)*lnx)
  445. + GSL_DBL_EPSILON * fabs(0.5*mu)*(fabs(lnxp1)+fabs(lnxm1));
  446. /* result = pre*|F|*|f| * cos(angle - tau * (log(x)+M_LN2))
  447. */
  448. gsl_sf_result cos_result;
  449. int stat_cos = gsl_sf_cos_e(angle + tau*(log(x) + M_LN2), &cos_result);
  450. int status = GSL_ERROR_SELECT_4(stat_cos, stat_gd, stat_gn, stat_F);
  451. if(cos_result.val == 0.0) {
  452. result->val = 0.0;
  453. result->err = 0.0;
  454. return status;
  455. }
  456. else {
  457. double lnFf_val = 0.5*log(reF*reF+imF*imF) + lgr_num.val - lgr_den.val;
  458. double lnFf_err = lgr_num.err + lgr_den.err + GSL_DBL_EPSILON * fabs(lnFf_val);
  459. double lnnoc_val = lnpre_const + lnpre_comm + lnFf_val;
  460. double lnnoc_err = lnpre_err + lnFf_err + GSL_DBL_EPSILON * fabs(lnnoc_val);
  461. int stat_e = gsl_sf_exp_mult_err_e(lnnoc_val, lnnoc_err,
  462. cos_result.val, cos_result.err,
  463. result);
  464. if(stat_e == GSL_SUCCESS) {
  465. *ln_multiplier = 0.0;
  466. }
  467. else {
  468. result->val = cos_result.val;
  469. result->err = cos_result.err;
  470. result->err += 2.0 * GSL_DBL_EPSILON * fabs(result->val);
  471. *ln_multiplier = lnnoc_val;
  472. }
  473. return status;
  474. }
  475. }
  476. /* P^{mu}_{-1/2 + I tau} first hypergeometric representation
  477. * -1 < x < 1
  478. * This is more effective for |x| small, however it will work w/o
  479. * reservation for any x < 0 because everything is positive
  480. * definite in that case.
  481. *
  482. * [Kolbig, (3)] (note typo in args of gamma functions)
  483. * [Bateman, (22)] (correct form)
  484. */
  485. static
  486. int
  487. conicalP_xlt1_hyperg_A(double mu, double tau, double x, gsl_sf_result * result)
  488. {
  489. double x2 = x*x;
  490. double err_amp = 1.0 + 1.0/(GSL_DBL_EPSILON + fabs(1.0-fabs(x)));
  491. double pre_val = M_SQRTPI / pow(0.5*sqrt(1-x2), mu);
  492. double pre_err = err_amp * GSL_DBL_EPSILON * (fabs(mu)+1.0) * fabs(pre_val) ;
  493. gsl_sf_result ln_g1, ln_g2, arg_g1, arg_g2;
  494. gsl_sf_result F1, F2;
  495. gsl_sf_result pre1, pre2;
  496. double t1_val, t1_err;
  497. double t2_val, t2_err;
  498. int stat_F1 = gsl_sf_hyperg_2F1_conj_e(0.25 - 0.5*mu, 0.5*tau, 0.5, x2, &F1);
  499. int stat_F2 = gsl_sf_hyperg_2F1_conj_e(0.75 - 0.5*mu, 0.5*tau, 1.5, x2, &F2);
  500. int status = GSL_ERROR_SELECT_2(stat_F1, stat_F2);
  501. gsl_sf_lngamma_complex_e(0.75 - 0.5*mu, -0.5*tau, &ln_g1, &arg_g1);
  502. gsl_sf_lngamma_complex_e(0.25 - 0.5*mu, -0.5*tau, &ln_g2, &arg_g2);
  503. gsl_sf_exp_err_e(-2.0*ln_g1.val, 2.0*ln_g1.err, &pre1);
  504. gsl_sf_exp_err_e(-2.0*ln_g2.val, 2.0*ln_g2.err, &pre2);
  505. pre2.val *= -2.0*x;
  506. pre2.err *= 2.0*fabs(x);
  507. pre2.err += GSL_DBL_EPSILON * fabs(pre2.val);
  508. t1_val = pre1.val * F1.val;
  509. t1_err = fabs(pre1.val) * F1.err + pre1.err * fabs(F1.val);
  510. t2_val = pre2.val * F2.val;
  511. t2_err = fabs(pre2.val) * F2.err + pre2.err * fabs(F2.val);
  512. result->val = pre_val * (t1_val + t2_val);
  513. result->err = pre_val * (t1_err + t2_err);
  514. result->err += pre_err * fabs(t1_val + t2_val);
  515. result->err += 2.0 * GSL_DBL_EPSILON * fabs(result->val);
  516. return status;
  517. }
  518. /* P^{mu}_{-1/2 + I tau}
  519. * defining hypergeometric representation
  520. * [Abramowitz+Stegun, 8.1.2]
  521. * 1 < x < 3
  522. * effective for x near 1
  523. *
  524. */
  525. #if 0
  526. static
  527. int
  528. conicalP_def_hyperg(double mu, double tau, double x, double * result)
  529. {
  530. double F;
  531. int stat_F = gsl_sf_hyperg_2F1_conj_renorm_e(0.5, tau, 1.0-mu, 0.5*(1.0-x), &F);
  532. *result = pow((x+1.0)/(x-1.0), 0.5*mu) * F;
  533. return stat_F;
  534. }
  535. #endif /* 0 */
  536. /* P^{mu}_{-1/2 + I tau} second hypergeometric representation
  537. * [Zhurina+Karmazina, (3.1)]
  538. * -1 < x < 3
  539. * effective for x near 1
  540. *
  541. */
  542. #if 0
  543. static
  544. int
  545. conicalP_xnear1_hyperg_C(double mu, double tau, double x, double * result)
  546. {
  547. double ln_pre, arg_pre;
  548. double ln_g1, arg_g1;
  549. double ln_g2, arg_g2;
  550. double F;
  551. int stat_F = gsl_sf_hyperg_2F1_conj_renorm_e(0.5+mu, tau, 1.0+mu, 0.5*(1.0-x), &F);
  552. gsl_sf_lngamma_complex_e(0.5+mu, tau, &ln_g1, &arg_g1);
  553. gsl_sf_lngamma_complex_e(0.5-mu, tau, &ln_g2, &arg_g2);
  554. ln_pre = mu*M_LN2 - 0.5*mu*log(fabs(x*x-1.0)) + ln_g1 - ln_g2;
  555. arg_pre = arg_g1 - arg_g2;
  556. *result = exp(ln_pre) * F;
  557. return stat_F;
  558. }
  559. #endif /* 0 */
  560. /* V0, V1 from Kolbig, m = 0
  561. */
  562. static
  563. int
  564. conicalP_0_V(const double t, const double f, const double tau, const double sgn,
  565. double * V0, double * V1)
  566. {
  567. double C[8];
  568. double T[8];
  569. double H[8];
  570. double V[12];
  571. int i;
  572. T[0] = 1.0;
  573. H[0] = 1.0;
  574. V[0] = 1.0;
  575. for(i=1; i<=7; i++) {
  576. T[i] = T[i-1] * t;
  577. H[i] = H[i-1] * (t*f);
  578. }
  579. for(i=1; i<=11; i++) {
  580. V[i] = V[i-1] * tau;
  581. }
  582. C[0] = 1.0;
  583. C[1] = (H[1]-1.0)/(8.0*T[1]);
  584. C[2] = (9.0*H[2] + 6.0*H[1] - 15.0 - sgn*8.0*T[2])/(128.0*T[2]);
  585. C[3] = 5.0*(15.0*H[3] + 27.0*H[2] + 21.0*H[1] - 63.0 - sgn*T[2]*(16.0*H[1]+24.0))/(1024.0*T[3]);
  586. C[4] = 7.0*(525.0*H[4] + 1500.0*H[3] + 2430.0*H[2] + 1980.0*H[1] - 6435.0
  587. + 192.0*T[4] - sgn*T[2]*(720.0*H[2]+1600.0*H[1]+2160.0)
  588. ) / (32768.0*T[4]);
  589. C[5] = 21.0*(2835.0*H[5] + 11025.0*H[4] + 24750.0*H[3] + 38610.0*H[2]
  590. + 32175.0*H[1] - 109395.0 + T[4]*(1984.0*H[1]+4032.0)
  591. - sgn*T[2]*(4800.0*H[3]+15120.0*H[2]+26400.0*H[1]+34320.0)
  592. ) / (262144.0*T[5]);
  593. C[6] = 11.0*(218295.0*H[6] + 1071630.0*H[5] + 3009825.0*H[4] + 6142500.0*H[3]
  594. + 9398025.0*H[2] + 7936110.0*H[1] - 27776385.0
  595. + T[4]*(254016.0*H[2]+749952.0*H[1]+1100736.0)
  596. - sgn*T[2]*(441000.0*H[4] + 1814400.0*H[3] + 4127760.0*H[2]
  597. + 6552000.0*H[1] + 8353800.0 + 31232.0*T[4]
  598. )
  599. ) / (4194304.0*T[6]);
  600. *V0 = C[0] + (-4.0*C[3]/T[1]+C[4])/V[4]
  601. + (-192.0*C[5]/T[3]+144.0*C[6]/T[2])/V[8]
  602. + sgn * (-C[2]/V[2]
  603. + (-24.0*C[4]/T[2]+12.0*C[5]/T[1]-C[6])/V[6]
  604. + (-1920.0*C[6]/T[4])/V[10]
  605. );
  606. *V1 = C[1]/V[1] + (8.0*(C[3]/T[2]-C[4]/T[1])+C[5])/V[5]
  607. + (384.0*C[5]/T[4] - 768.0*C[6]/T[3])/V[9]
  608. + sgn * ((2.0*C[2]/T[1]-C[3])/V[3]
  609. + (48.0*C[4]/T[3]-72.0*C[5]/T[2] + 18.0*C[6]/T[1])/V[7]
  610. + (3840.0*C[6]/T[5])/V[11]
  611. );
  612. return GSL_SUCCESS;
  613. }
  614. /* V0, V1 from Kolbig, m = 1
  615. */
  616. static
  617. int
  618. conicalP_1_V(const double t, const double f, const double tau, const double sgn,
  619. double * V0, double * V1)
  620. {
  621. double Cm1;
  622. double C[8];
  623. double T[8];
  624. double H[8];
  625. double V[12];
  626. int i;
  627. T[0] = 1.0;
  628. H[0] = 1.0;
  629. V[0] = 1.0;
  630. for(i=1; i<=7; i++) {
  631. T[i] = T[i-1] * t;
  632. H[i] = H[i-1] * (t*f);
  633. }
  634. for(i=1; i<=11; i++) {
  635. V[i] = V[i-1] * tau;
  636. }
  637. Cm1 = -1.0;
  638. C[0] = 3.0*(1.0-H[1])/(8.0*T[1]);
  639. C[1] = (-15.0*H[2]+6.0*H[1]+9.0+sgn*8.0*T[2])/(128.0*T[2]);
  640. C[2] = 3.0*(-35.0*H[3] - 15.0*H[2] + 15.0*H[1] + 35.0 + sgn*T[2]*(32.0*H[1]+8.0))/(1024.0*T[3]);
  641. C[3] = (-4725.0*H[4] - 6300.0*H[3] - 3150.0*H[2] + 3780.0*H[1] + 10395.0
  642. -1216.0*T[4] + sgn*T[2]*(6000.0*H[2]+5760.0*H[1]+1680.0)) / (32768.0*T[4]);
  643. C[4] = 7.0*(-10395.0*H[5] - 23625.0*H[4] - 28350.0*H[3] - 14850.0*H[2]
  644. +19305.0*H[1] + 57915.0 - T[4]*(6336.0*H[1]+6080.0)
  645. + sgn*T[2]*(16800.0*H[3] + 30000.0*H[2] + 25920.0*H[1] + 7920.0)
  646. ) / (262144.0*T[5]);
  647. C[5] = (-2837835.0*H[6] - 9168390.0*H[5] - 16372125.0*H[4] - 18918900*H[3]
  648. -10135125.0*H[2] + 13783770.0*H[1] + 43648605.0
  649. -T[4]*(3044160.0*H[2] + 5588352.0*H[1] + 4213440.0)
  650. +sgn*T[2]*(5556600.0*H[4] + 14817600.0*H[3] + 20790000.0*H[2]
  651. + 17297280.0*H[1] + 5405400.0 + 323072.0*T[4]
  652. )
  653. ) / (4194304.0*T[6]);
  654. C[6] = 0.0;
  655. *V0 = C[0] + (-4.0*C[3]/T[1]+C[4])/V[4]
  656. + (-192.0*C[5]/T[3]+144.0*C[6]/T[2])/V[8]
  657. + sgn * (-C[2]/V[2]
  658. + (-24.0*C[4]/T[2]+12.0*C[5]/T[1]-C[6])/V[6]
  659. );
  660. *V1 = C[1]/V[1] + (8.0*(C[3]/T[2]-C[4]/T[1])+C[5])/V[5]
  661. + (384.0*C[5]/T[4] - 768.0*C[6]/T[3])/V[9]
  662. + sgn * (Cm1*V[1] + (2.0*C[2]/T[1]-C[3])/V[3]
  663. + (48.0*C[4]/T[3]-72.0*C[5]/T[2] + 18.0*C[6]/T[1])/V[7]
  664. );
  665. return GSL_SUCCESS;
  666. }
  667. /*-*-*-*-*-*-*-*-*-*-*-* Functions with Error Codes *-*-*-*-*-*-*-*-*-*-*-*/
  668. /* P^0_{-1/2 + I lambda}
  669. */
  670. int
  671. gsl_sf_conicalP_0_e(const double lambda, const double x, gsl_sf_result * result)
  672. {
  673. /* CHECK_POINTER(result) */
  674. if(x <= -1.0) {
  675. DOMAIN_ERROR(result);
  676. }
  677. else if(x == 1.0) {
  678. result->val = 1.0;
  679. result->err = 0.0;
  680. return GSL_SUCCESS;
  681. }
  682. else if(lambda == 0.0) {
  683. gsl_sf_result K;
  684. int stat_K;
  685. if(x < 1.0) {
  686. const double th = acos(x);
  687. const double s = sin(0.5*th);
  688. stat_K = gsl_sf_ellint_Kcomp_e(s, GSL_MODE_DEFAULT, &K);
  689. result->val = 2.0/M_PI * K.val;
  690. result->err = 2.0/M_PI * K.err;
  691. result->err += 2.0 * GSL_DBL_EPSILON * fabs(result->val);
  692. return stat_K;
  693. }
  694. else {
  695. const double xi = acosh(x);
  696. const double c = cosh(0.5*xi);
  697. const double t = tanh(0.5*xi);
  698. stat_K = gsl_sf_ellint_Kcomp_e(t, GSL_MODE_DEFAULT, &K);
  699. result->val = 2.0/M_PI / c * K.val;
  700. result->err = 2.0/M_PI / c * K.err;
  701. result->err += 2.0 * GSL_DBL_EPSILON * fabs(result->val);
  702. return stat_K;
  703. }
  704. }
  705. else if( (x <= 0.0 && lambda < 1000.0)
  706. || (x < 0.1 && lambda < 17.0)
  707. || (x < 0.2 && lambda < 5.0 )
  708. ) {
  709. return conicalP_xlt1_hyperg_A(0.0, lambda, x, result);
  710. }
  711. else if( (x <= 0.2 && lambda < 17.0)
  712. || (x <= 1.5 && lambda < 20.0)
  713. ) {
  714. return gsl_sf_hyperg_2F1_conj_e(0.5, lambda, 1.0, (1.0-x)/2, result);
  715. }
  716. else if(1.5 < x && lambda < GSL_MAX(x,20.0)) {
  717. gsl_sf_result P;
  718. double lm;
  719. int stat_P = gsl_sf_conicalP_large_x_e(0.0, lambda, x,
  720. &P, &lm
  721. );
  722. int stat_e = gsl_sf_exp_mult_err_e(lm, 2.0*GSL_DBL_EPSILON * fabs(lm),
  723. P.val, P.err,
  724. result);
  725. return GSL_ERROR_SELECT_2(stat_e, stat_P);
  726. }
  727. else {
  728. double V0, V1;
  729. if(x < 1.0) {
  730. double th = acos(x);
  731. double sth = sqrt(1.0-x*x); /* sin(th) */
  732. gsl_sf_result I0, I1;
  733. int stat_I0 = gsl_sf_bessel_I0_scaled_e(th * lambda, &I0);
  734. int stat_I1 = gsl_sf_bessel_I1_scaled_e(th * lambda, &I1);
  735. int stat_I = GSL_ERROR_SELECT_2(stat_I0, stat_I1);
  736. int stat_V = conicalP_0_V(th, x/sth, lambda, -1.0, &V0, &V1);
  737. double bessterm = V0 * I0.val + V1 * I1.val;
  738. double besserr = fabs(V0) * I0.err + fabs(V1) * I1.err;
  739. double arg1 = th*lambda;
  740. double sqts = sqrt(th/sth);
  741. int stat_e = gsl_sf_exp_mult_err_e(arg1, 4.0 * GSL_DBL_EPSILON * fabs(arg1),
  742. sqts * bessterm, sqts * besserr,
  743. result);
  744. return GSL_ERROR_SELECT_3(stat_e, stat_V, stat_I);
  745. }
  746. else {
  747. double sh = sqrt(x-1.0)*sqrt(x+1.0); /* sinh(xi) */
  748. double xi = log(x + sh); /* xi = acosh(x) */
  749. gsl_sf_result J0, J1;
  750. int stat_J0 = gsl_sf_bessel_J0_e(xi * lambda, &J0);
  751. int stat_J1 = gsl_sf_bessel_J1_e(xi * lambda, &J1);
  752. int stat_J = GSL_ERROR_SELECT_2(stat_J0, stat_J1);
  753. int stat_V = conicalP_0_V(xi, x/sh, lambda, 1.0, &V0, &V1);
  754. double bessterm = V0 * J0.val + V1 * J1.val;
  755. double besserr = fabs(V0) * J0.err + fabs(V1) * J1.err;
  756. double pre_val = sqrt(xi/sh);
  757. double pre_err = 2.0 * fabs(pre_val);
  758. result->val = pre_val * bessterm;
  759. result->err = pre_val * besserr;
  760. result->err += pre_err * fabs(bessterm);
  761. result->err += 2.0 * GSL_DBL_EPSILON * fabs(result->val);
  762. return GSL_ERROR_SELECT_2(stat_V, stat_J);
  763. }
  764. }
  765. }
  766. /* P^1_{-1/2 + I lambda}
  767. */
  768. int
  769. gsl_sf_conicalP_1_e(const double lambda, const double x, gsl_sf_result * result)
  770. {
  771. /* CHECK_POINTER(result) */
  772. if(x <= -1.0) {
  773. DOMAIN_ERROR(result);
  774. }
  775. else if(lambda == 0.0) {
  776. gsl_sf_result K, E;
  777. int stat_K, stat_E;
  778. if(x == 1.0) {
  779. result->val = 0.0;
  780. result->err = 0.0;
  781. return GSL_SUCCESS;
  782. }
  783. else if(x < 1.0) {
  784. if(1.0-x < GSL_SQRT_DBL_EPSILON) {
  785. double err_amp = GSL_MAX_DBL(1.0, 1.0/(GSL_DBL_EPSILON + fabs(1.0-x)));
  786. result->val = 0.25/M_SQRT2 * sqrt(1.0-x) * (1.0 + 5.0/16.0 * (1.0-x));
  787. result->err = err_amp * 3.0 * GSL_DBL_EPSILON * fabs(result->val);
  788. return GSL_SUCCESS;
  789. }
  790. else {
  791. const double th = acos(x);
  792. const double s = sin(0.5*th);
  793. const double c2 = 1.0 - s*s;
  794. const double sth = sin(th);
  795. const double pre = 2.0/(M_PI*sth);
  796. stat_K = gsl_sf_ellint_Kcomp_e(s, GSL_MODE_DEFAULT, &K);
  797. stat_E = gsl_sf_ellint_Ecomp_e(s, GSL_MODE_DEFAULT, &E);
  798. result->val = pre * (E.val - c2 * K.val);
  799. result->err = pre * (E.err + fabs(c2) * K.err);
  800. result->err += 2.0 * GSL_DBL_EPSILON * fabs(result->val);
  801. return stat_K;
  802. }
  803. }
  804. else {
  805. if(x-1.0 < GSL_SQRT_DBL_EPSILON) {
  806. double err_amp = GSL_MAX_DBL(1.0, 1.0/(GSL_DBL_EPSILON + fabs(1.0-x)));
  807. result->val = -0.25/M_SQRT2 * sqrt(x-1.0) * (1.0 - 5.0/16.0 * (x-1.0));
  808. result->err = err_amp * 3.0 * GSL_DBL_EPSILON * fabs(result->val);
  809. return GSL_SUCCESS;
  810. }
  811. else {
  812. const double xi = acosh(x);
  813. const double c = cosh(0.5*xi);
  814. const double t = tanh(0.5*xi);
  815. const double sxi = sinh(xi);
  816. const double pre = 2.0/(M_PI*sxi) * c;
  817. stat_K = gsl_sf_ellint_Kcomp_e(t, GSL_MODE_DEFAULT, &K);
  818. stat_E = gsl_sf_ellint_Ecomp_e(t, GSL_MODE_DEFAULT, &E);
  819. result->val = pre * (E.val - K.val);
  820. result->err = pre * (E.err + K.err);
  821. result->err += 2.0 * GSL_DBL_EPSILON * fabs(result->val);
  822. return stat_K;
  823. }
  824. }
  825. }
  826. else if( (x <= 0.0 && lambda < 1000.0)
  827. || (x < 0.1 && lambda < 17.0)
  828. || (x < 0.2 && lambda < 5.0 )
  829. ) {
  830. return conicalP_xlt1_hyperg_A(1.0, lambda, x, result);
  831. }
  832. else if( (x <= 0.2 && lambda < 17.0)
  833. || (x < 1.5 && lambda < 20.0)
  834. ) {
  835. const double arg = fabs(x*x - 1.0);
  836. const double sgn = GSL_SIGN(1.0 - x);
  837. const double pre = 0.5*(lambda*lambda + 0.25) * sgn * sqrt(arg);
  838. gsl_sf_result F;
  839. int stat_F = gsl_sf_hyperg_2F1_conj_e(1.5, lambda, 2.0, (1.0-x)/2, &F);
  840. result->val = pre * F.val;
  841. result->err = fabs(pre) * F.err;
  842. result->err += 2.0 * GSL_DBL_EPSILON * fabs(result->val);
  843. return stat_F;
  844. }
  845. else if(1.5 <= x && lambda < GSL_MAX(x,20.0)) {
  846. gsl_sf_result P;
  847. double lm;
  848. int stat_P = gsl_sf_conicalP_large_x_e(1.0, lambda, x,
  849. &P, &lm
  850. );
  851. int stat_e = gsl_sf_exp_mult_err_e(lm, 2.0 * GSL_DBL_EPSILON * fabs(lm),
  852. P.val, P.err,
  853. result);
  854. return GSL_ERROR_SELECT_2(stat_e, stat_P);
  855. }
  856. else {
  857. double V0, V1;
  858. if(x < 1.0) {
  859. const double sqrt_1mx = sqrt(1.0 - x);
  860. const double sqrt_1px = sqrt(1.0 + x);
  861. const double th = acos(x);
  862. const double sth = sqrt_1mx * sqrt_1px; /* sin(th) */
  863. gsl_sf_result I0, I1;
  864. int stat_I0 = gsl_sf_bessel_I0_scaled_e(th * lambda, &I0);
  865. int stat_I1 = gsl_sf_bessel_I1_scaled_e(th * lambda, &I1);
  866. int stat_I = GSL_ERROR_SELECT_2(stat_I0, stat_I1);
  867. int stat_V = conicalP_1_V(th, x/sth, lambda, -1.0, &V0, &V1);
  868. double bessterm = V0 * I0.val + V1 * I1.val;
  869. double besserr = fabs(V0) * I0.err + fabs(V1) * I1.err
  870. + 2.0 * GSL_DBL_EPSILON * fabs(V0 * I0.val)
  871. + 2.0 * GSL_DBL_EPSILON * fabs(V1 * I1.val);
  872. double arg1 = th * lambda;
  873. double sqts = sqrt(th/sth);
  874. int stat_e = gsl_sf_exp_mult_err_e(arg1, 2.0 * GSL_DBL_EPSILON * fabs(arg1),
  875. sqts * bessterm, sqts * besserr,
  876. result);
  877. result->err *= 1.0/sqrt_1mx;
  878. return GSL_ERROR_SELECT_3(stat_e, stat_V, stat_I);
  879. }
  880. else {
  881. const double sqrt_xm1 = sqrt(x - 1.0);
  882. const double sqrt_xp1 = sqrt(x + 1.0);
  883. const double sh = sqrt_xm1 * sqrt_xp1; /* sinh(xi) */
  884. const double xi = log(x + sh); /* xi = acosh(x) */
  885. const double xi_lam = xi * lambda;
  886. gsl_sf_result J0, J1;
  887. const int stat_J0 = gsl_sf_bessel_J0_e(xi_lam, &J0);
  888. const int stat_J1 = gsl_sf_bessel_J1_e(xi_lam, &J1);
  889. const int stat_J = GSL_ERROR_SELECT_2(stat_J0, stat_J1);
  890. const int stat_V = conicalP_1_V(xi, x/sh, lambda, 1.0, &V0, &V1);
  891. const double bessterm = V0 * J0.val + V1 * J1.val;
  892. const double besserr = fabs(V0) * J0.err + fabs(V1) * J1.err
  893. + 512.0 * 2.0 * GSL_DBL_EPSILON * fabs(V0 * J0.val)
  894. + 512.0 * 2.0 * GSL_DBL_EPSILON * fabs(V1 * J1.val)
  895. + GSL_DBL_EPSILON * fabs(xi_lam * V0 * J1.val)
  896. + GSL_DBL_EPSILON * fabs(xi_lam * V1 * J0.val);
  897. const double pre = sqrt(xi/sh);
  898. result->val = pre * bessterm;
  899. result->err = pre * besserr * sqrt_xp1 / sqrt_xm1;
  900. result->err += 4.0 * GSL_DBL_EPSILON * fabs(result->val);
  901. return GSL_ERROR_SELECT_2(stat_V, stat_J);
  902. }
  903. }
  904. }
  905. /* P^{1/2}_{-1/2 + I lambda} (x)
  906. * [Abramowitz+Stegun 8.6.8, 8.6.12]
  907. * checked OK [GJ] Fri May 8 12:24:36 MDT 1998
  908. */
  909. int gsl_sf_conicalP_half_e(const double lambda, const double x,
  910. gsl_sf_result * result
  911. )
  912. {
  913. /* CHECK_POINTER(result) */
  914. if(x <= -1.0) {
  915. DOMAIN_ERROR(result);
  916. }
  917. else if(x < 1.0) {
  918. double err_amp = 1.0 + 1.0/(GSL_DBL_EPSILON + fabs(1.0-fabs(x)));
  919. double ac = acos(x);
  920. double den = sqrt(sqrt(1.0-x)*sqrt(1.0+x));
  921. result->val = Root_2OverPi_ / den * cosh(ac * lambda);
  922. result->err = err_amp * 3.0 * GSL_DBL_EPSILON * fabs(result->val);
  923. result->err *= fabs(ac * lambda) + 1.0;
  924. return GSL_SUCCESS;
  925. }
  926. else if(x == 1.0) {
  927. result->val = 0.0;
  928. result->err = 0.0;
  929. return GSL_SUCCESS;
  930. }
  931. else {
  932. /* x > 1 */
  933. double err_amp = 1.0 + 1.0/(GSL_DBL_EPSILON + fabs(1.0-fabs(x)));
  934. double sq_term = sqrt(x-1.0)*sqrt(x+1.0);
  935. double ln_term = log(x + sq_term);
  936. double den = sqrt(sq_term);
  937. double carg_val = lambda * ln_term;
  938. double carg_err = 2.0 * GSL_DBL_EPSILON * fabs(carg_val);
  939. gsl_sf_result cos_result;
  940. int stat_cos = gsl_sf_cos_err_e(carg_val, carg_err, &cos_result);
  941. result->val = Root_2OverPi_ / den * cos_result.val;
  942. result->err = err_amp * Root_2OverPi_ / den * cos_result.err;
  943. result->err += 4.0 * GSL_DBL_EPSILON * fabs(result->val);
  944. return stat_cos;
  945. }
  946. }
  947. /* P^{-1/2}_{-1/2 + I lambda} (x)
  948. * [Abramowitz+Stegun 8.6.9, 8.6.14]
  949. * checked OK [GJ] Fri May 8 12:24:43 MDT 1998
  950. */
  951. int gsl_sf_conicalP_mhalf_e(const double lambda, const double x, gsl_sf_result * result)
  952. {
  953. /* CHECK_POINTER(result) */
  954. if(x <= -1.0) {
  955. DOMAIN_ERROR(result);
  956. }
  957. else if(x < 1.0) {
  958. double ac = acos(x);
  959. double den = sqrt(sqrt(1.0-x)*sqrt(1.0+x));
  960. double arg = ac * lambda;
  961. double err_amp = 1.0 + 1.0/(GSL_DBL_EPSILON + fabs(1.0-fabs(x)));
  962. if(fabs(arg) < GSL_SQRT_DBL_EPSILON) {
  963. result->val = Root_2OverPi_ / den * ac;
  964. result->err = 2.0 * GSL_DBL_EPSILON * fabs(result->val);
  965. result->err *= err_amp;
  966. }
  967. else {
  968. result->val = Root_2OverPi_ / (den*lambda) * sinh(arg);
  969. result->err = GSL_DBL_EPSILON * (fabs(arg)+1.0) * fabs(result->val);
  970. result->err *= err_amp;
  971. result->err += 2.0 * GSL_DBL_EPSILON * fabs(result->val);
  972. }
  973. return GSL_SUCCESS;
  974. }
  975. else if(x == 1.0) {
  976. result->val = 0.0;
  977. result->err = 0.0;
  978. return GSL_SUCCESS;
  979. }
  980. else {
  981. /* x > 1 */
  982. double sq_term = sqrt(x-1.0)*sqrt(x+1.0);
  983. double ln_term = log(x + sq_term);
  984. double den = sqrt(sq_term);
  985. double arg_val = lambda * ln_term;
  986. double arg_err = 2.0 * GSL_DBL_EPSILON * fabs(arg_val);
  987. if(arg_val < GSL_SQRT_DBL_EPSILON) {
  988. result->val = Root_2OverPi_ / den * ln_term;
  989. result->err = 2.0 * GSL_DBL_EPSILON * fabs(result->val);
  990. return GSL_SUCCESS;
  991. }
  992. else {
  993. gsl_sf_result sin_result;
  994. int stat_sin = gsl_sf_sin_err_e(arg_val, arg_err, &sin_result);
  995. result->val = Root_2OverPi_ / (den*lambda) * sin_result.val;
  996. result->err = Root_2OverPi_ / fabs(den*lambda) * sin_result.err;
  997. result->err += 3.0 * GSL_DBL_EPSILON * fabs(result->val);
  998. return stat_sin;
  999. }
  1000. }
  1001. }
  1002. int gsl_sf_conicalP_sph_reg_e(const int l, const double lambda,
  1003. const double x,
  1004. gsl_sf_result * result
  1005. )
  1006. {
  1007. /* CHECK_POINTER(result) */
  1008. if(x <= -1.0 || l < -1) {
  1009. DOMAIN_ERROR(result);
  1010. }
  1011. else if(l == -1) {
  1012. return gsl_sf_conicalP_half_e(lambda, x, result);
  1013. }
  1014. else if(l == 0) {
  1015. return gsl_sf_conicalP_mhalf_e(lambda, x, result);
  1016. }
  1017. else if(x == 1.0) {
  1018. result->val = 0.0;
  1019. result->err = 0.0;
  1020. return GSL_SUCCESS;
  1021. }
  1022. else if(x < 0.0) {
  1023. double c = 1.0/sqrt(1.0-x*x);
  1024. gsl_sf_result r_Pellm1;
  1025. gsl_sf_result r_Pell;
  1026. int stat_0 = gsl_sf_conicalP_half_e(lambda, x, &r_Pellm1); /* P^( 1/2) */
  1027. int stat_1 = gsl_sf_conicalP_mhalf_e(lambda, x, &r_Pell); /* P^(-1/2) */
  1028. int stat_P = GSL_ERROR_SELECT_2(stat_0, stat_1);
  1029. double Pellm1 = r_Pellm1.val;
  1030. double Pell = r_Pell.val;
  1031. double Pellp1;
  1032. int ell;
  1033. for(ell=0; ell<l; ell++) {
  1034. double d = (ell+1.0)*(ell+1.0) + lambda*lambda;
  1035. Pellp1 = (Pellm1 - (2.0*ell+1.0)*c*x * Pell) / d;
  1036. Pellm1 = Pell;
  1037. Pell = Pellp1;
  1038. }
  1039. result->val = Pell;
  1040. result->err = (0.5*l + 1.0) * GSL_DBL_EPSILON * fabs(Pell);
  1041. result->err += GSL_DBL_EPSILON * l * fabs(result->val);
  1042. return stat_P;
  1043. }
  1044. else if(x < 1.0) {
  1045. const double xi = x/(sqrt(1.0-x)*sqrt(1.0+x));
  1046. gsl_sf_result rat;
  1047. gsl_sf_result Phf;
  1048. int stat_CF1 = conicalP_negmu_xlt1_CF1(0.5, l, lambda, x, &rat);
  1049. int stat_Phf = gsl_sf_conicalP_half_e(lambda, x, &Phf);
  1050. double Pellp1 = rat.val * GSL_SQRT_DBL_MIN;
  1051. double Pell = GSL_SQRT_DBL_MIN;
  1052. double Pellm1;
  1053. int ell;
  1054. for(ell=l; ell>=0; ell--) {
  1055. double d = (ell+1.0)*(ell+1.0) + lambda*lambda;
  1056. Pellm1 = (2.0*ell+1.0)*xi * Pell + d * Pellp1;
  1057. Pellp1 = Pell;
  1058. Pell = Pellm1;
  1059. }
  1060. result->val = GSL_SQRT_DBL_MIN * Phf.val / Pell;
  1061. result->err = GSL_SQRT_DBL_MIN * Phf.err / fabs(Pell);
  1062. result->err += fabs(rat.err/rat.val) * (l + 1.0) * fabs(result->val);
  1063. result->err += 2.0 * GSL_DBL_EPSILON * fabs(result->val);
  1064. return GSL_ERROR_SELECT_2(stat_Phf, stat_CF1);
  1065. }
  1066. else if(x == 1.0) {
  1067. result->val = 0.0;
  1068. result->err = 0.0;
  1069. return GSL_SUCCESS;
  1070. }
  1071. else {
  1072. /* x > 1.0 */
  1073. const double xi = x/sqrt((x-1.0)*(x+1.0));
  1074. gsl_sf_result rat;
  1075. int stat_CF1 = conicalP_negmu_xgt1_CF1(0.5, l, lambda, x, &rat);
  1076. int stat_P;
  1077. double Pellp1 = rat.val * GSL_SQRT_DBL_MIN;
  1078. double Pell = GSL_SQRT_DBL_MIN;
  1079. double Pellm1;
  1080. int ell;
  1081. for(ell=l; ell>=0; ell--) {
  1082. double d = (ell+1.0)*(ell+1.0) + lambda*lambda;
  1083. Pellm1 = (2.0*ell+1.0)*xi * Pell - d * Pellp1;
  1084. Pellp1 = Pell;
  1085. Pell = Pellm1;
  1086. }
  1087. if(fabs(Pell) > fabs(Pellp1)){
  1088. gsl_sf_result Phf;
  1089. stat_P = gsl_sf_conicalP_half_e(lambda, x, &Phf);
  1090. result->val = GSL_SQRT_DBL_MIN * Phf.val / Pell;
  1091. result->err = 2.0 * GSL_SQRT_DBL_MIN * Phf.err / fabs(Pell);
  1092. result->err += 2.0 * fabs(rat.err/rat.val) * (l + 1.0) * fabs(result->val);
  1093. result->err += 2.0 * GSL_DBL_EPSILON * fabs(result->val);
  1094. }
  1095. else {
  1096. gsl_sf_result Pmhf;
  1097. stat_P = gsl_sf_conicalP_mhalf_e(lambda, x, &Pmhf);
  1098. result->val = GSL_SQRT_DBL_MIN * Pmhf.val / Pellp1;
  1099. result->err = 2.0 * GSL_SQRT_DBL_MIN * Pmhf.err / fabs(Pellp1);
  1100. result->err += 2.0 * fabs(rat.err/rat.val) * (l + 1.0) * fabs(result->val);
  1101. result->err += 2.0 * GSL_DBL_EPSILON * fabs(result->val);
  1102. }
  1103. return GSL_ERROR_SELECT_2(stat_P, stat_CF1);
  1104. }
  1105. }
  1106. int gsl_sf_conicalP_cyl_reg_e(const int m, const double lambda,
  1107. const double x,
  1108. gsl_sf_result * result
  1109. )
  1110. {
  1111. /* CHECK_POINTER(result) */
  1112. if(x <= -1.0 || m < -1) {
  1113. DOMAIN_ERROR(result);
  1114. }
  1115. else if(m == -1) {
  1116. return gsl_sf_conicalP_1_e(lambda, x, result);
  1117. }
  1118. else if(m == 0) {
  1119. return gsl_sf_conicalP_0_e(lambda, x, result);
  1120. }
  1121. else if(x == 1.0) {
  1122. result->val = 0.0;
  1123. result->err = 0.0;
  1124. return GSL_SUCCESS;
  1125. }
  1126. else if(x < 0.0) {
  1127. double c = 1.0/sqrt(1.0-x*x);
  1128. gsl_sf_result r_Pkm1;
  1129. gsl_sf_result r_Pk;
  1130. int stat_0 = gsl_sf_conicalP_1_e(lambda, x, &r_Pkm1); /* P^1 */
  1131. int stat_1 = gsl_sf_conicalP_0_e(lambda, x, &r_Pk); /* P^0 */
  1132. int stat_P = GSL_ERROR_SELECT_2(stat_0, stat_1);
  1133. double Pkm1 = r_Pkm1.val;
  1134. double Pk = r_Pk.val;
  1135. double Pkp1;
  1136. int k;
  1137. for(k=0; k<m; k++) {
  1138. double d = (k+0.5)*(k+0.5) + lambda*lambda;
  1139. Pkp1 = (Pkm1 - 2.0*k*c*x * Pk) / d;
  1140. Pkm1 = Pk;
  1141. Pk = Pkp1;
  1142. }
  1143. result->val = Pk;
  1144. result->err = (m + 2.0) * GSL_DBL_EPSILON * fabs(Pk);
  1145. result->err += 2.0 * GSL_DBL_EPSILON * fabs(result->val);
  1146. return stat_P;
  1147. }
  1148. else if(x < 1.0) {
  1149. const double xi = x/(sqrt(1.0-x)*sqrt(1.0+x));
  1150. gsl_sf_result rat;
  1151. gsl_sf_result P0;
  1152. int stat_CF1 = conicalP_negmu_xlt1_CF1(0.0, m, lambda, x, &rat);
  1153. int stat_P0 = gsl_sf_conicalP_0_e(lambda, x, &P0);
  1154. double Pkp1 = rat.val * GSL_SQRT_DBL_MIN;
  1155. double Pk = GSL_SQRT_DBL_MIN;
  1156. double Pkm1;
  1157. int k;
  1158. for(k=m; k>0; k--) {
  1159. double d = (k+0.5)*(k+0.5) + lambda*lambda;
  1160. Pkm1 = 2.0*k*xi * Pk + d * Pkp1;
  1161. Pkp1 = Pk;
  1162. Pk = Pkm1;
  1163. }
  1164. result->val = GSL_SQRT_DBL_MIN * P0.val / Pk;
  1165. result->err = 2.0 * GSL_SQRT_DBL_MIN * P0.err / fabs(Pk);
  1166. result->err += 2.0 * fabs(rat.err/rat.val) * (m + 1.0) * fabs(result->val);
  1167. result->err += 2.0 * GSL_DBL_EPSILON * fabs(result->val);
  1168. return GSL_ERROR_SELECT_2(stat_P0, stat_CF1);
  1169. }
  1170. else if(x == 1.0) {
  1171. result->val = 0.0;
  1172. result->err = 0.0;
  1173. return GSL_SUCCESS;
  1174. }
  1175. else {
  1176. /* x > 1.0 */
  1177. const double xi = x/sqrt((x-1.0)*(x+1.0));
  1178. gsl_sf_result rat;
  1179. int stat_CF1 = conicalP_negmu_xgt1_CF1(0.0, m, lambda, x, &rat);
  1180. int stat_P;
  1181. double Pkp1 = rat.val * GSL_SQRT_DBL_MIN;
  1182. double Pk = GSL_SQRT_DBL_MIN;
  1183. double Pkm1;
  1184. int k;
  1185. for(k=m; k>-1; k--) {
  1186. double d = (k+0.5)*(k+0.5) + lambda*lambda;
  1187. Pkm1 = 2.0*k*xi * Pk - d * Pkp1;
  1188. Pkp1 = Pk;
  1189. Pk = Pkm1;
  1190. }
  1191. if(fabs(Pk) > fabs(Pkp1)){
  1192. gsl_sf_result P1;
  1193. stat_P = gsl_sf_conicalP_1_e(lambda, x, &P1);
  1194. result->val = GSL_SQRT_DBL_MIN * P1.val / Pk;
  1195. result->err = 2.0 * GSL_SQRT_DBL_MIN * P1.err / fabs(Pk);
  1196. result->err += 2.0 * fabs(rat.err/rat.val) * (m+2.0) * fabs(result->val);
  1197. result->err += 2.0 * GSL_DBL_EPSILON * fabs(result->val);
  1198. }
  1199. else {
  1200. gsl_sf_result P0;
  1201. stat_P = gsl_sf_conicalP_0_e(lambda, x, &P0);
  1202. result->val = GSL_SQRT_DBL_MIN * P0.val / Pkp1;
  1203. result->err = 2.0 * GSL_SQRT_DBL_MIN * P0.err / fabs(Pkp1);
  1204. result->err += 2.0 * fabs(rat.err/rat.val) * (m+2.0) * fabs(result->val);
  1205. result->err += 2.0 * GSL_DBL_EPSILON * fabs(result->val);
  1206. }
  1207. return GSL_ERROR_SELECT_2(stat_P, stat_CF1);
  1208. }
  1209. }
  1210. /*-*-*-*-*-*-*-*-*-* Functions w/ Natural Prototypes *-*-*-*-*-*-*-*-*-*-*/
  1211. #include "gsl_specfunc__eval.h"
  1212. double gsl_sf_conicalP_0(const double lambda, const double x)
  1213. {
  1214. EVAL_RESULT(gsl_sf_conicalP_0_e(lambda, x, &result));
  1215. }
  1216. double gsl_sf_conicalP_1(const double lambda, const double x)
  1217. {
  1218. EVAL_RESULT(gsl_sf_conicalP_1_e(lambda, x, &result));
  1219. }
  1220. double gsl_sf_conicalP_half(const double lambda, const double x)
  1221. {
  1222. EVAL_RESULT(gsl_sf_conicalP_half_e(lambda, x, &result));
  1223. }
  1224. double gsl_sf_conicalP_mhalf(const double lambda, const double x)
  1225. {
  1226. EVAL_RESULT(gsl_sf_conicalP_mhalf_e(lambda, x, &result));
  1227. }
  1228. double gsl_sf_conicalP_sph_reg(const int l, const double lambda, const double x)
  1229. {
  1230. EVAL_RESULT(gsl_sf_conicalP_sph_reg_e(l, lambda, x, &result));
  1231. }
  1232. double gsl_sf_conicalP_cyl_reg(const int m, const double lambda, const double x)
  1233. {
  1234. EVAL_RESULT(gsl_sf_conicalP_cyl_reg_e(m, lambda, x, &result));
  1235. }