pypeit.core.pydl module

pypeit.core.pydl.arange_ndim(dims, axis=0)[source]

Create an array where the values of the array are the index of the element along the selected axis.

Parameters:
  • dims (int, tuple, list) – Shape of the array to return along each dimension. The length of the tuple or list sets the number of dimensions in the returned array. If dims is an integer, this function is identical to numpy.arange.

  • axis (int, optional) – Axis along which to assign the index numbers.

Returns:

Integer array with the same shape as dims where each element equals its index number along axis.

Return type:

numpy.ndarray

Raises:

ValueError – Raised if axis is out of range for the provided dims.

Examples

>>> from pypeit.core.pydl import arange_ndim
>>> print(arange_ndim((4,3)))
[[0 0 0]
 [1 1 1]
 [2 2 2]
 [3 3 3]]
>>> print(arange_ndim((4,3), axis=1))
[[0 1 2]
 [0 1 2]
 [0 1 2]
 [0 1 2]]
class pypeit.core.pydl.chunks(ra, dec, minSize)[source]

Bases: object

chunks class

Functions for creating and manipulating spherical chunks are implemented as methods on this class.

assign(ra, dec, marginSize)[source]

Take the objects and the chunks (already defined in the constructor) and assign the objects to the appropriate chunks, with some leeway given by the parameter marginSize. Basically, at the end, each chunk should be associated with a list of the objects that belong to it.

chunkfriendsoffriends(ra, dec, chunkList, linkSep)[source]

Does friends-of-friends on the ra, dec that are defined by chunkList.

cosDecMin(i)[source]

Frequently used utility function.

friendsoffriends(ra, dec, linkSep)[source]

Friends-of-friends using chunked data.

get(ra, dec)[source]

Find the chunk to which a given point belongs.

getbounds(ra, dec, marginSize)[source]

Find the set of chunks a point (with margin) belongs to.

getraminmax(ra, raOffset)[source]

Utility function used by rarange.

rarange(ra, minSize)[source]

Finds the offset which yields the smallest raRange & returns both.

Notes

Warning

This is not (yet) well-defined for the case of only one point.

pypeit.core.pydl.djs_laxisnum(dims, iaxis=0)[source]

Returns an integer array where each element of the array is set equal to its index number in the specified axis.

Parameters:
  • dims (list) – Dimensions of the array to return.

  • iaxis (int, optional) – Index along this dimension.

Returns:

An array of indexes with dtype=int32.

Return type:

numpy.ndarray

Raises:

ValueError – If iaxis is greater than or equal to the number of dimensions, or if number of dimensions is greater than three.

Notes

For two or more dimensions, there is no difference between this routine and djs_laxisgen().

Examples

>>> from pypeit.core.pydl import djs_laxisnum
>>> print(djs_laxisnum((5,1), iaxis=0))
[[0]
 [1]
 [2]
 [3]
 [4]]
>>> print(djs_laxisnum(5))
[0 0 0 0 0]
>>> print(djs_laxisnum((5,1), iaxis=1))
[[0]
 [0]
 [0]
 [0]
 [0]]
>>> print(djs_laxisnum((1,5), iaxis=1))
[[0 1 2 3 4]]
>>> print(djs_laxisnum((4,3)))
[[0 0 0]
 [1 1 1]
 [2 2 2]
 [3 3 3]]
pypeit.core.pydl.djs_maskinterp(yval, mask, xval=None, axis=None, const=False)[source]

Interpolate over masked pixels in a vector, image or 3-D array.

Parameters:
  • yval (numpy.ndarray) – The input values

  • mask (numpy.ndarray) – The mask

  • xval (numpy.ndarray, optional) – If set, use these x values, otherwise use an array

  • axis (int, optional) – Must be set if yval has more than one dimension. If set to zero, interpolate along the first axis of the array, if set to one, interpolate along the second axis of the array, and so on.

  • const (bool, optional) – This value is passed to a helper function, djs_maskinterp1.

Returns:

The interpolated array.

Return type:

numpy.ndarray

pypeit.core.pydl.djs_maskinterp1(yval, mask, xval=None, const=False)[source]

Interpolate over a masked, 1-d array.

