bundles / scipy 1.17.1 / scipy / interpolate / _fitpack2 / RectSphereBivariateSpline
class
scipy.interpolate._fitpack2:RectSphereBivariateSpline
Signature
class RectSphereBivariateSpline ( u , v , r , s = 0.0 , pole_continuity = False , pole_values = None , pole_exact = False , pole_flat = False ) Members
Summary
Bivariate spline approximation over a rectangular mesh on a sphere.
Extended Summary
Can be used for smoothing data.
Parameters
u: array_like1-D array of colatitude coordinates in strictly ascending order. Coordinates must be given in radians and lie within the open interval
(0, pi).v: array_like1-D array of longitude coordinates in strictly ascending order. Coordinates must be given in radians. First element (
v[0]) must lie within the interval[-pi, pi). Last element (v[-1]) must satisfyv[-1] <= v[0] + 2*pi.r: array_like2-D array of data with shape
(u.size, v.size).s: float, optionalPositive smoothing factor defined for estimation condition (
s=0is for interpolation).pole_continuity: bool or (bool, bool), optionalOrder of continuity at the poles
u=0(pole_continuity[0]) andu=pi(pole_continuity[1]). The order of continuity at the pole will be 1 or 0 when this is True or False, respectively. Defaults to False.pole_values: float or (float, float), optionalData values at the poles
u=0andu=pi. Either the whole parameter or each individual element can be None. Defaults to None.pole_exact: bool or (bool, bool), optionalData value exactness at the poles
u=0andu=pi. If True, the value is considered to be the right function value, and it will be fitted exactly. If False, the value will be considered to be a data value just like the other data values. Defaults to False.pole_flat: bool or (bool, bool), optionalFor the poles at
u=0andu=pi, specify whether or not the approximation has vanishing derivatives. Defaults to False.
Notes
Currently, only the smoothing spline approximation (iopt[0] = 0 and iopt[0] = 1 in the FITPACK routine) is supported. The exact least-squares spline approximation is not implemented yet.
When actually performing the interpolation, the requested v values must lie within the same length 2pi interval that the original v values were chosen from.
For more information, see the FITPACK_ site about this function.
Array API Standard Support
RectSphereBivariateSpline 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 gridimport numpy as np lats = np.linspace(10, 170, 9) * np.pi / 180. lons = np.linspace(0, 350, 18) * np.pi / 180. data = np.dot(np.atleast_2d(90. - np.linspace(-80., 80., 18)).T, np.atleast_2d(180. - np.abs(np.linspace(0., 350., 9)))).T✓
new_lats = np.linspace(1, 180, 180) * np.pi / 180 new_lons = np.linspace(1, 360, 360) * np.pi / 180 new_lats, new_lons = np.meshgrid(new_lats, new_lons)✓
from scipy.interpolate import RectSphereBivariateSpline lut = RectSphereBivariateSpline(lats, lons, data)✓
data_interp = lut.ev(new_lats.ravel(), new_lons.ravel()).reshape((360, 180)).T✓
import matplotlib.pyplot as plt fig = plt.figure() ax1 = fig.add_subplot(211)✓
ax1.imshow(data, interpolation='nearest')
✗ax2 = fig.add_subplot(212)
✓ax2.imshow(data_interp, interpolation='nearest')
✗plt.show()
✓
fig2 = plt.figure() s = [3e9, 2e9, 1e9, 1e8]✓
for idx, sval in enumerate(s, 1): lut = RectSphereBivariateSpline(lats, lons, data, s=sval) data_interp = lut.ev(new_lats.ravel(), new_lons.ravel()).reshape((360, 180)).T ax = fig2.add_subplot(2, 2, idx) ax.imshow(data_interp, interpolation='nearest') ax.set_title(f"s = {sval:g}")✗
plt.show()
✓
See also
- BivariateSpline
a base class for bivariate splines.
- LSQBivariateSpline
a bivariate spline using weighted least-squares fitting
- LSQSphereBivariateSpline
a bivariate spline in spherical coordinates using weighted least-squares fitting
- RectBivariateSpline
a bivariate spline over a rectangular mesh.
- 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.RectSphereBivariateSpline