This repository has been archived by the owner on Jan 30, 2023. It is now read-only.
-
-
Notifications
You must be signed in to change notification settings - Fork 7
/
bessel.py
1580 lines (1238 loc) · 48.3 KB
/
bessel.py
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
r"""
Bessel Functions
This module provides symbolic Bessel Functions. These functions use the
`mpmath library`_ for numerical evaluation and Maxima, GiNaC, Pynac for
symbolics.
The main objects which are exported from this module are:
* ``bessel_J`` -- The Bessel J function
* ``bessel_Y`` -- The Bessel Y function
* ``bessel_I`` -- The Bessel I function
* ``bessel_K`` -- The Bessel K function
* ``Bessel`` -- A factory function for producing Bessel functions of
various kinds and orders
- Bessel functions, first defined by the Swiss mathematician
Daniel Bernoulli and named after Friedrich Bessel, are canonical
solutions y(x) of Bessel's differential equation:
.. math::
x^2 \frac{d^2 y}{dx^2} + x \frac{dy}{dx} + \left(x^2 - \nu^2\right)y =
0,
for an arbitrary complex number `\nu` (the order).
- In this module, `J_\nu` denotes the unique solution of Bessel's equation
which is non-singular at `x = 0`. This function is known as the Bessel
Function of the First Kind. This function also arises as a special case
of the hypergeometric function `{}_0F_1`:
.. math::
J_\nu(x) = \frac{x^n}{2^\nu \Gamma(\nu + 1)} {}_0F_1(\nu +
1, -\frac{x^2}{4}).
- The second linearly independent solution to Bessel's equation (which is
singular at `x=0`) is denoted by `Y_\nu` and is called the Bessel
Function of the Second Kind:
.. math::
Y_\nu(x) = \frac{ J_\nu(x) \cos(\pi \nu) -
J_{-\nu}(x)}{\sin(\pi \nu)}.
- There are also two commonly used combinations of the Bessel J and Y
Functions. The Bessel I Function, or the Modified Bessel Function of the
First Kind, is defined by:
.. math::
I_\nu(x) = i^{-\nu} J_\nu(ix).
The Bessel K Function, or the Modified Bessel Function of the Second Kind,
is defined by:
.. math::
K_\nu(x) = \frac{\pi}{2} \cdot \frac{I_{-\nu}(x) -
I_n(x)}{\sin(\pi \nu)}.
We should note here that the above formulas for Bessel Y and K functions
should be understood as limits when `\nu` is an integer.
- It follows from Bessel's differential equation that the derivative of
`J_n(x)` with respect to `x` is:
.. math::
\frac{d}{dx} J_n(x) = \frac{1}{x^n} \left(x^n J_{n-1}(x) - n x^{n-1}
J_n(z) \right)
- Another important formulation of the two linearly independent
solutions to Bessel's equation are the Hankel functions
`H_\nu^{(1)}(x)` and `H_\nu^{(2)}(x)`,
defined by:
.. math::
H_\nu^{(1)}(x) = J_\nu(x) + i Y_\nu(x)
.. math::
H_\nu^{(2)}(x) = J_\nu(x) - i Y_\nu(x)
where `i` is the imaginary unit (and `J_*` and
`Y_*` are the usual J- and Y-Bessel functions). These
linear combinations are also known as Bessel functions of the third
kind; they are also two linearly independent solutions of Bessel's
differential equation. They are named for Hermann Hankel.
EXAMPLES:
Evaluate the Bessel J function symbolically and numerically::
sage: bessel_J(0, x)
bessel_J(0, x)
sage: bessel_J(0, 0)
bessel_J(0, 0)
sage: bessel_J(0, x).diff(x)
-1/2*bessel_J(1, x) + 1/2*bessel_J(-1, x)
sage: N(bessel_J(0, 0), digits = 20)
1.0000000000000000000
sage: find_root(bessel_J(0,x), 0, 5)
2.404825557695773
Plot the Bessel J function::
sage: f(x) = Bessel(0)(x); f
x |--> bessel_J(0, x)
sage: plot(f, (x, 1, 10))
Visualize the Bessel Y function on the complex plane::
sage: complex_plot(bessel_Y(0, x), (-5, 5), (-5, 5))
Evaluate a combination of Bessel functions::
sage: f(x) = bessel_J(1, x) - bessel_Y(0, x)
sage: f(pi)
bessel_J(1, pi) - bessel_Y(0, pi)
sage: f(pi).n()
-0.0437509653365599
sage: f(pi).n(digits=50)
-0.043750965336559909054985168023342675387737118378169
Symbolically solve a second order differential equation with initial
conditions `y(1) = a` and `y'(1) = b` in terms of Bessel functions::
sage: y = function('y', x)
sage: a, b = var('a, b')
sage: diffeq = x^2*diff(y,x,x) + x*diff(y,x) + x^2*y == 0
sage: f = desolve(diffeq, y, [1, a, b]); f
(a*bessel_Y(1, 1) + b*bessel_Y(0, 1))*bessel_J(0, x)/(bessel_J(0,
1)*bessel_Y(1, 1) - bessel_J(1, 1)*bessel_Y(0, 1)) -
(a*bessel_J(1, 1) + b*bessel_J(0, 1))*bessel_Y(0, x)/(bessel_J(0,
1)*bessel_Y(1, 1) - bessel_J(1, 1)*bessel_Y(0, 1))
For more examples, see the docstring for :meth:`Bessel`.
AUTHORS:
- Benjamin Jones (2012-12-27): initial version
- Some of the documentation here has been adapted from David Joyner's
original documentation of Sage's special functions module (2006).
REFERENCES:
- Abramowitz and Stegun: Handbook of Mathematical Functions,
http://www.math.sfu.ca/~cbm/aands/
- http://en.wikipedia.org/wiki/Bessel_function
- mpmath Library `Bessel Functions`_
.. _`mpmath Library`: http://code.google.com/p/mpmath/
.. _`Bessel Functions`: http://mpmath.googlecode.com/svn/trunk/doc/build/functions/bessel.html
"""
#*****************************************************************************
# Copyright (C) 2013 Benjamin Jones <benjaminfjones@gmail.com>
#
# Distributed under the terms of the GNU General Public License (GPL)
#
# This code is distributed in the hope that it will be useful,
# but WITHOUT ANY WARRANTY; without even the implied warranty of
# MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
# General Public License for more details.
#
# The full text of the GPL is available at:
#
# http://www.gnu.org/licenses/
#*****************************************************************************
from sage.functions.other import sqrt
from sage.functions.log import exp
from sage.functions.hyperbolic import sinh, cosh
from sage.libs.mpmath import utils as mpmath_utils
from sage.misc.latex import latex
from sage.rings.all import RR, Integer
from sage.structure.coerce import parent
from sage.structure.element import get_coercion_model
from sage.symbolic.constants import pi
from sage.symbolic.function import BuiltinFunction, is_inexact
from sage.symbolic.expression import Expression
# remove after deprecation period
from sage.calculus.calculus import maxima
from sage.functions.other import real, imag
from sage.misc.sage_eval import sage_eval
from sage.rings.real_mpfr import RealField
from sage.plot.plot import plot
from sage.rings.all import ZZ
class Function_Bessel_J(BuiltinFunction):
r"""
The Bessel J Function, denoted by bessel_J(`\nu`, x) or `J_\nu(x)`.
As a Taylor series about `x=0` it is equal to:
.. math::
J_\nu(x) = \sum_{k=0}^\infty \frac{(-1)^k}{k! \Gamma(k+\nu+1)}
\left(\frac{x}{2}\right)^{2k+\nu}
The parameter `\nu` is called the order and may be any real or
complex number; however, integer and half-integer values are most
common. It is defined for all complex numbers `x` when `\nu`
is an integer or greater than zero and it diverges as `x \to 0`
for negative non-integer values of `\nu`.
For integer orders `\nu = n` there is an integral representation:
.. math::
J_n(x) = \frac{1}{\pi} \int_0^\pi \cos(n t - x \sin(t)) \; dt
This function also arises as a special case of the hypergeometric
function `{}_0F_1`:
.. math::
J_\nu(x) = \frac{x^n}{2^\nu \Gamma(\nu + 1)} {}_0F_1\left(\nu +
1, -\frac{x^2}{4}\right).
EXAMPLES::
sage: bessel_J(1.0, 1.0)
0.440050585744933
sage: bessel_J(2, I).n(digits=30)
-0.135747669767038281182852569995
sage: bessel_J(1, x)
bessel_J(1, x)
sage: n = var('n')
sage: bessel_J(n, x)
bessel_J(n, x)
Examples of symbolic manipulation::
sage: a = bessel_J(pi, bessel_J(1, I)); a
bessel_J(pi, bessel_J(1, I))
sage: N(a, digits=20)
0.00059023706363796717363 - 0.0026098820470081958110*I
sage: f = bessel_J(2, x)
sage: f.diff(x)
-1/2*bessel_J(3, x) + 1/2*bessel_J(1, x)
Comparison to a well-known integral representation of `J_1(1)`::
sage: A = numerical_integral(1/pi*cos(x - sin(x)), 0, pi)
sage: A[0] # abs tol 1e-14
0.44005058574493355
sage: bessel_J(1.0, 1.0) - A[0] < 1e-15
True
Currently, integration is not supported (directly) since we cannot
yet convert hypergeometric functions to and from Maxima::
sage: f = bessel_J(2, x)
sage: f.integrate(x)
Traceback (most recent call last):
...
TypeError: cannot coerce arguments: no canonical coercion from <type 'list'> to Symbolic Ring
sage: m = maxima(bessel_J(2, x))
sage: m.integrate(x)
hypergeometric([3/2],[5/2,3],-x^2/4)*x^3/24
Visualization::
sage: plot(bessel_J(1,x), (x,0,5), color='blue')
sage: complex_plot(bessel_J(1, x), (-5, 5), (-5, 5)) # long time
ALGORITHM:
Numerical evaluation is handled by the mpmath library. Symbolics are
handled by a combination of Maxima and Sage (Ginac/Pynac).
"""
def __init__(self):
"""
See the docstring for :meth:`Function_Bessel_J`.
EXAMPLES::
sage: sage.functions.bessel.Function_Bessel_J()
bessel_J
"""
BuiltinFunction.__init__(self, "bessel_J", nargs=2,
conversions=dict(mathematica='BesselJ',
maxima='bessel_j',
sympy='besselj'))
def _eval_(self, n, x):
"""
EXAMPLES::
sage: a, b = var('a, b')
sage: bessel_J(a, b)
bessel_J(a, b)
sage: bessel_J(1.0, 1.0)
0.440050585744933
"""
if (not isinstance(n, Expression) and
not isinstance(x, Expression) and
(is_inexact(n) or is_inexact(x))):
coercion_model = get_coercion_model()
n, x = coercion_model.canonical_coercion(n, x)
return self._evalf_(n, x, parent(n))
return None
def _evalf_(self, n, x, parent=None, algorithm=None):
"""
EXAMPLES::
sage: bessel_J(0.0, 1.0)
0.765197686557966
sage: bessel_J(0, 1).n(digits=20)
0.76519768655796655145
"""
import mpmath
return mpmath_utils.call(mpmath.besselj, n, x, parent=parent)
def _derivative_(self, n, x, diff_param):
"""
Return the derivative of the Bessel J function.
EXAMPLES::
sage: f(z) = bessel_J(10, z)
sage: derivative(f, z)
z |--> -1/2*bessel_J(11, z) + 1/2*bessel_J(9, z)
sage: nu = var('nu')
sage: bessel_J(nu, z).diff(nu)
Traceback (most recent call last):
...
NotImplementedError: derivative with respect to order
"""
if diff_param == 1:
return (bessel_J(n - 1, x) - bessel_J(n + 1, x)) / Integer(2)
else:
raise NotImplementedError('derivative with respect to order')
def _print_latex_(self, n, z):
"""
Custom _print_latex_ method.
EXAMPLES::
sage: latex(bessel_J(1, x))
\operatorname{J_{1}}(x)
"""
return r"\operatorname{J_{%s}}(%s)" % (latex(n), latex(z))
bessel_J = Function_Bessel_J()
class Function_Bessel_Y(BuiltinFunction):
r"""
The Bessel Y functions, also known as the Bessel functions of the second
kind, Weber functions, or Neumann functions.
`Y_\nu(z)` is a holomorphic function of `z` on the complex plane,
cut along the negative real axis. It is singular at `z = 0`. When `z`
is fixed, `Y_\nu(z)` is an entire function of the order `\nu`.
DEFINITION:
.. math::
Y_n(z) = \frac{J_\nu(z) \cos(\nu z) -
J_{-\nu}(z)}{\sin(\nu z)}
Its derivative with respect to `z` is:
.. math::
\frac{d}{dz} Y_n(z) = \frac{1}{z^n} \left(z^n Y_{n-1}(z) - n z^{n-1}
Y_n(z) \right)
EXAMPLES::
sage: bessel_Y(1, x)
bessel_Y(1, x)
sage: bessel_Y(1.0, 1.0)
-0.781212821300289
sage: n = var('n')
sage: bessel_Y(n, x)
bessel_Y(n, x)
sage: bessel_Y(2, I).n()
1.03440456978312 - 0.135747669767038*I
sage: bessel_Y(0, 0).n()
-infinity
Examples of symbolic manipulation::
sage: a = bessel_Y(pi, bessel_Y(1, I)); a
bessel_Y(pi, bessel_Y(1, I))
sage: N(a, digits=20)
4.2059146571791095708 + 21.307914215321993526*I
sage: f = bessel_Y(2, x)
sage: f.diff(x)
-1/2*bessel_Y(3, x) + 1/2*bessel_Y(1, x)
High precision and complex valued inputs (see :trac:`4230`)::
sage: bessel_Y(0, 1).n(128)
0.088256964215676957982926766023515162828
sage: bessel_Y(0, RealField(200)(1))
0.088256964215676957982926766023515162827817523090675546711044
sage: bessel_Y(0, ComplexField(200)(0.5+I))
0.077763160184438051408593468823822434235010300228009867784073 + 1.0142336049916069152644677682828326441579314239591288411739*I
Visualization::
sage: plot(bessel_Y(1,x), (x,0,5), color='blue')
sage: complex_plot(bessel_Y(1, x), (-5, 5), (-5, 5)) # long time
ALGORITHM:
Numerical evaluation is handled by the mpmath library. Symbolics are
handled by a combination of Maxima and Sage (Ginac/Pynac).
"""
def __init__(self):
"""
See the docstring for :meth:`Function_Bessel_Y`.
EXAMPLES::
sage: sage.functions.bessel.Function_Bessel_Y()(0, x)
bessel_Y(0, x)
"""
BuiltinFunction.__init__(self, "bessel_Y", nargs=2,
conversions=dict(mathematica='BesselY',
maxima='bessel_y',
sympy='bessely'))
def _eval_(self, n, x):
"""
EXAMPLES::
sage: a,b = var('a, b')
sage: bessel_Y(a, b)
bessel_Y(a, b)
sage: bessel_Y(0, 1).n(128)
0.088256964215676957982926766023515162828
"""
if (not isinstance(n, Expression) and not isinstance(x, Expression) and
(is_inexact(n) or is_inexact(x))):
coercion_model = get_coercion_model()
n, x = coercion_model.canonical_coercion(n, x)
return self._evalf_(n, x, parent(n))
return None # leaves the expression unevaluated
def _evalf_(self, n, x, parent=None, algorithm=None):
"""
EXAMPLES::
sage: bessel_Y(1.0+2*I, 3.0+4*I)
0.699410324467538 + 0.228917940896421*I
sage: bessel_Y(0, 1).n(256)
0.08825696421567695798292676602351516282781752309067554671104384761199978932351
"""
import mpmath
return mpmath_utils.call(mpmath.bessely, n, x, parent=parent)
def _derivative_(self, n, x, diff_param):
"""
Return the derivative of the Bessel Y function.
EXAMPLES::
sage: f(x) = bessel_Y(10, x)
sage: derivative(f, x)
x |--> -1/2*bessel_Y(11, x) + 1/2*bessel_Y(9, x)
sage: nu = var('nu')
sage: bessel_Y(nu, x).diff(nu)
Traceback (most recent call last):
...
NotImplementedError: derivative with respect to order
"""
if diff_param == 1:
return (bessel_Y(n - 1, x) - bessel_Y(n + 1, x)) / Integer(2)
else:
raise NotImplementedError('derivative with respect to order')
def _print_latex_(self, n, z):
"""
Custom _print_latex_ method.
EXAMPLES::
sage: latex(bessel_Y(1, x))
\operatorname{Y_{1}}(x)
"""
return r"\operatorname{Y_{%s}}(%s)" % (latex(n), latex(z))
bessel_Y = Function_Bessel_Y()
class Function_Bessel_I(BuiltinFunction):
r"""
The Bessel I function, or the Modified Bessel Function of the First Kind.
DEFINITION:
.. math::
I_\nu(x) = i^{-\nu} J_\nu(ix)
EXAMPLES::
sage: bessel_I(1, x)
bessel_I(1, x)
sage: bessel_I(1.0, 1.0)
0.565159103992485
sage: n = var('n')
sage: bessel_I(n, x)
bessel_I(n, x)
sage: bessel_I(2, I).n()
-0.114903484931900
Examples of symbolic manipulation::
sage: a = bessel_I(pi, bessel_I(1, I))
sage: N(a, digits=20)
0.00026073272117205890528 - 0.0011528954889080572266*I
sage: f = bessel_I(2, x)
sage: f.diff(x)
1/2*bessel_I(3, x) + 1/2*bessel_I(1, x)
Special identities that bessel_I satisfies::
sage: bessel_I(1/2, x)
sqrt(2)*sqrt(1/(pi*x))*sinh(x)
sage: eq = bessel_I(1/2, x) == bessel_I(0.5, x)
sage: eq.test_relation()
True
sage: bessel_I(-1/2, x)
sqrt(2)*sqrt(1/(pi*x))*cosh(x)
sage: eq = bessel_I(-1/2, x) == bessel_I(-0.5, x)
sage: eq.test_relation()
True
Examples of asymptotic behavior::
sage: limit(bessel_I(0, x), x=oo)
+Infinity
sage: limit(bessel_I(0, x), x=0)
1
High precision and complex valued inputs::
sage: bessel_I(0, 1).n(128)
1.2660658777520083355982446252147175376
sage: bessel_I(0, RealField(200)(1))
1.2660658777520083355982446252147175376076703113549622068081
sage: bessel_I(0, ComplexField(200)(0.5+I))
0.80644357583493619472428518415019222845373366024179916785502 + 0.22686958987911161141397453401487525043310874687430711021434*I
Visualization::
sage: plot(bessel_I(1,x), (x,0,5), color='blue')
sage: complex_plot(bessel_I(1, x), (-5, 5), (-5, 5)) # long time
ALGORITHM:
Numerical evaluation is handled by the mpmath library. Symbolics are
handled by a combination of Maxima and Sage (Ginac/Pynac).
"""
def __init__(self):
"""
See the docstring for :meth:`Function_Bessel_I`.
EXAMPLES::
sage: bessel_I(1,x)
bessel_I(1, x)
"""
BuiltinFunction.__init__(self, "bessel_I", nargs=2,
conversions=dict(mathematica='BesselI',
maxima='bessel_i',
sympy='besseli'))
def _eval_(self, n, x):
"""
EXAMPLES::
sage: y=var('y')
sage: bessel_I(y,x)
bessel_I(y, x)
sage: bessel_I(0.0, 1.0)
1.26606587775201
sage: bessel_I(1/2, 1)
sqrt(2)*sinh(1)/sqrt(pi)
sage: bessel_I(-1/2, pi)
sqrt(2)*cosh(pi)/pi
"""
if (not isinstance(n, Expression) and not isinstance(x, Expression) and
(is_inexact(n) or is_inexact(x))):
coercion_model = get_coercion_model()
n, x = coercion_model.canonical_coercion(n, x)
return self._evalf_(n, x, parent(n))
# special identities
if n == Integer(1) / Integer(2):
return sqrt(2 / (pi * x)) * sinh(x)
elif n == -Integer(1) / Integer(2):
return sqrt(2 / (pi * x)) * cosh(x)
return None # leaves the expression unevaluated
def _evalf_(self, n, x, parent=None, algorithm=None):
"""
EXAMPLES::
sage: bessel_I(1,3).n(digits=20)
3.9533702174026093965
"""
import mpmath
return mpmath_utils.call(mpmath.besseli, n, x, parent=parent)
def _derivative_(self, n, x, diff_param):
"""
Return the derivative of the Bessel I function `I_n(x)` with respect
to `x`.
EXAMPLES::
sage: f(z) = bessel_I(10, x)
sage: derivative(f, x)
z |--> 1/2*bessel_I(11, x) + 1/2*bessel_I(9, x)
sage: nu = var('nu')
sage: bessel_I(nu, x).diff(nu)
Traceback (most recent call last):
...
NotImplementedError: derivative with respect to order
"""
if diff_param == 1:
return (bessel_I(n - 1, x) + bessel_I(n + 1, x)) / Integer(2)
else:
raise NotImplementedError('derivative with respect to order')
def _print_latex_(self, n, z):
"""
Custom _print_latex_ method.
EXAMPLES::
sage: latex(bessel_I(1, x))
\operatorname{I_{1}}(x)
"""
return r"\operatorname{I_{%s}}(%s)" % (latex(n), latex(z))
bessel_I = Function_Bessel_I()
class Function_Bessel_K(BuiltinFunction):
r"""
The Bessel K function, or the modified Bessel function of the second kind.
DEFINITION:
.. math::
K_\nu(x) = \frac{\pi}{2} \frac{I_{-\nu}(x)-I_\nu(x)}{\sin(\nu \pi)}
EXAMPLES::
sage: bessel_K(1, x)
bessel_K(1, x)
sage: bessel_K(1.0, 1.0)
0.601907230197235
sage: n = var('n')
sage: bessel_K(n, x)
bessel_K(n, x)
sage: bessel_K(2, I).n()
-2.59288617549120 + 0.180489972066962*I
Examples of symbolic manipulation::
sage: a = bessel_K(pi, bessel_K(1, I)); a
bessel_K(pi, bessel_K(1, I))
sage: N(a, digits=20)
3.8507583115005220157 + 0.068528298579883425792*I
sage: f = bessel_K(2, x)
sage: f.diff(x)
1/2*bessel_K(3, x) + 1/2*bessel_K(1, x)
sage: bessel_K(1/2, x)
bessel_K(1/2, x)
sage: bessel_K(1/2, -1)
bessel_K(1/2, -1)
sage: bessel_K(1/2, 1)
sqrt(1/2)*sqrt(pi)*e^(-1)
Examples of asymptotic behavior::
sage: bessel_K(0, 0.0)
+infinity
sage: limit(bessel_K(0, x), x=0)
+Infinity
sage: limit(bessel_K(0, x), x=oo)
0
High precision and complex valued inputs::
sage: bessel_K(0, 1).n(128)
0.42102443824070833333562737921260903614
sage: bessel_K(0, RealField(200)(1))
0.42102443824070833333562737921260903613621974822666047229897
sage: bessel_K(0, ComplexField(200)(0.5+I))
0.058365979093103864080375311643360048144715516692187818271179 - 0.67645499731334483535184142196073004335768129348518210260256*I
Visualization::
sage: plot(bessel_K(1,x), (x,0,5), color='blue')
sage: complex_plot(bessel_K(1, x), (-5, 5), (-5, 5)) # long time
ALGORITHM:
Numerical evaluation is handled by the mpmath library. Symbolics are
handled by a combination of Maxima and Sage (Ginac/Pynac).
TESTS:
Verify that :trac:`3426` is fixed:
The Bessel K function can be evaluated numerically at complex orders::
sage: bessel_K(10 * I, 10).n()
9.82415743819925e-8
For a fixed imaginary order and increasing, real, second component the
value of Bessel K is exponentially decaying::
sage: for x in [10, 20, 50, 100, 200]: print bessel_K(5*I, x).n()
5.27812176514912e-6
3.11005908421801e-10
2.66182488515423e-23 - 8.59622057747552e-58*I
4.11189776828337e-45 - 1.01494840019482e-80*I
1.15159692553603e-88 - 6.75787862113718e-125*I
"""
def __init__(self):
"""
See the docstring for :meth:`Function_Bessel_K`.
EXAMPLES::
sage: sage.functions.bessel.Function_Bessel_K()
bessel_K
"""
BuiltinFunction.__init__(self, "bessel_K", nargs=2,
conversions=dict(mathematica='BesselK',
maxima='bessel_k',
sympy='besselk'))
def _eval_(self, n, x):
"""
EXAMPLES::
sage: bessel_K(1,0)
bessel_K(1, 0)
sage: bessel_K(1.0, 0.0)
+infinity
sage: bessel_K(-1, 1).n(128)
0.60190723019723457473754000153561733926
"""
if (not isinstance(n, Expression) and not isinstance(x, Expression) and
(is_inexact(n) or is_inexact(x))):
coercion_model = get_coercion_model()
n, x = coercion_model.canonical_coercion(n, x)
return self._evalf_(n, x, parent(n))
# special identity
if n == Integer(1) / Integer(2) and x > 0:
return sqrt(pi / 2) * exp(-x) * x ** (-Integer(1) / Integer(2))
return None # leaves the expression unevaluated
def _evalf_(self, n, x, parent=None, algorithm=None):
"""
EXAMPLES::
sage: bessel_K(0.0, 1.0)
0.421024438240708
sage: bessel_K(0, RealField(128)(1))
0.42102443824070833333562737921260903614
"""
import mpmath
return mpmath_utils.call(mpmath.besselk, n, x, parent=parent)
def _derivative_(self, n, x, diff_param):
"""
Return the derivative of the Bessel K function.
EXAMPLES::
sage: f(x) = bessel_K(10, x)
sage: derivative(f, x)
x |--> 1/2*bessel_K(11, x) + 1/2*bessel_K(9, x)
sage: nu = var('nu')
sage: bessel_K(nu, x).diff(nu)
Traceback (most recent call last):
...
NotImplementedError: derivative with respect to order
"""
if diff_param == 1:
return (bessel_K(n - 1, x) + bessel_K(n + 1, x)) / Integer(2)
else:
raise NotImplementedError('derivative with respect to order')
def _print_latex_(self, n, z):
"""
Custom _print_latex_ method.
EXAMPLES::
sage: latex(bessel_K(1, x))
\operatorname{K_{1}}(x)
"""
return r"\operatorname{K_{%s}}(%s)" % (latex(n), latex(z))
bessel_K = Function_Bessel_K()
# dictionary used in Bessel
bessel_type_dict = {'I': bessel_I, 'J': bessel_J, 'K': bessel_K, 'Y': bessel_Y}
def Bessel(*args, **kwds):
"""
A function factory that produces symbolic I, J, K, and Y Bessel functions.
There are several ways to call this function:
- ``Bessel(order, type)``
- ``Bessel(order)`` -- type defaults to 'J'
- ``Bessel(order, typ=T)``
- ``Bessel(typ=T)`` -- order is unspecified, this is a 2-parameter
function
- ``Bessel()`` -- order is unspecified, type is 'J'
where ``order`` can be any integer and T must be one of the strings 'I',
'J', 'K', or 'Y'.
See the EXAMPLES below.
EXAMPLES:
Construction of Bessel functions with various orders and types::
sage: Bessel()
bessel_J
sage: Bessel(1)(x)
bessel_J(1, x)
sage: Bessel(1, 'Y')(x)
bessel_Y(1, x)
sage: Bessel(-2, 'Y')(x)
bessel_Y(-2, x)
sage: Bessel(typ='K')
bessel_K
sage: Bessel(0, typ='I')(x)
bessel_I(0, x)
Evaluation::
sage: f = Bessel(1)
sage: f(3.0)
0.339058958525936
sage: f(3)
bessel_J(1, 3)
sage: f(3).n(digits=50)
0.33905895852593645892551459720647889697308041819801
sage: g = Bessel(typ='J')
sage: g(1,3)
bessel_J(1, 3)
sage: g(2, 3+I).n()
0.634160370148554 + 0.0253384000032695*I
sage: abs(numerical_integral(1/pi*cos(3*sin(x)), 0.0, pi)[0] - Bessel(0, 'J')(3.0)) < 1e-15
True
Symbolic calculus::
sage: f(x) = Bessel(0, 'J')(x)
sage: derivative(f, x)
x |--> -1/2*bessel_J(1, x) + 1/2*bessel_J(-1, x)
sage: derivative(f, x, x)
x |--> 1/4*bessel_J(2, x) - 1/2*bessel_J(0, x) + 1/4*bessel_J(-2, x)
Verify that `J_0` satisfies Bessel's differential equation numerically
using the ``test_relation()`` method::
sage: y = bessel_J(0, x)
sage: diffeq = x^2*derivative(y,x,x) + x*derivative(y,x) + x^2*y == 0
sage: diffeq.test_relation(proof=False)
True
Conversion to other systems::
sage: x,y = var('x,y')
sage: f = maxima(Bessel(typ='K')(x,y))
sage: f.derivative('x')
%pi*csc(%pi*x)*('diff(bessel_i(-x,y),x,1)-'diff(bessel_i(x,y),x,1))/2-%pi*bessel_k(x,y)*cot(%pi*x)
sage: f.derivative('y')
-(bessel_k(x+1,y)+bessel_k(x-1,y))/2
Compute the particular solution to Bessel's Differential Equation that
satisfies `y(1) = 1` and `y'(1) = 1`, then verify the initial conditions
and plot it::
sage: y = function('y', x)
sage: diffeq = x^2*diff(y,x,x) + x*diff(y,x) + x^2*y == 0
sage: f = desolve(diffeq, y, [1, 1, 1]); f
(bessel_Y(1, 1) + bessel_Y(0, 1))*bessel_J(0, x)/(bessel_J(0,
1)*bessel_Y(1, 1) - bessel_J(1, 1)*bessel_Y(0, 1)) - (bessel_J(1,
1) + bessel_J(0, 1))*bessel_Y(0, x)/(bessel_J(0, 1)*bessel_Y(1, 1)
- bessel_J(1, 1)*bessel_Y(0, 1))
sage: f.subs(x=1).n() # numerical verification
1.00000000000000
sage: fp = f.diff(x)
sage: fp.subs(x=1).n()
1.00000000000000
sage: f.subs(x=1).simplify_full() # symbolic verification
1
sage: fp = f.diff(x)
sage: fp.subs(x=1).simplify_full()
1
sage: plot(f, (x,0,5))
Plotting::
sage: f(x) = Bessel(0)(x); f
x |--> bessel_J(0, x)
sage: plot(f, (x, 1, 10))
sage: plot([ Bessel(i, 'J') for i in range(5) ], 2, 10)
sage: G = Graphics()
sage: G += sum([ plot(Bessel(i), 0, 4*pi, rgbcolor=hue(sin(pi*i/10))) for i in range(5) ])
sage: show(G)
A recreation of Abramowitz and Stegun Figure 9.1::
sage: G = plot(Bessel(0, 'J'), 0, 15, color='black')
sage: G += plot(Bessel(0, 'Y'), 0, 15, color='black')
sage: G += plot(Bessel(1, 'J'), 0, 15, color='black', linestyle='dotted')
sage: G += plot(Bessel(1, 'Y'), 0, 15, color='black', linestyle='dotted')
sage: show(G, ymin=-1, ymax=1)
"""
# Determine the order and type of function from the arguments and keywords.
# These are recored in local variables: _type, _order, _system, _nargs.
_type = None
if len(args) == 0: # no order specified
_order = None
_nargs = 2
elif len(args) == 1: # order is specified
_order = args[0]
_nargs = 1
elif len(args) == 2: # both order and type are positional arguments
_order = args[0]
_type = args[1]
_nargs = 1
else:
from sage.misc.superseded import deprecation
deprecation(4102, 'precision argument is deprecated; algorithm '
'argument is currently deprecated, but will be '
'available as a named keyword in the future')
return _Bessel(*args, **kwds)
# check for type inconsistency
if _type is not None and 'typ' in kwds and _type != kwds['typ']:
raise ValueError("inconsistent types given")
# record the function type
if _type is None:
if 'typ' in kwds:
_type = kwds['typ']
else:
_type = 'J'
if not (_type in ['I', 'J', 'K', 'Y']):
raise ValueError("type must be one of I, J, K, Y")
# record the numerical evaluation system
if 'algorithm' in kwds:
_system = kwds['algorithm']
else: