scipy.special.jve#

scipy.special.jve(v, z, out=None) = <ufunc 'jve'>#

Exponentially scaled Bessel function of the first kind of order v.

Defined as:

jve(v, z) = jv(v, z) * exp(-abs(z.imag))
Parameters:
varray_like

Order (float).

zarray_like

Argument (float or complex).

outndarray, optional

Optional output array for the function values

Returns:
Jscalar or ndarray

Value of the exponentially scaled Bessel function.

See also

jv

Unscaled Bessel function of the first kind

Notes

For positive v values, the computation is carried out using the AMOS [1] zbesj routine, which exploits the connection to the modified Bessel function \(I_v\),

\[ \begin{align}\begin{aligned}J_v(z) = \exp(v\pi\imath/2) I_v(-\imath z)\qquad (\Im z > 0)\\J_v(z) = \exp(-v\pi\imath/2) I_v(\imath z)\qquad (\Im z < 0)\end{aligned}\end{align} \]

For negative v values the formula,

\[J_{-v}(z) = J_v(z) \cos(\pi v) - Y_v(z) \sin(\pi v)\]

is used, where \(Y_v(z)\) is the Bessel function of the second kind, computed using the AMOS routine zbesy. Note that the second term is exactly zero for integer v; to improve accuracy the second term is explicitly omitted for v values such that v = floor(v).

Exponentially scaled Bessel functions are useful for large arguments z: for these, the unscaled Bessel functions can easily under-or overflow.

References

[1]

Donald E. Amos, “AMOS, A Portable Package for Bessel Functions of a Complex Argument and Nonnegative Order”, http://netlib.org/amos/

Examples

Compare the output of jv and jve for large complex arguments for z by computing their values for order v=1 at z=1000j. We see that jv overflows but jve returns a finite number:

>>> import numpy as np
>>> from scipy.special import jv, jve
>>> v = 1
>>> z = 1000j
>>> jv(v, z), jve(v, z)
((inf+infj), (7.721967686709077e-19+0.012610930256928629j))

For real arguments for z, jve returns the same as jv.

>>> v, z = 1, 1000
>>> jv(v, z), jve(v, z)
(0.004728311907089523, 0.004728311907089523)

The function can be evaluated for several orders at the same time by providing a list or NumPy array for v:

>>> jve([1, 3, 5], 1j)
array([1.27304208e-17+2.07910415e-01j, -4.99352086e-19-8.15530777e-03j,
       6.11480940e-21+9.98657141e-05j])

In the same way, the function can be evaluated at several points in one call by providing a list or NumPy array for z:

>>> jve(1, np.array([1j, 2j, 3j]))
array([1.27308412e-17+0.20791042j, 1.31814423e-17+0.21526929j,
       1.20521602e-17+0.19682671j])

It is also possible to evaluate several orders at several points at the same time by providing arrays for v and z with compatible shapes for broadcasting. Compute jve for two different orders v and three points z resulting in a 2x3 array.

>>> v = np.array([[1], [3]])
>>> z = np.array([1j, 2j, 3j])
>>> v.shape, z.shape
((2, 1), (3,))
>>> jve(v, z)
array([[1.27304208e-17+0.20791042j,  1.31810070e-17+0.21526929j,
        1.20517622e-17+0.19682671j],
       [-4.99352086e-19-0.00815531j, -1.76289571e-18-0.02879122j,
        -2.92578784e-18-0.04778332j]])