pypeit.bspline.bspline module

Implements the bspline class

class pypeit.bspline.bspline.bspline(x, fullbkpt=None, nord=4, npoly=1, bkpt=None, bkspread=1.0, verbose=False, from_dict=None, **kwargs)[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

Fit coefficients

funcname

str

Function of fit

icoeff

numpy.ndarray

numpy.floating

??

mask

numpy.ndarray

numpy.bool

Mask

nord

int

Order of the bspline fit

npoly

int

Order of the bspline polynomial

xmax

float

Normalization for input data

xmin

float

Normalization for input data

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) – The data.

  • nord (int, optional) – To be documented.

  • npoly (int, optional) – To be documented.

  • bkpt (numpy.ndarray, optional) – To be documented.

  • bkspread (float, optional) – To be documented.

  • verbose (bool, optional.) – If True print extra information.

breakpoints

Breakpoints for bspline, spacing for these breakpoints are determined by keywords inputs;

nord

Order of bspline; [default=4]

npoly

Polynomial order to fit over 2nd variable (when specified as x2): [default=1]

mask

Output mask, set =1 for good points, =0 for bad points;

coeff

Output coefficient of the bspline;

icoeff

Cholesky band matrix used to solve for the bspline coefficients;

xmin

Normalization minimum for x2; [default max(xdata)]

xmax

Normalization maximum for x2; [default min(xdata)]

funcname

Function for the second variable; [default ‘legendre’]

from_dict

If not None, create a bspline from a dictionary created by to_dict(). [default ‘None’] It is possible to instantiate a bspline from a dict without the x data: new_bspline = bspline(None, from_dict=dictionary)

_bundle()[source]

Overload for the HDU name

Return type:

list

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': 'Fit coefficients', 'otype': <class 'numpy.ndarray'>}, 'funcname': {'descr': 'Function of fit', 'otype': <class 'str'>}, 'icoeff': {'atype': <class 'numpy.floating'>, 'descr': '??', 'otype': <class 'numpy.ndarray'>}, 'mask': {'atype': <class 'numpy.bool_'>, 'descr': 'Mask', 'otype': <class 'numpy.ndarray'>}, 'nord': {'descr': 'Order of the bspline fit', 'otype': <class 'int'>}, 'npoly': {'descr': 'Order of the bspline polynomial', 'otype': <class 'int'>}, 'xmax': {'descr': 'Normalization for input data', 'otype': <class 'float'>}, 'xmin': {'descr': 'Normalization for input data', 'otype': <class 'float'>}}

Provides the class data model. This is undefined in the abstract class and should be overwritten in the derived classes.

The format of the datamodel needed for each implementation of a DataContainer derived class is as follows.

The datamodel itself is a class attribute (i.e., it is a member of the class, not just of an instance of the class). The datamodel is a dictionary of dictionaries: Each key of the datamodel dictionary provides the name of a given datamodel element, and the associated item (dictionary) for the datamodel element provides the type and description information for that datamodel element. For each datamodel element, the dictionary item must provide:

  • otype: This is the type of the object for this datamodel item. E.g., for a float or a numpy.ndarray, you would set otype=float and otype=np.ndarray, respectively.

  • descr: This provides a text description of the datamodel element. This is used to construct the datamodel tables in the pypeit documentation.

If the object type is a numpy.ndarray, you should also provide the atype keyword that sets the type of the data contained within the array. E.g., for a floating point array containing an image, your datamodel could be simply:

datamodel = {'image' : dict(otype=np.ndarray, atype=float, descr='My image')}

More advanced examples are given in the top-level module documentation.

Currently, datamodel components are restricted to have otype that are tuple, int, float, numpy.integer, numpy.floating, numpy.ndarray, or astropy.table.Table objects. E.g., datamodel values for otype cannot be dict.

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

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]
to_dict()[source]

Write bspline attributes to a dict.

Attributes returned are: breakpoints, nord, npoly, mask, coeff, icoeff, xmin, xmax, and funcname.

Note

numpy.ndarray objects are converted to lists in the dictionary to make it JSON compatible.

Returns:

A dictionary with the above keys and items.

Return type:

dict

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'

Provides the string representation of the class version.

This is currently put to minimal use so far, but will used for I/O verification in the future.

Each derived class should provide a version to guard against data model changes during development.

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]