NAME
Math::GSL::SF - Special Functions
SYNOPSIS
use Math::GSL::SF qw/:all/;
DESCRIPTION
This module contains a data structure named gsl_sf_result. To create a new one use
$r = Math::GSL::SF::gsl_sf_result_struct->new;
You can then access the elements of the structure in this way :
my $val = $r->{val};
my $error = $r->{err};
Here is a list of all included functions:
gsl_sf_airy_Ai_e($x, $mode)
gsl_sf_airy_Ai($x, $mode, $result)
-
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.
gsl_sf_airy_Bi_e($x, $mode, $result)
gsl_sf_airy_Bi($x, $mode)
-
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.
gsl_sf_airy_Ai_scaled_e($x, $mode, $result)
gsl_sf_airy_Ai_scaled($x, $mode)
-
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.
gsl_sf_airy_Bi_scaled_e($x, $mode, $result)
gsl_sf_airy_Bi_scaled($x, $mode)
-
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.
gsl_sf_airy_Ai_deriv_e($x, $mode, $result)
gsl_sf_airy_Ai_deriv($x, $mode)
-
These routines compute the Airy function derivative Ai'($x) with an accuracy specified by $mode. $result is a gsl_sf_result structure.
gsl_sf_airy_Bi_deriv_e($x, $mode, $result)
gsl_sf_airy_Bi_deriv($x, $mode)
-
These routines compute the Airy function derivative Bi'($x) with an accuracy specified by $mode. $result is a gsl_sf_result structure.
gsl_sf_airy_Ai_deriv_scaled_e($x, $mode, $result)
gsl_sf_airy_Ai_deriv_scaled($x, $mode)
-
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.
gsl_sf_airy_Bi_deriv_scaled_e($x, $mode, $result)
gsl_sf_airy_Bi_deriv_scaled($x, $mode)
-
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.
gsl_sf_airy_zero_Ai_e($s, $result)
gsl_sf_airy_zero_Ai($s)
-
These routines compute the location of the s-th zero of the Airy function Ai($x). $result is a gsl_sf_result structure.
gsl_sf_airy_zero_Bi_e($s, $result)
gsl_sf_airy_zero_Bi($s)
-
These routines compute the location of the s-th zero of the Airy function Bi($x). $result is a gsl_sf_result structure.
gsl_sf_airy_zero_Ai_deriv_e($s, $result)
gsl_sf_airy_zero_Ai_deriv($s)
-
These routines compute the location of the s-th zero of the Airy function derivative Ai'(x). $result is a gsl_sf_result structure.
gsl_sf_airy_zero_Bi_deriv_e($s, $result)
gsl_sf_airy_zero_Bi_deriv($s)
-
These routines compute the location of the s-th zero of the Airy function derivative Bi'(x). $result is a gsl_sf_result structure.
gsl_sf_bessel_J0_e($x, $result)
gsl_sf_bessel_J0($x)
-
These routines compute the regular cylindrical Bessel function of zeroth order, J_0($x). $result is a gsl_sf_result structure.
gsl_sf_bessel_J1_e($x, $result)
gsl_sf_bessel_J1($x)
-
These routines compute the regular cylindrical Bessel function of first order, J_1($x). $result is a gsl_sf_result structure.
gsl_sf_bessel_Jn_e($n, $x, $result)
gsl_sf_bessel_Jn($n, $x)
-
These routines compute the regular cylindrical Bessel function of order n, J_n($x).
gsl_sf_bessel_Jn_array($nmin, $nmax, $x)
-
This routine computes the values of the regular cylindrical Bessel functions J_n($x) for n from $nmin to $nmax inclusive, returning the results in an array reference. The values are computed using recurrence relations for efficiency, and therefore may differ slightly from the exact values.
gsl_sf_bessel_Y0_e($x, $result)
gsl_sf_bessel_Y0($x)
-
These routines compute the irregular spherical Bessel function of zeroth order, y_0(x) = -\cos(x)/x.
gsl_sf_bessel_Y1_e($x, $result)
gsl_sf_bessel_Y1($x)
-
These routines compute the irregular spherical Bessel function of first order, y_1(x) = -(\cos(x)/x + \sin(x))/x.
gsl_sf_bessel_Yn_e
($n, $x, $result)gsl_sf_bessel_Yn($n, $x)
-
These routines compute the irregular cylindrical Bessel function of order $n, Y_n(x), for x>0.
gsl_sf_bessel_Yn_array($nmin, $nmax, $x)
-
This routine computes the values of the irregular cylindrical Bessel functions Y_n(x) for n from $nmin to $nmax inclusive, returning the results in an array reference. The domain of the function is $x>0. The values are computed using recurrence relations for efficiency, and therefore may differ slightly from the exact values.
gsl_sf_bessel_I0_e($x, $result)
gsl_sf_bessel_I0($x)
-
These routines compute the regular modified cylindrical Bessel function of zeroth order, I_0(x).
gsl_sf_bessel_In_e($n, $x, $result)
gsl_sf_bessel_In($n, $x)
-
These routines compute the regular modified cylindrical Bessel function of order $n, I_n(x).
gsl_sf_bessel_In_array($nmin, $nmax, $x)
-
This routine computes the values of the regular modified cylindrical Bessel functions I_n(x) for n from $nmin to $nmax inclusive, returning the results in an array reference. The start of the range nmin must be positive or zero. The values are computed using recurrence relations for efficiency, and therefore may differ slightly from the exact values.
gsl_sf_bessel_I0_scaled_e($x, $result)
gsl_sf_bessel_I0_scaled($x)
-
These routines compute the scaled regular modified cylindrical Bessel function of zeroth order \exp(-|x|) I_0(x).
gsl_sf_bessel_I1_scaled_e($x, $result)
gsl_sf_bessel_I1_scaled($x)
-
These routines compute the scaled regular modified cylindrical Bessel function of first order \exp(-|x|) I_1(x).
gsl_sf_bessel_In_scaled_e($n, $x, $result)
gsl_sf_bessel_In_scaled($n, $x)
-
These routines compute the scaled regular modified cylindrical Bessel function of order $n, \exp(-|x|) I_n(x)
gsl_sf_bessel_In_scaled_array($nmin, $nmax, $x)
-
This routine computes the values of the scaled regular cylindrical Bessel functions exp(-|$x|) I_n($x) for n from $nmin to $nmax inclusive, returning the results in an array reference. The start of the range nmin must be positive or zero. The values are computed using recurrence relations for efficiency, and therefore may differ slightly from the exact values.
gsl_sf_bessel_K0_e($x, $result)
gsl_sf_bessel_K0($x)
-
These routines compute the irregular modified cylindrical Bessel function of zeroth order, K_0(x), for x > 0.
gsl_sf_bessel_K1_e($x, $result)
gsl_sf_bessel_K1($x)
-
These routines compute the irregular modified cylindrical Bessel function of first order, K_1(x), for x > 0.
gsl_sf_bessel_Kn_e($n, $x, $result)
gsl_sf_bessel_Kn($n, $x)
-
These routines compute the irregular modified cylindrical Bessel function of order $n, K_n(x), for x > 0.
gsl_sf_bessel_Kn_array($nmin, $nmax, $x)
-
This routine computes the values of the Irregular Modified Cylindrical Bessel Functions K_n($x) for n from $nmin to $nmax inclusive, returning the results in an array reference. The values are computed using recurrence relations for efficiency, and therefore may differ slightly from the exact values. This function is only defined for positive $x and with throw an exception otherwise.
gsl_sf_bessel_K0_scaled_e($x, $result)
gsl_sf_bessel_K0_scaled($x)
-
These routines compute the scaled irregular modified cylindrical Bessel function of zeroth order \exp(x) K_0(x) for x>0.
gsl_sf_bessel_Kn_scaled_array($nmin, $max, $x)
-
This routine computes the values of the scaled irregular cylindrical Bessel functions exp(x) K_n(x) for n from nmin to nmax inclusive, storing the results in the array result_array. The start of the range nmin must be positive or zero. The domain of the function is x>0. The values are computed using recurrence relations for efficiency, and therefore may differ slightly from the exact values.
gsl_sf_bessel_jl_array($lmax, $x)
-
This routine computes the values of the regular spherical Bessel functions j_l(x) for l from 0 to
$lmax
inclusive for$lmax
>= 0 and $x >= 0, returning the results in an array reference. The values are computed using recurrence relations for efficiency, and therefore may differ slightly from the exact values.
gsl_sf_bessel_jl_steed_array($lmax, $x)
-
This routine uses Steed’s method to compute the values of the regular spherical Bessel functions j_l(x) for l from 0 to $lmax inclusive for $lmax >= 0 and $x >= 0, storing the results in the array result_array. The Steed/Barnett algorithm is described in Comp. Phys. Comm. 21, 297 (1981). Steed’s method is more stable than the recurrence used in the other functions but is also slower.
gsl_sf_bessel_yl_array($lmax, $x)
-
This routine computes the values of the irregular spherical Bessel functions y_l(x) for l from 0 to $lmax inclusive for lmax >= 0, returning the results in an array reference. The values are computed using recurrence relations for efficiency, and therefore may differ slightly from the exact values.
gsl_sf_bessel_il_scaled_array($lmax, $x)
-
This routine computes the values of the scaled regular modified spherical Bessel functions exp(-|x|) i_l(x) for l from 0 to $lmax inclusive for $lmax >= 0, returning the results in an array reference. The values are computed using recurrence relations for efficiency, and therefore may differ slightly from the exact values.
gsl_sf_bessel_kl_scaled_array($lmax, $x)
-
This routine computes the values of the scaled irregular modified spherical Bessel functions exp($x) k_l($x) for l from 0 to lmax inclusive for $lmax >= 0 and $x>0, returning the results in an array reference. The values are computed using recurrence relations for efficiency, and therefore may differ slightly from the exact values.
gsl_sf_bessel_sequence_Jnu_e
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.gsl_sf_coulomb_wave_F_array
-gsl_sf_coulomb_wave_FG_array
-gsl_sf_coulomb_wave_FGp_array
-gsl_sf_coulomb_wave_sphF_array
-gsl_sf_coulomb_CL_e($L, $eta, $result)
- This function computes the Coulomb wave function normalization constant C_L($eta) for $L > -1.gsl_sf_coulomb_CL_arrayi
-
gsl_sf_coupling_3j_e($two_ja, $two_jb, $two_jc, $two_ma, $two_mb, $two_mc, $result)
gsl_sf_coupling_3j($two_ja, $two_jb, $two_jc, $two_ma, $two_mb, $two_mc)
-
These routines compute the Wigner 3-j coefficient, (ja jb jc ma mb mc) where the arguments are given in half-integer units, ja = $two_ja/2, ma = $two_ma/2, etc.
gsl_sf_coupling_6j_e($two_ja, $two_jb, $two_jc, $two_jd, $two_je, $two_jf, $result)
gsl_sf_coupling_6j($two_ja, $two_jb, $two_jc, $two_jd, $two_je, $two_jf)
-
These routines compute the Wigner 6-j coefficient, {ja jb jc jd je jf} where the arguments are given in half-integer units, ja = $two_ja/2, ma = $two_ma/2, etc.
gsl_sf_coupling_RacahW_e
gsl_sf_coupling_RacahW
gsl_sf_coupling_9j_e($two_ja, $two_jb, $two_jc, $two_jd, $two_je, $two_jf, $two_jg, $two_jh, $two_ji, $result)
gsl_sf_coupling_9j($two_ja, $two_jb, $two_jc, $two_jd, $two_je, $two_jf, $two_jg, $two_jh, $two_ji)
-
These routines compute the Wigner 9-j coefficient,
{ja jb jc jd je jf jg jh ji} where the arguments are given in half-integer units, ja = two_ja/2, ma = two_ma/2, etc.
gsl_sf_dawson_e($x, $result)
gsl_sf_dawson($x)
-
These routines compute the value of Dawson's integral for $x.
gsl_sf_debye_1_e($x, $result)
gsl_sf_debye_1($x)
-
These routines compute the first-order Debye function D_1(x) = (1/x) \int_0^x dt (t/(e^t - 1)).
gsl_sf_debye_2_e($x, $result)
gsl_sf_debye_2($x)
-
These routines compute the second-order Debye function D_2(x) = (2/x^2) \int_0^x dt (t^2/(e^t - 1)).
gsl_sf_debye_3_e($x, $result)
gsl_sf_debye_3($x)
-
These routines compute the third-order Debye function D_3(x) = (3/x^3) \int_0^x dt (t^3/(e^t - 1)).
gsl_sf_debye_4_e($x, $result)
gsl_sf_debye_4($x)
-
These routines compute the fourth-order Debye function D_4(x) = (4/x^4) \int_0^x dt (t^4/(e^t - 1)).
gsl_sf_debye_5_e($x, $result)
gsl_sf_debye_5($x)
-
These routines compute the fifth-order Debye function D_5(x) = (5/x^5) \int_0^x dt (t^5/(e^t - 1)).
gsl_sf_debye_6_e($x, $result)
gsl_sf_debye_6($x)
-
These routines compute the sixth-order Debye function D_6(x) = (6/x^6) \int_0^x dt (t^6/(e^t - 1)).
gsl_sf_dilog_e ($x, $result)
gsl_sf_dilog($x)
-
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).
gsl_sf_complex_dilog_xy_e
-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.gsl_sf_complex_spence_xy_e
-
gsl_sf_multiply
gsl_sf_multiply_e($x, $y, $result)
- This function multiplies $x and $y storing the product and its associated error in $result.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.
gsl_sf_ellint_Kcomp_e($k, $mode, $result)
gsl_sf_ellint_Kcomp($k, $mode)
-
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.
gsl_sf_ellint_Dcomp_e
gsl_sf_ellint_Dcomp
gsl_sf_ellint_F_e($phi, $k, $mode, $result)
gsl_sf_ellint_F($phi, $k, $mode)
-
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.
gsl_sf_ellint_E_e($phi, $k, $mode, $result)
gsl_sf_ellint_E($phi, $k, $mode)
-
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.
gsl_sf_ellint_P_e($phi, $k, $n, $mode, $result)
gsl_sf_ellint_P($phi, $k, $n, $mode)
-
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.
gsl_sf_ellint_D_e($phi, $k, $n, $mode, $result)
gsl_sf_ellint_D($phi, $k, $n, $mode)
-
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.
gsl_sf_ellint_RC_e($x, $y, $mode, $result)
gsl_sf_ellint_RC($x, $y, $mode)
-
These routines compute the incomplete elliptic integral RC($x,$y) to the accuracy specified by the mode variable $mode.
gsl_sf_ellint_RD_e($x, $y, $z, $mode, $result)
gsl_sf_ellint_RD($x, $y, $z, $mode)
-
These routines compute the incomplete elliptic integral RD($x,$y,$z) to the accuracy specified by the mode variable $mode.
gsl_sf_ellint_RF_e($x, $y, $z, $mode, $result)
gsl_sf_ellint_RF($x, $y, $z, $mode)
-
These routines compute the incomplete elliptic integral RF($x,$y,$z) to the accuracy specified by the mode variable $mode.
gsl_sf_ellint_RJ_e($x, $y, $z, $p, $mode, $result)
gsl_sf_ellint_RJ($x, $y, $z, $p, $mode)
-
These routines compute the incomplete elliptic integral RJ($x,$y,$z,$p) to the accuracy specified by the mode variable $mode.
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.gsl_sf_erfc_e($x, $result)
gsl_sf_erfc($x)
-
These routines compute the complementary error function erfc(x) = 1 - erf(x) = (2/\sqrt(\pi)) \int_x^\infty \exp(-t^2).
gsl_sf_log_erfc_e($x, $result)
gsl_sf_log_erfc($x)
-
These routines compute the logarithm of the complementary error function \log(\erfc(x)).
gsl_sf_erf_e($x, $result)
gsl_sf_erf($x)
-
These routines compute the error function erf(x), where erf(x) = (2/\sqrt(\pi)) \int_0^x dt \exp(-t^2).
gsl_sf_erf_Z_e($x, $result)
gsl_sf_erf_Z($x)
-
These routines compute the Gaussian probability density function Z(x) = (1/\sqrt{2\pi}) \exp(-x^2/2).
gsl_sf_erf_Q_e($x, $result)
gsl_sf_erf_Q($x)
-
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.
gsl_sf_hazard_e($x, $result)
gsl_sf_hazard($x)
-
These routines compute the hazard function for the normal distribution.
gsl_sf_exp_e($x, $result)
gsl_sf_exp($x)
-
These routines provide an exponential function \exp(x) using GSL semantics and error checking.
gsl_sf_exp_mult_e
gsl_sf_exp_mult
gsl_sf_expm1_e($x, $result)
gsl_sf_expm1($x)
-
These routines compute the quantity \exp(x)-1 using an algorithm that is accurate for small x.
gsl_sf_exprel_e($x, $result)
gsl_sf_exprel($x)
-
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.
gsl_sf_exprel_2_e($x, $result)
gsl_sf_exprel_2($x)
-
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.
gsl_sf_exprel_n_e($x, $result)
gsl_sf_exprel_n($x)
-
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, exprel_N(x) = N!/x^N (\exp(x) - \sum_{k=0}^{N-1} x^k/k!) = 1 + x/(N+1) + x^2/((N+1)(N+2)) + ... = 1F1 (1,1+N,x)
gsl_sf_exp_err_e($x, $dx, $result)
- This function exponentiates $x with an associated absolute error $dx.gsl_sf_exp_err_e10_e
-gsl_sf_exp_mult_err_e($x, $dx, $y, $dy, $result)
-gsl_sf_exp_mult_err_e10_e
-
gsl_sf_expint_E1_e($x, $result)
gsl_sf_expint_E1($x)
-
These routines compute the exponential integral E_1(x), E_1(x) := \Re \int_1^\infty dt \exp(-xt)/t.
gsl_sf_expint_E2_e($x, $result)
gsl_sf_expint_E2($x)
-
These routines compute the second-order exponential integral E_2(x), E_2(x) := \Re \int_1^\infty dt \exp(-xt)/t^2.
gsl_sf_expint_En_e($n, $x, $result)
gsl_sf_expint_En($n, $x)
-
These routines compute the exponential integral E_n(x) of order n, E_n(x) := \Re \int_1^\infty dt \exp(-xt)/t^n.
gsl_sf_expint_E1_scaled_e
gsl_sf_expint_E1_scaled
gsl_sf_expint_E2_scaled_e
gsl_sf_expint_E2_scaled
gsl_sf_expint_En_scaled_e
gsl_sf_expint_En_scaled
gsl_sf_expint_Ei_e($x, $result)
gsl_sf_expint_Ei($x)
-
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.
gsl_sf_expint_Ei_scaled_e
gsl_sf_expint_Ei_scaled
gsl_sf_Shi_e($x, $result)
gsl_sf_Shi($x)
-
These routines compute the integral Shi(x) = \int_0^x dt \sinh(t)/t.
gsl_sf_Chi_e($x, $result)
gsl_sf_Chi($x)
-
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).
gsl_sf_expint_3_e($x, $result)
gsl_sf_expint_3($x)
-
These routines compute the third-order exponential integral Ei_3(x) = \int_0^xdt \exp(-t^3) for x >= 0.
gsl_sf_Si_e($x, $result)
gsl_sf_Si($x)
-
These routines compute the Sine integral Si(x) = \int_0^x dt \sin(t)/t.
gsl_sf_Ci_e($x, $result)
gsl_sf_Ci($x)
-
These routines compute the Cosine integral Ci(x) = -\int_x^\infty dt \cos(t)/t for x > 0.
gsl_sf_fermi_dirac_m1_e($x, $result)
gsl_sf_fermi_dirac_m1($x)
-
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).
gsl_sf_fermi_dirac_0_e($x, $result)
gsl_sf_fermi_dirac_0($x)
-
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).
gsl_sf_fermi_dirac_1_e($x, $result)
gsl_sf_fermi_dirac_1($x)
-
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)).
gsl_sf_fermi_dirac_2_e($x, $result)
gsl_sf_fermi_dirac_2($x)
-
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)).
gsl_sf_fermi_dirac_int_e($j, $x, $result)
gsl_sf_fermi_dirac_int($j, $x)
-
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)).
gsl_sf_fermi_dirac_mhalf_e($x, $result)
gsl_sf_fermi_dirac_mhalf($x)
-
These routines compute the complete Fermi-Dirac integral F_{-1/2}(x).
gsl_sf_fermi_dirac_half_e($x, $result)
gsl_sf_fermi_dirac_half($x)
-
These routines compute the complete Fermi-Dirac integral F_{1/2}(x).
gsl_sf_fermi_dirac_3half_e($x, $result)
gsl_sf_fermi_dirac_3half($x)
-
These routines compute the complete Fermi-Dirac integral F_{3/2}(x).
gsl_sf_fermi_dirac_inc_0_e($x, $b, $result)
gsl_sf_fermi_dirac_inc_0($x, $b, $result)
-
These routines compute the incomplete Fermi-Dirac integral with an index of zero, F_0(x,b) = \ln(1 + e^{b-x}) - (b-x).
gsl_sf_legendre_Pl_e($l, $x, $result)
gsl_sf_legendre_Pl($l, $x)
-
These functions evaluate the Legendre polynomial P_l(x) for a specific value of l, x subject to l >= 0, |x| <= 1
gsl_sf_legendre_Pl_array
gsl_sf_legendre_Pl_deriv_array
gsl_sf_legendre_P1_e($x, $result)
gsl_sf_legendre_P2_e($x, $result)
gsl_sf_legendre_P3_e($x, $result)
gsl_sf_legendre_P1($x)
gsl_sf_legendre_P2($x)
gsl_sf_legendre_P3($x)
-
These functions evaluate the Legendre polynomials P_l(x) using explicit representations for l=1, 2, 3.
gsl_sf_legendre_Q0_e($x, $result)
gsl_sf_legendre_Q0($x)
-
These routines compute the Legendre function Q_0(x) for x > -1, x != 1.
gsl_sf_legendre_Q1_e($x, $result)
gsl_sf_legendre_Q1($x)
-
These routines compute the Legendre function Q_1(x) for x > -1, x != 1.
gsl_sf_legendre_Ql_e($l, $x, $result)
gsl_sf_legendre_Ql($l, $x)
-
These routines compute the Legendre function Q_l(x) for x > -1, x != 1 and l >= 0.
gsl_sf_legendre_Plm_e($l, $m, $x, $result)
gsl_sf_legendre_Plm($l, $m, $x)
-
These routines compute the associated Legendre polynomial P_l^m(x) for m >= 0, l >= m, |x| <= 1.
gsl_sf_legendre_Plm_array
gsl_sf_legendre_Plm_deriv_array
gsl_sf_legendre_sphPlm_e($l, $m, $x, $result)
gsl_sf_legendre_sphPlm($l, $m, $x)
-
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).
gsl_sf_legendre_sphPlm_array
gsl_sf_legendre_sphPlm_deriv_array
gsl_sf_lngamma_e($x, $result)
gsl_sf_lngamma($x)
-
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.
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).
gsl_sf_gamma_e
gsl_sf_gamma
gsl_sf_gammastar_e
gsl_sf_gammastar
gsl_sf_gammainv_e
gsl_sf_gammainv
gsl_sf_lngamma_complex_e
gsl_sf_gamma_inc_Q_e
gsl_sf_gamma_inc_Q
gsl_sf_gamma_inc_P_e
gsl_sf_gamma_inc_P
gsl_sf_gamma_inc_e
gsl_sf_gamma_inc
gsl_sf_taylorcoeff_e
gsl_sf_taylorcoeff
gsl_sf_fact_e
gsl_sf_fact
gsl_sf_doublefact_e
gsl_sf_doublefact
gsl_sf_lnfact_e
gsl_sf_lnfact
gsl_sf_lndoublefact_e
gsl_sf_lndoublefact
gsl_sf_lnchoose_e
gsl_sf_lnchoose
gsl_sf_choose_e
gsl_sf_choose
gsl_sf_lnpoch_e
gsl_sf_lnpoch
gsl_sf_lnpoch_sgn_e
gsl_sf_poch_e
gsl_sf_poch
gsl_sf_pochrel_e
gsl_sf_pochrel
gsl_sf_lnbeta_e
gsl_sf_lnbeta
gsl_sf_lnbeta_sgn_e
gsl_sf_beta_e
gsl_sf_beta
gsl_sf_beta_inc_e
gsl_sf_beta_inc
gsl_sf_gegenpoly_1_e
gsl_sf_gegenpoly_2_e
gsl_sf_gegenpoly_3_e
gsl_sf_gegenpoly_1
gsl_sf_gegenpoly_2
gsl_sf_gegenpoly_3
gsl_sf_gegenpoly_n_e
gsl_sf_gegenpoly_n
gsl_sf_gegenpoly_array
gsl_sf_hyperg_0F1_e
gsl_sf_hyperg_0F1
gsl_sf_hyperg_1F1_int_e
gsl_sf_hyperg_1F1_int
gsl_sf_hyperg_1F1_e
gsl_sf_hyperg_1F1
gsl_sf_hyperg_U_int_e
gsl_sf_hyperg_U_int
gsl_sf_hyperg_U_int_e10_e
gsl_sf_hyperg_U_e
gsl_sf_hyperg_U
gsl_sf_hyperg_U_e10_e
gsl_sf_hyperg_2F1_e
gsl_sf_hyperg_2F1
gsl_sf_hyperg_2F1_conj_e
gsl_sf_hyperg_2F1_conj
gsl_sf_hyperg_2F1_renorm_e
gsl_sf_hyperg_2F1_renorm
gsl_sf_hyperg_2F1_conj_renorm_e
gsl_sf_hyperg_2F1_conj_renorm
gsl_sf_hyperg_2F0_e
gsl_sf_hyperg_2F0
gsl_sf_laguerre_1_e
gsl_sf_laguerre_2_e
gsl_sf_laguerre_3_e
gsl_sf_laguerre_1
gsl_sf_laguerre_2
gsl_sf_laguerre_3
gsl_sf_laguerre_n_e
gsl_sf_laguerre_n
gsl_sf_lambert_W0_e
gsl_sf_lambert_W0
gsl_sf_lambert_Wm1_e
gsl_sf_lambert_Wm1
gsl_sf_conicalP_half_e
gsl_sf_conicalP_half
gsl_sf_conicalP_mhalf_e
gsl_sf_conicalP_mhalf
gsl_sf_conicalP_0_e
gsl_sf_conicalP_0
gsl_sf_conicalP_1_e
gsl_sf_conicalP_1
gsl_sf_conicalP_sph_reg_e
gsl_sf_conicalP_sph_reg
gsl_sf_conicalP_cyl_reg_e
gsl_sf_conicalP_cyl_reg
gsl_sf_legendre_H3d_0_e
gsl_sf_legendre_H3d_0
gsl_sf_legendre_H3d_1_e
gsl_sf_legendre_H3d_1
gsl_sf_legendre_H3d_e
gsl_sf_legendre_H3d
gsl_sf_legendre_H3d_array
gsl_sf_log_e
gsl_sf_log
gsl_sf_log_abs_e
gsl_sf_log_abs
gsl_sf_complex_log_e
gsl_sf_log_1plusx_e
gsl_sf_log_1plusx
gsl_sf_log_1plusx_mx_e
gsl_sf_log_1plusx_mx
gsl_sf_mathieu_a_array
gsl_sf_mathieu_b_array
gsl_sf_mathieu_a
gsl_sf_mathieu_b
gsl_sf_mathieu_a_coeff
gsl_sf_mathieu_b_coeff
gsl_sf_mathieu_alloc
gsl_sf_mathieu_free
gsl_sf_mathieu_ce
gsl_sf_mathieu_se
gsl_sf_mathieu_ce_array
gsl_sf_mathieu_se_array
gsl_sf_mathieu_Mc
gsl_sf_mathieu_Ms
gsl_sf_mathieu_Mc_array
gsl_sf_mathieu_Ms_array
gsl_sf_pow_int_e
gsl_sf_pow_int
gsl_sf_psi_int_e
gsl_sf_psi_int
gsl_sf_psi_e
gsl_sf_psi
gsl_sf_psi_1piy_e
gsl_sf_psi_1piy
gsl_sf_complex_psi_e
gsl_sf_psi_1_int_e
gsl_sf_psi_1_int
gsl_sf_psi_1_e
gsl_sf_psi_1
gsl_sf_psi_n_e
gsl_sf_psi_n
gsl_sf_result_smash_e
gsl_sf_synchrotron_1_e
gsl_sf_synchrotron_1
gsl_sf_synchrotron_2_e
gsl_sf_synchrotron_2
gsl_sf_transport_2_e
gsl_sf_transport_2
gsl_sf_transport_3_e
gsl_sf_transport_3
gsl_sf_transport_4_e
gsl_sf_transport_4
gsl_sf_transport_5_e
gsl_sf_transport_5
gsl_sf_sin_e
gsl_sf_sin
gsl_sf_cos_e
gsl_sf_cos
gsl_sf_hypot_e
gsl_sf_hypot
gsl_sf_complex_sin_e
gsl_sf_complex_cos_e
gsl_sf_complex_logsin_e
gsl_sf_sinc_e
gsl_sf_sinc
gsl_sf_lnsinh_e
gsl_sf_lnsinh
gsl_sf_lncosh_e
gsl_sf_lncosh
gsl_sf_polar_to_rect
gsl_sf_rect_to_polar
gsl_sf_sin_err_e
gsl_sf_cos_err_e
gsl_sf_angle_restrict_symm_e
gsl_sf_angle_restrict_symm
gsl_sf_angle_restrict_pos_e
gsl_sf_angle_restrict_pos
gsl_sf_angle_restrict_symm_err_e
gsl_sf_angle_restrict_pos_err_e
-
gsl_sf_atanint_e
gsl_sf_atanint
-
These routines compute the Arctangent integral, which is defined as AtanInt(x) = \int_0^x dt \arctan(t)/t.
gsl_sf_zeta_int_e
gsl_sf_zeta_int
gsl_sf_zeta_e gsl_sf_zeta
gsl_sf_zetam1_e
gsl_sf_zetam1
gsl_sf_zetam1_int_e
gsl_sf_zetam1_int
gsl_sf_hzeta_e
gsl_sf_hzeta
gsl_sf_eta_int_e
gsl_sf_eta_int
gsl_sf_eta_e
gsl_sf_eta
This module also contains the following constants used as mode in various of those functions :
GSL_PREC_DOUBLE - Double-precision, a relative accuracy of approximately 2 * 10^-16.
GSL_PREC_SINGLE - Single-precision, a relative accuracy of approximately 10^-7.
GSL_PREC_APPROX - Approximate values, a relative accuracy of approximately 5 * 10^-4.
You can import the functions that you want to use by giving a space separated
list to Math::GSL::SF when you use the package. You can also write
use Math::GSL::SF qw/:all/
to use all available functions of the module. Note that
the tag names begin with a colon. Other tags are also available, here is a
complete list of all tags for this module :
airy
bessel
clausen
hydrogenic
coulumb
coupling
dawson
debye
dilog
factorial
misc
elliptic
error
hypergeometric
laguerre
legendre
gamma
transport
trig
zeta
eta
vars
For more informations on the functions, we refer you to the GSL official documentation: http://www.gnu.org/software/gsl/manual/html_node/
EXAMPLES
This example computes the dilogarithm of 1/10 :
use Math::GSL::SF qw/dilog/;
my $x = gsl_sf_dilog(0.1);
print "gsl_sf_dilog(0.1) = $x\n";
An example using Math::GSL::SF and gnuplot is in the examples/sf folder of the source code.
AUTHORS
Jonathan "Duke" Leto <jonathan@leto.net> and Thierry Moisan <thierry.moisan@gmail.com>
COPYRIGHT AND LICENSE
Copyright (C) 2008-2023 Jonathan "Duke" Leto and Thierry Moisan
This program is free software; you can redistribute it and/or modify it under the same terms as Perl itself.