qnm.angular

Solve the angular Teukolsky equation via spectral decomposition.

TODO Documentation.

Functions

C_and_sep_const_closest(A0, s, c, m, l_max) Get a single eigenvalue and eigenvector of decomposition matrix, where the eigenvalue is closest to some guess A0.
M_matrix(s, c, m, l_max) Spherical-spheroidal decomposition matrix truncated at l_max.
M_matrix_elem(s, c, m, l, lprime) The (l, lprime) matrix element from the spherical-spheroidal decomposition matrix from Eq.
calA(s, l, m) Eq.
calB(s, l, m) Eq.
calC(s, l, m) Eq.
calD(s, l, m) Eq.
calE(s, l, m) Eq.
calF(s, l, m) Eq.
calG(s, l, m) Eq.
calH(s, l, m) Eq.
give_M_matrix_elem_ufunc(s, c, m) Gives ufunc that implements matrix elements of the spherical-spheroidal decomposition matrix.
l_min(s, m) Minimum allowed value of l for a given s, m.
sep_const_closest(A0, s, c, m, l_max) Gives the separation constant that is closest to A0.
sep_consts(s, c, m, l_max) Finds eigenvalues of decomposition matrix, i.e.
swsphericalh_A(s, l, m) Angular separation constant at a=0.
qnm.angular.C_and_sep_const_closest(A0, s, c, m, l_max)[source]

Get a single eigenvalue and eigenvector of decomposition matrix, where the eigenvalue is closest to some guess A0.

Parameters:
A0: complex

Value close to the desired separation constant.

s: int

Spin-weight of interest

c: complex

Oblateness of spheroidal harmonic

m: int

Magnetic quantum number

l_max: int

Maximum angular quantum number

Returns:
complex, complex ndarray

The first element of the tuple is the eigenvalue that is closest in value to A0. The second element of the tuple is the corresponding eigenvector.

qnm.angular.M_matrix(s, c, m, l_max)[source]

Spherical-spheroidal decomposition matrix truncated at l_max.

Parameters:
s: int

Spin-weight of interest

c: complex

Oblateness of the spheroidal harmonic

m: int

Magnetic quantum number

l_max: int

Maximum angular quantum number

Returns:
complex ndarray

Decomposition matrix

qnm.angular.M_matrix_elem(s, c, m, l, lprime)[source]

The (l, lprime) matrix element from the spherical-spheroidal decomposition matrix from Eq. (55).

Parameters:
s: int

Spin-weight of interest

c: complex

Oblateness of the spheroidal harmonic

m: int

Magnetic quantum number

l: int

Angular quantum number of interest

lprime: int

Primed quantum number of interest

Returns:
complex

Matrix element M_{l, lprime}

qnm.angular.calA(s, l, m)[source]

Eq. (53a)

qnm.angular.calB(s, l, m)[source]

Eq. (53c)

qnm.angular.calC(s, l, m)[source]

Eq. (53e)

qnm.angular.calD(s, l, m)[source]

Eq. (53b)

qnm.angular.calE(s, l, m)[source]

Eq. (53d)

qnm.angular.calF(s, l, m)[source]

Eq. (52b)

qnm.angular.calG(s, l, m)[source]

Eq. (52c)

qnm.angular.calH(s, l, m)[source]

Eq. (52d)

qnm.angular.give_M_matrix_elem_ufunc(s, c, m)[source]

Gives ufunc that implements matrix elements of the spherical-spheroidal decomposition matrix.

Parameters:
s: int

Spin-weight of interest

c: complex

Oblateness of the spheroidal harmonic

m: int

Magnetic quantum number

Returns:
ufunc

Implements elements of M matrix

qnm.angular.l_min(s, m)[source]

Minimum allowed value of l for a given s, m.

The formula is l_min = max(|m|,|s|).

Parameters:
s: int

Spin-weight of interest

m: int

Magnetic quantum number

Returns:
int

l_min

qnm.angular.sep_const_closest(A0, s, c, m, l_max)[source]

Gives the separation constant that is closest to A0.

Parameters:
A0: complex

Value close to the desired separation constant.

s: int

Spin-weight of interest

c: complex

Oblateness of spheroidal harmonic

m: int

Magnetic quantum number

l_max: int

Maximum angular quantum number

Returns:
complex

Separation constant that is the closest to A0.

qnm.angular.sep_consts(s, c, m, l_max)[source]

Finds eigenvalues of decomposition matrix, i.e. the separation constants, As.

Parameters:
s: int

Spin-weight of interest

c: complex

Oblateness of spheroidal harmonic

m: int

Magnetic quantum number

l_max: int

Maximum angular quantum number

Returns:
complex ndarray

Eigenvalues of spherical-spheroidal decomposition matrix

qnm.angular.swsphericalh_A(s, l, m)[source]

Angular separation constant at a=0.

Eq. (50). Has no dependence on m. The formula is
A_0 = l(l+1) - s(s+1)
Parameters:
s: int

Spin-weight of interest

l: int

Angular quantum number of interest

m: int

Magnetic quantum number, ignored

Returns:
int

Value of A(a=0) = l(l+1) - s(s+1)