{ } Raw JSON

bundles / scipy latest / scipy / interpolate / _fitpack2 / LSQSphereBivariateSpline

class

scipy.interpolate._fitpack2:LSQSphereBivariateSpline

source: /scipy/interpolate/_fitpack2.py :1988

Signature

class   LSQSphereBivariateSpline ( theta phi r tt tp w = None eps = 1e-16 )

Members

Summary

Weighted least-squares bivariate spline approximation in spherical coordinates.

Extended Summary

Determines a smoothing bicubic spline according to a given set of knots in the theta and phi directions.

Parameters

theta, phi, r : array_like

1-D sequences of data points (order is not important). Coordinates must be given in radians. Theta must lie within the interval [0, pi], and phi must lie within the interval [0, 2pi].

tt, tp : array_like

Strictly ordered 1-D sequences of knots coordinates. Coordinates must satisfy 0 < tt[i] < pi, 0 < tp[i] < 2*pi.

w : array_like, optional

Positive 1-D sequence of weights, of the same length as theta, phi and r.

eps : float, optional

A threshold for determining the effective rank of an over-determined linear system of equations. eps should have a value within the open interval (0, 1), the default is 1e-16.

Notes

For more information, see the FITPACK_ site about this function.

Array API Standard Support

LSQSphereBivariateSpline is not in-scope for support of Python Array API Standard compatible backends other than NumPy.

See dev-arrayapi for more information.

Examples

Suppose we have global data on a coarse grid (the input data does not have to be on a grid):
from scipy.interpolate import LSQSphereBivariateSpline
import numpy as np
import matplotlib.pyplot as plt
theta = np.linspace(0, np.pi, num=7)
phi = np.linspace(0, 2*np.pi, num=9)
data = np.empty((theta.shape[0], phi.shape[0]))
data[:,0], data[0,:], data[-1,:] = 0., 0., 0.
data[1:-1,1], data[1:-1,-1] = 1., 1.
data[1,1:-1], data[-2,1:-1] = 1., 1.
data[2:-2,2], data[2:-2,-2] = 2., 2.
data[2,2:-2], data[-3,2:-2] = 2., 2.
data[3,3:-2] = 3.
data = np.roll(data, 4, 1)
We need to set up the interpolator object. Here, we must also specify the coordinates of the knots to use.
lats, lons = np.meshgrid(theta, phi)
knotst, knotsp = theta.copy(), phi.copy()
knotst[0] += .0001
knotst[-1] -= .0001
knotsp[0] += .0001
knotsp[-1] -= .0001
lut = LSQSphereBivariateSpline(lats.ravel(), lons.ravel(),
                               data.T.ravel(), knotst, knotsp)
As a first test, we'll see what the algorithm returns when run on the input coordinates
data_orig = lut(theta, phi)
Finally we interpolate the data to a finer grid
fine_lats = np.linspace(0., np.pi, 70)
fine_lons = np.linspace(0., 2*np.pi, 90)
data_lsq = lut(fine_lats, fine_lons)
fig = plt.figure()
ax1 = fig.add_subplot(131)
ax1.imshow(data, interpolation='nearest')
ax2 = fig.add_subplot(132)
ax2.imshow(data_orig, interpolation='nearest')
ax3 = fig.add_subplot(133)
ax3.imshow(data_lsq, interpolation='nearest')
plt.show()
fig-696c38780b0598dd.png

See also

BivariateSpline

a base class for bivariate splines.

LSQBivariateSpline

a bivariate spline using weighted least-squares fitting

RectBivariateSpline

a bivariate spline over a rectangular mesh.

RectSphereBivariateSpline

a bivariate spline over a rectangular mesh on a sphere

SmoothBivariateSpline

a smoothing bivariate spline through the given points

SmoothSphereBivariateSpline

a smoothing bivariate spline in spherical coordinates

UnivariateSpline

a smooth univariate spline to fit a given set of data points.

bisplev

a function to evaluate a bivariate B-spline and its derivatives

bisplrep

a function to find a bivariate B-spline representation of a surface

Aliases

  • scipy.interpolate.LSQSphereBivariateSpline