Special Functions¶
AUTHORS:
- David Joyner (2006-13-06): initial version
- David Joyner (2006-30-10): bug fixes to pari wrappers of Bessel functions, hypergeometric_U
- William Stein (2008-02): Impose some sanity checks.
- David Joyner (2008-04-23): addition of elliptic integrals
This module provides easy access to many of Maxima and PARI’s special functions.
Maxima’s special functions package (which includes spherical harmonic functions, spherical Bessel functions (of the 1st and 2nd kind), and spherical Hankel functions (of the 1st and 2nd kind)) was written by Barton Willis of the University of Nebraska at Kearney. It is released under the terms of the General Public License (GPL).
Support for elliptic functions and integrals was written by Raymond Toy. It is placed under the terms of the General Public License (GPL) that governs the distribution of Maxima.
Next, we summarize some of the properties of the functions implemented here.
Spherical harmonics: Laplace’s equation in spherical coordinates is:
Note that the spherical coordinates
and
are defined here as follows:
is the colatitude or polar angle, ranging from
and
the azimuth or longitude, ranging from
.
The general solution which remains finite towards infinity is a linear combination of functions of the form
and
where
are the associated Legendre polynomials, and with integer parameters
and
from
to
. Put in another way, the solutions with integer parameters
and
, can be written as linear combinations of:
where the functions
are the spherical harmonic functions with parameters
,
, which can be written as:
The spherical harmonics obey the normalisation condition
When solving for separable solutions of Laplace’s equation in spherical coordinates, the radial equation has the form:
The spherical Bessel functions
and
, are two linearly independent solutions to this equation. They are related to the ordinary Bessel functions
and
by:
For
, the confluent hypergeometric function
is defined to be the solution to Kummer’s differential equation
which satisfies
, as
. (There is a linearly independent solution, called Kummer’s function
, which is not implemented.)
The incomplete elliptic integrals (of the first kind, etc.) are:
and the complete ones are obtained by taking
.
REFERENCES:
- Abramowitz and Stegun: Handbook of Mathematical Functions, http://www.math.sfu.ca/~cbm/aands/
- http://en.wikipedia.org/wiki/Spherical_harmonics
- http://en.wikipedia.org/wiki/Helmholtz_equation
- Online Encyclopedia of Special Function http://algo.inria.fr/esf/index.html
TODO: Resolve weird bug in commented out code in hypergeometric_U below.
AUTHORS:
- David Joyner and William Stein
Added 16-02-2008 (wdj): optional calls to scipy and replace all ‘#random’ by ‘...’ (both at the request of William Stein)
Warning
SciPy’s versions are poorly documented and seem less accurate than the Maxima and PARI versions; typically they are limited by hardware floats precision.
-
class
sage.functions.special.
EllipticE
¶ Bases:
sage.functions.special.MaximaFunction
This returns the value of the “incomplete elliptic integral of the second kind,”
i.e.,
integrate(sqrt(1 - m*sin(x)^2), x, 0, phi)
. Takinggives
elliptic_ec
.EXAMPLES:
sage: z = var("z") sage: # this is still wrong: must be abs(sin(z)) + 2*round(z/pi) sage: elliptic_e(z, 1) 2*round(z/pi) + sin(z) sage: elliptic_e(z, 0) z sage: elliptic_e(0.5, 0.1) # abs tol 2e-15 0.498011394498832
-
class
sage.functions.special.
EllipticEC
¶ Bases:
sage.functions.special.MaximaFunction
This returns the value of the “complete elliptic integral of the second kind,”
EXAMPLES:
sage: elliptic_ec(0.1) 1.53075763689776 sage: elliptic_ec(x).diff() 1/2*(elliptic_ec(x) - elliptic_kc(x))/x sage: loads(dumps(elliptic_ec)) elliptic_ec
-
class
sage.functions.special.
EllipticEU
¶ Bases:
sage.functions.special.MaximaFunction
Return the value of the “incomplete elliptic integral of the second kind,”
where
.
EXAMPLES:
sage: elliptic_eu (0.5, 0.1) 0.496054551286597
-
class
sage.functions.special.
EllipticF
¶ Bases:
sage.functions.special.MaximaFunction
This returns the value of the “incomplete elliptic integral of the first kind,”
i.e.,
integrate(1/sqrt(1 - m*sin(x)^2), x, 0, phi)
. Takinggives
elliptic_kc
.EXAMPLES:
sage: z = var("z") sage: elliptic_f (z, 0) z sage: elliptic_f (z, 1) log(tan(1/4*pi + 1/2*z)) sage: elliptic_f (0.2, 0.1) 0.200132506747543
-
class
sage.functions.special.
EllipticKC
¶ Bases:
sage.functions.special.MaximaFunction
This returns the value of the “complete elliptic integral of the first kind,”
EXAMPLES:
sage: elliptic_kc(0.5) 1.85407467730137 sage: elliptic_f(RR(pi/2), 0.5) 1.85407467730137
-
class
sage.functions.special.
EllipticPi
¶ Bases:
sage.functions.special.MaximaFunction
This returns the value of the “incomplete elliptic integral of the third kind,”
INPUT:
n
– a real number, called the “characteristic”t
– a real number, called the “amplitude”m
– a real number, called the “parameter”
EXAMPLES:
sage: N(elliptic_pi(1, pi/4, 1)) 1.14779357469632
Compare the value computed by Maxima to the definition as a definite integral (using GSL):
sage: elliptic_pi(0.1, 0.2, 0.3) 0.200665068220979 sage: numerical_integral(1/(1-0.1*sin(x)^2)/sqrt(1-0.3*sin(x)^2), 0.0, 0.2) (0.2006650682209791, 2.227829789769088e-15)
ALGORITHM:
Numerical evaluation and symbolic manipulation are provided by Maxima.
REFERENCES:
- Abramowitz and Stegun: Handbook of Mathematical Functions, section 17.7 http://www.math.sfu.ca/~cbm/aands/
- Elliptic Functions in Maxima
-
class
sage.functions.special.
MaximaFunction
(name, nargs=2, conversions={})¶ Bases:
sage.symbolic.function.BuiltinFunction
EXAMPLES:
sage: from sage.functions.special import MaximaFunction sage: f = MaximaFunction("jacobi_sn") sage: f(1,1) tanh(1) sage: f(1/2,1/2).n() 0.470750473655657
-
class
sage.functions.special.
SphericalHarmonic
¶ Bases:
sage.symbolic.function.BuiltinFunction
Returns the spherical harmonic function
.
For integers
,
, simplification is done automatically. Numeric evaluation is supported for complex
and
.
EXAMPLES:
sage: x, y = var('x, y') sage: spherical_harmonic(3, 2, x, y) 15/4*sqrt(7/30)*cos(x)*e^(2*I*y)*sin(x)^2/sqrt(pi) sage: spherical_harmonic(3, 2, 1, 2) 15/4*sqrt(7/30)*cos(1)*e^(4*I)*sin(1)^2/sqrt(pi) sage: spherical_harmonic(3 + I, 2., 1, 2) -0.351154337307488 - 0.415562233975369*I sage: latex(spherical_harmonic(3, 2, x, y, hold=True)) Y_{3}^{2}\left(x, y\right) sage: spherical_harmonic(1, 2, x, y) 0
-
sage.functions.special.
elliptic_j
(z)¶ Returns the elliptic modular
-function evaluated at
.
INPUT:
z
(complex) – a complex number with positive imaginary part.
OUTPUT:
(complex) The value of
.
ALGORITHM:
Calls the
pari
functionellj()
.AUTHOR:
John Cremona
EXAMPLES:
sage: elliptic_j(CC(i)) 1728.00000000000 sage: elliptic_j(sqrt(-2.0)) 8000.00000000000 sage: z = ComplexField(100)(1,sqrt(11))/2 sage: elliptic_j(z) -32768.000... sage: elliptic_j(z).real().round() -32768 :: sage: tau = (1 + sqrt(-163))/2 sage: (-elliptic_j(tau.n(100)).real().round())^(1/3) 640320
-
sage.functions.special.
error_fcn
(t)¶ The complementary error function
(t belongs to RR). This function is currently always evaluated immediately.
EXAMPLES:
sage: error_fcn(6) 2.15197367124989e-17 sage: error_fcn(RealField(100)(1/2)) 0.47950012218695346231725334611
Note this is literally equal to
:
sage: 1 - error_fcn(0.5) 0.520499877813047 sage: erf(0.5) 0.520499877813047
-
sage.functions.special.
hypergeometric_U
(alpha, beta, x, algorithm='pari', prec=53)¶ Default is a wrap of PARI’s hyperu(alpha,beta,x) function. Optionally, algorithm = “scipy” can be used.
The confluent hypergeometric function
is defined to be the solution to Kummer’s differential equation
This satisfies
, as
, and is sometimes denoted
x^{-a}2_F_0(a,1+a-b,-1/x)
. This is not the same as Kummer’s-hypergeometric function, denoted sometimes as
_1F_1(alpha,beta,x)
, though it satisfies the same DE thatdoes.
Warning
In the literature, both are called “Kummer confluent hypergeometric” functions.
EXAMPLES:
sage: hypergeometric_U(1,1,1,"scipy") 0.596347362323... sage: hypergeometric_U(1,1,1) 0.59634736232319... sage: hypergeometric_U(1,1,1,"pari",70) 0.59634736232319407434...
-
sage.functions.special.
maxima_function
(name)¶ Returns a function which is evaluated both symbolically and numerically via Maxima. In particular, it returns an instance of
MaximaFunction
.Note
This function is cached so that duplicate copies of the same function are not created.
EXAMPLES:
sage: spherical_hankel2(2,i) -e
-
sage.functions.special.
meval
(x)¶ Return
x
evaluated in Maxima, then returned to Sage.This is used to evaluate several of these special functions.
TEST:
sage: from sage.functions.special import spherical_bessel_J sage: spherical_bessel_J(2.,3.) # rel tol 1e-10 0.2986374970757335
-
sage.functions.special.
spherical_bessel_J
(n, var, algorithm='maxima')¶ Returns the spherical Bessel function of the first kind for integers n >= 1.
Reference: AS 10.1.8 page 437 and AS 10.1.15 page 439.
EXAMPLES:
sage: spherical_bessel_J(2,x) ((3/x^2 - 1)*sin(x) - 3*cos(x)/x)/x sage: spherical_bessel_J(1, 5.2, algorithm='scipy') -0.12277149950007... sage: spherical_bessel_J(1, 3, algorithm='scipy') 0.345677499762355...
-
sage.functions.special.
spherical_bessel_Y
(n, var, algorithm='maxima')¶ Returns the spherical Bessel function of the second kind for integers n -1.
Reference: AS 10.1.9 page 437 and AS 10.1.15 page 439.
EXAMPLES:
sage: x = PolynomialRing(QQ, 'x').gen() sage: spherical_bessel_Y(2,x) -((3/x^2 - 1)*cos(x) + 3*sin(x)/x)/x
-
sage.functions.special.
spherical_hankel1
(n, var)¶ Returns the spherical Hankel function of the first kind for integers
, written as a string. Reference: AS 10.1.36 page 439.
EXAMPLES:
sage: spherical_hankel1(2, x) (I*x^2 - 3*x - 3*I)*e^(I*x)/x^3
-
sage.functions.special.
spherical_hankel2
(n, x)¶ Returns the spherical Hankel function of the second kind for integers
, written as a string. Reference: AS 10.1.17 page 439.
EXAMPLES:
sage: spherical_hankel2(2, x) (-I*x^2 - 3*x + 3*I)*e^(-I*x)/x^3
Here I = sqrt(-1).