pypeit.bspline.bspline module

Implements the bspline class

class pypeit.bspline.bspline.bspline(x=None, fullbkpt=None, nord=4, npoly=1, bkpt=None, bkspread=1.0, bkspace=None, nbkpts=None, everyn=None, funcname='legendre')[source]

Bases: DataContainer

Bspline class.

Functions in the bspline library are implemented as methods on this class.

The datamodel attributes are:

Version: 1.0.0

Attribute

Type

Array Type

Description

breakpoints

numpy.ndarray

numpy.floating

Breakpoint locations

coeff

numpy.ndarray

numpy.floating

Output fit coefficients

funcname

str

Function type for the 2nd variable (when x2 is specified)

icoeff

numpy.ndarray

numpy.floating

Cholesky band matrix used to solve for the bspline coefficients

mask

numpy.ndarray

numpy.bool

Output mask

nord

int

Order of the bspline fit

npoly

int

Order of polynomial to fit over 2nd variable (when x2 is specified)

xmax

float

Normalization maximum for x2

xmin

float

Normalization minimum for x2

When written to an output-file HDU, all numpy.ndarray elements are bundled into an astropy.io.fits.BinTableHDU, and the other elements are written as header keywords. Any datamodel elements that are None are not included in the output.

Parameters:
  • x (numpy.ndarray, optional) – Independent variable for the definition of the b-spline. If None, fullbkpt must be provided.

  • fullbkpt (numpy.ndarray, optional) – The full set of breakpoints. The input vector is sorted and cast as a float. If the length of the vector is less than twice nord, it is also padded with nord-1 values, as needed for the construction of the b-spline. If None, x must be provided.

  • nord (int, optional) – Order of the b-spline.

  • npoly (int, optional) – Polynomial order to fit over 2nd variable (when specified using x2; see fit()).

  • bkpt (numpy.ndarray, optional) – A precalculated set of breakpoints within the range of x. The input vector is sorted and any points beyond the range of x are omitted. If only one breakpoint is provided (or not omitted), one breakpoint is set at each end of x. If the breakpoints do not cover the full range, the first and last breakpoints are moved such that they do. If None, the breakpoints are determined using the keywords below.

  • bkspread (float, optional) – Scale factor for the separation between breakpoints.

  • bkspace (float, optional) – Defines the separation between breakpoints. If provided, nbkpts and everyn are ignored.

  • nbkpts (int, optional) – Defines the number of breakpoints used to span the full range of x. Only used if bkspace is None. If provided, everyn is ignored.

  • everyn (int, float, optional) – Places a breakpoint at every Nth value of x. Only used if bkspace and nbkpts are both None.

  • funcname (str, optional) – Function for the second variable (when specified using x2; see fit()).

_bundle()[source]

Overload the base class method (see _bundle()) to set the HDU name explicitly to BSPLINE.

static _fill_bkpt(bkpt, nord, bkspread)[source]

Helper function used to pad the breakpoint vector according to the order of the b-spline.

Parameters:
  • bkpt (numpy.ndarray) – The current set of breakpoints.

  • nord (int) – Order of the b-spline.

  • bkspread (float) – Scale factor for the separation between breakpoints.

Returns:

The padded set of breakpoints (typically fullbkpt as used by the class).

Return type:

numpy.ndarray

action(x, x2=None)[source]

Construct banded bspline matrix, with dimensions [ndata, bandwidth].

Parameters:
Returns:

A 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.

Return type:

tuple

bsplvn(x, ileft)[source]

To be documented.

Parameters:
Returns:

vnikx – To be documented.

Return type:

numpy.ndarray

copy()[source]

Return a copied instance of the object.

datamodel = {'breakpoints': {'atype': <class 'numpy.floating'>, 'descr': 'Breakpoint locations', 'otype': <class 'numpy.ndarray'>}, 'coeff': {'atype': <class 'numpy.floating'>, 'descr': 'Output fit coefficients', 'otype': <class 'numpy.ndarray'>}, 'funcname': {'descr': 'Function type for the 2nd variable (when x2 is specified)', 'otype': <class 'str'>}, 'icoeff': {'atype': <class 'numpy.floating'>, 'descr': 'Cholesky band matrix used to solve for the bspline coefficients', 'otype': <class 'numpy.ndarray'>}, 'mask': {'atype': <class 'numpy.bool'>, 'descr': 'Output mask', 'otype': <class 'numpy.ndarray'>}, 'nord': {'descr': 'Order of the bspline fit', 'otype': <class 'int'>}, 'npoly': {'descr': 'Order of polynomial to fit over 2nd variable (when x2 is specified)', 'otype': <class 'int'>}, 'xmax': {'descr': 'Normalization maximum for x2', 'otype': <class 'float'>}, 'xmin': {'descr': 'Normalization minimum for x2', 'otype': <class 'float'>}}

Datamodel components.

fit(xdata, ydata, invvar, x2=None)[source]

