An example code fragment using the ZOOM SpecialFunctions wrappers around the GNU Scientific Library (aka GSL) looks like
. . #include "SpecialFunctions/Bessel.h" . . // Compute J3(2.5) two different ways and compare double J3One = Bessel::Jn( 3, 2.5 ); Dpair J3Pair = Bessel::JnErr( 3, -2.5 ); double J3Two = J3Pair.first; double J3Err = J3Pair.second; if( abs(J3One-J3Two) > J3Err ) cerr << "Bad news, Bubba!" << endl; . .The non-native return types shown here and in the tables below are all defined through typedef declarations in the relevant header files. For the record, they are also given here,
- typedef std::vector<double> Dvector;
- typedef std::pair<double,double> Dpair;
- typedef std::complex<double> Dcomplex;
- typedef std::pair<Dcomplex,Dcomplex> Cpair;
The GSL library makes heavy use of a two-member struct gsl_sf_result such that the first member is the numerical result and the second is the estimated error for that result, when requested. We have uniformly replaced that usage in the following way. When no error estimate is requested, the return value is either a simple double or a Dcomplex, as appropriate. When an error estimate is requested, we return a Dpair or, in the case of a complex result, a Cpair where the first member of the pair is the result and the second is the estimated error.
Note that the tables do not show return types for any of the GSL functions. In standard C fashion, the GSL functions all return an integer status code with a zero value to indicate success and something non-zero to indicate one of a number of possible error conditions as described in the GSL Reference Manual. The SpecialFunctions wrappers intercept the GSL error returns and turn them into ZMExceptions of varying severity levels. The user can select any of the several ZMException handlers or supply one of his/her own to deal with these error returns in the same way as any other ZMexception and does not need to build up special code to deal with the specifics of the GSL design. Note also that, except for the Chebyshev class, all our wrappers have been implemented in the form of static methods. In fact, rather than define classes, we have used structs instead. Consequently, there are no constructors, no destructors and no member data. There is nothing to instantiate and the methods are uniformly invoked as shown in the tables without reference to any object.
With this list of correspondences, the user can refer to the existing GNU Scientific Library Reference Manual (a PostScript file in the doc subdirectory of the ZOOM SpecialFunctions package) for most of the operational details. In large measure, this relieves us of the need to duplicate existing documentation and relieves the user of having to cope with two manuals that inevitably get out of synch with each other.
It is only fair to emphasize that the whole GSL is much larger than the portion we have wrapped here. The goal of this package was limited to making the set of special functions included in the GSL available to the C++ programmer. As it stands now, this package is based on what the GSL developers call release 0.5++. It is clear that GSL is a "work-in-progress." For instance, there are functions in various implementation files with no corresponding prototype in the associated header file and vice versa. All such anomalies have been excluded from the SpecialFunctions package at the cost of a small loss of functionality.
Some portions of this library, notably Exponential, IntegerPowers and Trigonometry are redundant in that the standard library is supposed to provide all that. They are present in the GSL partly to cover those cases where the standard library is deficient and, as noted in the GSL manual, to provide uniform semantics across platforms. For the most part, we have retained all that.
| SpecialFunctions library | GSL library |
| double Airy::AiryAi( const double x ) | gsl_sf_airy_Ai_impl( x, GSL_MODE_DEFAULT, &_result ) |
| Dpair Airy::AiryAiErr( const double x ) | gsl_sf_airy_Ai_e( x, GSL_MODE_DEFAULT, &_result ) |
| double Airy::AiryBi( const double x ) | gsl_sf_airy_Bi_impl( x, GSL_MODE_DEFAULT, &_result ) |
| Dpair Airy::AiryBiErr( const double x ) | gsl_sf_airy_Bi_e( x, GSL_MODE_DEFAULT, &_result ) |
| double Airy::AiryAiScaled( const double x ) | gsl_sf_airy_Ai_scaled_impl( x, GSL_MODE_DEFAULT, &_result ) |
| Dpair Airy::AiryAiScaledErr( const double x ) | gsl_sf_airy_Ai_scaled_e( x, GSL_MODE_DEFAULT, &_result ) |
| double Airy::AiryBiScaled( const double x ) | gsl_sf_airy_Bi_scaled_impl( x, GSL_MODE_DEFAULT, &_result ) |
| Dpair Airy::AiryBiScaledErr( const double x ) | gsl_sf_airy_Bi_scaled_e( x, GSL_MODE_DEFAULT, &_result ) |
| double Airy::AiryAiPrime( const double x ) | gsl_sf_airy_Ai_deriv_impl( x, GSL_MODE_DEFAULT, &_result ) |
| Dpair Airy::AiryAiPrimeErr( const double x ) | gsl_sf_airy_Ai_deriv_e( x, GSL_MODE_DEFAULT, &_result ) |
| double Airy::AiryBiPrime( const double x ) | gsl_sf_airy_Bi_deriv_impl( x, GSL_MODE_DEFAULT, &_result ) |
| Dpair Airy::AiryBiPrimeErr( const double x ) | gsl_sf_airy_Bi_deriv_e( x, GSL_MODE_DEFAULT, &_result ) |
| double Airy::AiryAiPrimeScaled( const double x ) | gsl_sf_airy_Ai_deriv_scaled_impl( x, GSL_MODE_DEFAULT, &_result ) |
| Dpair Airy::AiryAiPrimeScaledErr( const double x ) | gsl_sf_airy_Ai_deriv_scaled_e( x, GSL_MODE_DEFAULT, &_result ) |
| double Airy::AiryBiPrimeScaled( const double x ) | gsl_sf_airy_Bi_deriv_scaled_impl( x, GSL_MODE_DEFAULT, &_result |
| Dpair Airy::AiryBiPrimeScaledErr( const double x ) | gsl_sf_airy_Bi_deriv_scaled_e( x, GSL_MODE_DEFAULT, &_result ) |
| double Airy::AiryAiZero( const int n ) | gsl_sf_airy_zero_Ai_impl( n, &_result ) |
| Dpair Airy::AiryAiZeroErr( const int n ) | gsl_sf_airy_zero_Ai_e( n, &_result ) |
| double Airy::AiryBiZero( const int n ) | gsl_sf_airy_zero_Bi_impl( n, &_result ) |
| Dpair Airy::AiryBiZeroErr( const int n ) | gsl_sf_airy_zero_Bi_e( n, &_result ) |
| double Airy::AiryAiPrimeZero( const int n ) | gsl_sf_airy_zero_Ai_deriv_impl( n, &_result ) |
| Dpair Airy::AiryAiPrimeZeroErr( const int n ) | gsl_sf_airy_zero_Ai_deriv_e( n, &_result ) |
| double Airy::AiryBiPrimeZero( const int n ) | gsl_sf_airy_zero_Bi_deriv_impl( n, &_result ) |
| Dpair Airy::AiryBiPrimeZeroErr( const int n ) | gsl_sf_airy_zero_Bi_deriv_e( n, &_result ) |
| SpecialFunctions library | GSL library |
| double Bessel::J0( const double x ) | gsl_sf_bessel_J0_impl( x, &_result ) |
| Dpair Bessel::J0Err( const double x ) | gsl_sf_bessel_J0_e( x, &_result ) |
| double Bessel::J1( const double x ) | gsl_sf_bessel_J1_impl( x, &_result ) |
| Dpair Bessel::J1Err( const double x ) | gsl_sf_bessel_J1_e( x, &_result ) |
| double Bessel::Jn( const int n, const double x ) | gsl_sf_bessel_Jn_impl( n, x, &_result ) |
| Dpair Bessel::JnErr( const int n, const double x ) | gsl_sf_bessel_Jn_e( n, x, &_result ) |
| Dvector Bessel::JnArray( const int nmin, const int nmax, const double x ) | gsl_sf_bessel_Jn_array_impl( nmin, nmax, x, temp ) |
| Dvector Bessel::JnErrArray( const int nmin, const int nmax, const double x ) | gsl_sf_bessel_Jn_array_e( nmin, nmax, x, temp ) |
| double Bessel::Y0( const double x ) | gsl_sf_bessel_Y0_impl( x, &_result ) |
| Dpair Bessel::Y0Err( const double x ) | gsl_sf_bessel_Y0_e( x, &_result ) |
| double Bessel::Y1( const double x ) | gsl_sf_bessel_Y1_impl( x, &_result ) |
| Dpair Bessel::Y1Err( const double x ) | gsl_sf_bessel_Y1_e( x, &_result ) |
| double Bessel::Yn( const int n, const double x ) | gsl_sf_bessel_Yn_impl( n, x, &_result ) |
| Dpair Bessel::YnErr( const int n, const double x ) | gsl_sf_bessel_Yn_e( n, x, &_result ) |
| Dvector Bessel::YnArray( const int nmin, const int nmax, const double x ) | gsl_sf_bessel_Yn_array_impl( nmin, nmax, x, temp ) |
| Dvector Bessel::YnErrArray( const int nmin, const int nmax, const double x ) | gsl_sf_bessel_Yn_array_e( nmin, nmax, x, temp ) |
| double Bessel::I0( const double x ) | gsl_sf_bessel_I0_impl( x, &_result ) |
| Dpair Bessel::I0Err( const double x ) | gsl_sf_bessel_I0_e( x, &_result ) |
| double Bessel::I1( const double x ) | gsl_sf_bessel_I1_impl( x, &_result ) |
| Dpair Bessel::I1Err( const double x ) | gsl_sf_bessel_I1_e( x, &_result ) |
| double Bessel::In( const int n, const double x ) | gsl_sf_bessel_In_impl( n, x, &_result ) |
| Dpair Bessel::InErr( const int n, const double x ) | gsl_sf_bessel_In_e( n, x, &_result ) |
| Dvector Bessel::InArray( const int nmin, const int nmax, const double x ) | gsl_sf_bessel_In_array_impl( nmin, nmax, x, temp ) |
| Dvector Bessel::InErrArray( const int nmin, const int nmax, const double x ) | gsl_sf_bessel_In_array_e( nmin, nmax, x, temp ) |
| double Bessel::I0Scaled( const double x ) | gsl_sf_bessel_I0_scaled_impl( x, &_result ) |
| Dpair Bessel::I0ErrScaled( const double x ) | gsl_sf_bessel_I0_scaled_e( x, &_result ) |
| double Bessel::I1Scaled( const double x ) | gsl_sf_bessel_I1_scaled_impl( x, &_result ) |
| Dpair Bessel::I1ErrScaled( const double x ) | gsl_sf_bessel_I1_scaled_e( x, &_result ) |
| double Bessel::InScaled( const int n, const double x ) | gsl_sf_bessel_In_scaled_impl( n, x, &_result ) |
| Dpair Bessel::InErrScaled( const int n, const double x ) | gsl_sf_bessel_In_scaled_e( n, x, &_result ) |
| Dvector Bessel::InArrayScaled( const int nmin, const int nmax, const double x ) | gsl_sf_bessel_In_scaled_array_impl( nmin, nmax, x, temp ) |
| Dvector Bessel::InErrArrayScaled( const int nmin, const int nmax, const double x ) | gsl_sf_bessel_In_scaled_array_e( nmin, nmax, x, temp ) |
| double Bessel::K0( const double x ) | gsl_sf_bessel_K0_impl( x, &_result ) |
| Dpair Bessel::K0Err( const double x ) | gsl_sf_bessel_K0_e( x, &_result ) |
| double Bessel::K1( const double x ) | gsl_sf_bessel_K1_impl( x, &_result ) |
| Dpair Bessel::K1Err( const double x ) | gsl_sf_bessel_K1_e( x, &_result ) |
| double Bessel::Kn( const int n, const double x ) | gsl_sf_bessel_Kn_impl( n, x, &_result ) |
| Dpair Bessel::KnErr( const int n, const double x ) | gsl_sf_bessel_Kn_e( n, x, &_result ) |
| Dvector Bessel::KnArray( const int nmin, const int nmax, const double x ) | gsl_sf_bessel_Kn_array_impl( nmin, nmax, x, temp ) |
| Dvector Bessel::KnErrArray( const int nmin, const int nmax, const double x ) | gsl_sf_bessel_Kn_array_e( nmin, nmax, x, temp ) |
| double Bessel::K0Scaled( const double x ) | gsl_sf_bessel_K0_scaled_impl( x, &_result ) |
| Dpair Bessel::K0ErrScaled( const double x ) | gsl_sf_bessel_K0_scaled_e( x, &_result ) |
| double Bessel::K1Scaled( const double x ) | gsl_sf_bessel_K1_scaled_impl( x, &_result ) |
| Dpair Bessel::K1ErrScaled( const double x ) | gsl_sf_bessel_K1_scaled_e( x, &_result ) |
| double Bessel::KnScaled( const int n, const double x ) | gsl_sf_bessel_Kn_scaled_impl( n, x, &_result ) |
| Dpair Bessel::KnErrScaled( const int n, const double x ) | gsl_sf_bessel_Kn_scaled_e( n, x, &_result ) |
| Dvector Bessel::KnArrayScaled( const int nmin, const int nmax, const double x ) | gsl_sf_bessel_Kn_scaled_array_impl( nmin, nmax, x, temp ) |
| Dvector Bessel::KnErrArrayScaled( const int nmin, const int nmax, const double x ) | gsl_sf_bessel_Kn_scaled_array_e( nmin, nmax, x, temp ) |
| double Bessel::RegularSphericalJ0( const double x ) | gsl_sf_bessel_j0_impl( x, &_result ) |
| Dpair Bessel::RegularSphericalJ0Err( const double x ) | gsl_sf_bessel_j0_e( x, &_result ) |
| double Bessel::RegularSphericalJ1( const double x ) | gsl_sf_bessel_j1_impl( x, &_result ) |
| Dpair Bessel::RegularSphericalJ1Err( const double x ) | gsl_sf_bessel_j1_e( x, &_result ) |
| double Bessel::RegularSphericalJ2( const double x ) | gsl_sf_bessel_j2_impl( x, &_result ) |
| Dpair Bessel::RegularSphericalJ2Err( const double x ) | gsl_sf_bessel_j2_e( x, &_result ) |
| double Bessel::RegularSphericalJn( const int n, const double x ) | gsl_sf_bessel_jl_impl( n, x, &_result ) |
| Dpair Bessel::RegularSphericalJnErr( const int n, const double x ) | gsl_sf_bessel_jl_e( n, x, &_result ) |
| Dvector Bessel::RegularSphericalJnArray( const int nmax, const double x ) | gsl_sf_bessel_jl_array_impl( nmax, x, temp ) |
| Dvector Bessel::RegularSphericalJnSteedArray( const int nmax, const double x ) | gsl_sf_bessel_jl_array_impl( nmax, x, temp ) |
| Dvector Bessel::RegularSphericalJnErrArray( const int nmax, const double x ) | gsl_sf_bessel_jl_array_e( nmax, x, temp ) |
| double Bessel::IrregularSphericalY0( const double x ) | gsl_sf_bessel_y0_impl( x, &_result ) |
| Dpair Bessel::IrregularSphericalY0Err( const double x ) | gsl_sf_bessel_y0_e( x, &_result ) |
| double Bessel::IrregularSphericalY1( const double x ) | gsl_sf_bessel_y1_impl( x, &_result ) |
| Dpair Bessel::IrregularSphericalY1Err( const double x ) | gsl_sf_bessel_y1_e( x, &_result ) |
| double Bessel::IrregularSphericalY2( const double x ) | gsl_sf_bessel_y2_impl( x, &_result ) |
| Dpair Bessel::IrregularSphericalY2Err( const double x ) | gsl_sf_bessel_y2_e( x, &_result ) |
| double Bessel::IrregularSphericalYn( const int n, const double x ) | gsl_sf_bessel_yl_impl( n, x, &_result ) |
| Dpair Bessel::IrregularSphericalYnErr( const int n, const double x ) | gsl_sf_bessel_yl_e( n, x, &_result ) |
| Dvector Bessel::IrregularSphericalYnArray( const int nmax, const double x ) | gsl_sf_bessel_yl_array_impl( nmax, x, temp ) |
| Dvector Bessel::IrregularSphericalYnErrArray( const int nmax, const double x ) | gsl_sf_bessel_yl_array_e( nmax, x, temp ) |
| double Bessel::RegularModifiedSphericalJ0( const double x ) | gsl_sf_bessel_i0_scaled_impl( x, &_result ) |
| Dpair Bessel::RegularModifiedSphericalJ0Err( const double x ) | gsl_sf_bessel_i0_scaled_e( x, &_result ) |
| double Bessel::RegularModifiedSphericalJ1( const double x ) | gsl_sf_bessel_i1_scaled_impl( x, &_result ) |
| Dpair Bessel::RegularModifiedSphericalJ1Err( const double x ) | gsl_sf_bessel_i1_scaled_e( x, &_result ) |
| double Bessel::RegularModifiedSphericalJ2( const double x ) | gsl_sf_bessel_i2_scaled_impl( x, &_result ) |
| Dpair Bessel::RegularModifiedSphericalJ2Err( const double x ) | gsl_sf_bessel_i2_scaled_e( x, &_result ) |
| double Bessel::RegularModifiedSphericalJn( const int n, const double x ) | gsl_sf_bessel_il_scaled_impl( n, x, &_result ) |
| Dpair Bessel::RegularModifiedSphericalJnErr( const int n, const double x ) | gsl_sf_bessel_il_scaled_e( n, x, &_result ) |
| Dvector Bessel::RegularModifiedSphericalJnArray( const int nmax, const double x ) | gsl_sf_bessel_il_scaled_array_impl( nmax, x, temp ) |
| Dvector Bessel::RegularModifiedSphericalJnErrArray( const int nmax, const double x ) | gsl_sf_bessel_il_scaled_array_e( nmax, x, temp ) |
| double Bessel::IrregularModifiedSphericalY0( const double x ) | gsl_sf_bessel_k0_scaled_impl( x, &_result ) |
| Dpair Bessel::IrregularModifiedSphericalY0Err( const double x ) | gsl_sf_bessel_k0_scaled_e( x, &_result ) |
| double Bessel::IrregularModifiedSphericalY1( const double x ) | gsl_sf_bessel_k1_scaled_impl( x, &_result ) |
| Dpair Bessel::IrregularModifiedSphericalY1Err( const double x ) | gsl_sf_bessel_k1_scaled_e( x, &_result ) |
| double Bessel::IrregularModifiedSphericalY2( const double x ) | gsl_sf_bessel_k2_scaled_impl( x, &_result ) |
| Dpair Bessel::IrregularModifiedSphericalY2Err( const double x ) | gsl_sf_bessel_k2_scaled_e( x, &_result ) |
| double Bessel::IrregularModifiedSphericalYn( const int n, const double x ) | gsl_sf_bessel_kl_scaled_impl( n, x, &_result ) |
| Dpair Bessel::IrregularModifiedSphericalYnErr( const int n, const double x ) | gsl_sf_bessel_kl_scaled_e( n, x, &_result ) |
| Dvector Bessel::IrregularModifiedSphericalYnArray( const int nmax, const double x ) | gsl_sf_bessel_kl_scaled_array_impl( nmax, x, temp ) |
| Dvector Bessel::IrregularModifiedSphericalYnErrArray( const int nmax, const double x ) | gsl_sf_bessel_kl_scaled_array_e( nmax, x, temp ) |
| double Bessel::Jnu( const double nu, const double x ) | gsl_sf_bessel_Jnu_impl( nu, x, &_result ) |
| Dpair Bessel::JnuErr( const double nu, const double x ) | gsl_sf_bessel_Jnu_e( nu, x, &_result ) |
| Dvector Bessel::JnuArray( const double nu, const Dvector x ) | gsl_sf_bessel_sequence_Jnu_impl( nu, GSL_MODE_DEFAULT, x.size(), temp ) |
| Dvector Bessel::JnuErrArray( const double nu, const Dvector x ) | gsl_sf_bessel_sequence_Jnu_e( nu, GSL_MODE_DEFAULT, x.size(), temp ) |
| double Bessel::Ynu( const double nu, const double x ) | gsl_sf_bessel_Ynu_impl( nu, x, &_result ) |
| Dpair Bessel::YnuErr( const double nu, const double x ) | gsl_sf_bessel_Ynu_e( nu, x, &_result ) |
| double Bessel::Inu( const double nu, const double x ) | gsl_sf_bessel_Inu_impl( nu, x, &_result ) |
| Dpair Bessel::InuErr( const double nu, const double x ) | gsl_sf_bessel_Inu_e( nu, x, &_result ) |
| double Bessel::InuScaled( const double nu, const double x ) | gsl_sf_bessel_Inu_scaled_impl( nu, x, &_result ) |
| Dpair Bessel::InuErrScaled( const double nu, const double x ) | gsl_sf_bessel_Inu_scaled_e( nu, x, &_result ) |
| double Bessel::Knu( const double nu, const double x ) | gsl_sf_bessel_Knu_impl( nu, x, &_result ) |
| Dpair Bessel::KnuErr( const double nu, const double x ) | gsl_sf_bessel_Knu_e( nu, x, &_result ) |
| double Bessel::lnKnu( const double nu, const double x ) | gsl_sf_bessel_lnKnu_impl( nu, x, &_result ) |
| Dpair Bessel::lnKnuErr( const double nu, const double x ) | gsl_sf_bessel_lnKnu_e( nu, x, &_result ) |
| double Bessel::KnuScaled( const double nu, const double x ) | gsl_sf_bessel_Knu_scaled_impl( nu, x, &_result ) |
| Dpair Bessel::KnuErrScaled( const double nu, const double x ) | gsl_sf_bessel_Knu_scaled_e( nu, x, &_result ) |
| double Bessel::J0Zero( const int n ) | gsl_sf_bessel_zero_J0_impl( n, &_result ) |
| Dpair Bessel::J0ZeroErr( const int n ) | gsl_sf_bessel_zero_J0_e( n, &_result ) |
| double Bessel::J1Zero( const int n ) | gsl_sf_bessel_zero_J1_impl( n, &_result ) |
| double Bessel::JnuZero( const double nu, const int n ) | gsl_sf_bessel_zero_Jnu_impl( nu, n, &_result ) |
| Dpair Bessel::JnuZeroErr( const double nu, const int n ) | gsl_sf_bessel_zero_Jnu_e( nu, n, &_result ) |
| SpecialFunctions library | GSL library |
| double ClausenIntegral::Clausen( const double x ) | gsl_sf_clausen_impl( x, &_result ) |
| Dpair ClausenIntegral::ClausenErr( const double x ) | gsl_sf_clausen_e( x, &_result ) |
| SpecialFunctions library | GSL library |
| double Coulomb::HydrogenicR1( const double Z, const double r ) | gsl_sf_hydrogenicR_1_impl( Z, r, &_result ) |
| Dpair Coulomb::HydrogenicR1Err( const double Z, const double r ) | gsl_sf_hydrogenicR_1_e( Z, r, &_result ) |
| double Coulomb::HydrogenicRn( const int n, const int l, const double Z, const double r ) | gsl_sf_hydrogenicR_impl( n, l, Z, r, &_result ) |
| Dpair Coulomb::HydrogenicRnErr( const int n, const int l, const double Z, const double r ) | gsl_sf_hydrogenicR_e( n, l, Z, r, &_result ) |
| bool Coulomb::CoulombWaveFG( const double eta, const double x, const double lambdaF, const int klambdaG, double* F, double* Fp, double* G, double* Gp, double* Fscale, double* Gscale ) | gsl_sf_coulomb_wave_FG_impl( eta, x, lambdaF, klambdaG, &_resultF, &_resultFp, &_resultG, &_resultGp, Fscale, Gscale ) |
| SpecialFunctions library | GSL library |
| double DawsonIntegral::Dawson( const double x ) | gsl_sf_dawson_impl( x, &_result ) |
| Dpair DawsonIntegral::DawsonErr( const double x ) | gsl_sf_dawson_e( x, &_result ) |
| SpecialFunctions library | GSL library |
| double DebyeIntegral::Debye1( const double x ) | gsl_sf_debye_1_impl( x, &_result ) |
| Dpair DebyeIntegral::Debye1Err( const double x ) | gsl_sf_debye_1_e( x, &_result ) |
| double DebyeIntegral::Debye2( const double x ) | gsl_sf_debye_2_impl( x, &_result ) |
| Dpair DebyeIntegral::Debye2Err( const double x ) | gsl_sf_debye_2_e( x, &_result ) |
| double DebyeIntegral::Debye3( const double x ) | gsl_sf_debye_3_impl( x, &_result ) |
| Dpair DebyeIntegral::Debye3Err( const double x ) | gsl_sf_debye_3_e( x, &_result ) |
| double DebyeIntegral::Debye4( const double x ) | gsl_sf_debye_4_impl( x, &_result ) |
| Dpair DebyeIntegral::Debye4Err( const double x ) | gsl_sf_debye_4_e( x, &_result ) |
| SpecialFunctions library | GSL library |
| double DiLogarithm::DiLogReal( const double x ) | gsl_sf_dilog_impl( x, &_result ) |
| Dpair DiLogarithm::DiLogRealErr( const double x ) | gsl_sf_dilog_e( x, &_result ) |
| Dcomplex DiLogarithm::DiLog( const Dcomplex z ) | gsl_sf_complex_dilog_impl( abs(z), arg(z), &_result, &_result2 ) |
| Cpair DiLogarithm::DiLogErr( const Dcomplex z ) | gsl_sf_complex_dilog_e( abs(z), arg(z), &_result, &_result2 ) |
| SpecialFunctions library | GSL library |
| double EllipticIntegrals::CompleteK( const double k ) | gsl_sf_ellint_Kcomp_impl( k, GSL_MODE_DEFAULT, &_result ) |
| Dpair EllipticIntegrals::CompleteKErr( const double k ) | gsl_sf_ellint_Kcomp_e( k, GSL_MODE_DEFAULT, &_result ) |
| double EllipticIntegrals::CompleteE( const double k ) | gsl_sf_ellint_Ecomp_impl( k, GSL_MODE_DEFAULT, &_result ) |
| Dpair EllipticIntegrals::CompleteEErr( const double k ) | gsl_sf_ellint_Ecomp_e( k, GSL_MODE_DEFAULT, &_result ) |
| double EllipticIntegrals::IncompleteF( const double phi, const double k ) | gsl_sf_ellint_F_impl( phi, k, GSL_MODE_DEFAULT, &_result ) |
| Dpair EllipticIntegrals::IncompleteFErr( const double phi, const double k ) | gsl_sf_ellint_F_e( phi, k, GSL_MODE_DEFAULT, &_result ) |
| double EllipticIntegrals::IncompleteE( const double phi, const double k ) | gsl_sf_ellint_E_impl( phi, k, GSL_MODE_DEFAULT, &_result ) |
| Dpair EllipticIntegrals::IncompleteEErr( const double phi, const double k ) | gsl_sf_ellint_E_e( phi, k, GSL_MODE_DEFAULT, &_result ) |
| double EllipticIntegrals::IncompleteP( const double phi, const double k, const double n ) | gsl_sf_ellint_P_impl( phi, k, n, GSL_MODE_DEFAULT, &_result ) |
| Dpair EllipticIntegrals::IncompletePErr( const double phi, const double k, const double n ) | gsl_sf_ellint_P_e( phi, k, n, GSL_MODE_DEFAULT, &_result ) |
| double EllipticIntegrals::IncompleteD( const double phi, const double k, const double n ) | gsl_sf_ellint_D_impl( phi, k, n, GSL_MODE_DEFAULT, &_result ) |
| Dpair EllipticIntegrals::IncompleteDErr( const double phi, const double k, const double n ) | gsl_sf_ellint_D_e( phi, k, n, GSL_MODE_DEFAULT, &_result ) |
| double EllipticIntegrals::CarlsonRC( const double x, const double y ) | gsl_sf_ellint_RC_impl( x, y, GSL_MODE_DEFAULT, &_result ) |
| Dpair EllipticIntegrals::CarlsonRCErr( const double x, const double y ) | gsl_sf_ellint_RC_e( x, y, GSL_MODE_DEFAULT, &_result ) |
| double EllipticIntegrals::CarlsonRD( const double x, const double y, const double z ) | gsl_sf_ellint_RD_impl( x, y, z, GSL_MODE_DEFAULT, &_result ) |
| Dpair EllipticIntegrals::CarlsonRDErr( const double x, const double y, const double z ) | gsl_sf_ellint_RD_e( x, y, z, GSL_MODE_DEFAULT, &_result ) |
| double EllipticIntegrals::CarlsonRF( const double x, const double y, const double z ) | gsl_sf_ellint_RF_impl( x, y, z, GSL_MODE_DEFAULT, &_result ) |
| Dpair EllipticIntegrals::CarlsonRFErr( const double x, const double y, const double z ) | gsl_sf_ellint_RF_e( x, y, z, GSL_MODE_DEFAULT, &_result ) |
| double EllipticIntegrals::CarlsonRJ( const double x, const double y, const double z, const double p ) | gsl_sf_ellint_RJ_impl( x, y, z, p, GSL_MODE_DEFAULT, &_result ) |
| Dpair EllipticIntegrals::CarlsonRJErr( const double x, const double y, const double z, const double p ) | gsl_sf_ellint_RJ_e( x, y, z, p, GSL_MODE_DEFAULT, &_result ) |
| SpecialFunctions library | GSL library |
| double ErrorFunction::Erfc( const double x ) | gsl_sf_erfc_impl( x, &_result ) |
| Dpair ErrorFunction::ErfcErr( const double x ) | gsl_sf_erfc_e( x, &_result ) |
| double ErrorFunction::lnErfc( const double x ) | gsl_sf_log_erfc_impl( x, &_result ) |
| Dpair ErrorFunction::lnErfcErr( const double x ) | gsl_sf_log_erfc_e( x, &_result ) |
| double ErrorFunction::Erf( const double x ) | gsl_sf_erf_impl( x, &_result ) |
| Dpair ErrorFunction::ErfErr( const double x ) | gsl_sf_erf_e( x, &_result ) |
| double ErrorFunction::ProbZ( const double x ) | gsl_sf_erf_Z_impl( x, &_result ) |
| Dpair ErrorFunction::ProbZErr( const double x ) | gsl_sf_erf_Z_e( x, &_result ) |
| double ErrorFunction::ProbQ( const double x ) | gsl_sf_erf_Q_impl( x, &_result ) |
| Dpair ErrorFunction::ProbQErr( const double x ) | gsl_sf_erf_Q_e( x, &_result ) |
| SpecialFunctions library | GSL library |
| double ExpIntegrals::E1( const double x ) | gsl_sf_expint_E1_impl( x, &_result ) | Dpair ExpIntegrals::E1Err( const double x ) | gsl_sf_expint_E1_e( x, &_result ) |
| double ExpIntegrals::E2( const double x ) | gsl_sf_expint_E2_impl( x, &_result ) |
| Dpair ExpIntegrals::E2Err( const double x ) | gsl_sf_expint_E2_e( x, &_result ) |
| double ExpIntegrals::Ei( const double x ) | gsl_sf_expint_Ei_impl( x, &_result ) |
| Dpair ExpIntegrals::EiErr( const double x ) | gsl_sf_expint_Ei_e( x, &_result ) |
| double ExpIntegrals::Shi( const double x ) | gsl_sf_Shi_impl( x, &_result ) |
| Dpair ExpIntegrals::ShiErr( const double x ) | gsl_sf_Shi_e( x, &_result ) |
| double ExpIntegrals::Chi( const double x ) | gsl_sf_Chi_impl( x, &_result ) |
| Dpair ExpIntegrals::ChiErr( const double x ) | gsl_sf_Chi_e( x, &_result ) |
| double ExpIntegrals::Ei3( const double x ) | gsl_sf_expint_3_impl( x, &_result ) |
| Dpair ExpIntegrals::Ei3Err( const double x ) | gsl_sf_expint_3_e( x, &_result ) |
| double ExpIntegrals::Si( const double x ) | gsl_sf_Si_impl( x, &_result ) |
| Dpair ExpIntegrals::SiErr( const double x ) | gsl_sf_Si_e( x, &_result ) |
| double ExpIntegrals::Ci( const double x ) | gsl_sf_Ci_impl( x, &_result ) |
| Dpair ExpIntegrals::CiErr( const double x ) | gsl_sf_Ci_e( x, &_result ) |
| double ExpIntegrals::AtanInt( const double x ) | gsl_sf_atanint_impl( x, &_result ) |
| Dpair ExpIntegrals::AtanIntErr( const double x ) | gsl_sf_atanint_e( x, &_result ) |
| SpecialFunctions library | GSL library |
| double Exponential::Exp( const double x ) | gsl_sf_exp_impl( x, &_result ) |
| Dpair Exponential::ExpErr( const double x ) | gsl_sf_exp_e( x, &_result ) |
| double Exponential::yExp( const double x, double y ) | gsl_sf_exp_mult_impl( x, y, &_result ) |
| Dpair Exponential::yExpErr( const double x, double y ) | gsl_sf_exp_mult_e( x, y, &_result ) |
| double Exponential::ExpM1( const double x ) | gsl_sf_expm1_impl( x, &_result ) |
| Dpair Exponential::ExpM1Err( const double x ) | gsl_sf_expm1_e( x, &_result ) |
| double Exponential::RelExp( const double x ) | gsl_sf_exprel_impl( x, &_result ) |
| Dpair Exponential::RelExpErr( const double x ) | gsl_sf_exprel_e( x, &_result ) |
| double Exponential::RelExp2( const double x ) | gsl_sf_exprel_2_impl( x, &_result ) |
| Dpair Exponential::RelExp2Err( const double x ) | gsl_sf_exprel_2_e( x, &_result ) |
| double Exponential::RelExpN( const int n, const double x ) | gsl_sf_exprel_n_impl( n, x, &_result ) |
| Dpair Exponential::RelExpNErr( const int n, const double x ) | gsl_sf_exprel_n_e( n, x, &_result ) |
| Dpair Exponential::ExpxDx( const double x, const double dx ) | gsl_sf_exp_err_impl( x, dx, &_result ) |
| Dpair Exponential::yExpxDxDy( const double x, const double dx, const double y, const double dy ) | gsl_sf_exp_mult_err_impl( x, dx, y, dy, &_result ) |
| SpecialFunctions library | GSL library |
| double FermiDirac::CompleteFDm1( const double x ) | gsl_sf_fermi_dirac_m1_impl( x, &_result ) |
| Dpair FermiDirac::CompleteFDm1Err( const double x ) | gsl_sf_fermi_dirac_m1_e( x, &_result ) |
| double FermiDirac::CompleteFD0( const double x ) | gsl_sf_fermi_dirac_0_impl( x, &_result ) |
| Dpair FermiDirac::CompleteFD0Err( const double x ) | gsl_sf_fermi_dirac_0_e( x, &_result ) |
| double FermiDirac::CompleteFD1( const double x ) | gsl_sf_fermi_dirac_1_impl( x, &_result ) |
| Dpair FermiDirac::CompleteFD1Err( const double x ) | gsl_sf_fermi_dirac_1_e( x, &_result ) |
| double FermiDirac::CompleteFD2( const double x ) | gsl_sf_fermi_dirac_2_impl( x, &_result ) |
| Dpair FermiDirac::CompleteFD2Err( const double x ) | gsl_sf_fermi_dirac_2_e( x, &_result ) |
| double FermiDirac::CompleteFDj( const int j, const double x ) | gsl_sf_fermi_dirac_int_impl( j, x, &_result ) |
| Dpair FermiDirac::CompleteFDjErr( const int j, const double x ) | gsl_sf_fermi_dirac_int_e( j, x, &_result ) |
| double FermiDirac::CompleteFDmHalf( const double x ) | gsl_sf_fermi_dirac_mhalf_impl( x, &_result ) |
| Dpair FermiDirac::CompleteFDmHalfErr( const double x ) | gsl_sf_fermi_dirac_mhalf_e( x, &_result ) |
| double FermiDirac::CompleteFDHalf( const double x ) | gsl_sf_fermi_dirac_half_impl( x, &_result ) |
| Dpair FermiDirac::CompleteFDHalfErr( const double x ) | gsl_sf_fermi_dirac_half_e( x, &_result ) |
| double FermiDirac::CompleteFD3Half( const double x ) | gsl_sf_fermi_dirac_3half_impl( x, &_result ) |
| Dpair FermiDirac::CompleteFD3HalfErr( const double x ) | gsl_sf_fermi_dirac_3half_e( x, &_result ) |
| double FermiDirac::IncompleteFD0( const double x, const double b ) | gsl_sf_fermi_dirac_inc_0_impl( x, b, &_result ) |
| Dpair FermiDirac::IncompleteFD0Err( const double x, const double b ) | gsl_sf_fermi_dirac_inc_0_e( x, b, &_result ) |
| SpecialFunctions library | GSL library |
| double Gamma::lnGammaAbs( const double x ) | gsl_sf_lngamma_impl( x, &_result ) |
| Dpair Gamma::lnGammaAbsErr( const double x ) | gsl_sf_lngamma_e( x, &_result ) |
| Dpair Gamma::lnGamma( const double x ) | gsl_sf_lngamma_sgn_impl( x, &_result, &sgn ) |
| Dpair Gamma::lnGammaErr( const double x ) | gsl_sf_lngamma_sgn_e( x, &_result, &sgn ) |
| double Gamma::GammaFunction( const double x ) | gsl_sf_gamma_impl( x, &_result ) |
| Dpair Gamma::GammaFunctionErr( const double x ) | gsl_sf_gamma_e( x, &_result ) |
| double Gamma::GammaStar( const double x ) | gsl_sf_gammastar_impl( x, &_result ) |
| Dpair Gamma::GammaStarErr( const double x ) | gsl_sf_gammastar_e( x, &_result ) |
| double Gamma::ReciprocalGamma( const double x ) | gsl_sf_gammainv_impl( x, &_result ) |
| Dpair Gamma::ReciprocalGammaErr( const double x ) | gsl_sf_gammainv_e( x, &_result ) |
| Dpair Gamma::lnGammaZ( const Dcomplex z ) | gsl_sf_lngamma_complex_impl( real(z), imag(z), &_result, &_result2 ) |
| Dpair Gamma::lnGammaZErr( const Dcomplex z ) | gsl_sf_lngamma_complex_e( real(z), imag(z), &_result, &_result2 ) |
| double Gamma::TaylorN( const int n, const double x ) | gsl_sf_taylorcoeff_impl( n, x, &_result ) |
| Dpair Gamma::TaylorNErr( const int n, const double x ) | gsl_sf_taylorcoeff_e( n, x, &_result ) |
| double Gamma::Factorial( const unsigned int n ) | gsl_sf_fact_impl( n, &_result ) |
| Dpair Gamma::FactorialErr( const unsigned int n ) | gsl_sf_fact_e( n, &_result ) |
| double Gamma::DoubleFactorial( const unsigned int n ) | gsl_sf_doublefact_impl( n, &_result ) |
| Dpair Gamma::DoubleFactorialErr( const unsigned int n ) | gsl_sf_doublefact_e( n, &_result ) |
| double Gamma::lnFactorial( const unsigned int n ) | gsl_sf_lnfact_impl( n, &_result ) |
| Dpair Gamma::lnFactorialErr( const unsigned int n ) | gsl_sf_lnfact_e( n, &_result ) |
| double Gamma::lnDoubleFactorial( const unsigned int n ) | gsl_sf_lndoublefact_impl( n, &_result ) |
| Dpair Gamma::lnDoubleFactorialErr( const unsigned int n ) | gsl_sf_lndoublefact_e( n, &_result ) |
| double Gamma::lnBinomialCoefficient( const unsigned int n, const unsigned int m ) | gsl_sf_lnchoose_impl( n, m, &_result ) |
| Dpair Gamma::lnBinomialCoefficientErr( const unsigned int n, const unsigned int m ) | gsl_sf_lnchoose_e( n, m, &_result ) |
| double Gamma::BinomialCoefficient( const unsigned int n, const unsigned int m ) | gsl_sf_choose_impl( n, m, &_result ) |
| Dpair Gamma::BinomialCoefficientErr( const unsigned int n, const unsigned int m ) | gsl_sf_choose_e( n, m, &_result ) |
| double Gamma::lnPochammer( const double a, const double x ) | gsl_sf_lnpoch_impl( a, x, &_result ) |
| Dpair Gamma::lnPochammerErr( const double a, const double x ) | gsl_sf_lnpoch_e( a, x, &_result ) |
| Dpair Gamma::lnPochammerSign( const double a, const double x ) | gsl_sf_lnpoch_sgn_impl( a, x, &_result, &sgn ) |
| Dpair Gamma::lnPochammerSignErr( const double a, const double x ) | gsl_sf_lnpoch_sgn_e( a, x, &_result, &sgn ) |
| double Gamma::Pochammer( const double a, const double x ) | gsl_sf_poch_impl( a, x, &_result ) |
| Dpair Gamma::PochammerErr( const double a, const double x ) | gsl_sf_poch_e( a, x, &_result ) |
| double Gamma::RelativePochammer( const double a, const double x ) | gsl_sf_pochrel_impl( a, x, &_result ) |
| Dpair Gamma::RelativePochammerErr( const double a, const double x ) | gsl_sf_pochrel_e( a, x, &_result ) |
| double Gamma::IncompleteGamma( const double a, const double x ) | gsl_sf_gamma_inc_Q_impl( a, x, &_result ) |
| Dpair Gamma::IncompleteGammaErr( const double a, const double x ) | gsl_sf_gamma_inc_Q_e( a, x, &_result ) |
| double Gamma::CompIncompleteGamma( const double a, const double x ) | gsl_sf_gamma_inc_P_impl( a, x, &_result ) |
| Dpair Gamma::CompIncompleteGammaErr( const double a, const double x ) | gsl_sf_gamma_inc_P_e( a, x, &_result ) |
| double Gamma::lnBeta( const double a, const double b ) | gsl_sf_lnbeta_impl( a, b, &_result ) |
| Dpair Gamma::lnBetaErr( const double a, const double b ) | gsl_sf_lnbeta_e( a, b, &_result ) |
| double Gamma::Beta( const double a, const double b ) | gsl_sf_beta_impl( a, b, &_result ) |
| Dpair Gamma::BetaErr( const double a, const double b ) | gsl_sf_beta_e( a, b, &_result ) |
| SpecialFunctions library | GSL library |
| double Gegenbauer::GegenPoly1( const double lambda, const double x ) | gsl_sf_gegenpoly_1_impl(lambda, x, &_result ) |
| Dpair Gegenbauer::GegenPoly1Err( const double lambda, const double x ) | gsl_sf_gegenpoly_1_e(lambda, x, &_result ) |
| double Gegenbauer::GegenPoly2( const double lambda, const double x ) | gsl_sf_gegenpoly_2_impl(lambda, x, &_result ) |
| Dpair Gegenbauer::GegenPoly2Err( const double lambda, const double x ) | gsl_sf_gegenpoly_2_e(lambda, x, &_result ) |
| double Gegenbauer::GegenPoly3( const double lambda, const double x ) | gsl_sf_gegenpoly_3_impl(lambda, x, &_result ) |
| Dpair Gegenbauer::GegenPoly3Err( const double lambda, const double x ) | gsl_sf_gegenpoly_3_e(lambda, x, &_result ) |
| double Gegenbauer::GegenPolyN(const int n, const double lambda, const double x ) | gsl_sf_gegenpoly_n_impl( n, lambda, x, &_result ) |
| Dpair Gegenbauer::GegenPolyNErr( const int n, const double lambda, const double x ) | gsl_sf_gegenpoly_n_e( n, lambda, x, &_result ) |
| Dvector Gegenbauer::GegenPolyArray( const int nmax, const double lambda, const double x ) | gsl_sf_gegenpoly_array_impl( nmax, lambda, x, temp ) |
| Dvector Gegenbauer::GegenPolyArrayErr( const int nmax, const double lambda, const double x ) | gsl_sf_gegenpoly_array_e( nmax, lambda, x, temp ) |
| SpecialFunctions library | GSL library |
| double HyperGeometric::HyperG0F1( const double c, const double x ) | gsl_sf_hyperg_0F1_impl( c, x, &_result ) |
| Dpair HyperGeometric::HyperG0F1Err( const double c, const double x ) | gsl_sf_hyperg_0F1_e( c, x, &_result ) |
| double HyperGeometric::HyperG1F1( const int m, const int n, const double x ) | gsl_sf_hyperg_1F1_int_impl( m, n, x, &_result ) |
| Dpair HyperGeometric::HyperG1F1Err( const int m, const int n, const double x ) | gsl_sf_hyperg_1F1_int_e( m, n, x, &_result ) |
| double HyperGeometric::HyperG1F1( const double a, const double b, const double x ) | gsl_sf_hyperg_1F1_impl( a, b, x, &_result ) |
| Dpair HyperGeometric::HyperG1F1Err( const double a, const double b, const double x ) | gsl_sf_hyperg_1F1_e( a, b, x, &_result ) |
| double HyperGeometric::HyperGU( const int m, const int n, const double x ) | gsl_sf_hyperg_U_int_impl( m, n, x, &_result ) |
| Dpair HyperGeometric::HyperGUErr( const int m, const int n, const double x ) | gsl_sf_hyperg_U_int_e( m, n, x, &_result ) |
| double HyperGeometric::HyperGU( const double a, const double b, const double x ) | gsl_sf_hyperg_U_impl( a, b, x, &_result ) |
| Dpair HyperGeometric::HyperGUErr( const double a, const double b, const double x ) | gsl_sf_hyperg_U_e( a, b, x, &_result ) |
| double HyperGeometric::HyperGauss( const double a, const double b, const double c, const double x ) | gsl_sf_hyperg_2F1_impl( a, b, c, x, &_result ) |
| Dpair HyperGeometric::HyperGaussErr( const double a, const double b, const double c, const double x ) | gsl_sf_hyperg_2F1_e( a, b, c, x, &_result ) |
| double HyperGeometric::HyperGaussConjugate( const Dcomplex a, const double c, const double x ) | gsl_sf_hyperg_2F1_conj_impl( real(a), imag(a), c, x, &_result ) |
| Dpair HyperGeometric::HyperGaussConjugateErr( const Dcomplex a, const double c, const double x ) | gsl_sf_hyperg_2F1_conj_e( real(a), imag(a), c, x, &_result ) |
| double HyperGeometric::NormedHyperGauss( const double a, const double b, const double c, const double x ) | gsl_sf_hyperg_2F1_renorm_impl( a, b, c, x, &_result ) |
| Dpair HyperGeometric::NormedHyperGaussErr( const double a, const double b, const double c, const double x ) | gsl_sf_hyperg_2F1_renorm_e( a, b, c, x, &_result ) |
| double HyperGeometric::NormedHyperGaussConjugate( const Dcomplex a, const double c, const double x ) | gsl_sf_hyperg_2F1_conj_renorm_impl( real(a), imag(a), c, x, &_result ) |
| Dpair HyperGeometric::NormedHyperGaussConjugateErr( const Dcomplex a, const double c, const double x ) | gsl_sf_hyperg_2F1_conj_renorm_e( real(a), imag(a), c, x, &_result ) |
| double HyperGeometric::HyperG2F0( const double a, const double b, const double x ) | gsl_sf_hyperg_2F0_impl( a, b, x, &_result ) |
| Dpair HyperGeometric::HyperG2F0Err( const double a, const double b, const double x ) | gsl_sf_hyperg_2F0_e( a, b, x, &_result ) |
| SpecialFunctions library | GSL library |
| double IntegerPowers::Pow2( const double x ) | gsl_sf_pow_2( x ) |
| double IntegerPowers::Pow3( const double x ) | gsl_sf_pow_3( x ) |
| double IntegerPowers::Pow4( const double x ) | gsl_sf_pow_4( x ) |
| double IntegerPowers::Pow5( const double x ) | gsl_sf_pow_5( x ) |
| double IntegerPowers::Pow6( const double x ) | gsl_sf_pow_6( x ) |
| double IntegerPowers::Pow7( const double x ) | gsl_sf_pow_7( x ) |
| double IntegerPowers::Pow8( const double x ) | gsl_sf_pow_8( x ) |
| double IntegerPowers::Pow9( const double x ) | gsl_sf_pow_9( x ) |
| double IntegerPowers::PowN( const double x, const int n ) | gsl_sf_pow_int_impl( x, n, &_result ) |
| Dpair IntegerPowers::PowNErr( const double x, const int n ) | gsl_sf_pow_int_e( x, n, &_result ) |
| SpecialFunctions library | GSL library | |
| void JacobiElliptics::Jacobi( const double u, const double m, double* sn, double* cn, double* dn ) | gsl_sf_elljac_impl( u, m, sn, cn, dn ) | |
| void JacobiElliptics::JacobiErr( const double u, const double m, double* sn, double* cn, double* dn ) | gsl_sf_elljac_impl( u, m, sn, cn, dn ) |
| SpecialFunctions library | GSL library |
| double Laguerre::Laguerre1( const double a, const double x ) | gsl_sf_laguerre_1_impl( a, x, &_result ) |
| Dpair Laguerre::Laguerre1Err( const double a, const double x ) | gsl_sf_laguerre_1_e( a, x, &_result ) |
| double Laguerre::Laguerre2( const double a, const double x ) | gsl_sf_laguerre_2_impl( a, x, &_result ) |
| Dpair Laguerre::Laguerre2Err( const double a, const double x ) | gsl_sf_laguerre_2_e( a, x, &_result ) |
| double Laguerre::Laguerre3( const double a, const double x ) | gsl_sf_laguerre_3_impl( a, x, &_result ) |
| Dpair Laguerre::Laguerre3Err( const double a, const double x ) | gsl_sf_laguerre_3_e( a, x, &_result ) |
| double Laguerre::LaguerreN(const int n, const double a, const double x ) | gsl_sf_laguerre_n_impl( n, a, x, &_result ) |
| Dpair Laguerre::LaguerreNErr( const int n, const double a, const double x ) | gsl_sf_laguerre_n_e( n, a, x, &_result ) |
| SpecialFunctions library | GSL library |
| double Legendre::P1( const double x ) | gsl_sf_legendre_P1_impl( x, &_result ) |
| Dpair Legendre::P1Err( const double x ) | gsl_sf_legendre_P1_e( x, &_result ) |
| double Legendre::P2( const double x ) | gsl_sf_legendre_P2_impl( x, &_result ) |
| Dpair Legendre::P2Err( const double x ) | gsl_sf_legendre_P2_e( x, &_result ) |
| double Legendre::P3( const double x ) | gsl_sf_legendre_P3_impl( x, &_result ) |
| Dpair Legendre::P3Err( const double x ) | gsl_sf_legendre_P3_e( x, &_result ) |
| double Legendre::Pl( const int l, const double x ) | gsl_sf_legendre_Pl_impl( l, x, &_result ) |
| Dpair Legendre::PlErr( const int l, const double x ) | gsl_sf_legendre_Pl_e( l, x, &_result ) |
| Dvector Legendre::PlArray( const int lmax, const double x ) | gsl_sf_legendre_Pl_array_impl( lmax, x, temp ) |
| Dvector Legendre::PlErrArray( const int lmax, const double x ) | gsl_sf_legendre_Pl_array_e( lmax, x, temp ) |
| double Legendre::Q0( const double x ) | gsl_sf_legendre_Q0_impl( x, &_result ) |
| Dpair Legendre::Q0Err( const double x ) | gsl_sf_legendre_Q0_e( x, &_result ) |
| double Legendre::Q1( const double x ) | gsl_sf_legendre_Q1_impl( x, &_result ) |
| Dpair Legendre::Q1Err( const double x ) | gsl_sf_legendre_Q1_e( x, &_result ) |
| double Legendre::Ql( const int l, const double x ) | gsl_sf_legendre_Ql_impl( l, x, &_result ) |
| Dpair Legendre::QlErr( const int l, const double x ) | gsl_sf_legendre_Ql_e( l, x, &_result ) |
| double Legendre::Plm( const int l, const int m, const double x ) | gsl_sf_legendre_Plm_impl( l, m, x, &_result ) |
| Dpair Legendre::PlmErr( const int l, const int m, const double x ) | gsl_sf_legendre_Plm_e( l, m, x, &_result ) |
| Dvector Legendre::PlmArray( const int lmax, const int m, const double x ) | gsl_sf_legendre_Plm_array_impl( lmax, m, x, temp ) |
| double Legendre::SphericalPlm( const int l, const int m, const double x ) | gsl_sf_legendre_sphPlm_impl( l, m, x, &_result ) |
| Dpair Legendre::SphericalPlmErr( const int l, const int m, const double x ) | gsl_sf_legendre_sphPlm_e( l, m, x, &_result ) |
| Dvector Legendre::SphericalPlmArray( const int lmax, const int m, const double x ) | gsl_sf_legendre_sphPlm_array_impl( lmax, m, x, temp ) |
| Dvector Legendre::SphericalPlmErrArray( const int lmax, const int m, const double x ) | gsl_sf_legendre_sphPlm_array_e( lmax, m, x, temp ) |
| double Legendre::RegularSphericalConical( const double lambda, const double x ) | gsl_sf_conicalP_mhalf_impl( lambda, x, &_result ) |
| Dpair Legendre::RegularSphericalConicalErr( const double lambda, const double x ) | gsl_sf_conicalP_mhalf_e( lambda, x, &_result ) |
| double Legendre::IrregularSphericalConical( const double lambda, const double x ) | gsl_sf_conicalP_half_impl( lambda, x, &_result ) |
| Dpair Legendre::IrregularSphericalConicalErr( const double lambda, const double x ) | gsl_sf_conicalP_half_e( lambda, x, &_result ) |
| double Legendre::ConicalP0( const double lambda, const double x ) | gsl_sf_conicalP_0_impl( lambda, x, &_result ) |
| Dpair Legendre::ConicalP0Err( const double lambda, const double x ) | gsl_sf_conicalP_0_e( lambda, x, &_result ) |
| double Legendre::ConicalP1( const double lambda, const double x ) | gsl_sf_conicalP_1_impl( lambda, x, &_result ) |
| Dpair Legendre::ConicalP1Err( const double lambda, const double x ) | gsl_sf_conicalP_1_e( lambda, x, &_result ) |
| double Legendre::RegularSphericalConicalPl( const int l, const double lambda, const double x ) | gsl_sf_conicalP_sph_reg_impl( l, lambda, x, &_result ) |
| Dpair Legendre::RegularSphericalConicalPlErr( const int l, const double lambda, const double x ) | gsl_sf_conicalP_sph_reg_e( l, lambda, x, &_result ) |
| double Legendre::RegularCylindricalConicalPl( const int l, const double lambda, const double x ) | gsl_sf_conicalP_cyl_reg_impl( l, lambda, x, &_result ) |
| Dpair Legendre::RegularCylindricalConicalPlErr( const int l, const double lambda, const double x ) | gsl_sf_conicalP_cyl_reg_e( l, lambda, x, &_result ) |
| double Legendre::H3d0( const double lambda, const double eta ) | gsl_sf_legendre_H3d_0_impl( lambda, eta, &_result ) |
| Dpair Legendre::H3d0Err( const double lambda, const double eta ) | gsl_sf_legendre_H3d_0_e( lambda, eta, &_result ) |
| double Legendre::H3d1( const double lambda, const double eta ) | gsl_sf_legendre_H3d_1_impl( lambda, eta, &_result ) |
| Dpair Legendre::H3d1Err( const double lambda, const double eta ) | gsl_sf_legendre_H3d_1_e( lambda, eta, &_result ) |
| double Legendre::H3dl( const int l, const double lambda, const double eta ) | gsl_sf_legendre_H3d_impl( l, lambda, eta, &_result ) |
| Dpair Legendre::H3dlErr( const int l, const double lambda, const double eta ) | gsl_sf_legendre_H3d_e( l, lambda, eta, &_result ) |
| Dvector Legendre::H3dlArray( const int lmax, const double lambda, const double eta ) | gsl_sf_legendre_H3d_array_impl( lmax, lambda, eta, temp ) |
| Dvector Legendre::H3dlErrArray( const int lmax, const double lambda, const double eta ) | gsl_sf_legendre_H3d_array_e( lmax, lambda, eta, temp ) |
| SpecialFunctions library | GSL library |
| double Logarithm::Log( const double x ) | gsl_sf_log_impl( x, &_result ) |
| Dpair Logarithm::LogErr( const double x ) | gsl_sf_log_e( x, &_result ) |
| double Logarithm::LogAbs( const double x ) | gsl_sf_log_abs_impl( x, &_result ) |
| Dpair Logarithm::LogAbsErr( const double x ) | gsl_sf_log_abs_e( x, &_result ) |
| Dcomplex Logarithm::Log( const Dcomplex z ) | gsl_sf_complex_log_impl( real(z), imag(z), &_result, &_result2 ) |
| Cpair Logarithm::LogErr( const Dcomplex z ) | gsl_sf_complex_log_e( real(z), imag(z), &_result, &_result2 ) |
| double Logarithm::Log1p( const double x ) | gsl_sf_log_1plusx_impl( x, &_result ) |
| Dpair Logarithm::Log1pErr( const double x ) | gsl_sf_log_1plusx_e( x, &_result ) |
| double Logarithm::Log1pm( const double x ) | gsl_sf_log_1plusx_mx_impl( x, &_result ) |
| Dpair Logarithm::Log1pmErr( const double x ) | gsl_sf_log_1plusx_mx_e( x, &_result ) |
| SpecialFunctions library | GSL library |
| double PolyGamma::DiGamma( const int n ) | gsl_sf_psi_int_impl( n, &_result ) |
| Dpair PolyGamma::DiGammaErr( const int n ) | gsl_sf_psi_e( n, &_result ) |
| double PolyGamma::DiGamma( const double x ) | gsl_sf_psi_impl( x, &_result ) |
| Dpair PolyGamma::DiGammaErr( const double x ) | gsl_sf_psi_e( x, &_result ) |
| double PolyGamma::DiGammaI( const double x ) | gsl_sf_psi_1piy_impl( x, &_result ) |
| Dpair PolyGamma::DiGammaIErr( const double x ) | gsl_sf_psi_1piy_e( x, &_result ) |
| double PolyGamma::TriGamma( const int n ) | gsl_sf_psi_1_int_impl( n, &_result ) |
| Dpair PolyGamma::TriGammaErr( const int n ) | gsl_sf_psi_1_int_e( n, &_result ) |
| double PolyGamma::PolyGammaN( const int n, const double x ) | gsl_sf_psi_n_impl( n, x, &_result ) |
| Dpair PolyGamma::PolyGammaNErr( const int n, const double x ) | gsl_sf_psi_n_e( n, x, &_result ) |
| SpecialFunctions library | GSL library |
| double Synchrotron::Synchrotron1( const double x ) | gsl_sf_synchrotron_1_impl( x, &_result ) |
| Dpair Synchrotron::Synchrotron1Err( const double x ) | gsl_sf_synchrotron_1_e( x, &_result ) |
| double Synchrotron::Synchrotron2( const double x ) | gsl_sf_synchrotron_2_impl( x, &_result ) |
| Dpair Synchrotron::Synchrotron2Err( const double x ) | gsl_sf_synchrotron_2_e( x, &_result ) |
| SpecialFunctions library | GSL library |
| double Transport::Transport2( const double x ) | gsl_sf_transport_2_impl( x, &_result ) |
| Dpair Transport::Transport2Err( const double x ) | gsl_sf_transport_2_e( x, &_result ) |
| double Transport::Transport3( const double x ) | gsl_sf_transport_3_impl( x, &_result ) |
| Dpair Transport::Transport3Err( const double x ) | gsl_sf_transport_3_e( x, &_result ) |
| double Transport::Transport4( const double x ) | gsl_sf_transport_4_impl( x, &_result ) |
| Dpair Transport::Transport4Err( const double x ) | gsl_sf_transport_4_e( x, &_result ) |
| double Transport::Transport5( const double x ) | gsl_sf_transport_5_impl( x, &_result ) |
| Dpair Transport::Transport5Err( const double x ) | gsl_sf_transport_5_e( x, &_result ) |
| SpecialFunctions library | GSL library |
| double Trigonometry::Sin( const double x ) | gsl_sf_sin_impl( x, &_result ) |
| Dpair Trigonometry::SinErr( const double x ) | gsl_sf_sin_e( x, &_result ) |
| double Trigonometry::Cos( const double x ) | gsl_sf_cos_impl( x, &_result ) |
| Dpair Trigonometry::CosErr( const double x ) | gsl_sf_cos_e( x, &_result ) |
| double Trigonometry::Hypot( const double x, double y ) | gsl_sf_hypot_impl( x, y, &_result ) |
| Dpair Trigonometry::HypotErr( const double x, double y ) | gsl_sf_hypot_e( x, y, &_result ) |
| Dcomplex Trigonometry::Sin( const Dcomplex z ) | gsl_sf_complex_sin_impl( real(z), imag(z), &_result, &_result2 ) |
| Cpair Trigonometry::SinErr( const Dcomplex z ) | gsl_sf_complex_sin_e( real(z), imag(z), &_result, &_result2 ) |
| Dcomplex Trigonometry::Cos( const Dcomplex z ) | gsl_sf_complex_cos_impl( real(z), imag(z), &_result, &_result2 ) |
| Cpair Trigonometry::CosErr( const Dcomplex z ) | gsl_sf_complex_cos_e( real(z), imag(z), &_result, &_result2 ) |
| Dcomplex Trigonometry::lnSin( const Dcomplex z ) | gsl_sf_complex_logsin_impl( real(z), imag(z), &_result, &_result2 ) |
| Cpair Trigonometry::lnSinErr( const Dcomplex z ) | gsl_sf_complex_logsin_e( real(z), imag(z), &_result, &_result2 ) |
| double Trigonometry::Sinc( const double x ) | gsl_sf_sinc_impl( x, &_result ) |
| Dpair Trigonometry::SincErr( const double x ) | gsl_sf_sinc_e( x, &_result ) |
| double Trigonometry::lnSinh( const double x ) | gsl_sf_lnsinh_impl( x, &_result ) |
| Dpair Trigonometry::lnSinhErr( const double x ) | gsl_sf_lnsinh_e( x, &_result ) |
| double Trigonometry::lnCosh( const double x ) | gsl_sf_lncosh_impl( x, &_result ) |
| Dpair Trigonometry::lnCoshErr( const double x ) | gsl_sf_lncosh_e( x, &_result ) |
| Dpair Trigonometry::PolarToRect( Dpair rtheta ) | gsl_sf_polar_to_rect_impl( rtheta.first, rtheta.second, &_result, &_result2 ) |
| Dpair Trigonometry::PolarToRectErr( Dpair rtheta, double& dx, double& dy ) | gsl_sf_polar_to_rect_e( rtheta.first, rtheta.second, &_result, &_result2 ) |
| Dpair Trigonometry::RectToPolar( Dpair xy ) | gsl_sf_rect_to_polar_impl( xy.first, xy.second, &_result, &_result2 ) |
| Dpair Trigonometry::RectToPolarErr( Dpair xy, double& dr, double& dtheta ) | gsl_sf_rect_to_polar_e( xy.first, xy.second, &_result, &_result2 ) |
| Dpair Trigonometry::Sin( const double x, const double dx ) | gsl_sf_sin_err_impl( x, dx, &_result ) |
| Dpair Trigonometry::Cos( const double x, const double dx ) | gsl_sf_cos_err_impl( x, dx, &_result ) |
| double Trigonometry::StandardizeSymm( const double theta ) | gsl_sf_angle_restrict_symm_impl( &th ) |
| Dpair Trigonometry::StandardizeSymmErr( const double theta ) | gsl_sf_angle_restrict_symm_err_impl( theta, &_result ) |
| double Trigonometry::StandardizePos( const double theta ) | gsl_sf_angle_restrict_pos_impl( &th ) |
| Dpair Trigonometry::StandardizePosErr( const double theta ) | gsl_sf_angle_restrict_pos_err_impl( theta, &_result ) |
| SpecialFunctions library | GSL library |
| double Wigner::ThreeJ( const int two_ja, const int two_jb, const int two_jc, const int two_ma, const int two_mb, const int two_mc ) | gsl_sf_coupling_3j_impl( two_ja, two_jb, two_jc, two_ma, two_mb, two_mc, &_result ) |
| Dpair Wigner::ThreeJErr( const int two_ja, const int two_jb, const int two_jc, const int two_ma, const int two_mb, const int two_mc ) | gsl_sf_coupling_3j_e( two_ja, two_jb, two_jc, two_ma, two_mb, two_mc, &_result ) |
| double Wigner::SixJ( const int two_ja, const int two_jb, const int two_jc, const int two_jd, const int two_je, const int two_jf ) | gsl_sf_coupling_6j_impl( two_ja, two_jb, two_jc, two_jd, two_je, two_jf, &_result ) |
| Dpair Wigner::SixJErr( const int two_ja, const int two_jb, const int two_jc, const int two_jd, const int two_je, const int two_jf ) | gsl_sf_coupling_6j_e( two_ja, two_jb, two_jc, two_jd, two_je, two_jf, &_result ) |
| double Wigner::NineJ( const int two_ja, const int two_jb, const int two_jc, const int two_jd, const int two_je, const int two_jf, const int two_jg, const int two_jh, const int two_ji ) | gsl_sf_coupling_9j_impl( two_ja, two_jb, two_jc, two_jd, two_je, two_jf, two_jg, two_jh, two_ji, &_result ) |
| Dpair Wigner::NineJErr( const int two_ja, const int two_jb, const int two_jc, const int two_jd, const int two_je, const int two_jf, const int two_jg, const int two_jh, const int two_ji ) | gsl_sf_coupling_9j_impl( two_ja, two_jb, two_jc, two_jd, two_je, two_jf, two_jg, two_jh, two_ji, &_result ) |
| SpecialFunctions library | GSL library |
| double Zeta::RZeta( const int n ) | gsl_sf_zeta_int_impl( n, &_result ) |
| Dpair Zeta::RZetaErr( const int n ) | gsl_sf_zeta_int_e( n, &_result ) |
| double Zeta::RZeta( const double s ) | gsl_sf_zeta_impl( s, &_result ) |
| Dpair Zeta::RZetaErr( const double s ) | gsl_sf_zeta_e( s, &_result ) |
| double Zeta::HZeta( const double s, const double q ) | gsl_sf_hzeta_impl( s, q, &_result ) |
| Dpair Zeta::HZetaErr( const double s, const double q ) | gsl_sf_hzeta_e( s, q, &_result ) |
| double Zeta::Eta( const int n ) | gsl_sf_eta_int_impl( n, &_result ) |
| Dpair Zeta::EtaErr( const int n ) | gsl_sf_eta_int_e( n, &_result ) |
| double Zeta::Eta( const double s ) | gsl_sf_eta_impl( s, &_result ) |
| Dpair Zeta::EtaErr( const double s ) | gsl_sf_eta_e( s, &_result ) |