interpolated_coordinates.utils.spline

scipy Splines, with Units

A module for scipy splines classes with units support.

scipy [scipy], [Dierckx] splines do not support Quantity because they do not understand Unit. A standard workaround solution when one needs to interpolate is to strip quantities of their units, apply the interpolation, then add units back.

As an example:

>>> import numpy as np
>>> import astropy.units as u
>>> x = np.linspace(-3, 3, 50) * u.s
>>> y = 8 * u.m / (x.value**2 + 4)
>>> xs = np.linspace(-2, 2, 10) * u.s  # for evaluating spline

>>> from scipy.interpolate import InterpolatedUnivariateSpline
>>> spl = InterpolatedUnivariateSpline(x.to_value(u.s), y.to_value(u.m))
>>> spl(xs.to_value(u.s)) * u.m  # evaluate, adding back units
<Quantity [1.00000009, 1.24615404, 1.52830261, 1.79999996, 1.97560874,
           1.97560874, 1.79999996, 1.52830261, 1.24615404, 1.00000009] m>

This is fine, but a bit of a hassle. Instead, we can wrap the unit stripping / adding process into a unit-aware version of the spline interpolation classes.

The same example as above, but with the new class:

>>> from interpolated_coordinates.utils import InterpolatedUnivariateSplinewithUnits
>>> spl = InterpolatedUnivariateSplinewithUnits(x, y)
>>> spl(xs)
<Quantity [1.00000009, 1.24615404, 1.52830261, 1.79999996, 1.97560874,
           1.97560874, 1.79999996, 1.52830261, 1.24615404, 1.00000009] m>

Using InterpolatedUnivariateSplinewithUnits, interpolation with ndarray AND Quantity inputs just work.

Plotting this example:

(Source code, png, hires.png, pdf)

example spline plot.

References

[Dierckx]

Paul Dierckx, Curve and Surface Fitting with Splines, Oxford University Press, 1993

[scipy]

Virtanen, P., Gommers, R., Oliphant, M., Reddy, T., Cournapeau, E., Peterson, P., Weckesser, J., Walt, M., Wilson, J., Millman, N., Nelson, A., Jones, R., Larson, E., Carey, ., Feng, Y., Moore, J., Laxalde, D., Perktold, R., Henriksen, I., Quintero, C., Archibald, A., Pedregosa, P., & SciPy 1.0 Contributors (2020). SciPy 1.0: Fundamental Algorithms for Scientific Computing in Python. Nature Methods, 17, 261-272.

Reference/API

interpolated_coordinates.utils.spline Module

Classes

UnivariateSplinewithUnits(x, y[, w, bbox, ...])

1-D smoothing spline fit to a given set of data points.

InterpolatedUnivariateSplinewithUnits(x, y)

1-D interpolating spline for a given set of data points, with units.

LSQUnivariateSplinewithUnits(x, y, t[, w, ...])

1-D spline with explicit internal knots.

Class Inheritance Diagram

digraph inheritance70f3cd3701 { bgcolor=transparent; rankdir=LR; size="8.0, 12.0"; "InterpolatedUnivariateSpline" [URL="https://docs.scipy.org/doc/scipy/reference/generated/scipy.interpolate.InterpolatedUnivariateSpline.html#scipy.interpolate.InterpolatedUnivariateSpline",fillcolor=white,fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",fontsize=10,height=0.25,shape=box,style="setlinewidth(0.5),filled",target="_top",tooltip="1-D interpolating spline for a given set of data points."]; "UnivariateSpline" -> "InterpolatedUnivariateSpline" [arrowsize=0.5,style="setlinewidth(0.5)"]; "InterpolatedUnivariateSplinewithUnits" [URL="../api/interpolated_coordinates.utils.spline.InterpolatedUnivariateSplinewithUnits.html#interpolated_coordinates.utils.spline.InterpolatedUnivariateSplinewithUnits",fillcolor=white,fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",fontsize=10,height=0.25,shape=box,style="setlinewidth(0.5),filled",target="_top",tooltip="1-D interpolating spline for a given set of data points, with units."]; "UnivariateSplinewithUnits" -> "InterpolatedUnivariateSplinewithUnits" [arrowsize=0.5,style="setlinewidth(0.5)"]; "InterpolatedUnivariateSpline" -> "InterpolatedUnivariateSplinewithUnits" [arrowsize=0.5,style="setlinewidth(0.5)"]; "LSQUnivariateSpline" [URL="https://docs.scipy.org/doc/scipy/reference/generated/scipy.interpolate.LSQUnivariateSpline.html#scipy.interpolate.LSQUnivariateSpline",fillcolor=white,fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",fontsize=10,height=0.25,shape=box,style="setlinewidth(0.5),filled",target="_top",tooltip="1-D spline with explicit internal knots."]; "UnivariateSpline" -> "LSQUnivariateSpline" [arrowsize=0.5,style="setlinewidth(0.5)"]; "LSQUnivariateSplinewithUnits" [URL="../api/interpolated_coordinates.utils.spline.LSQUnivariateSplinewithUnits.html#interpolated_coordinates.utils.spline.LSQUnivariateSplinewithUnits",fillcolor=white,fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",fontsize=10,height=0.25,shape=box,style="setlinewidth(0.5),filled",target="_top",tooltip="1-D spline with explicit internal knots."]; "UnivariateSplinewithUnits" -> "LSQUnivariateSplinewithUnits" [arrowsize=0.5,style="setlinewidth(0.5)"]; "LSQUnivariateSpline" -> "LSQUnivariateSplinewithUnits" [arrowsize=0.5,style="setlinewidth(0.5)"]; "UnivariateSpline" [URL="https://docs.scipy.org/doc/scipy/reference/generated/scipy.interpolate.UnivariateSpline.html#scipy.interpolate.UnivariateSpline",fillcolor=white,fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",fontsize=10,height=0.25,shape=box,style="setlinewidth(0.5),filled",target="_top",tooltip="1-D smoothing spline fit to a given set of data points."]; "UnivariateSplinewithUnits" [URL="../api/interpolated_coordinates.utils.spline.UnivariateSplinewithUnits.html#interpolated_coordinates.utils.spline.UnivariateSplinewithUnits",fillcolor=white,fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",fontsize=10,height=0.25,shape=box,style="setlinewidth(0.5),filled",target="_top",tooltip="1-D smoothing spline fit to a given set of data points."]; "UnivariateSpline" -> "UnivariateSplinewithUnits" [arrowsize=0.5,style="setlinewidth(0.5)"]; }