Skip to content

Commit 9230fa8

Browse files
committed
Use new scipy.special.assoc_legendre_p_all for getting Associated Legendre functions in SCFPotential
1 parent 4fc5db9 commit 9230fa8

File tree

1 file changed

+43
-13
lines changed

1 file changed

+43
-13
lines changed

galpy/potential/SCFPotential.py

+43-13
Original file line numberDiff line numberDiff line change
@@ -1,9 +1,17 @@
11
import hashlib
22

33
import numpy
4+
import scipy
45
from numpy.polynomial.legendre import leggauss
6+
from packaging.version import parse as parse_version
57
from scipy import integrate
6-
from scipy.special import gamma, gammaln, lpmn
8+
from scipy.special import gamma, gammaln
9+
10+
_SCIPY_VERSION = parse_version(scipy.__version__)
11+
if _SCIPY_VERSION < parse_version("1.15"): # pragma: no cover
12+
from scipy.special import lpmn
13+
else:
14+
from scipy.special import assoc_legendre_p_all
715

816
from ..util import conversion, coords
917
from ..util._optional_deps import _APY_LOADED
@@ -412,7 +420,17 @@ def _compute(self, funcTilde, R, z, phi):
412420
N, L, M = Acos.shape
413421
r, theta, phi = coords.cyl_to_spher(R, z, phi)
414422

415-
PP = lpmn(M - 1, L - 1, numpy.cos(theta))[0].T ##Get the Legendre polynomials
423+
## Get the Legendre polynomials
424+
if _SCIPY_VERSION < parse_version("1.15"): # pragma: no cover
425+
PP = lpmn(M - 1, L - 1, numpy.cos(theta))[0].T
426+
else:
427+
PP = numpy.swapaxes(
428+
assoc_legendre_p_all(L - 1, M - 1, numpy.cos(theta), branch_cut=2)[
429+
0, :, :M
430+
],
431+
0,
432+
1,
433+
).T
416434
func_tilde = funcTilde(r, N, L) ## Tilde of the function of interest
417435

418436
func = numpy.zeros(
@@ -515,9 +533,15 @@ def _computeforce(self, R, z, phi=0, t=0):
515533
dPhi_dphi = self._cached_dPhi_dphi
516534

517535
else:
518-
PP, dPP = lpmn(
519-
M - 1, L - 1, numpy.cos(theta)
520-
) ##Get the Legendre polynomials
536+
## Get the Legendre polynomials
537+
if _SCIPY_VERSION < parse_version("1.15"): # pragma: no cover
538+
PP, dPP = lpmn(M - 1, L - 1, numpy.cos(theta))
539+
else:
540+
PP, dPP = assoc_legendre_p_all(
541+
L - 1, M - 1, numpy.cos(theta), branch_cut=2, diff_n=1
542+
)
543+
PP = numpy.swapaxes(PP[:, :M], 0, 1)
544+
dPP = numpy.swapaxes(dPP[:, :M], 0, 1)
521545
PP = PP.T[None, :, :]
522546
dPP = dPP.T[None, :, :]
523547
phi_tilde = self._phiTilde(r, N, L)[:, :, numpy.newaxis]
@@ -924,14 +948,13 @@ def integrand(xi, costheta):
924948
r = _xiToR(xi, a)
925949
R = r * numpy.sqrt(1 - costheta**2.0)
926950
z = r * costheta
927-
Legandre = lpmn(0, L - 1, costheta)[0].T[numpy.newaxis, :, 0]
951+
if _SCIPY_VERSION < parse_version("1.15"): # pragma: no cover
952+
PP = lpmn(0, L - 1, costheta)[0].T[numpy.newaxis, :, 0]
953+
else:
954+
PP = assoc_legendre_p_all(L - 1, 0, costheta, branch_cut=2)[0].T
928955
dV = (1.0 + xi) ** 2.0 * numpy.power(1.0 - xi, -4.0)
929956
phi_nl = (
930-
a**3
931-
* (1.0 + xi) ** l
932-
* (1.0 - xi) ** (l + 1.0)
933-
* _C(xi, N, L)[:, :]
934-
* Legandre
957+
a**3 * (1.0 + xi) ** l * (1.0 - xi) ** (l + 1.0) * _C(xi, N, L)[:, :] * PP
935958
)
936959
param[0] = R
937960
param[1] = z
@@ -1089,15 +1112,22 @@ def integrand(xi, costheta, phi):
10891112
r = _xiToR(xi, a)
10901113
R = r * numpy.sqrt(1 - costheta**2.0)
10911114
z = r * costheta
1092-
Legandre = lpmn(L - 1, L - 1, costheta)[0].T[numpy.newaxis, :, :]
1115+
if _SCIPY_VERSION < parse_version("1.15"): # pragma: no cover
1116+
PP = lpmn(L - 1, L - 1, costheta)[0].T[numpy.newaxis, :, :]
1117+
else:
1118+
PP = numpy.swapaxes(
1119+
assoc_legendre_p_all(L - 1, L - 1, costheta, branch_cut=2)[0][:, :L],
1120+
0,
1121+
1,
1122+
).T[numpy.newaxis, :, :]
10931123
dV = (1.0 + xi) ** 2.0 * numpy.power(1.0 - xi, -4.0)
10941124

10951125
phi_nl = (
10961126
-(a**3)
10971127
* (1.0 + xi) ** l
10981128
* (1.0 - xi) ** (l + 1.0)
10991129
* _C(xi, N, L)[:, :, numpy.newaxis]
1100-
* Legandre
1130+
* PP
11011131
)
11021132

11031133
return (

0 commit comments

Comments
 (0)