bspline¶
- class pydl.pydlutils.bspline.bspline(x, nord=4, npoly=1, bkpt=None, bkspread=1.0, placed=None, bkspace=None, nbkpts=None, everyn=None)[source]¶
Bases:
objectB-spline class.
Functions in the idlutils bspline library are implemented as methods on this class.
- Parameters:
- x
numpy.ndarray The data.
- nord
int, optional The order of the B-spline. Default is 4, which is cubic.
- npoly
int, optional Polynomial order to fit over 2nd variable, if supplied. If not supplied the order is 1.
- bkpt
numpy.ndarray, optional An initial breakpoint vector. If omitted, breakpoints will be computed from other options.
- bkspread
float, optional Applies a scale factor to the difference between successive value of
bkpt.- placed
numpy.ndarray, optional Precalculated breakpoint positions.
- bkspace
float, optional Spacing of breakpoints in units of
x.- nbkpts
int, optional Number of breakpoints to span
xrange; minimum is 2 (the endpoints).- everyn
int, optional Spacing of breakpoints in good pixels.
- x
Notes
Although this code has been tested (in the sense that it reproduces the original idlutils code) and appears to work for a subset 1-dimensional cases, it should not be used for 2 dimensions or more exotic 1-dimensional cases.
- Attributes:
- breakpoints
The processed breakpoints for the fit, based on the initialization options.
- nord
The order of the B-spline. Default is 4, which is cubic.
- npoly
Polynomial order to fit over 2nd variable, if supplied. If not supplied the order is 1.
- mask
A mask for
breakpoints.- coeff
Internal attribute used by
fit().- icoeff
Internal attribute used by
fit().- xmin
Normalization minimum for the second variable in 2-dimensional fits.
- xmax
Normalization maximum for the second variable in 2-dimensional fits.
- funcname
For 2-dimensional fits, this is the function for the second variable. The default is ‘legendre’.
Init creates an object whose attributes are similar to the structure returned by the
create_bsplineset()function.Methods Summary
action(x[, x2])Construct banded B-spline matrix, with dimensions [ndata, bandwidth].
bsplvn(x, ileft)Calculates the value of all possibly nonzero B-splines at
xof a certain order.fit(xdata, ydata, invvar[, x2])Calculate a B-spline in the least-squares sense.
intrv(x)Find the segment between breakpoints which contain each value in the array
x.maskpoints(err)Perform simple logic of which breakpoints to mask.
value(x[, x2, action, lower, upper])Evaluate a B-spline at specified values.
Methods Documentation
- action(x, x2=None)[source]¶
Construct banded B-spline matrix, with dimensions [ndata, bandwidth].
- Parameters:
- x
numpy.ndarray Independent variable.
- x2
numpy.ndarray, optional Orthogonal dependent variable for 2d fits.
- x
- Returns:
tupleA tuple containing the B-spline action matrix; the ‘lower’ parameter, a list of pixel positions, each corresponding to the first occurence of position greater than breakpoint indx; and ‘upper’, Same as lower, but denotes the upper pixel positions.
- bsplvn(x, ileft)[source]¶
Calculates the value of all possibly nonzero B-splines at
xof a certain order.- Parameters:
- x
numpy.ndarray Independent variable.
- ileft
int Breakpoint segements that contain
x.
- x
- Returns:
numpy.ndarrayB-spline values.
- fit(xdata, ydata, invvar, x2=None)[source]¶
Calculate a B-spline in the least-squares sense.
Fit is based on two variables:
xdatawhich is sorted and spans a large range where breakpoints are requiredydatawhich can be described with a low order polynomial.- Parameters:
- xdata
numpy.ndarray Independent variable.
- ydata
numpy.ndarray Dependent variable.
- invvar
numpy.ndarray Inverse variance of
ydata.- x2
numpy.ndarray, optional Orthogonal dependent variable for 2d fits.
- xdata
- Returns:
tupleA tuple containing an integer error code, and the evaluation of the b-spline at the input values. An error code of -2 is a failure, -1 indicates dropped breakpoints, 0 is success, and positive integers indicate ill-conditioned breakpoints.
- intrv(x)[source]¶
Find the segment between breakpoints which contain each value in the array
x.The minimum breakpoint is
nbkptord - 1, and the maximum isnbkpt - nbkptord - 1.- Parameters:
- x
numpy.ndarray Data values, assumed to be monotonically increasing.
- x
- Returns:
numpy.ndarrayPosition of array elements with respect to breakpoints.
- maskpoints(err)[source]¶
Perform simple logic of which breakpoints to mask.
- Parameters:
- err
numpy.ndarray The list of indexes returned by the cholesky routines.
- err
- Returns:
intAn integer indicating the results of the masking. -1 indicates that the error points were successfully masked. -2 indicates failure; the calculation should be aborted.
Notes
The mask attribute is modified, assuming it is possible to create the mask.
- value(x, x2=None, action=None, lower=None, upper=None)[source]¶
Evaluate a B-spline at specified values.
- Parameters:
- x
numpy.ndarray Independent variable.
- x2
numpy.ndarray, optional Orthogonal dependent variable for 2d fits.
- action
numpy.ndarray, optional Action matrix to use. If not supplied it is calculated.
- lower
numpy.ndarray, optional If the action parameter is supplied, this parameter must also be supplied.
- upper
numpy.ndarray, optional If the action parameter is supplied, this parameter must also be supplied.
- x
- Returns:
tupleA tuple containing the results of the bspline evaluation and a mask indicating where the evaluation was good.