f15.3.10 {hypergeo}R Documentation

Transformations of the hypergeometric function

Description

Transformations of the hypergeometric function detailed in AMS-55, page 559-560.

Usage

f15.3.10       (A, B,    z, tol = 0, maxiter = 2000, method = "a")
f15.3.10_a     (A, B,    z, tol = 0, maxiter = 2000              )
f15.3.10_b     (A, B,    z, tol = 0, maxiter = 2000              )
f15.3.11       (A, B, m, z, tol = 0, maxiter = 2000, method = "a")
f15.3.11_bit1  (A, B, m, z, tol = 0                              )
f15.3.11_bit2_a(A, B, m, z, tol = 0, maxiter = 2000              )
f15.3.11_bit2_b(A, B, m, z, tol = 0, maxiter = 2000              )
f15.3.12       (A, B, m, z, tol = 0, maxiter = 2000, method = "a")
f15.3.12_bit1  (A, B, m, z, tol = 0                              )
f15.3.12_bit2_a(A, B, m, z, tol = 0, maxiter = 2000              )
f15.3.12_bit2_b(A, B, m, z, tol = 0, maxiter = 2000              )
f15.3.13       (A, C,    z, tol = 0, maxiter = 2000, method = "a")
f15.3.13_a     (A, C,    z, tol = 0, maxiter = 2000              )
f15.3.13_b     (A, C,    z, tol = 0, maxiter = 2000              )
f15.3.14       (A, C, m, z, tol = 0, maxiter = 2000, method = "a")
f15.3.14_bit1_a(A, C, m, z, tol = 0, maxiter = 2000              )
f15.3.14_bit1_b(A, C, m, z, tol = 0, maxiter = 2000              )
f15.3.14_bit2  (A, C, m, z, tol = 0                              )
f15.3.13_14    (A, C, m, z, tol = 0, maxiter = 2000, method = "a")
f15.3.10_11_12 (A, B, m, z, tol = 0, maxiter = 2000, method = "a")
f15.1.1        (A, B, C, z, tol = 0, maxiter = 2000              )

Arguments

A, B, C

Parameters of the hypergeometric function

m

Integer linking A, B, C as set out in AMS-55, page 559,560

z

primary complex argument

tol, maxiter

numerical parameters

method

Length 1 character vector specifying the method. See details

Details

Naming scheme (functions and arguments) follows AMS-55, pages 559-560.

The method argument to (eg) f15.3.14() specifies whether to use psigamma() directly (method “a”), or the recurrence 6.3.5 (method “b”). Press et al recommend method “b”, presumably on the grounds of execution speed. I'm not so sure (method “a” seems to be more accurate in the sense that it returns values closer to those of Maple).

Method “c” means to use a totally dull, slow, direct (but clearly correct) summation, for the purposes of debugging. This is only used for the functions documented under wolfram.Rd

Functions f15.3.13_14() and f15.3.10_11_12() are convenience wrappers. For example, function f15.3.13_14() dispatches to either f15.3.13() or f15.3.14() depending on the value of m.

Note

These functions are not really designed to be called by the user: use hypergeo() instead, or hypergeo_cover[123]() for specific cases.

Author(s)

Robin K. S. Hankin

References

M. Abramowitz and I. A. Stegun 1965. Handbook of mathematical functions. New York: Dover

See Also

hypergeo,wolfram,hypergeo_cover1

Examples


f15.3.10_11_12(A=1.1, B=pi, m= +3, z=.1+.1i)
f15.3.10_11_12(A=1.1, B=pi, m= -3, z=.1+.1i)


[Package hypergeo version 1.2-13 Index]