8 gsl_sf_airy_Ai_scaled_e
10 gsl_sf_airy_Bi_scaled_e
12 gsl_sf_airy_Ai_deriv_e
14 gsl_sf_airy_Bi_deriv_e
16 gsl_sf_airy_Ai_deriv_scaled_e
17 gsl_sf_airy_Ai_deriv_scaled
18 gsl_sf_airy_Bi_deriv_scaled_e
19 gsl_sf_airy_Bi_deriv_scaled
24 gsl_sf_airy_zero_Ai_deriv_e
25 gsl_sf_airy_zero_Ai_deriv
26 gsl_sf_airy_zero_Bi_deriv_e
27 gsl_sf_airy_zero_Bi_deriv
36 gsl_sf_bessel_Jn_array
43 gsl_sf_bessel_Yn_array
50 gsl_sf_bessel_In_array
51 gsl_sf_bessel_I0_scaled_e
52 gsl_sf_bessel_I0_scaled
53 gsl_sf_bessel_I1_scaled_e
54 gsl_sf_bessel_I1_scaled
55 gsl_sf_bessel_In_scaled_e
56 gsl_sf_bessel_In_scaled
57 gsl_sf_bessel_In_scaled_array
64 gsl_sf_bessel_Kn_array
65 gsl_sf_bessel_K0_scaled_e
66 gsl_sf_bessel_K0_scaled
67 gsl_sf_bessel_K1_scaled_e
68 gsl_sf_bessel_K1_scaled
69 gsl_sf_bessel_Kn_scaled_e
70 gsl_sf_bessel_Kn_scaled
71 gsl_sf_bessel_Kn_scaled_array
80 gsl_sf_bessel_jl_array
81 gsl_sf_bessel_jl_steed_array
90 gsl_sf_bessel_yl_array
91 gsl_sf_bessel_i0_scaled_e
92 gsl_sf_bessel_i0_scaled
93 gsl_sf_bessel_i1_scaled_e
94 gsl_sf_bessel_i1_scaled
95 gsl_sf_bessel_i2_scaled_e
96 gsl_sf_bessel_i2_scaled
97 gsl_sf_bessel_il_scaled_e
98 gsl_sf_bessel_il_scaled
99 gsl_sf_bessel_il_scaled_array
100 gsl_sf_bessel_k0_scaled_e
101 gsl_sf_bessel_k0_scaled
102 gsl_sf_bessel_k1_scaled_e
103 gsl_sf_bessel_k1_scaled
104 gsl_sf_bessel_k2_scaled_e
105 gsl_sf_bessel_k2_scaled
106 gsl_sf_bessel_kl_scaled_e
107 gsl_sf_bessel_kl_scaled
108 gsl_sf_bessel_kl_scaled_array
113 gsl_sf_bessel_sequence_Jnu_e
114 gsl_sf_bessel_Inu_scaled_e
115 gsl_sf_bessel_Inu_scaled
118 gsl_sf_bessel_Knu_scaled_e
119 gsl_sf_bessel_Knu_scaled
122 gsl_sf_bessel_lnKnu_e
124 gsl_sf_bessel_zero_J0_e
125 gsl_sf_bessel_zero_J0
126 gsl_sf_bessel_zero_J1_e
127 gsl_sf_bessel_zero_J1
128 gsl_sf_bessel_zero_Jnu_e
129 gsl_sf_bessel_zero_Jnu
131 @EXPORT_clausen = qw/
135 @EXPORT_hydrogenic = qw/
136 gsl_sf_hydrogenicR_1_e
141 @EXPORT_coulumb = qw/
142 gsl_sf_coulomb_wave_FG_e
143 gsl_sf_coulomb_wave_F_array
144 gsl_sf_coulomb_wave_FG_array
145 gsl_sf_coulomb_wave_FGp_array
146 gsl_sf_coulomb_wave_sphF_array
148 gsl_sf_coulomb_CL_array
150 @EXPORT_coupling = qw/
155 gsl_sf_coupling_RacahW_e
156 gsl_sf_coupling_RacahW
159 gsl_sf_coupling_6j_INCORRECT_e
160 gsl_sf_coupling_6j_INCORRECT
183 gsl_sf_complex_dilog_xy_e
184 gsl_sf_complex_dilog_e
188 gsl_sf_complex_spence_xy_e
191 gsl_sf_multiply_err_e
193 @EXPORT_elliptic = qw/
194 gsl_sf_ellint_Kcomp_e
196 gsl_sf_ellint_Ecomp_e
198 gsl_sf_ellint_Pcomp_e
200 gsl_sf_ellint_Dcomp_e
234 push @EXPORT_misc, qw/
240 gsl_sf_exp_mult_e10_e
251 gsl_sf_exp_mult_err_e
252 gsl_sf_exp_mult_err_e10_e
259 gsl_sf_expint_E1_scaled_e
260 gsl_sf_expint_E1_scaled
261 gsl_sf_expint_E2_scaled_e
262 gsl_sf_expint_E2_scaled
263 gsl_sf_expint_En_scaled_e
264 gsl_sf_expint_En_scaled
267 gsl_sf_expint_Ei_scaled_e
268 gsl_sf_expint_Ei_scaled
280 @EXPORT_fermi_dirac = qw/
281 gsl_sf_fermi_dirac_m1_e
282 gsl_sf_fermi_dirac_m1
283 gsl_sf_fermi_dirac_0_e
285 gsl_sf_fermi_dirac_1_e
287 gsl_sf_fermi_dirac_2_e
289 gsl_sf_fermi_dirac_int_e
290 gsl_sf_fermi_dirac_int
291 gsl_sf_fermi_dirac_mhalf_e
292 gsl_sf_fermi_dirac_mhalf
293 gsl_sf_fermi_dirac_half_e
294 gsl_sf_fermi_dirac_half
295 gsl_sf_fermi_dirac_3half_e
296 gsl_sf_fermi_dirac_3half
297 gsl_sf_fermi_dirac_inc_0_e
298 gsl_sf_fermi_dirac_inc_0
300 @EXPORT_legendre = qw/
303 gsl_sf_legendre_Pl_array
304 gsl_sf_legendre_Pl_deriv_array
317 gsl_sf_legendre_Plm_e
319 gsl_sf_legendre_Plm_array
320 gsl_sf_legendre_Plm_deriv_array
321 gsl_sf_legendre_sphPlm_e
322 gsl_sf_legendre_sphPlm
323 gsl_sf_legendre_sphPlm_array
324 gsl_sf_legendre_sphPlm_deriv_array
325 gsl_sf_legendre_array_size
326 gsl_sf_legendre_H3d_0_e
327 gsl_sf_legendre_H3d_0
328 gsl_sf_legendre_H3d_1_e
329 gsl_sf_legendre_H3d_1
330 gsl_sf_legendre_H3d_e
332 gsl_sf_legendre_H3d_array
344 gsl_sf_lngamma_complex_e
352 @EXPORT_factorial = qw/
359 gsl_sf_lndoublefact_e
362 @EXPORT_hypergeometric = qw/
365 gsl_sf_hyperg_1F1_int_e
366 gsl_sf_hyperg_1F1_int
369 gsl_sf_hyperg_U_int_e
371 gsl_sf_hyperg_U_int_e10_e
374 gsl_sf_hyperg_U_e10_e
377 gsl_sf_hyperg_2F1_conj_e
378 gsl_sf_hyperg_2F1_conj
379 gsl_sf_hyperg_2F1_renorm_e
380 gsl_sf_hyperg_2F1_renorm
381 gsl_sf_hyperg_2F1_conj_renorm_e
382 gsl_sf_hyperg_2F1_conj_renorm
386 @EXPORT_laguerre = qw/
396 push @EXPORT_misc, qw/
425 gsl_sf_gegenpoly_array
430 gsl_sf_conicalP_half_e
432 gsl_sf_conicalP_mhalf_e
433 gsl_sf_conicalP_mhalf
438 gsl_sf_conicalP_sph_reg_e
439 gsl_sf_conicalP_sph_reg
440 gsl_sf_conicalP_cyl_reg_e
441 gsl_sf_conicalP_cyl_reg
449 gsl_sf_log_1plusx_mx_e
466 gsl_sf_result_smash_e
467 gsl_sf_synchrotron_1_e
469 gsl_sf_synchrotron_2_e
472 @EXPORT_mathieu = qw/
473 gsl_sf_mathieu_a_array
474 gsl_sf_mathieu_b_array
477 gsl_sf_mathieu_a_coeff
478 gsl_sf_mathieu_b_coeff
483 gsl_sf_mathieu_ce_array
484 gsl_sf_mathieu_se_array
487 gsl_sf_mathieu_Mc_array
488 gsl_sf_mathieu_Ms_array
490 @EXPORT_transport = qw/
511 gsl_sf_complex_logsin_e
522 gsl_sf_angle_restrict_symm_e
523 gsl_sf_angle_restrict_symm
524 gsl_sf_angle_restrict_pos_e
525 gsl_sf_angle_restrict_pos
526 gsl_sf_angle_restrict_symm_err_e
527 gsl_sf_angle_restrict_pos_err_e
552 GSL_SF_DOUBLEFACT_NMAX
557 @EXPORT_airy, @EXPORT_bessel, @EXPORT_clausen, @EXPORT_hydrogenic,
558 @EXPORT_coulumb, @EXPORT_coupling, @EXPORT_dawson, @EXPORT_debye,
559 @EXPORT_dilog, @EXPORT_misc, @EXPORT_elliptic, @EXPORT_error, @EXPORT_legendre,
560 @EXPORT_gamma, @EXPORT_transport, @EXPORT_trig, @EXPORT_zeta, @EXPORT_eta,
561 @EXPORT_vars, @EXPORT_mathieu, @EXPORT_hypergeometric
565 all => [ @EXPORT_OK ],
566 airy => [ @EXPORT_airy ],
567 bessel => [ @EXPORT_bessel ],
568 clausen => [ @EXPORT_clausen ],
569 coulumb => [ @EXPORT_coulumb ],
570 coupling => [ @EXPORT_coupling ],
571 dawson => [ @EXPORT_dawson ],
572 debye => [ @EXPORT_debye ],
573 dilog => [ @EXPORT_dilog ],
574 eta => [ @EXPORT_eta ],
575 elliptic => [ @EXPORT_elliptic ],
576 error => [ @EXPORT_error ],
577 factorial => [ @EXPORT_factorial ],
578 gamma => [ @EXPORT_gamma ],
579 hydrogenic => [ @EXPORT_hydrogenic ],
580 hypergeometric => [ @EXPORT_hypergeometric ],
581 laguerre => [ @EXPORT_laguerre ],
582 legendre => [ @EXPORT_legendre ],
583 mathieu => [ @EXPORT_mathieu ],
584 misc => [ @EXPORT_misc ],
585 transport => [ @EXPORT_transport ],
586 trig => [ @EXPORT_trig ],
587 vars => [ @EXPORT_vars ],
588 zeta => [ @EXPORT_zeta ],
595 Math::GSL::SF - Special Functions
599 use Math::GSL::SF qw /:all/;
603 This module contains a data structure named gsl_sf_result. To create a new one use
605 $r = Math::GSL::SF::gsl_sf_result_struct->new;
607 You can then access the elements of the structure in this way :
611 my $error = $r->{err};
613 Here is a list of all included functions:
617 =item C<gsl_sf_airy_Ai_e($x, $mode)>
619 =item C<gsl_sf_airy_Ai($x, $mode, $result)>
621 - These routines compute the Airy function Ai($x) with an accuracy specified by $mode. $mode should be $GSL_PREC_DOUBLE, $GSL_PREC_SINGLE or $GSL_PREC_APPROX. $result is a gsl_sf_result structure.
627 =item C<gsl_sf_airy_Bi_e($x, $mode, $result)>
629 =item C<gsl_sf_airy_Bi($x, $mode)>
631 - These routines compute the Airy function Bi($x) with an accuracy specified by $mode. $mode should be $GSL_PREC_DOUBLE, $GSL_PREC_SINGLE or $GSL_PREC_APPROX. $result is a gsl_sf_result structure.
637 =item C<gsl_sf_airy_Ai_scaled_e($x, $mode, $result)>
639 =item C<gsl_sf_airy_Ai_scaled($x, $mode)>
641 - These routines compute a scaled version of the Airy function S_A($x) Ai($x). For $x>0 the scaling factor S_A($x) is \exp(+(2/3) $x**(3/2)), and is 1 for $x<0. $result is a gsl_sf_result structure.
647 =item C<gsl_sf_airy_Bi_scaled_e($x, $mode, $result)>
649 =item C<gsl_sf_airy_Bi_scaled($x, $mode)>
651 - These routines compute a scaled version of the Airy function S_B($x) Bi($x). For $x>0 the scaling factor S_B($x) is exp(-(2/3) $x**(3/2)), and is 1 for $x<0. $result is a gsl_sf_result structure.
657 =item C<gsl_sf_airy_Ai_deriv_e($x, $mode, $result)>
659 =item C<gsl_sf_airy_Ai_deriv($x, $mode)>
661 - These routines compute the Airy function derivative Ai'($x) with an accuracy specified by $mode. $result is a gsl_sf_result structure.
667 =item C<gsl_sf_airy_Bi_deriv_e($x, $mode, $result)>
669 =item C<gsl_sf_airy_Bi_deriv($x, $mode)>
671 -These routines compute the Airy function derivative Bi'($x) with an accuracy specified by $mode. $result is a gsl_sf_result structure.
677 =item C<gsl_sf_airy_Ai_deriv_scaled_e($x, $mode, $result)>
679 =item C<gsl_sf_airy_Ai_deriv_scaled($x, $mode)>
681 -These routines compute the scaled Airy function derivative S_A(x) Ai'(x). For x>0 the scaling factor S_A(x) is \exp(+(2/3) x^(3/2)), and is 1 for x<0. $result is a gsl_sf_result structure.
687 =item C<gsl_sf_airy_Bi_deriv_scaled_e($x, $mode, $result)>
689 =item C<gsl_sf_airy_Bi_deriv_scaled($x, $mode)>
691 -These routines compute the scaled Airy function derivative S_B(x) Bi'(x). For x>0 the scaling factor S_B(x) is exp(-(2/3) x^(3/2)), and is 1 for x<0. $result is a gsl_sf_result structure.
697 =item C<gsl_sf_airy_zero_Ai_e($s, $result)>
699 =item C<gsl_sf_airy_zero_Ai($s)>
701 -These routines compute the location of the s-th zero of the Airy function Ai($x). $result is a gsl_sf_result structure.
707 =item C<gsl_sf_airy_zero_Bi_e($s, $result)>
709 =item C<gsl_sf_airy_zero_Bi($s)>
711 -These routines compute the location of the s-th zero of the Airy function Bi($x). $result is a gsl_sf_result structure.
717 =item C<gsl_sf_airy_zero_Ai_deriv_e($s, $result)>
719 =item C<gsl_sf_airy_zero_Ai_deriv($s)>
721 -These routines compute the location of the s-th zero of the Airy function derivative Ai'(x). $result is a gsl_sf_result structure.
727 =item C<gsl_sf_airy_zero_Bi_deriv_e($s, $result)>
729 =item C<gsl_sf_airy_zero_Bi_deriv($s)>
731 - These routines compute the location of the s-th zero of the Airy function derivative Bi'(x). $result is a gsl_sf_result structure.
737 =item C<gsl_sf_bessel_J0_e($x, $result)>
739 =item C<gsl_sf_bessel_J0($x)>
741 -These routines compute the regular cylindrical Bessel function of zeroth order, J_0($x). $result is a gsl_sf_result structure.
747 =item C<gsl_sf_bessel_J1_e($x, $result)>
749 =item C<gsl_sf_bessel_J1($x)>
751 - These routines compute the regular cylindrical Bessel function of first order, J_1($x). $result is a gsl_sf_result structure.
757 =item C<gsl_sf_bessel_Jn_e($n, $x, $result)>
759 =item C<gsl_sf_bessel_Jn($n, $x)>
761 -These routines compute the regular cylindrical Bessel function of order n, J_n($x).
767 =item C<gsl_sf_bessel_Jn_array($nmin, $nmax, $x, $result_array)> - This routine computes the values of the regular cylindrical Bessel functions J_n($x) for n from $nmin to $nmax inclusive, storing the results in the array $result_array. The values are computed using recurrence relations for efficiency, and therefore may differ slightly from the exact values.
773 =item C<gsl_sf_bessel_Y0_e($x, $result)>
775 =item C<gsl_sf_bessel_Y0($x)>
777 - These routines compute the irregular spherical Bessel function of zeroth order, y_0(x) = -\cos(x)/x.
783 =item C<gsl_sf_bessel_Y1_e($x, $result)>
785 =item C<gsl_sf_bessel_Y1($x)>
787 -These routines compute the irregular spherical Bessel function of first order, y_1(x) = -(\cos(x)/x + \sin(x))/x.
793 =item C<gsl_sf_bessel_Yn_e>($n, $x, $result)
795 =item C<gsl_sf_bessel_Yn($n, $x)>
797 -These routines compute the irregular cylindrical Bessel function of order $n, Y_n(x), for x>0.
803 =item C<gsl_sf_bessel_Yn_array>
811 =item C<gsl_sf_bessel_I0_e($x, $result)>
813 =item C<gsl_sf_bessel_I0($x)>
815 -These routines compute the regular modified cylindrical Bessel function of zeroth order, I_0(x).
821 =item C<gsl_sf_bessel_I1_e($x, $result)>
823 =item C<gsl_sf_bessel_I1($x)>
825 -These routines compute the regular modified cylindrical Bessel function of first order, I_1(x).
831 =item C<gsl_sf_bessel_In_e($n, $x, $result)>
833 =item C<gsl_sf_bessel_In($n, $x)>
835 -These routines compute the regular modified cylindrical Bessel function of order $n, I_n(x).
841 =item C<gsl_sf_bessel_In_array>
849 =item C<gsl_sf_bessel_I0_scaled_e($x, $result)>
851 =item C<gsl_sf_bessel_I0_scaled($x)>
853 -These routines compute the scaled regular modified cylindrical Bessel function of zeroth order \exp(-|x|) I_0(x).
859 =item C<gsl_sf_bessel_I1_scaled_e($x, $result)>
861 =item C<gsl_sf_bessel_I1_scaled($x)>
863 -These routines compute the scaled regular modified cylindrical Bessel function of first order \exp(-|x|) I_1(x).
869 =item C<gsl_sf_bessel_In_scaled_e($n, $x, $result)>
871 =item C<gsl_sf_bessel_In_scaled($n, $x)>
873 -These routines compute the scaled regular modified cylindrical Bessel function of order $n, \exp(-|x|) I_n(x)
879 =item C<gsl_sf_bessel_In_scaled_array>
887 =item C<gsl_sf_bessel_K0_e($x, $result)>
889 =item C<gsl_sf_bessel_K0($x)>
891 -These routines compute the irregular modified cylindrical Bessel function of zeroth order, K_0(x), for x > 0.
897 =item C<gsl_sf_bessel_K1_e($x, $result)>
899 =item C<gsl_sf_bessel_K1($x)>
901 -These routines compute the irregular modified cylindrical Bessel function of first order, K_1(x), for x > 0.
907 =item C<gsl_sf_bessel_Kn_e($n, $x, $result)>
909 =item C<gsl_sf_bessel_Kn($n, $x)>
911 -These routines compute the irregular modified cylindrical Bessel function of order $n, K_n(x), for x > 0.
917 =item C<gsl_sf_bessel_Kn_array>
925 =item C<gsl_sf_bessel_K0_scaled_e($x, $result)>
927 =item C<gsl_sf_bessel_K0_scaled($x)>
929 -These routines compute the scaled irregular modified cylindrical Bessel function of zeroth order \exp(x) K_0(x) for x>0.
935 =item C<gsl_sf_bessel_K1_scaled_e($x, $result)>
937 =item C<gsl_sf_bessel_K1_scaled($x)>
945 =item C<gsl_sf_bessel_Kn_scaled_e($n, $x, $result)>
947 =item C<gsl_sf_bessel_Kn_scaled($n, $x)>
955 =item C<gsl_sf_bessel_Kn_scaled_array >
963 =item C<gsl_sf_bessel_j0_e($x, $result)>
965 =item C<gsl_sf_bessel_j0($x)>
973 =item C<gsl_sf_bessel_j1_e($x, $result)>
975 =item C<gsl_sf_bessel_j1($x)>
983 =item C<gsl_sf_bessel_j2_e($x, $result)>
985 =item C<gsl_sf_bessel_j2($x)>
993 =item C<gsl_sf_bessel_jl_e($l, $x, $result)>
995 =item C<gsl_sf_bessel_jl($l, $x)>
1003 =item C<gsl_sf_bessel_jl_array>
1011 =item C<gsl_sf_bessel_jl_steed_array>
1019 =item C<gsl_sf_bessel_y0_e($x, $result)>
1021 =item C<gsl_sf_bessel_y0($x)>
1029 =item C<gsl_sf_bessel_y1_e($x, $result)>
1031 =item C<gsl_sf_bessel_y1($x)>
1039 =item C<gsl_sf_bessel_y2_e($x, $result)>
1041 =item C<gsl_sf_bessel_y2($x)>
1049 =item C<gsl_sf_bessel_yl_e($l, $x, $result)>
1051 =item C<gsl_sf_bessel_yl($l, $x)>
1059 =item C<gsl_sf_bessel_yl_array>
1067 =item C<gsl_sf_bessel_i0_scaled_e($x, $result)>
1069 =item C<gsl_sf_bessel_i0_scaled($x)>
1077 =item C<gsl_sf_bessel_i1_scaled_e($x, $result)>
1079 =item C<gsl_sf_bessel_i1_scaled($x)>
1087 =item C<gsl_sf_bessel_i2_scaled_e($x, $result)>
1089 =item C<gsl_sf_bessel_i2_scaled($x)>
1097 =item C<gsl_sf_bessel_il_scaled_e($l, $x, $result)>
1099 =item C<gsl_sf_bessel_il_scaled($x)>
1107 =item C<gsl_sf_bessel_il_scaled_array>
1115 =item C<gsl_sf_bessel_k0_scaled_e($x, $result)>
1117 =item C<gsl_sf_bessel_k0_scale($x)>
1125 =item C<gsl_sf_bessel_k1_scaled_e($x, $result)>
1127 =item C<gsl_sf_bessel_k1_scaled($x)>
1135 =item C<gsl_sf_bessel_k2_scaled_e($x, $result) >
1137 =item C<gsl_sf_bessel_k2_scaled($x)>
1145 =item C<gsl_sf_bessel_kl_scaled_e($l, $x, $result)>
1147 =item C<gsl_sf_bessel_kl_scaled($l, $x)>
1155 =item C<gsl_sf_bessel_kl_scaled_array>
1163 =item C<gsl_sf_bessel_Jnu_e($nu, $x, $result)>
1165 =item C<gsl_sf_bessel_Jnu($nu, $x)>
1173 =item C<gsl_sf_bessel_sequence_Jnu_e >
1181 =item C<gsl_sf_bessel_Ynu_e($nu, $x, $result)>
1183 =item C<gsl_sf_bessel_Ynu($nu, $x)>
1191 =item C<gsl_sf_bessel_Inu_scaled_e($nu, $x, $result)>
1193 =item C<gsl_sf_bessel_Inu_scaled($nu, $x)>
1201 =item C<gsl_sf_bessel_Inu_e($nu, $x, $result)>
1203 =item C<gsl_sf_bessel_Inu($nu, $x)>
1211 =item C<gsl_sf_bessel_Knu_scaled_e($nu, $x, $result)>
1213 =item C<gsl_sf_bessel_Knu_scaled($nu, $x)>
1221 =item C<gsl_sf_bessel_Knu_e($nu, $x, $result)>
1223 =item C<gsl_sf_bessel_Knu($nu, $x)>
1231 =item C<gsl_sf_bessel_lnKnu_e($nu, $x, $result)>
1233 =item C<gsl_sf_bessel_lnKnu($nu, $x)>
1241 =item C<gsl_sf_bessel_zero_J0_e($s, $result)>
1243 =item C<gsl_sf_bessel_zero_J0($s)>
1251 =item C<gsl_sf_bessel_zero_J1_e($s, $result)>
1253 =item C<gsl_sf_bessel_zero_J1($s)>
1261 =item C<gsl_sf_bessel_zero_Jnu_e($nu, $s, $result)>
1263 =item C<gsl_sf_bessel_zero_Jnu($nu, $s)>
1271 =item C<gsl_sf_clausen_e($x, $result)>
1273 =item C<gsl_sf_clausen($x)>
1281 =item C<gsl_sf_hydrogenicR_1_e($Z, $r, $result)>
1283 =item C<gsl_sf_hydrogenicR_1($Z, $r)>
1291 =item C<gsl_sf_hydrogenicR_e($n, $l, $Z, $r, $result)>
1293 =item C<gsl_sf_hydrogenicR($n, $l, $Z, $r)>
1301 =item C<gsl_sf_coulomb_wave_FG_e($eta, $x, $L_F, $k, $F, gsl_sf_result * Fp, gsl_sf_result * G, $Gp)> - This function computes the Coulomb wave functions F_L(\eta,x), G_{L-k}(\eta,x) and their derivatives F'_L(\eta,x), G'_{L-k}(\eta,x) with respect to $x. The parameters are restricted to L, L-k > -1/2, x > 0 and integer $k. Note that L itself is not restricted to being an integer. The results are stored in the parameters $F, $G for the function values and $Fp, $Gp for the derivative values. $F, $G, $Fp, $Gp are all gsl_result structs. If an overflow occurs, $GSL_EOVRFLW is returned and scaling exponents are returned as second and third values.
1303 =item C<gsl_sf_coulomb_wave_F_array > -
1305 =item C<gsl_sf_coulomb_wave_FG_array> -
1307 =item C<gsl_sf_coulomb_wave_FGp_array> -
1309 =item C<gsl_sf_coulomb_wave_sphF_array> -
1311 =item C<gsl_sf_coulomb_CL_e($L, $eta, $result)> - This function computes the Coulomb wave function normalization constant C_L($eta) for $L > -1.
1313 =item C<gsl_sf_coulomb_CL_arrayi> -
1319 =item C<gsl_sf_coupling_3j_e($two_ja, $two_jb, $two_jc, $two_ma, $two_mb, $two_mc, $result)>
1321 =item C<gsl_sf_coupling_3j($two_ja, $two_jb, $two_jc, $two_ma, $two_mb, $two_mc)>
1323 - These routines compute the Wigner 3-j coefficient,
1326 where the arguments are given in half-integer units, ja = $two_ja/2, ma = $two_ma/2, etc.
1332 =item C<gsl_sf_coupling_6j_e($two_ja, $two_jb, $two_jc, $two_jd, $two_je, $two_jf, $result)>
1334 =item C<gsl_sf_coupling_6j($two_ja, $two_jb, $two_jc, $two_jd, $two_je, $two_jf)>
1336 - These routines compute the Wigner 6-j coefficient,
1339 where the arguments are given in half-integer units, ja = $two_ja/2, ma = $two_ma/2, etc.
1345 =item C<gsl_sf_coupling_RacahW_e>
1347 =item C<gsl_sf_coupling_RacahW>
1355 =item C<gsl_sf_coupling_9j_e($two_ja, $two_jb, $two_jc, $two_jd, $two_je, $two_jf, $two_jg, $two_jh, $two_ji, $result)>
1357 =item C<gsl_sf_coupling_9j($two_ja, $two_jb, $two_jc, $two_jd, $two_je, $two_jf, $two_jg, $two_jh, $two_ji)>
1359 -These routines compute the Wigner 9-j coefficient,
1364 where the arguments are given in half-integer units, ja = two_ja/2, ma = two_ma/2, etc.
1370 =item C<gsl_sf_dawson_e($x, $result)>
1372 =item C<gsl_sf_dawson($x)>
1374 -These routines compute the value of Dawson's integral for $x.
1380 =item C<gsl_sf_debye_1_e($x, $result)>
1382 =item C<gsl_sf_debye_1($x)>
1384 -These routines compute the first-order Debye function D_1(x) = (1/x) \int_0^x dt (t/(e^t - 1)).
1390 =item C<gsl_sf_debye_2_e($x, $result)>
1392 =item C<gsl_sf_debye_2($x)>
1394 -These routines compute the second-order Debye function D_2(x) = (2/x^2) \int_0^x dt (t^2/(e^t - 1)).
1400 =item C<gsl_sf_debye_3_e($x, $result)>
1402 =item C<gsl_sf_debye_3($x)>
1404 -These routines compute the third-order Debye function D_3(x) = (3/x^3) \int_0^x dt (t^3/(e^t - 1)).
1410 =item C<gsl_sf_debye_4_e($x, $result)>
1412 =item C<gsl_sf_debye_4($x)>
1414 -These routines compute the fourth-order Debye function D_4(x) = (4/x^4) \int_0^x dt (t^4/(e^t - 1)).
1420 =item C<gsl_sf_debye_5_e($x, $result)>
1422 =item C<gsl_sf_debye_5($x)>
1424 -These routines compute the fifth-order Debye function D_5(x) = (5/x^5) \int_0^x dt (t^5/(e^t - 1)).
1430 =item C<gsl_sf_debye_6_e($x, $result)>
1432 =item C<gsl_sf_debye_6($x)>
1434 -These routines compute the sixth-order Debye function D_6(x) = (6/x^6) \int_0^x dt (t^6/(e^t - 1)).
1440 =item C<gsl_sf_dilog_e ($x, $result)>
1442 =item C<gsl_sf_dilog($x)>
1444 - These routines compute the dilogarithm for a real argument. In Lewin's notation this is Li_2(x), the real part of the dilogarithm of a real x. It is defined by the integral representation Li_2(x) = - \Re \int_0^x ds \log(1-s) / s. Note that \Im(Li_2(x)) = 0 for x <= 1, and -\pi\log(x) for x > 1. Note that Abramowitz & Stegun refer to the Spence integral S(x)=Li_2(1-x) as the dilogarithm rather than Li_2(x).
1450 =item C<gsl_sf_complex_dilog_xy_e> -
1452 =item C<gsl_sf_complex_dilog_e($r, $theta, $result_re, $result_im)> - This function computes the full complex-valued dilogarithm for the complex argument z = r \exp(i \theta). The real and imaginary parts of the result are returned in the $result_re and $result_im gsl_result structs.
1454 =item C<gsl_sf_complex_spence_xy_e> -
1460 =item C<gsl_sf_multiply>
1462 =item C<gsl_sf_multiply_e($x, $y, $result)> - This function multiplies $x and $y storing the product and its associated error in $result.
1464 =item C<gsl_sf_multiply_err_e($x, $dx, $y, $dy, $result)> - This function multiplies $x and $y with associated absolute errors $dx and $dy. The product xy +/- xy \sqrt((dx/x)^2 +(dy/y)^2) is stored in $result.
1473 =item C<gsl_sf_ellint_Kcomp_e($k, $mode, $result)>
1475 =item C<gsl_sf_ellint_Kcomp($k, $mode)>
1477 -These routines compute the complete elliptic integral K($k) to the accuracy specified by the mode variable mode. Note that Abramowitz & Stegun define this function in terms of the parameter m = k^2.
1483 =item C<gsl_sf_ellint_Ecomp_e($k, $mode, $result)>
1485 =item C<gsl_sf_ellint_Ecomp($k, $mode)>
1493 =item C<gsl_sf_ellint_Pcomp_e($k, $n, $mode, $result)>
1495 =item C<gsl_sf_ellint_Pcomp($k, $n, $mode)>
1503 =item C<gsl_sf_ellint_Dcomp_e>
1505 =item C<gsl_sf_ellint_Dcomp >
1513 =item C<gsl_sf_ellint_F_e($phi, $k, $mode, $result)>
1515 =item C<gsl_sf_ellint_F($phi, $k, $mode)>
1517 -These routines compute the incomplete elliptic integral F($phi,$k) to the accuracy specified by the mode variable mode. Note that Abramowitz & Stegun define this function in terms of the parameter m = k^2.
1523 =item C<gsl_sf_ellint_E_e($phi, $k, $mode, $result)>
1525 =item C<gsl_sf_ellint_E($phi, $k, $mode)>
1527 -These routines compute the incomplete elliptic integral E($phi,$k) to the accuracy specified by the mode variable mode. Note that Abramowitz & Stegun define this function in terms of the parameter m = k^2.
1533 =item C<gsl_sf_ellint_P_e($phi, $k, $n, $mode, $result)>
1535 =item C<gsl_sf_ellint_P($phi, $k, $n, $mode)>
1537 -These routines compute the incomplete elliptic integral \Pi(\phi,k,n) to the accuracy specified by the mode variable mode. Note that Abramowitz & Stegun define this function in terms of the parameters m = k^2 and \sin^2(\alpha) = k^2, with the change of sign n \to -n.
1543 =item C<gsl_sf_ellint_D_e($phi, $k, $n, $mode, $result)>
1545 =item C<gsl_sf_ellint_D($phi, $k, $n, $mode)>
1547 -These functions compute the incomplete elliptic integral D(\phi,k) which is defined through the Carlson form RD(x,y,z) by the following relation, D(\phi,k,n) = (1/3)(\sin(\phi))^3 RD (1-\sin^2(\phi), 1-k^2 \sin^2(\phi), 1). The argument $n is not used and will be removed in a future release.
1553 =item C<gsl_sf_ellint_RC_e($x, $y, $mode, $result)>
1555 =item C<gsl_sf_ellint_RC($x, $y, $mode)>
1557 - These routines compute the incomplete elliptic integral RC($x,$y) to the accuracy specified by the mode variable $mode.
1563 =item C<gsl_sf_ellint_RD_e($x, $y, $z, $mode, $result)>
1565 =item C<gsl_sf_ellint_RD($x, $y, $z, $mode)>
1567 - These routines compute the incomplete elliptic integral RD($x,$y,$z) to the accuracy specified by the mode variable $mode.
1573 =item C<gsl_sf_ellint_RF_e($x, $y, $z, $mode, $result)>
1575 =item C<gsl_sf_ellint_RF($x, $y, $z, $mode)>
1577 - These routines compute the incomplete elliptic integral RF($x,$y,$z) to the accuracy specified by the mode variable $mode.
1583 =item C<gsl_sf_ellint_RJ_e($x, $y, $z, $p, $mode, $result)>
1585 =item C<gsl_sf_ellint_RJ($x, $y, $z, $p, $mode)>
1587 - These routines compute the incomplete elliptic integral RJ($x,$y,$z,$p) to the accuracy specified by the mode variable $mode.
1593 =item C<gsl_sf_elljac_e($u, $m)> - This function computes the Jacobian elliptic functions sn(u|m), cn(u|m), dn(u|m) by descending Landen transformations. The function returns 0 if the operation succeded, 1 otherwise and then returns the result of sn, cn and dn in this order.
1595 =item C<gsl_sf_erfc_e($x, $result)>
1597 =item C<gsl_sf_erfc($x)>
1599 -These routines compute the complementary error function erfc(x) = 1 - erf(x) = (2/\sqrt(\pi)) \int_x^\infty \exp(-t^2).
1605 =item C<gsl_sf_log_erfc_e($x, $result)>
1607 =item C<gsl_sf_log_erfc($x)>
1609 -These routines compute the logarithm of the complementary error function \log(\erfc(x)).
1615 =item C<gsl_sf_erf_e($x, $result)>
1617 =item C<gsl_sf_erf($x)>
1619 -These routines compute the error function erf(x), where erf(x) = (2/\sqrt(\pi)) \int_0^x dt \exp(-t^2).
1625 =item C<gsl_sf_erf_Z_e($x, $result)>
1627 =item C<gsl_sf_erf_Z($x)>
1629 -These routines compute the Gaussian probability density function Z(x) = (1/\sqrt{2\pi}) \exp(-x^2/2).
1635 =item C<gsl_sf_erf_Q_e($x, $result)>
1637 =item C<gsl_sf_erf_Q($x)>
1639 - These routines compute the upper tail of the Gaussian probability function Q(x) = (1/\sqrt{2\pi}) \int_x^\infty dt \exp(-t^2/2). The hazard function for the normal distribution, also known as the inverse Mill's ratio, is defined as, h(x) = Z(x)/Q(x) = \sqrt{2/\pi} \exp(-x^2 / 2) / \erfc(x/\sqrt 2) It decreases rapidly as x approaches -\infty and asymptotes to h(x) \sim x as x approaches +\infty.
1645 =item C<gsl_sf_hazard_e($x, $result)>
1647 =item C<gsl_sf_hazard($x)>
1649 - These routines compute the hazard function for the normal distribution.
1655 =item C<gsl_sf_exp_e($x, $result)>
1657 =item C<gsl_sf_exp($x)>
1659 - These routines provide an exponential function \exp(x) using GSL semantics and error checking.
1665 =item C<gsl_sf_exp_e10_e> -
1671 =item C<gsl_sf_exp_mult_e >
1673 =item C<gsl_sf_exp_mult>
1681 =item C<gsl_sf_exp_mult_e10_e> -
1687 =item C<gsl_sf_expm1_e($x, $result)>
1689 =item C<gsl_sf_expm1($x)>
1691 -These routines compute the quantity \exp(x)-1 using an algorithm that is accurate for small x.
1697 =item C<gsl_sf_exprel_e($x, $result)>
1699 =item C<gsl_sf_exprel($x)>
1701 -These routines compute the quantity (\exp(x)-1)/x using an algorithm that is accurate for small x. For small x the algorithm is based on the expansion (\exp(x)-1)/x = 1 + x/2 + x^2/(2*3) + x^3/(2*3*4) + \dots.
1707 =item C<gsl_sf_exprel_2_e($x, $result)>
1709 =item C<gsl_sf_exprel_2($x)>
1711 -These routines compute the quantity 2(\exp(x)-1-x)/x^2 using an algorithm that is accurate for small x. For small x the algorithm is based on the expansion 2(\exp(x)-1-x)/x^2 = 1 + x/3 + x^2/(3*4) + x^3/(3*4*5) + \dots.
1717 =item C<gsl_sf_exprel_n_e($x, $result)>
1719 =item C<gsl_sf_exprel_n($x)>
1721 -These routines compute the N-relative exponential, which is the n-th generalization of the functions gsl_sf_exprel and gsl_sf_exprel2. The N-relative exponential is given by,
1722 exprel_N(x) = N!/x^N (\exp(x) - \sum_{k=0}^{N-1} x^k/k!)
1723 = 1 + x/(N+1) + x^2/((N+1)(N+2)) + ...
1730 =item C<gsl_sf_exp_err_e($x, $dx, $result)> - This function exponentiates $x with an associated absolute error $dx.
1732 =item C<gsl_sf_exp_err_e10_e> -
1734 =item C<gsl_sf_exp_mult_err_e($x, $dx, $y, $dy, $result)> -
1736 =item C<gsl_sf_exp_mult_err_e10_e> -
1742 =item C<gsl_sf_expint_E1_e($x, $result)>
1744 =item C<gsl_sf_expint_E1($x)>
1746 -These routines compute the exponential integral E_1(x), E_1(x) := \Re \int_1^\infty dt \exp(-xt)/t.
1752 =item C<gsl_sf_expint_E2_e($x, $result)>
1754 =item C<gsl_sf_expint_E2($x)>
1756 -These routines compute the second-order exponential integral E_2(x),
1757 E_2(x) := \Re \int_1^\infty dt \exp(-xt)/t^2.
1764 =item C<gsl_sf_expint_En_e($n, $x, $result)>
1766 =item C<gsl_sf_expint_En($n, $x)>
1768 -These routines compute the exponential integral E_n(x) of order n,
1769 E_n(x) := \Re \int_1^\infty dt \exp(-xt)/t^n.
1776 =item C<gsl_sf_expint_E1_scaled_e >
1778 =item C<gsl_sf_expint_E1_scaled>
1786 =item C<gsl_sf_expint_E2_scaled_e>
1788 =item C<gsl_sf_expint_E2_scaled >
1796 =item C<gsl_sf_expint_En_scaled_e>
1798 =item C<gsl_sf_expint_En_scaled>
1806 =item C<gsl_sf_expint_Ei_e($x, $result)>
1808 =item C<gsl_sf_expint_Ei($x)>
1810 -These routines compute the exponential integral Ei(x), Ei(x) := - PV(\int_{-x}^\infty dt \exp(-t)/t) where PV denotes the principal value of the integral.
1816 =item C<gsl_sf_expint_Ei_scaled_e>
1818 =item C<gsl_sf_expint_Ei_scaled >
1826 =item C<gsl_sf_Shi_e($x, $result)>
1828 =item C<gsl_sf_Shi($x)>
1830 -These routines compute the integral Shi(x) = \int_0^x dt \sinh(t)/t.
1836 =item C<gsl_sf_Chi_e($x, $result)>
1838 =item C<gsl_sf_Chi($x)>
1840 -These routines compute the integral Chi(x) := \Re[ \gamma_E + \log(x) + \int_0^x dt (\cosh[t]-1)/t] , where \gamma_E is the Euler constant (available as $M_EULER from the Math::GSL::Const module).
1846 =item C<gsl_sf_expint_3_e($x, $result)>
1848 =item C<gsl_sf_expint_3($x)>
1850 -These routines compute the third-order exponential integral Ei_3(x) = \int_0^xdt \exp(-t^3) for x >= 0.
1856 =item C<gsl_sf_Si_e($x, $result)>
1858 =item C<gsl_sf_Si($x)>
1860 -These routines compute the Sine integral Si(x) = \int_0^x dt \sin(t)/t.
1866 =item C<gsl_sf_Ci_e($x, $result)>
1868 =item C<gsl_sf_Ci($x)>
1870 -These routines compute the Cosine integral Ci(x) = -\int_x^\infty dt \cos(t)/t for x > 0.
1876 =item C<gsl_sf_fermi_dirac_m1_e($x, $result)>
1878 =item C<gsl_sf_fermi_dirac_m1($x)>
1880 -These routines compute the complete Fermi-Dirac integral with an index of -1. This integral is given by F_{-1}(x) = e^x / (1 + e^x).
1886 =item C<gsl_sf_fermi_dirac_0_e($x, $result)>
1888 =item C<gsl_sf_fermi_dirac_0($x)>
1890 -These routines compute the complete Fermi-Dirac integral with an index of 0. This integral is given by F_0(x) = \ln(1 + e^x).
1896 =item C<gsl_sf_fermi_dirac_1_e($x, $result)>
1898 =item C<gsl_sf_fermi_dirac_1($x)>
1900 -These routines compute the complete Fermi-Dirac integral with an index of 1, F_1(x) = \int_0^\infty dt (t /(\exp(t-x)+1)).
1906 =item C<gsl_sf_fermi_dirac_2_e($x, $result)>
1908 =item C<gsl_sf_fermi_dirac_2($x)>
1910 -These routines compute the complete Fermi-Dirac integral with an index of 2, F_2(x) = (1/2) \int_0^\infty dt (t^2 /(\exp(t-x)+1)).
1916 =item C<gsl_sf_fermi_dirac_int_e($j, $x, $result)>
1918 =item C<gsl_sf_fermi_dirac_int($j, $x)>
1920 -These routines compute the complete Fermi-Dirac integral with an integer index of j, F_j(x) = (1/\Gamma(j+1)) \int_0^\infty dt (t^j /(\exp(t-x)+1)).
1926 =item C<gsl_sf_fermi_dirac_mhalf_e($x, $result)>
1928 =item C<gsl_sf_fermi_dirac_mhalf($x)>
1930 -These routines compute the complete Fermi-Dirac integral F_{-1/2}(x).
1936 =item C<gsl_sf_fermi_dirac_half_e($x, $result)>
1938 =item C<gsl_sf_fermi_dirac_half($x)>
1940 -These routines compute the complete Fermi-Dirac integral F_{1/2}(x).
1946 =item C<gsl_sf_fermi_dirac_3half_e($x, $result)>
1948 =item C<gsl_sf_fermi_dirac_3half($x)>
1950 -These routines compute the complete Fermi-Dirac integral F_{3/2}(x).
1956 =item C<gsl_sf_fermi_dirac_inc_0_e($x, $b, $result)>
1958 =item C<gsl_sf_fermi_dirac_inc_0($x, $b, $result)>
1960 -These routines compute the incomplete Fermi-Dirac integral with an index of zero, F_0(x,b) = \ln(1 + e^{b-x}) - (b-x).
1966 =item C<gsl_sf_legendre_Pl_e($l, $x, $result)>
1968 =item C<gsl_sf_legendre_Pl($l, $x)>
1970 -These functions evaluate the Legendre polynomial P_l(x) for a specific value of l, x subject to l >= 0, |x| <= 1
1976 =item C<gsl_sf_legendre_Pl_array>
1978 =item C<gsl_sf_legendre_Pl_deriv_array>
1986 =item C<gsl_sf_legendre_P1_e($x, $result)>
1988 =item C<gsl_sf_legendre_P2_e($x, $result)>
1990 =item C<gsl_sf_legendre_P3_e($x, $result)>
1992 =item C<gsl_sf_legendre_P1($x)>
1994 =item C<gsl_sf_legendre_P2($x)>
1996 =item C<gsl_sf_legendre_P3($x)>
1998 -These functions evaluate the Legendre polynomials P_l(x) using explicit representations for l=1, 2, 3.
2004 =item C<gsl_sf_legendre_Q0_e($x, $result)>
2006 =item C<gsl_sf_legendre_Q0($x)>
2008 -These routines compute the Legendre function Q_0(x) for x > -1, x != 1.
2014 =item C<gsl_sf_legendre_Q1_e($x, $result)>
2016 =item C<gsl_sf_legendre_Q1($x)>
2018 -These routines compute the Legendre function Q_1(x) for x > -1, x != 1.
2024 =item C<gsl_sf_legendre_Ql_e($l, $x, $result)>
2026 =item C<gsl_sf_legendre_Ql($l, $x)>
2028 -These routines compute the Legendre function Q_l(x) for x > -1, x != 1 and l >= 0.
2034 =item C<gsl_sf_legendre_Plm_e($l, $m, $x, $result)>
2036 =item C<gsl_sf_legendre_Plm($l, $m, $x)>
2038 -These routines compute the associated Legendre polynomial P_l^m(x) for m >= 0, l >= m, |x| <= 1.
2044 =item C<gsl_sf_legendre_Plm_array>
2046 =item C<gsl_sf_legendre_Plm_deriv_array >
2054 =item C<gsl_sf_legendre_sphPlm_e($l, $m, $x, $result)>
2056 =item C<gsl_sf_legendre_sphPlm($l, $m, $x)>
2058 -These routines compute the normalized associated Legendre polynomial $\sqrt{(2l+1)/(4\pi)} \sqrt{(l-m)!/(l+m)!} P_l^m(x)$ suitable for use in spherical harmonics. The parameters must satisfy m >= 0, l >= m, |x| <= 1. Theses routines avoid the overflows that occur for the standard normalization of P_l^m(x).
2064 =item C<gsl_sf_legendre_sphPlm_array >
2066 =item C<gsl_sf_legendre_sphPlm_deriv_array>
2074 =item C<gsl_sf_legendre_array_size> -
2080 =item C<gsl_sf_lngamma_e($x, $result)>
2082 =item C<gsl_sf_lngamma($x)>
2084 -These routines compute the logarithm of the Gamma function, \log(\Gamma(x)), subject to x not being a negative integer or zero. For x<0 the real part of \log(\Gamma(x)) is returned, which is equivalent to \log(|\Gamma(x)|). The function is computed using the real Lanczos method.
2090 =item C<gsl_sf_lngamma_sgn_e($x, $result_lg)> - This routine returns the sign of the gamma function and the logarithm of its magnitude into this order, subject to $x not being a negative integer or zero. The function is computed using the real Lanczos method. The value of the gamma function can be reconstructed using the relation \Gamma(x) = sgn * \exp(resultlg).
2096 =item C<gsl_sf_gamma_e >
2098 =item C<gsl_sf_gamma>
2106 =item C<gsl_sf_gammastar_e>
2108 =item C<gsl_sf_gammastar >
2116 =item C<gsl_sf_gammainv_e>
2118 =item C<gsl_sf_gammainv>
2126 =item C<gsl_sf_lngamma_complex_e >
2134 =item C<gsl_sf_gamma_inc_Q_e>
2136 =item C<gsl_sf_gamma_inc_Q>
2144 =item C<gsl_sf_gamma_inc_P_e >
2146 =item C<gsl_sf_gamma_inc_P>
2154 =item C<gsl_sf_gamma_inc_e>
2156 =item C<gsl_sf_gamma_inc >
2164 =item C<gsl_sf_taylorcoeff_e>
2166 =item C<gsl_sf_taylorcoeff>
2174 =item C<gsl_sf_fact_e >
2176 =item C<gsl_sf_fact>
2184 =item C<gsl_sf_doublefact_e>
2186 =item C<gsl_sf_doublefact >
2194 =item C<gsl_sf_lnfact_e>
2196 =item C<gsl_sf_lnfact>
2204 =item C<gsl_sf_lndoublefact_e >
2206 =item C<gsl_sf_lndoublefact>
2214 =item C<gsl_sf_lnchoose_e>
2216 =item C<gsl_sf_lnchoose >
2224 =item C<gsl_sf_choose_e>
2226 =item C<gsl_sf_choose>
2234 =item C<gsl_sf_lnpoch_e >
2236 =item C<gsl_sf_lnpoch>
2244 =item C<gsl_sf_lnpoch_sgn_e>
2252 =item C<gsl_sf_poch_e >
2254 =item C<gsl_sf_poch>
2262 =item C<gsl_sf_pochrel_e>
2264 =item C<gsl_sf_pochrel >
2272 =item C<gsl_sf_lnbeta_e>
2274 =item C<gsl_sf_lnbeta>
2282 =item C<gsl_sf_lnbeta_sgn_e >
2290 =item C<gsl_sf_beta_e>
2292 =item C<gsl_sf_beta>
2300 =item C<gsl_sf_beta_inc_e >
2302 =item C<gsl_sf_beta_inc>
2310 =item C<gsl_sf_gegenpoly_1_e>
2312 =item C<gsl_sf_gegenpoly_2_e >
2314 =item C<gsl_sf_gegenpoly_3_e>
2316 =item C<gsl_sf_gegenpoly_1>
2318 =item C<gsl_sf_gegenpoly_2 >
2320 =item C<gsl_sf_gegenpoly_3>
2328 =item C<gsl_sf_gegenpoly_n_e>
2330 =item C<gsl_sf_gegenpoly_n >
2338 =item C<gsl_sf_gegenpoly_array>
2340 =item C<gsl_sf_hyperg_0F1_e>
2342 =item C<gsl_sf_hyperg_0F1 >
2350 =item C<gsl_sf_hyperg_1F1_int_e>
2352 =item C<gsl_sf_hyperg_1F1_int>
2360 =item C<gsl_sf_hyperg_1F1_e >
2362 =item C<gsl_sf_hyperg_1F1>
2370 =item C<gsl_sf_hyperg_U_int_e>
2372 =item C<gsl_sf_hyperg_U_int >
2380 =item C<gsl_sf_hyperg_U_int_e10_e>
2388 =item C<gsl_sf_hyperg_U_e>
2390 =item C<gsl_sf_hyperg_U >
2398 =item C<gsl_sf_hyperg_U_e10_e>
2406 =item C<gsl_sf_hyperg_2F1_e>
2408 =item C<gsl_sf_hyperg_2F1 >
2416 =item C<gsl_sf_hyperg_2F1_conj_e>
2418 =item C<gsl_sf_hyperg_2F1_conj>
2426 =item C<gsl_sf_hyperg_2F1_renorm_e >
2428 =item C<gsl_sf_hyperg_2F1_renorm>
2436 =item C<gsl_sf_hyperg_2F1_conj_renorm_e>
2438 =item C<gsl_sf_hyperg_2F1_conj_renorm >
2446 =item C<gsl_sf_hyperg_2F0_e>
2448 =item C<gsl_sf_hyperg_2F0>
2456 =item C<gsl_sf_laguerre_1_e >
2458 =item C<gsl_sf_laguerre_2_e>
2460 =item C<gsl_sf_laguerre_3_e>
2462 =item C<gsl_sf_laguerre_1 >
2464 =item C<gsl_sf_laguerre_2>
2466 =item C<gsl_sf_laguerre_3>
2474 =item C<gsl_sf_laguerre_n_e >
2476 =item C<gsl_sf_laguerre_n>
2484 =item C<gsl_sf_lambert_W0_e>
2486 =item C<gsl_sf_lambert_W0 >
2494 =item C<gsl_sf_lambert_Wm1_e>
2496 =item C<gsl_sf_lambert_Wm1>
2504 =item C<gsl_sf_conicalP_half_e >
2506 =item C<gsl_sf_conicalP_half>
2514 =item C<gsl_sf_conicalP_mhalf_e>
2516 =item C<gsl_sf_conicalP_mhalf >
2524 =item C<gsl_sf_conicalP_0_e>
2526 =item C<gsl_sf_conicalP_0>
2534 =item C<gsl_sf_conicalP_1_e >
2536 =item C<gsl_sf_conicalP_1>
2544 =item C<gsl_sf_conicalP_sph_reg_e>
2546 =item C<gsl_sf_conicalP_sph_reg >
2554 =item C<gsl_sf_conicalP_cyl_reg_e>
2556 =item C<gsl_sf_conicalP_cyl_reg>
2564 =item C<gsl_sf_legendre_H3d_0_e >
2566 =item C<gsl_sf_legendre_H3d_0>
2574 =item C<gsl_sf_legendre_H3d_1_e>
2576 =item C<gsl_sf_legendre_H3d_1 >
2584 =item C<gsl_sf_legendre_H3d_e>
2586 =item C<gsl_sf_legendre_H3d>
2594 =item C<gsl_sf_legendre_H3d_array >
2602 =item C<gsl_sf_log_e>
2612 =item C<gsl_sf_log_abs_e >
2614 =item C<gsl_sf_log_abs>
2622 =item C<gsl_sf_complex_log_e>
2630 =item C<gsl_sf_log_1plusx_e >
2632 =item C<gsl_sf_log_1plusx>
2640 =item C<gsl_sf_log_1plusx_mx_e>
2642 =item C<gsl_sf_log_1plusx_mx >
2650 =item C<gsl_sf_mathieu_a_array>
2652 =item C<gsl_sf_mathieu_b_array>
2660 =item C<gsl_sf_mathieu_a >
2662 =item C<gsl_sf_mathieu_b>
2670 =item C<gsl_sf_mathieu_a_coeff>
2672 =item C<gsl_sf_mathieu_b_coeff >
2680 =item C<gsl_sf_mathieu_alloc>
2688 =item C<gsl_sf_mathieu_free>
2696 =item C<gsl_sf_mathieu_ce >
2698 =item C<gsl_sf_mathieu_se>
2706 =item C<gsl_sf_mathieu_ce_array>
2708 =item C<gsl_sf_mathieu_se_array >
2716 =item C<gsl_sf_mathieu_Mc>
2718 =item C<gsl_sf_mathieu_Ms>
2726 =item C<gsl_sf_mathieu_Mc_array >
2728 =item C<gsl_sf_mathieu_Ms_array>
2736 =item C<gsl_sf_pow_int_e>
2738 =item C<gsl_sf_pow_int >
2746 =item C<gsl_sf_psi_int_e>
2748 =item C<gsl_sf_psi_int>
2756 =item C<gsl_sf_psi_e >
2766 =item C<gsl_sf_psi_1piy_e>
2768 =item C<gsl_sf_psi_1piy >
2776 =item C<gsl_sf_complex_psi_e>
2784 =item C<gsl_sf_psi_1_int_e>
2786 =item C<gsl_sf_psi_1_int >
2794 =item C<gsl_sf_psi_1_e >
2796 =item C<gsl_sf_psi_1>
2804 =item C<gsl_sf_psi_n_e >
2806 =item C<gsl_sf_psi_n>
2814 =item C<gsl_sf_result_smash_e>
2822 =item C<gsl_sf_synchrotron_1_e >
2824 =item C<gsl_sf_synchrotron_1>
2832 =item C<gsl_sf_synchrotron_2_e>
2834 =item C<gsl_sf_synchrotron_2 >
2842 =item C<gsl_sf_transport_2_e>
2844 =item C<gsl_sf_transport_2>
2852 =item C<gsl_sf_transport_3_e >
2854 =item C<gsl_sf_transport_3>
2862 =item C<gsl_sf_transport_4_e>
2864 =item C<gsl_sf_transport_4 >
2872 =item C<gsl_sf_transport_5_e>
2874 =item C<gsl_sf_transport_5>
2882 =item C<gsl_sf_sin_e >
2892 =item C<gsl_sf_cos_e>
2894 =item C<gsl_sf_cos >
2902 =item C<gsl_sf_hypot_e>
2904 =item C<gsl_sf_hypot>
2912 =item C<gsl_sf_complex_sin_e >
2920 =item C<gsl_sf_complex_cos_e>
2928 =item C<gsl_sf_complex_logsin_e>
2936 =item C<gsl_sf_sinc_e >
2938 =item C<gsl_sf_sinc>
2946 =item C<gsl_sf_lnsinh_e>
2948 =item C<gsl_sf_lnsinh >
2956 =item C<gsl_sf_lncosh_e>
2958 =item C<gsl_sf_lncosh>
2966 =item C<gsl_sf_polar_to_rect >
2974 =item C<gsl_sf_rect_to_polar>
2982 =item C<gsl_sf_sin_err_e>
2984 =item C<gsl_sf_cos_err_e >
2992 =item C<gsl_sf_angle_restrict_symm_e>
2994 =item C<gsl_sf_angle_restrict_symm>
3002 =item C<gsl_sf_angle_restrict_pos_e >
3004 =item C<gsl_sf_angle_restrict_pos>
3012 =item C<gsl_sf_angle_restrict_symm_err_e>
3014 =item C<gsl_sf_angle_restrict_pos_err_e >
3018 =item C<gsl_sf_atanint_e>
3020 =item C<gsl_sf_atanint>
3022 -These routines compute the Arctangent integral, which is defined as AtanInt(x) = \int_0^x dt \arctan(t)/t.
3026 =item C<gsl_sf_zeta_int_e >
3028 =item C<gsl_sf_zeta_int>
3030 =item C<gsl_sf_zeta_e gsl_sf_zeta >
3032 =item C<gsl_sf_zetam1_e>
3034 =item C<gsl_sf_zetam1>
3036 =item C<gsl_sf_zetam1_int_e >
3038 =item C<gsl_sf_zetam1_int>
3040 =item C<gsl_sf_hzeta_e>
3042 =item C<gsl_sf_hzeta >
3044 =item C<gsl_sf_eta_int_e>
3046 =item C<gsl_sf_eta_int>
3048 =item C<gsl_sf_eta_e>
3050 =item C<gsl_sf_eta >
3054 This module also contains the following constants used as mode in various of those functions :
3058 =item * GSL_PREC_DOUBLE - Double-precision, a relative accuracy of approximately 2 * 10^-16.
3060 =item * GSL_PREC_SINGLE - Single-precision, a relative accuracy of approximately 10^-7.
3062 =item * GSL_PREC_APPROX - Approximate values, a relative accuracy of approximately 5 * 10^-4.
3066 You can import the functions that you want to use by giving a space separated
3067 list to Math::GSL::SF when you use the package. You can also write
3068 use Math::GSL::SF qw/:all/
3069 to use all avaible functions of the module. Note that
3070 the tag names begin with a colon. Other tags are also available, here is a
3071 complete list of all tags for this module :
3101 =item C<hypergeometric>
3121 For more informations on the functions, we refer you to the GSL offcial
3122 documentation: L<http://www.gnu.org/software/gsl/manual/html_node/>
3124 Tip : search on google: site:http://www.gnu.org/software/gsl/manual/html_node/name_of_the_function_you_want
3128 This example computes the dilogarithm of 1/10 :
3130 use Math::GSL::SF qw/dilog/;
3131 my $x = gsl_sf_dilog(0.1);
3132 print "gsl_sf_dilog(0.1) = $x\n";
3134 An example using Math::GSL::SF and gnuplot is in the B<examples/sf> folder of the source code.
3138 Jonathan Leto <jonathan@leto.net> and Thierry Moisan <thierry.moisan@gmail.com>
3140 =head1 COPYRIGHT AND LICENSE
3142 Copyright (C) 2008-2009 Jonathan Leto and Thierry Moisan
3144 This program is free software; you can redistribute it and/or modify it
3145 under the same terms as Perl itself.