Calculate a B-spline in the least-squares sense.

Fit is based on two variables: x which is sorted and spans a large range where bkpts are required y which can be described with a low order polynomial.

Parameters:
Returns:

A 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.

Return type:

tuple

static get_breakpoints(x=None, bkpt=None, fullbkpt=None, nord=4, bkspread=1.0, bkspace=None, nbkpts=None, everyn=None)[source]

Generate the set of breakpoints for the b-spline.

Parameters:
  • x (numpy.ndarray, optional) – Independent variable for the definition of the b-spline. If None, fullbkpt must be provided.

  • bkpt (numpy.ndarray, optional) – A precalculated set of breakpoints within the range of x. The input vector is sorted and any points beyond the range of x are omitted. If only one breakpoint is provided (or not omitted), one breakpoint is set at each end of x. If the breakpoints do not cover the full range, the first and last breakpoints are moved such that they do. If None, the breakpoints are determined using the keywords below.

  • fullbkpt (numpy.ndarray, optional) – The full set of breakpoints. The input vector is sorted and cast as a float. If the length of the vector is less than twice nord, it is also padded with nord-1 values, as needed for the construction of the b-spline. If None, x must be provided.

  • nord (int, optional) – Order of the b-spline.

  • bkspread (float, optional) – Scale factor for the separation between breakpoints.

  • bkspace (float, optional) – Defines the separation between breakpoints. If provided, nbkpts and everyn are ignored.

  • nbkpts (int, optional) – Defines the number of breakpoints used to span the full range of x. Only used if bkspace is None. If provided, everyn is ignored.

  • everyn (int, float, optional) – Places a breakpoint at every Nth value of x. Only used if bkspace and nbkpts are both None.

Returns:

Vector with the breakpoints

Return type:

numpy.ndarray

Raises:

ValueError – Raised if neither fullbkpt nor x are provided.

maskpoints(err)[source]

Perform simple logic of which breakpoints to mask.

Parameters:

err (numpy.ndarray, int) – The list of indexes returned by the cholesky routines. This is indexed to the set of currently good breakpoints (i.e. self.mask=True) and the first nord are skipped.

Returns:

An integer indicating the results of the masking. -1 indicates that the error points were successfully masked. -2 indicates failure; the calculation should be aborted.

Return type:

int

Notes

The mask attribute is modified, assuming it is possible to create the mask.

reinit_coeff()[source]
value(x, x2=None, action=None, lower=None, upper=None)[source]

Evaluate a bspline 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.

Returns:

  • yfit (numpy.ndarray) – Results of the bspline evaluation

  • mask (numpy.ndarray) – Mask indicating where the evaluation was good (i.e., True is good).

version = '1.0.0'

Datamodel version

workit(xdata, ydata, invvar, action, lower, upper)[source]

An internal routine for bspline_extract and bspline_radial which solve a general banded correlation matrix which is represented by the variable “action”. This routine only solves the linear system once, and stores the coefficients in sset. A non-zero return value signifies a failed inversion

Parameters:
  • xdata (numpy.ndarray) – Independent variable.

  • ydata (numpy.ndarray) – Dependent variable.

  • invvar (numpy.ndarray) – Inverse variance of ydata.

  • action (numpy.ndarray) – Banded correlation matrix

  • lower (numpy.ndarray) – A list of pixel positions, each corresponding to the first occurence of position greater than breakpoint indx

  • upper (numpy.ndarray) – Same as lower, but denotes the upper pixel positions

Returns:

  • success (int) – Method error code: 0 is good; -1 is dropped breakpoints, try again; -2 is failure, should abort.

  • yfit (numpy.ndarray) – Evaluation of the b-spline yfit at the input values.

pypeit.bspline.bspline.uniq(x, index=None)[source]

Return the indices of the last occurrence of the unique elements in a sorted array.

The input vector must be sorted before being passed to this function. This can be done by sorting x directly or by passing the array that sorts x (index).

Replicates the IDL UNIQ() function.

Parameters:
  • x (array-like) – Search this array for unique items.

  • index (array-like, optional) –

    This array provides the array subscripts that sort x. That is:

    index = np.argsort(x)
    

Returns:

result – The indices of the last occurence in x of its unique values.

Return type:

numpy.ndarray

Notes

Given a sorted array, and assuming that there is a set of adjacent identical items, uniq() will return the subscript of the last unique item. This charming feature is retained for reproducibility.

References

http://www.harrisgeospatial.com/docs/uniq.html

Speed improvement thanks to discussion here: https://stackoverflow.com/questions/47495510/numpy-in-a-sorted-list-find-the-first-and-the-last-index-for-each-unique-value

Examples

>>> import numpy as np
>>> from pypeit.core.pydl import uniq
>>> data = np.array([ 1, 2, 3, 1, 5, 6, 1, 7, 3, 2, 5, 9, 11, 1 ])
>>> print(uniq(np.sort(data)))
[ 3  5  7  9 10 11 12 13]