Book a Demo!
CoCalc Logo Icon
StoreFeaturesDocsShareSupportNewsAboutPoliciesSign UpSign In
sagemath
GitHub Repository: sagemath/sagesmc
Path: blob/master/src/sage/functions/special.py
8815 views
1
r"""
2
Special Functions
3
4
AUTHORS:
5
6
- David Joyner (2006-13-06): initial version
7
8
- David Joyner (2006-30-10): bug fixes to pari wrappers of Bessel
9
functions, hypergeometric_U
10
11
- William Stein (2008-02): Impose some sanity checks.
12
13
- David Joyner (2008-04-23): addition of elliptic integrals
14
15
This module provides easy access to many of Maxima and PARI's
16
special functions.
17
18
Maxima's special functions package (which includes spherical
19
harmonic functions, spherical Bessel functions (of the 1st and 2nd
20
kind), and spherical Hankel functions (of the 1st and 2nd kind))
21
was written by Barton Willis of the University of Nebraska at
22
Kearney. It is released under the terms of the General Public
23
License (GPL).
24
25
Support for elliptic functions and integrals was written by Raymond
26
Toy. It is placed under the terms of the General Public License
27
(GPL) that governs the distribution of Maxima.
28
29
Next, we summarize some of the properties of the functions
30
implemented here.
31
32
33
- Airy function The function `Ai(x)` and the related
34
function `Bi(x)`, which is also called an Airy function,
35
are solutions to the differential equation
36
37
38
.. math::
39
40
y'' - xy = 0,
41
42
known as the Airy equation. They belong to the class of 'Bessel functions of
43
fractional order'. The initial conditions
44
`Ai(0) = (\Gamma(2/3)3^{2/3})^{-1}`,
45
`Ai'(0) = -(\Gamma(1/3)3^{1/3})^{-1}` define
46
`Ai(x)`. The initial conditions
47
`Bi(0) = 3^{1/2}Ai(0)`, `Bi'(0) = -3^{1/2}Ai'(0)`
48
define `Bi(x)`.
49
50
They are named after the British astronomer George Biddell Airy.
51
52
- Spherical harmonics: Laplace's equation in spherical coordinates
53
is:
54
55
.. math::
56
57
{\frac{1}{r^2}}{\frac{\partial}{\partial r}} \left(r^2 {\frac{\partial f}{\partial r}}\right) + {\frac{1}{r^2}\sin\theta}{\frac{\partial}{\partial \theta}} \left(\sin\theta {\frac{\partial f}{\partial \theta}}\right) + {\frac{1}{r^2\sin^2\theta}}{\frac{\partial^2 f}{\partial \varphi^2}} = 0.
58
59
60
Note that the spherical coordinates `\theta` and
61
`\varphi` are defined here as follows: `\theta` is
62
the colatitude or polar angle, ranging from
63
`0\leq\theta\leq\pi` and `\varphi` the azimuth or
64
longitude, ranging from `0\leq\varphi<2\pi`.
65
66
The general solution which remains finite towards infinity is a
67
linear combination of functions of the form
68
69
.. math::
70
71
r^{-1-\ell} \cos (m \varphi) P_\ell^m (\cos{\theta} )
72
73
74
and
75
76
.. math::
77
78
r^{-1-\ell} \sin (m \varphi) P_\ell^m (\cos{\theta} )
79
80
81
where `P_\ell^m` are the associated Legendre polynomials,
82
and with integer parameters `\ell \ge 0` and `m`
83
from `0` to `\ell`. Put in another way, the
84
solutions with integer parameters `\ell \ge 0` and
85
`- \ell\leq m\leq \ell`, can be written as linear
86
combinations of:
87
88
.. math::
89
90
U_{\ell,m}(r,\theta , \varphi ) = r^{-1-\ell} Y_\ell^m( \theta , \varphi )
91
92
93
where the functions `Y` are the spherical harmonic
94
functions with parameters `\ell`, `m`, which can be
95
written as:
96
97
.. math::
98
99
Y_\ell^m( \theta , \varphi ) = \sqrt{{\frac{(2\ell+1)}{4\pi}}{\frac{(\ell-m)!}{(\ell+m)!}}} \cdot e^{i m \varphi } \cdot P_\ell^m ( \cos{\theta} ) .
100
101
102
103
The spherical harmonics obey the normalisation condition
104
105
106
.. math::
107
108
\int_{\theta=0}^\pi\int_{\varphi=0}^{2\pi} Y_\ell^mY_{\ell'}^{m'*}\,d\Omega =\delta_{\ell\ell'}\delta_{mm'}\quad\quad d\Omega =\sin\theta\,d\varphi\,d\theta .
109
110
111
112
- When solving for separable solutions of Laplace's equation in
113
spherical coordinates, the radial equation has the form:
114
115
.. math::
116
117
x^2 \frac{d^2 y}{dx^2} + 2x \frac{dy}{dx} + [x^2 - n(n+1)]y = 0.
118
119
120
The spherical Bessel functions `j_n` and `y_n`,
121
are two linearly independent solutions to this equation. They are
122
related to the ordinary Bessel functions `J_n` and
123
`Y_n` by:
124
125
.. math::
126
127
j_n(x) = \sqrt{\frac{\pi}{2x}} J_{n+1/2}(x),
128
129
130
131
.. math::
132
133
y_n(x) = \sqrt{\frac{\pi}{2x}} Y_{n+1/2}(x) = (-1)^{n+1} \sqrt{\frac{\pi}{2x}} J_{-n-1/2}(x).
134
135
136
137
- For `x>0`, the confluent hypergeometric function
138
`y = U(a,b,x)` is defined to be the solution to Kummer's
139
differential equation
140
141
142
.. math::
143
144
xy'' + (b-x)y' - ay = 0,
145
146
which satisfies `U(a,b,x) \sim x^{-a}`, as
147
`x\rightarrow \infty`. (There is a linearly independent
148
solution, called Kummer's function `M(a,b,x)`, which is not
149
implemented.)
150
151
- Jacobi elliptic functions can be thought of as generalizations
152
of both ordinary and hyperbolic trig functions. There are twelve
153
Jacobian elliptic functions. Each of the twelve corresponds to an
154
arrow drawn from one corner of a rectangle to another.
155
156
::
157
158
n ------------------- d
159
| |
160
| |
161
| |
162
s ------------------- c
163
164
Each of the corners of the rectangle are labeled, by convention, s,
165
c, d and n. The rectangle is understood to be lying on the complex
166
plane, so that s is at the origin, c is on the real axis, and n is
167
on the imaginary axis. The twelve Jacobian elliptic functions are
168
then pq(x), where p and q are one of the letters s,c,d,n.
169
170
The Jacobian elliptic functions are then the unique
171
doubly-periodic, meromorphic functions satisfying the following
172
three properties:
173
174
175
#. There is a simple zero at the corner p, and a simple pole at the
176
corner q.
177
178
#. The step from p to q is equal to half the period of the function
179
pq(x); that is, the function pq(x) is periodic in the direction pq,
180
with the period being twice the distance from p to q. Also, pq(x)
181
is also periodic in the other two directions as well, with a period
182
such that the distance from p to one of the other corners is a
183
quarter period.
184
185
#. If the function pq(x) is expanded in terms of x at one of the
186
corners, the leading term in the expansion has a coefficient of 1.
187
In other words, the leading term of the expansion of pq(x) at the
188
corner p is x; the leading term of the expansion at the corner q is
189
1/x, and the leading term of an expansion at the other two corners
190
is 1.
191
192
193
We can write
194
195
.. math::
196
197
pq(x)=\frac{pr(x)}{qr(x)}
198
199
200
where `p`, `q`, and `r` are any of the
201
letters `s`, `c`, `d`, `n`, with
202
the understanding that `ss=cc=dd=nn=1`.
203
204
Let
205
206
.. math::
207
208
u=\int_0^\phi \frac{d\theta} {\sqrt {1-m \sin^2 \theta}}
209
210
211
Then the *Jacobi elliptic function* `sn(u)` is given by
212
213
.. math::
214
215
{sn}\; u = \sin \phi
216
217
and `cn(u)` is given by
218
219
.. math::
220
221
{cn}\; u = \cos \phi
222
223
and
224
225
.. math::
226
227
{dn}\; u = \sqrt {1-m\sin^2 \phi}.
228
229
To emphasize the dependence on `m`, one can write
230
`sn(u,m)` for example (and similarly for `cn` and
231
`dn`). This is the notation used below.
232
233
For a given `k` with `0 < k < 1` they therefore are
234
solutions to the following nonlinear ordinary differential
235
equations:
236
237
238
- `\mathrm{sn}\,(x;k)` solves the differential equations
239
240
.. math::
241
242
\frac{\mathrm{d}^2 y}{\mathrm{d}x^2} + (1+k^2) y - 2 k^2 y^3 = 0,
243
244
and
245
246
.. math::
247
248
\left(\frac{\mathrm{d} y}{\mathrm{d}x}\right)^2 = (1-y^2) (1-k^2 y^2).
249
250
- `\mathrm{cn}\,(x;k)` solves the differential equations
251
252
253
.. math::
254
255
\frac{\mathrm{d}^2 y}{\mathrm{d}x^2} + (1-2k^2) y + 2 k^2 y^3 = 0,
256
257
and `\left(\frac{\mathrm{d} y}{\mathrm{d}x}\right)^2 = (1-y^2) (1-k^2 + k^2 y^2)`.
258
259
- `\mathrm{dn}\,(x;k)` solves the differential equations
260
261
.. math::
262
263
\frac{\mathrm{d}^2 y}{\mathrm{d}x^2} - (2 - k^2) y + 2 y^3 = 0,
264
265
and `\left(\frac{\mathrm{d} y}{\mathrm{d}x}\right)^2= y^2 (1 - k^2 - y^2)`.
266
267
If `K(m)` denotes the complete elliptic integral of the
268
first kind (denoted ``elliptic_kc``), the elliptic functions
269
`sn (x,m)` and `cn (x,m)` have real periods
270
`4K(m)`, whereas `dn (x,m)` has a period
271
`2K(m)`. The limit `m\rightarrow 0` gives
272
`K(0) = \pi/2` and trigonometric functions:
273
`sn(x, 0) = \sin x`, `cn(x, 0) = \cos x`,
274
`dn(x, 0) = 1`. The limit `m \rightarrow 1` gives
275
`K(1) \rightarrow \infty` and hyperbolic functions:
276
`sn(x, 1) = \tanh x`,
277
`cn(x, 1) = \mbox{\rm sech} x`,
278
`dn(x, 1) = \mbox{\rm sech} x`.
279
280
- The incomplete elliptic integrals (of the first kind, etc.) are:
281
282
.. math::
283
284
\begin{array}{c} \displaystyle\int_0^\phi \frac{1}{\sqrt{1 - m\sin(x)^2}}\, dx,\\ \displaystyle\int_0^\phi \sqrt{1 - m\sin(x)^2}\, dx,\\ \displaystyle\int_0^\phi \frac{\sqrt{1-mt^2}}{\sqrt(1 - t^2)}\, dx,\\ \displaystyle\int_0^\phi \frac{1}{\sqrt{1 - m\sin(x)^2\sqrt{1 - n\sin(x)^2}}}\, dx, \end{array}
285
286
and the complete ones are obtained by taking `\phi =\pi/2`.
287
288
289
REFERENCES:
290
291
- Abramowitz and Stegun: Handbook of Mathematical Functions,
292
http://www.math.sfu.ca/~cbm/aands/
293
294
- http://en.wikipedia.org/wiki/Airy_function
295
296
- http://en.wikipedia.org/wiki/Spherical_harmonics
297
298
- http://en.wikipedia.org/wiki/Helmholtz_equation
299
300
- http://en.wikipedia.org/wiki/Jacobi's_elliptic_functions
301
302
- A. Khare, U. Sukhatme, 'Cyclic Identities Involving
303
Jacobi Elliptic Functions', Math ArXiv, math-ph/0201004
304
305
- Online Encyclopedia of Special Function
306
http://algo.inria.fr/esf/index.html
307
308
TODO: Resolve weird bug in commented out code in hypergeometric_U
309
below.
310
311
AUTHORS:
312
313
- David Joyner and William Stein
314
315
Added 16-02-2008 (wdj): optional calls to scipy and replace all
316
'#random' by '...' (both at the request of William Stein)
317
318
.. warning::
319
320
SciPy's versions are poorly documented and seem less
321
accurate than the Maxima and PARI versions.
322
"""
323
324
#*****************************************************************************
325
# Copyright (C) 2006 William Stein <[email protected]>
326
# 2006 David Joyner <[email protected]>
327
#
328
# Distributed under the terms of the GNU General Public License (GPL)
329
#
330
# This code is distributed in the hope that it will be useful,
331
# but WITHOUT ANY WARRANTY; without even the implied warranty of
332
# MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
333
# General Public License for more details.
334
#
335
# The full text of the GPL is available at:
336
#
337
# http://www.gnu.org/licenses/
338
#*****************************************************************************
339
340
from sage.plot.plot import plot
341
from sage.rings.real_mpfr import RealField
342
from sage.rings.complex_field import ComplexField
343
from sage.misc.sage_eval import sage_eval
344
from sage.rings.all import ZZ, RR, RDF
345
from sage.functions.other import real, imag, log_gamma
346
from sage.symbolic.function import BuiltinFunction
347
from sage.calculus.calculus import maxima
348
349
_done = False
350
def _init():
351
"""
352
Internal function which checks if Maxima has loaded the
353
"orthopoly" package. All functions using this in this
354
file should call this function first.
355
356
TEST:
357
358
The global starts ``False``::
359
360
sage: sage.functions.special._done
361
False
362
363
Then after using one of these functions, it changes::
364
365
sage: from sage.functions.special import airy_ai
366
sage: airy_ai(1.0)
367
0.135292416313
368
sage: sage.functions.special._done
369
True
370
"""
371
global _done
372
if _done:
373
return
374
maxima.eval('load("orthopoly");')
375
maxima.eval('orthopoly_returns_intervals:false;')
376
_done = True
377
378
def meval(x):
379
"""
380
Returns ``x`` evaluated in Maxima, then returned to Sage.
381
This is used to evaluate several of these special functions.
382
383
TEST::
384
385
sage: from sage.functions.special import airy_ai
386
sage: airy_bi(1.0)
387
1.20742359495
388
"""
389
return maxima(x).sage()
390
391
392
class MaximaFunction(BuiltinFunction):
393
"""
394
EXAMPLES::
395
396
sage: from sage.functions.special import MaximaFunction
397
sage: f = MaximaFunction("jacobi_sn")
398
sage: f(1,1)
399
tanh(1)
400
sage: f(1/2,1/2).n()
401
0.470750473655657
402
"""
403
def __init__(self, name, nargs=2, conversions={}):
404
"""
405
EXAMPLES::
406
407
sage: from sage.functions.special import MaximaFunction
408
sage: f = MaximaFunction("jacobi_sn")
409
sage: f(1,1)
410
tanh(1)
411
sage: f(1/2,1/2).n()
412
0.470750473655657
413
"""
414
c = dict(maxima=name)
415
c.update(conversions)
416
BuiltinFunction.__init__(self, name=name, nargs=nargs,
417
conversions=c)
418
419
def _maxima_init_evaled_(self, *args):
420
"""
421
Returns a string which represents this function evaluated at
422
*args* in Maxima.
423
424
EXAMPLES::
425
426
sage: from sage.functions.special import MaximaFunction
427
sage: f = MaximaFunction("jacobi_sn")
428
sage: f._maxima_init_evaled_(1/2, 1/2)
429
'jacobi_sn(1/2, 1/2)'
430
431
TESTS:
432
433
Check if complex numbers in the arguments are converted to maxima
434
correctly :trac:`7557`::
435
436
sage: t = jacobi('sn',1.2+2*I*elliptic_kc(1-.5),.5)
437
sage: t._maxima_init_(maxima)
438
'0.88771548861927...*%i'
439
sage: t.n() # abs tol 1e-13
440
0.887715488619280 - 1.79195288804672e-15*I
441
442
"""
443
args_maxima = []
444
for a in args:
445
if isinstance(a, str):
446
args_maxima.append(a)
447
elif hasattr(a, '_maxima_init_'):
448
args_maxima.append(a._maxima_init_())
449
else:
450
args_maxima.append(str(a))
451
return "%s(%s)"%(self.name(), ', '.join(args_maxima))
452
453
def _evalf_(self, *args, **kwds):
454
"""
455
Returns a numerical approximation of this function using
456
Maxima. Currently, this is limited to 53 bits of precision.
457
458
EXAMPLES::
459
460
sage: from sage.functions.special import MaximaFunction
461
sage: f = MaximaFunction("jacobi_sn")
462
sage: f(1/2,1/2)
463
jacobi_sn(1/2, 1/2)
464
sage: f(1/2,1/2).n()
465
0.470750473655657
466
467
TESTS::
468
469
sage: f(1/2,1/2).n(150)
470
Traceback (most recent call last):
471
...
472
NotImplementedError: jacobi_sn not implemented for precision > 53
473
"""
474
parent = kwds['parent']
475
if hasattr(parent, 'prec') and parent.prec() > 53:
476
raise NotImplementedError, "%s not implemented for precision > 53"%self.name()
477
_init()
478
return parent(maxima("%s, numer"%self._maxima_init_evaled_(*args)))
479
480
def _eval_(self, *args):
481
"""
482
Returns a string which represents this function evaluated at
483
*args* in Maxima.
484
485
EXAMPLES::
486
487
sage: from sage.functions.special import MaximaFunction
488
sage: f = MaximaFunction("jacobi_sn")
489
sage: f(1,1)
490
tanh(1)
491
492
sage: f._eval_(1,1)
493
tanh(1)
494
495
Here arccoth doesn't have 1 in its domain, so we just hold the expression:
496
497
sage: elliptic_e(arccoth(1), x^2*e)
498
elliptic_e(arccoth(1), x^2*e)
499
"""
500
_init()
501
try:
502
s = maxima(self._maxima_init_evaled_(*args))
503
except TypeError:
504
return None
505
if self.name() in repr(s):
506
return None
507
else:
508
return s.sage()
509
510
from sage.misc.cachefunc import cached_function
511
512
@cached_function
513
def maxima_function(name):
514
"""
515
Returns a function which is evaluated both symbolically and
516
numerically via Maxima. In particular, it returns an instance
517
of :class:`MaximaFunction`.
518
519
.. note::
520
521
This function is cached so that duplicate copies of the same
522
function are not created.
523
524
EXAMPLES::
525
526
sage: spherical_hankel2(2,i)
527
-e
528
"""
529
# The superclass of MaximaFunction, BuiltinFunction, assumes that there
530
# will be only one symbolic function with the same name and class.
531
# We create a new class for each Maxima function wrapped.
532
class NewMaximaFunction(MaximaFunction):
533
def __init__(self):
534
"""
535
Constructs an object that wraps a Maxima function.
536
537
TESTS::
538
539
sage: spherical_hankel2(2,x)
540
(-I*x^2 - 3*x + 3*I)*e^(-I*x)/x^3
541
"""
542
MaximaFunction.__init__(self, name)
543
544
return NewMaximaFunction()
545
546
547
def airy_ai(x):
548
r"""
549
The function `Ai(x)` and the related function `Bi(x)`,
550
which is also called an *Airy function*, are
551
solutions to the differential equation
552
553
.. math::
554
555
y'' - xy = 0,
556
557
known as the *Airy equation*. The initial conditions
558
`Ai(0) = (\Gamma(2/3)3^{2/3})^{-1}`,
559
`Ai'(0) = -(\Gamma(1/3)3^{1/3})^{-1}` define `Ai(x)`.
560
The initial conditions `Bi(0) = 3^{1/2}Ai(0)`,
561
`Bi'(0) = -3^{1/2}Ai'(0)` define `Bi(x)`.
562
563
They are named after the British astronomer George Biddell Airy.
564
They belong to the class of "Bessel functions of fractional order".
565
566
EXAMPLES::
567
568
sage: airy_ai(1.0) # last few digits are random
569
0.135292416312881400
570
sage: airy_bi(1.0) # last few digits are random
571
1.20742359495287099
572
573
REFERENCE:
574
575
- Abramowitz and Stegun: Handbook of Mathematical Functions,
576
http://www.math.sfu.ca/~cbm/aands/
577
578
- http://en.wikipedia.org/wiki/Airy_function
579
"""
580
_init()
581
return RDF(meval("airy_ai(%s)"%RDF(x)))
582
583
def airy_bi(x):
584
r"""
585
The function `Ai(x)` and the related function `Bi(x)`,
586
which is also called an *Airy function*, are
587
solutions to the differential equation
588
589
.. math::
590
591
y'' - xy = 0,
592
593
known as the *Airy equation*. The initial conditions
594
`Ai(0) = (\Gamma(2/3)3^{2/3})^{-1}`,
595
`Ai'(0) = -(\Gamma(1/3)3^{1/3})^{-1}` define `Ai(x)`.
596
The initial conditions `Bi(0) = 3^{1/2}Ai(0)`,
597
`Bi'(0) = -3^{1/2}Ai'(0)` define `Bi(x)`.
598
599
They are named after the British astronomer George Biddell Airy.
600
They belong to the class of "Bessel functions of fractional order".
601
602
EXAMPLES::
603
604
sage: airy_ai(1) # last few digits are random
605
0.135292416312881400
606
sage: airy_bi(1) # last few digits are random
607
1.20742359495287099
608
609
REFERENCE:
610
611
- Abramowitz and Stegun: Handbook of Mathematical Functions,
612
http://www.math.sfu.ca/~cbm/aands/
613
614
- http://en.wikipedia.org/wiki/Airy_function
615
"""
616
_init()
617
return RDF(meval("airy_bi(%s)"%RDF(x)))
618
619
620
def hypergeometric_U(alpha,beta,x,algorithm="pari",prec=53):
621
r"""
622
Default is a wrap of PARI's hyperu(alpha,beta,x) function.
623
Optionally, algorithm = "scipy" can be used.
624
625
The confluent hypergeometric function `y = U(a,b,x)` is
626
defined to be the solution to Kummer's differential equation
627
628
.. math::
629
630
xy'' + (b-x)y' - ay = 0.
631
632
This satisfies `U(a,b,x) \sim x^{-a}`, as
633
`x\rightarrow \infty`, and is sometimes denoted
634
``x^{-a}2_F_0(a,1+a-b,-1/x)``. This is not the same as Kummer's
635
`M`-hypergeometric function, denoted sometimes as
636
``_1F_1(alpha,beta,x)``, though it satisfies the same DE that
637
`U` does.
638
639
.. warning::
640
641
In the literature, both are called "Kummer confluent
642
hypergeometric" functions.
643
644
EXAMPLES::
645
646
sage: hypergeometric_U(1,1,1,"scipy")
647
0.596347362323...
648
sage: hypergeometric_U(1,1,1)
649
0.59634736232319...
650
sage: hypergeometric_U(1,1,1,"pari",70)
651
0.59634736232319407434...
652
"""
653
if algorithm=="scipy":
654
if prec != 53:
655
raise ValueError, "for the scipy algorithm the precision must be 53"
656
import scipy.special
657
ans = str(scipy.special.hyperu(float(alpha),float(beta),float(x)))
658
ans = ans.replace("(","")
659
ans = ans.replace(")","")
660
ans = ans.replace("j","*I")
661
return sage_eval(ans)
662
elif algorithm=='pari':
663
from sage.libs.pari.all import pari
664
R = RealField(prec)
665
return R(pari(R(alpha)).hyperu(R(beta), R(x), precision=prec))
666
else:
667
raise ValueError, "unknown algorithm '%s'"%algorithm
668
669
def spherical_bessel_J(n, var, algorithm="maxima"):
670
r"""
671
Returns the spherical Bessel function of the first kind for
672
integers n >= 1.
673
674
Reference: AS 10.1.8 page 437 and AS 10.1.15 page 439.
675
676
EXAMPLES::
677
678
sage: spherical_bessel_J(2,x)
679
((3/x^2 - 1)*sin(x) - 3*cos(x)/x)/x
680
sage: spherical_bessel_J(1, 5.2, algorithm='scipy')
681
-0.12277149950007...
682
sage: spherical_bessel_J(1, 3, algorithm='scipy')
683
0.345677499762355...
684
"""
685
if algorithm=="scipy":
686
from scipy.special.specfun import sphj
687
return sphj(int(n), float(var))[1][-1]
688
elif algorithm == 'maxima':
689
_init()
690
return meval("spherical_bessel_j(%s,%s)"%(ZZ(n),var))
691
else:
692
raise ValueError, "unknown algorithm '%s'"%algorithm
693
694
def spherical_bessel_Y(n,var, algorithm="maxima"):
695
r"""
696
Returns the spherical Bessel function of the second kind for
697
integers n -1.
698
699
Reference: AS 10.1.9 page 437 and AS 10.1.15 page 439.
700
701
EXAMPLES::
702
703
sage: x = PolynomialRing(QQ, 'x').gen()
704
sage: spherical_bessel_Y(2,x)
705
-((3/x^2 - 1)*cos(x) + 3*sin(x)/x)/x
706
"""
707
if algorithm=="scipy":
708
import scipy.special
709
ans = str(scipy.special.sph_yn(int(n),float(var)))
710
ans = ans.replace("(","")
711
ans = ans.replace(")","")
712
ans = ans.replace("j","*I")
713
return sage_eval(ans)
714
elif algorithm == 'maxima':
715
_init()
716
return meval("spherical_bessel_y(%s,%s)"%(ZZ(n),var))
717
else:
718
raise ValueError, "unknown algorithm '%s'"%algorithm
719
720
def spherical_hankel1(n, var):
721
r"""
722
Returns the spherical Hankel function of the first kind for
723
integers `n > -1`, written as a string. Reference: AS
724
10.1.36 page 439.
725
726
EXAMPLES::
727
728
sage: spherical_hankel1(2, x)
729
(I*x^2 - 3*x - 3*I)*e^(I*x)/x^3
730
"""
731
return maxima_function("spherical_hankel1")(ZZ(n), var)
732
733
def spherical_hankel2(n,x):
734
r"""
735
Returns the spherical Hankel function of the second kind for
736
integers `n > -1`, written as a string. Reference: AS 10.1.17 page
737
439.
738
739
EXAMPLES::
740
741
sage: spherical_hankel2(2, x)
742
(-I*x^2 - 3*x + 3*I)*e^(-I*x)/x^3
743
744
Here I = sqrt(-1).
745
"""
746
return maxima_function("spherical_hankel2")(ZZ(n), x)
747
748
def spherical_harmonic(m,n,x,y):
749
r"""
750
Returns the spherical Harmonic function of the second kind for
751
integers `n > -1`, `|m|\leq n`. Reference:
752
Merzbacher 9.64.
753
754
EXAMPLES::
755
756
sage: x,y = var('x,y')
757
sage: spherical_harmonic(3,2,x,y)
758
15/4*sqrt(7/30)*cos(x)*e^(2*I*y)*sin(x)^2/sqrt(pi)
759
sage: spherical_harmonic(3,2,1,2)
760
15/4*sqrt(7/30)*cos(1)*e^(4*I)*sin(1)^2/sqrt(pi)
761
"""
762
_init()
763
return meval("spherical_harmonic(%s,%s,%s,%s)"%(ZZ(m),ZZ(n),x,y))
764
765
####### elliptic functions and integrals
766
767
def elliptic_j(z):
768
r"""
769
Returns the elliptic modular `j`-function evaluated at `z`.
770
771
INPUT:
772
773
- ``z`` (complex) -- a complex number with positive imaginary part.
774
775
OUTPUT:
776
777
(complex) The value of `j(z)`.
778
779
ALGORITHM:
780
781
Calls the ``pari`` function ``ellj()``.
782
783
AUTHOR:
784
785
John Cremona
786
787
EXAMPLES::
788
789
sage: elliptic_j(CC(i))
790
1728.00000000000
791
sage: elliptic_j(sqrt(-2.0))
792
8000.00000000000
793
sage: z = ComplexField(100)(1,sqrt(11))/2
794
sage: elliptic_j(z)
795
-32768.000...
796
sage: elliptic_j(z).real().round()
797
-32768
798
799
"""
800
CC = z.parent()
801
from sage.rings.complex_field import is_ComplexField
802
if not is_ComplexField(CC):
803
CC = ComplexField()
804
try:
805
z = CC(z)
806
except ValueError:
807
raise ValueError, "elliptic_j only defined for complex arguments."
808
from sage.libs.all import pari
809
return CC(pari(z).ellj())
810
811
def jacobi(sym,x,m):
812
r"""
813
Here sym = "pq", where p,q in c,d,n,s. This returns the value of
814
the Jacobi function pq(x,m), as described in the documentation for
815
Sage's "special" module. There are a total of 12 functions
816
described by this.
817
818
EXAMPLES::
819
820
sage: jacobi("sn",1,1)
821
tanh(1)
822
sage: jacobi("cd",1,1/2)
823
jacobi_cd(1, 1/2)
824
sage: RDF(jacobi("cd",1,1/2))
825
0.724009721659
826
sage: RDF(jacobi("cn",1,1/2)); RDF(jacobi("dn",1,1/2)); RDF(jacobi("cn",1,1/2)/jacobi("dn",1,1/2))
827
0.595976567672
828
0.823161001632
829
0.724009721659
830
sage: jsn = jacobi("sn",x,1)
831
sage: P = plot(jsn,0,1)
832
833
To view this, type P.show().
834
"""
835
return maxima_function("jacobi_%s"%sym)(x,m)
836
837
def inverse_jacobi(sym,x,m):
838
"""
839
Here sym = "pq", where p,q in c,d,n,s. This returns the value of
840
the inverse Jacobi function `pq^{-1}(x,m)`. There are a
841
total of 12 functions described by this.
842
843
EXAMPLES::
844
845
sage: jacobi("sn",1/2,1/2)
846
jacobi_sn(1/2, 1/2)
847
sage: float(jacobi("sn",1/2,1/2))
848
0.4707504736556572
849
sage: float(inverse_jacobi("sn",0.47,1/2))
850
0.4990982313222196
851
sage: float(inverse_jacobi("sn",0.4707504,0.5))
852
0.4999999114665548
853
sage: P = plot(inverse_jacobi('sn', x, 0.5), 0, 1, plot_points=20)
854
855
Now to view this, just type show(P).
856
"""
857
return maxima_function("inverse_jacobi_%s"%sym)(x,m)
858
859
#### elliptic integrals
860
861
class EllipticE(MaximaFunction):
862
r"""
863
This returns the value of the "incomplete elliptic integral of the
864
second kind,"
865
866
.. math::
867
868
\int_0^\phi \sqrt{1 - m\sin(x)^2}\, dx,
869
870
i.e., ``integrate(sqrt(1 - m*sin(x)^2), x, 0, phi)``. Taking `\phi
871
= \pi/2` gives ``elliptic_ec``.
872
873
EXAMPLES::
874
875
sage: z = var("z")
876
sage: # this is still wrong: must be abs(sin(z)) + 2*round(z/pi)
877
sage: elliptic_e(z, 1)
878
2*round(z/pi) + sin(z)
879
sage: elliptic_e(z, 0)
880
z
881
sage: elliptic_e(0.5, 0.1)
882
0.498011394498832
883
884
sage: loads(dumps(elliptic_e))
885
elliptic_e
886
"""
887
def __init__(self):
888
"""
889
EXAMPLES::
890
891
sage: elliptic_e(0.5, 0.1)
892
0.498011394498832
893
"""
894
MaximaFunction.__init__(self, "elliptic_e")
895
896
elliptic_e = EllipticE()
897
898
class EllipticEC(MaximaFunction):
899
"""
900
This returns the value of the "complete elliptic integral of the
901
second kind,"
902
903
.. math::
904
905
\int_0^{\pi/2} \sqrt{1 - m\sin(x)^2}\, dx.
906
907
EXAMPLES::
908
909
sage: elliptic_ec(0.1)
910
1.53075763689776
911
sage: elliptic_ec(x).diff()
912
1/2*(elliptic_ec(x) - elliptic_kc(x))/x
913
914
sage: loads(dumps(elliptic_ec))
915
elliptic_ec
916
"""
917
def __init__(self):
918
"""
919
EXAMPLES::
920
921
sage: elliptic_ec(0.1)
922
1.53075763689776
923
"""
924
MaximaFunction.__init__(self, "elliptic_ec", nargs=1)
925
926
def _derivative_(self, *args, **kwds):
927
"""
928
EXAMPLES::
929
930
sage: elliptic_ec(x).diff()
931
1/2*(elliptic_ec(x) - elliptic_kc(x))/x
932
"""
933
diff_param = kwds['diff_param']
934
assert diff_param == 0
935
x = args[diff_param]
936
return (elliptic_ec(x) - elliptic_kc(x))/(2*x)
937
938
elliptic_ec = EllipticEC()
939
940
class EllipticEU(MaximaFunction):
941
r"""
942
This returns the value of the "incomplete elliptic integral of the
943
second kind,"
944
945
.. math::
946
947
\int_0^u \mathrm{dn}(x,m)^2\, dx = \int_0^\tau
948
{\sqrt{1-m x^2}\over\sqrt{1-x^2}}\, dx.
949
950
where `\tau = \mathrm{sn}(u, m)`.
951
952
EXAMPLES::
953
954
sage: elliptic_eu (0.5, 0.1)
955
0.496054551286597
956
"""
957
def __init__(self):
958
r"""
959
EXAMPLES::
960
961
sage: elliptic_eu (0.5, 0.1)
962
0.496054551286597
963
"""
964
MaximaFunction.__init__(self, "elliptic_eu")
965
966
elliptic_eu = EllipticEU()
967
968
class EllipticF(MaximaFunction):
969
r"""
970
This returns the value of the "incomplete elliptic integral of the
971
first kind,"
972
973
.. math::
974
975
\int_0^\phi \frac{dx}{\sqrt{1 - m\sin(x)^2}},
976
977
i.e., ``integrate(1/sqrt(1 - m*sin(x)^2), x, 0, phi)``. Taking
978
`\phi = \pi/2` gives ``elliptic_kc``.
979
980
EXAMPLES::
981
982
sage: z = var("z")
983
sage: elliptic_f (z, 0)
984
z
985
sage: elliptic_f (z, 1)
986
log(tan(1/4*pi + 1/2*z))
987
sage: elliptic_f (0.2, 0.1)
988
0.200132506747543
989
"""
990
def __init__(self):
991
r"""
992
EXAMPLES::
993
994
sage: elliptic_f (0.2, 0.1)
995
0.200132506747543
996
"""
997
MaximaFunction.__init__(self, "elliptic_f")
998
999
elliptic_f = EllipticF()
1000
1001
class EllipticKC(MaximaFunction):
1002
r"""
1003
This returns the value of the "complete elliptic integral of the
1004
first kind,"
1005
1006
.. math::
1007
1008
\int_0^{\pi/2} \frac{dx}{\sqrt{1 - m\sin(x)^2}}.
1009
1010
EXAMPLES::
1011
1012
sage: elliptic_kc(0.5)
1013
1.85407467730137
1014
sage: elliptic_f(RR(pi/2), 0.5)
1015
1.85407467730137
1016
"""
1017
def __init__(self):
1018
r"""
1019
EXAMPLES::
1020
1021
sage: elliptic_kc(0.5)
1022
1.85407467730137
1023
sage: elliptic_f(RR(pi/2), 0.5)
1024
1.85407467730137
1025
"""
1026
MaximaFunction.__init__(self, "elliptic_kc", nargs=1)
1027
1028
elliptic_kc = EllipticKC()
1029
1030
class EllipticPi(MaximaFunction):
1031
r"""
1032
This returns the value of the "incomplete elliptic integral of the
1033
third kind,"
1034
1035
.. math::
1036
1037
\text{elliptic\_pi}(n, t, m) = \int_0^t \frac{dx}{(1 - n \sin(x)^2)
1038
\sqrt{1 - m \sin(x)^2}}.
1039
1040
INPUT:
1041
1042
- ``n`` -- a real number, called the "characteristic"
1043
1044
- ``t`` -- a real number, called the "amplitude"
1045
1046
- ``m`` -- a real number, called the "parameter"
1047
1048
EXAMPLES::
1049
1050
sage: N(elliptic_pi(1, pi/4, 1))
1051
1.14779357469632
1052
1053
Compare the value computed by Maxima to the definition as a definite integral
1054
(using GSL)::
1055
1056
sage: elliptic_pi(0.1, 0.2, 0.3)
1057
0.200665068220979
1058
sage: numerical_integral(1/(1-0.1*sin(x)^2)/sqrt(1-0.3*sin(x)^2), 0.0, 0.2)
1059
(0.2006650682209791, 2.227829789769088e-15)
1060
1061
ALGORITHM:
1062
1063
Numerical evaluation and symbolic manipulation are provided by `Maxima`_.
1064
1065
REFERENCES:
1066
1067
- Abramowitz and Stegun: Handbook of Mathematical Functions, section 17.7
1068
http://www.math.sfu.ca/~cbm/aands/
1069
- Elliptic Functions in `Maxima`_
1070
1071
.. _`Maxima`: http://maxima.sourceforge.net/docs/manual/en/maxima_16.html#SEC91
1072
"""
1073
def __init__(self):
1074
r"""
1075
EXAMPLES::
1076
1077
sage: elliptic_pi(0.1, 0.2, 0.3)
1078
0.200665068220979
1079
"""
1080
MaximaFunction.__init__(self, "elliptic_pi", nargs=3)
1081
1082
elliptic_pi = EllipticPi()
1083
1084
1085
def lngamma(t):
1086
r"""
1087
This method is deprecated, please use
1088
:meth:`~sage.functions.other.log_gamma` instead.
1089
1090
See the :meth:`~sage.functions.other.log_gamma` function for '
1091
documentation and examples.
1092
1093
EXAMPLES::
1094
1095
sage: lngamma(RR(6))
1096
doctest:...: DeprecationWarning: The method lngamma() is deprecated. Use log_gamma() instead.
1097
See http://trac.sagemath.org/6992 for details.
1098
4.78749174278205
1099
"""
1100
from sage.misc.superseded import deprecation
1101
deprecation(6992, "The method lngamma() is deprecated. Use log_gamma() instead.")
1102
return log_gamma(t)
1103
1104
def error_fcn(t):
1105
r"""
1106
The complementary error function
1107
`\frac{2}{\sqrt{\pi}}\int_t^\infty e^{-x^2} dx` (t belongs
1108
to RR). This function is currently always
1109
evaluated immediately.
1110
1111
EXAMPLES::
1112
1113
sage: error_fcn(6)
1114
2.15197367124989e-17
1115
sage: error_fcn(RealField(100)(1/2))
1116
0.47950012218695346231725334611
1117
1118
Note this is literally equal to `1 - erf(t)`::
1119
1120
sage: 1 - error_fcn(0.5)
1121
0.520499877813047
1122
sage: erf(0.5)
1123
0.520499877813047
1124
"""
1125
try:
1126
return t.erfc()
1127
except AttributeError:
1128
from sage.rings.real_mpfr import RR
1129
try:
1130
return RR(t).erfc()
1131
except Exception:
1132
raise NotImplementedError
1133
1134
1135
1136
1137