Skip to content

nschloe/scipyx

Repository files navigation

scipyx

PyPi Version PyPI pyversions GitHub stars PyPi downloads

gh-actions codecov LGTM Code style: black

SciPy is large library used everywhere in scientific computing. That's why breaking backwards-compatibility comes as a significant cost and is almost always avoided, even if the API of some methods is arguably lacking. This package provides drop-in wrappers "fixing" those.

npx does the same for NumPy.

If you have a fix for a SciPy method that can't go upstream for some reason, feel free to PR here.

Krylov methods

import numpy as np
import scipy.sparse
import scipyx as spx

# create tridiagonal (-1, 2, -1) matrix
n = 100
data = -np.ones((3, n))
data[1] = 2.0
A = scipy.sparse.spdiags(data, [-1, 0, 1], n, n)
A = A.tocsr()
b = np.ones(n)


sol, info = spx.cg(A, b, tol=1.0e-10)
sol, info = spx.minres(A, b, tol=1.0e-10)
sol, info = spx.gmres(A, b, tol=1.0e-10)
sol, info = spx.bicg(A, b, tol=1.0e-10)
sol, info = spx.bicgstab(A, b, tol=1.0e-10)
sol, info = spx.cgs(A, b, tol=1.0e-10)
sol, info = spx.qmr(A, b, tol=1.0e-10)

sol is the solution of the linear system A @ x = b (or None if no convergence), and info contains some useful data, e.g., info.resnorms. The solution sol and all callback x have the shape of x0/b. The methods are wrappers around SciPy's iterative solvers.

Relevant issues:

Optimization

import scipyx as spx


def f(x):
    return (x ** 2 - 2) ** 2


x0 = 1.5
out = spx.minimize(f, x0)
print(out.x)

x0 = -3.2
x, _ = spx.leastsq(f, x0)
print(x)

In scipyx, all intermediate values x and the result from a minimization out.x will have the same shape as x0. (In SciPy, they always have shape (n,), no matter the input vector.)

Relevant issues:

Rolling Lagrange interpolation

import numpy as np
import scipyx as spx


x = np.linspace(0.0, 1.0, 11)
y = np.sin(7.0 * x)

poly = spx.interp_rolling_lagrange(x, y, order=3)

Given an array of coordinates x and an array of values y, you can use scipyx to compute a piecewise polynomial Lagrange approximation. The order + 1 closest coordinates x/y are considered for each interval.

Order 0 Order 1 Order 2

Jacobi elliptic functions with complex argument

SciPy supports Jacobi elliptic functions as ellipj. Unfortunately, only real-valued argument u and parameter m are allowed. scipyx expands support to complex-valued argument u.

import scipyx as spx

u = 1.0 + 2.0j
m = 0.8
# sn, cn, dn, ph = scipy.special.ellipj(x, m)  # not working
sn, cn, dn, ph = spx.ellipj(u, m)

Relevant bug reports:

License

This software is published under the BSD-3-Clause license.