Parameters:
  • yval (numpy.ndarray) – The input values.

  • mask (numpy.ndarray) – The mask.

  • xval (numpy.ndarray, optional) – If set, use these x values, otherwise use an array.

  • const (bool, optional) – If set to True, bad values around the edges of the array will be set to a constant value. Because of the default behavior of numpy.interp(), this value actually makes no difference in the output.

Returns:

The yval array with masked values replaced by interpolated values.

Return type:

numpy.ndarray

pypeit.core.pydl.djs_reject(data, model, outmask=None, inmask=None, invvar=None, lower=None, upper=None, percentile=False, maxdev=None, maxrej=None, groupdim=None, groupsize=None, groupbadpix=False, grow=0, sticky=False, use_mad=False)[source]

Routine to reject points when doing an iterative fit to data.

Parameters:
  • data (numpy.ndarray) – The data

  • model (numpy.ndarray) – The model, must have the same number of dimensions as data.

  • outmask (numpy.ndarray, optional) – Output mask, generated by a previous call to djs_reject. If sticky=True, then bad points accumulate in this mask between calls. Otherwise, this mask is only used to determine if the rejection iterations are complete (e.g. to set qdone). Although this parameter is technically optional, it will almost always be set. If not supplied, this mask will be initialized to a mask that masks nothing and qdone will always be returned as True.

  • inmask (numpy.ndarray, optional) – Input mask. Bad points are marked with a value that evaluates to False. Must have the same number of dimensions as data. Points masked as bad “False” in the inmask will also always evaluate to “False” in the outmask

  • invvar (numpy.ndarray, optional) – Inverse variance of the data, used to reject points based on the values of upper and lower.

  • lower (int or float, optional) – If set, reject points with data < model - lower * sigm, where sigma = 1.0/sqrt(invvar)

  • upper (int or float, optional) – If set, reject points with data > model + upper * sigma, where sigma = 1.0/sqrt(invvar)

  • maxdev (int or float, optional) – If set, reject points with abs(data-model) > maxdev. It is permitted to set all three of lower, upper and maxdev.

  • maxrej (int or numpy.ndarray, optional) – Maximum number of points to reject in this iteration. If groupsize or groupdim are set to arrays, this should be an array as well.

  • groupdim (class: int) – Dimension along which to group the data; set to 1 to group along the 1st dimension, 2 for the 2nd dimension, etc. If data has shape [100,200], then setting GROUPDIM=2 is equivalent to grouping the data with groupsize=100. In either case, there are 200 groups, specified by [*,i]. NOT WELL TESTED IN PYTHON!

  • groupsize (class: int) – If this and maxrej are set, then reject a maximum of maxrej points per group of groupsize points, where the grouping is performed in the along the dimension of the data vector. (For use in curve fitting, one probably wants to make sure that data is sorted according to the indpeendent variable. For multi-dimensional arrays where one desires this grouping along each dimension, then groupdim should be set. If groupdim is also set, then this specifies sub-groups within that.

  • groupbadpix (bool, optional) – If set to True, consecutive sets of bad pixels are considered groups, overriding the values of groupsize.

  • grow (int, optional, default = 0) – If set to a non-zero integer, N, the N nearest neighbors of rejected pixels will also be rejected.

  • sticky (bool, optional) – If set to True then points rejected in outmask from a previous call to djs_reject are kept rejected. If set to False, if a fit (model) changes between iterations, points can alternate from being rejected to not rejected.

  • use_mad – It set to True, compute the median of the maximum absolute deviation between the data and use this for the rejection instead of the default which is to compute the standard deviation of the yarray - modelfit. Note that it is not possible to specify use_mad=True and also pass in values invvar, and the code will return an error if this is done.

Returns:

  • outmask (np.ndarray, boolean) – mask where rejected data values are False

  • qdone (boolean) – a value set to “True” if djs_reject believes there is no further rejection to be done. This will be set to “False” if the points marked as rejected in the outmask have changed. It will be set to “True” when the same points are rejected in outmask as from a previous call. It will also be set to “False” if model is set to None. Recall that outmask is also an optional input parameter. If it is not set, then qdone will simply return true, so outmask needs to be input from the previous iteration for the routine to do something meaningful.

Raises:

ValueError – If dimensions of various inputs do not match.

pypeit.core.pydl.func_fit(x, y, ncoeff, invvar=None, function_name='legendre', ia=None, inputans=None, inputfunc=None)[source]

Fit x, y positions to a functional form.

Parameters:
  • x (array-like) – X values (independent variable).

  • y (array-like) – Y values (dependent variable).

  • ncoeff (int) – Number of coefficients to fit.

  • invvar (array-like, optional) – Weight values; inverse variance.

  • function_name (str, optional) – Function name, default ‘legendre’.

  • ia (array-like, optional) – An array of bool of length ncoeff specifying free (True) and fixed (False) parameters.

  • inputans (array-like, optional) – An array of values of length ncoeff specifying the values of the fixed parameters.

  • inputfunc (array-like, optional) – Multiply the function fit by these values.

Returns:

Fit coefficients, length ncoeff; fitted values.

Return type:

tuple() of array-like

Raises:
  • KeyError – If an invalid function type is selected.

  • ValueError – If input dimensions do not agree.

pypeit.core.pydl.gcirc(ra1, dec1, ra2, dec2, units=2)[source]

Computes rigorous great circle arc distances.

Parameters:
  • ra1 (float or array-like) – RA and Dec of two points.

  • dec1 (float or array-like) – RA and Dec of two points.

  • ra2 (float or array-like) – RA and Dec of two points.

  • dec2 (float or array-like) – RA and Dec of two points.

  • units ({ 0, 1, 2 }, optional) –

    • units = 0: everything is already in radians

    • units = 1: RA in hours, dec in degrees, distance in arcsec.

    • units = 2: RA, dec in degrees, distance in arcsec (default)

Returns:

The angular distance. Units of the value returned depend on the input value of units.

Return type:

float or array-like

Notes

The formula below is the one best suited to handling small angular separations. See: http://en.wikipedia.org/wiki/Great-circle_distance

class pypeit.core.pydl.groups(coordinates, distance, separation='euclid')[source]

Bases: object

Group a set of objects (a list of coordinates in some space) based on a friends-of-friends algorithm

static euclid(x1, x2)[source]

Pythagorean theorem in Euclidean space with arbitrary number of dimensions.

static sphereradec(x1, x2)[source]

Separation of two points on a 2D-sphere, assuming they are in longitude-latitude or right ascension-declination form. Assumes everything is already in radians.

pypeit.core.pydl.spheregroup(ra, dec, linklength, chunksize=None)[source]

Perform friends-of-friends grouping given ra/dec coordinates.

Parameters:
  • ra (numpy.ndarray) – Arrays of coordinates to group in decimal degrees.

  • dec (numpy.ndarray) – Arrays of coordinates to group in decimal degrees.

  • linklength (float) – Linking length for the groups in decimal degrees.

  • chunksize (float, optional) – Break up the sphere into chunks of this size in decimal degrees.

Returns:

A tuple containing the group number of each object, the multiplicity of each group, the first member of each group, and the next member of the group for each object.

Return type:

tuple()

Raises:

raise PypeItError – If the array of coordinates only contains one point.

Notes

It is important that chunksize >= 4 * linklength. This is enforced.

Warning

Behavior at the poles is not well tested.

pypeit.core.pydl.spherematch(ra1, dec1, ra2, dec2, matchlength, chunksize=None, maxmatch=1)[source]

Match points on a sphere.

Parameters:
  • ra1 (numpy.ndarray) – The sets of coordinates to match. Assumed to be in decimal degrees

  • dec1 (numpy.ndarray) – The sets of coordinates to match. Assumed to be in decimal degrees

  • ra2 (numpy.ndarray) – The sets of coordinates to match. Assumed to be in decimal degrees

  • dec2 (numpy.ndarray) – The sets of coordinates to match. Assumed to be in decimal degrees

  • matchlength (float) – Two points closer than this separation are matched. Assumed to be in decimal degrees.

  • chunksize (float, optional) – Value to pass to chunk assignment.

  • maxmatch (int, optional) – Allow up to maxmatch matches per coordinate. Default 1. If set to zero, All possible matches will be returned.

Returns:

A tuple containing the indices into the first set of points, the indices into the second set of points and the match distance in decimal degrees.

Return type:

tuple()

Notes

If you have sets of coordinates that differ in size, call this function with the larger list first. This exploits the inherent asymmetry in the underlying code to reduce memory use.

Warning

Behavior at the poles is not well tested.