TessTargetPixelFile

class lightkurve.targetpixelfile.TessTargetPixelFile(path, quality_bitmask='default', **kwargs)

Bases: lightkurve.targetpixelfile.TargetPixelFile

Represents pixel data products created by NASA’s TESS pipeline.

This class enables extraction of custom light curves and centroid positions.

Parameters
pathstr

Path to a Kepler Target Pixel (FITS) File.

quality_bitmaskstr or int

Bitmask specifying quality flags of cadences that should be ignored.

kwargsdict

Keyword arguments passed to astropy.io.fits.open().

Attributes Summary

astropy_time

Returns an AstroPy Time object for all good-quality cadences.

background_mask

Returns the background mask used by the TESS pipeline.

cadenceno

Return the cadence number for all good-quality cadences.

camera

TESS Camera number (‘CAMERA’ header keyword).

ccd

TESS CCD number (‘CCD’ header keyword).

column

CCD pixel column number (‘1CRV5P’ header keyword).

dec

Declination of target (‘DEC_OBJ’ header keyword).

flux

Returns the flux for all good-quality cadences.

flux_bkg

Returns the background flux for all good-quality cadences.

flux_bkg_err

flux_err

Returns the flux uncertainty for all good-quality cadences.

hdu

header

Returns the header of the primary extension.

mission

nan_time_mask

Returns a boolean mask flagging cadences whose time is nan.

pipeline_mask

Returns the optimal aperture mask used by the pipeline.

pos_corr1

Returns the column position correction.

pos_corr2

Returns the row position correction.

quality

Returns the quality flag integer of every good cadence.

ra

Right Ascension of target (‘RA_OBJ’ header keyword).

row

CCD pixel row number (‘2CRV5P’ header keyword).

sector

TESS Sector number (‘SECTOR’ header keyword).

shape

Return the cube dimension shape.

time

Returns the time for all good-quality cadences.

wcs

Returns an astropy.wcs.WCS object with the World Coordinate System solution for the target pixel file.

Methods Summary

create_threshold_mask(self[, threshold, …])

Returns an aperture mask creating using the thresholding method.

cutout(self[, center, size])

Cut a rectangle out of the Target Pixel File.

estimate_centroids(self[, aperture_mask, method])

Returns the flux center of an object inside aperture_mask.

extract_aperture_photometry(self[, …])

Performs aperture photometry.

get_bkg_lightcurve(self[, aperture_mask])

get_coordinates(self[, cadence])

Returns two 3D arrays of RA and Dec values in decimal degrees.

get_keyword(self, keyword[, hdu, default])

Returns a header keyword value.

interact(self[, notebook_url, max_cadences, …])

Display an interactive Jupyter Notebook widget to inspect the pixel data.

interact_sky(self[, notebook_url, …])

Display a Jupyter Notebook widget showing Gaia DR2 positions on top of the pixels.

plot(self[, ax, frame, cadenceno, bkg, …])

Plot the pixel data for a single frame (i.e.

show_properties(self)

Prints a description of all non-callable attributes.

to_corrector(self[, method])

Returns a Corrector instance to remove systematics.

to_fits(self[, output_fn, overwrite])

Writes the TPF to a FITS file on disk.

to_lightcurve(self[, method])

Performs photometry on the pixel data and returns a LightCurve object.

Attributes Documentation

astropy_time

Returns an AstroPy Time object for all good-quality cadences.

background_mask

Returns the background mask used by the TESS pipeline.

cadenceno

Return the cadence number for all good-quality cadences.

camera

TESS Camera number (‘CAMERA’ header keyword).

ccd

TESS CCD number (‘CCD’ header keyword).

column

CCD pixel column number (‘1CRV5P’ header keyword).

dec

Declination of target (‘DEC_OBJ’ header keyword).

flux

Returns the flux for all good-quality cadences.

flux_bkg

Returns the background flux for all good-quality cadences.

flux_bkg_err
flux_err

Returns the flux uncertainty for all good-quality cadences.

hdu
header

Returns the header of the primary extension.

mission
nan_time_mask

Returns a boolean mask flagging cadences whose time is nan.

pipeline_mask

Returns the optimal aperture mask used by the pipeline.

pos_corr1

Returns the column position correction.

pos_corr2

Returns the row position correction.

quality

Returns the quality flag integer of every good cadence.

ra

Right Ascension of target (‘RA_OBJ’ header keyword).

row

CCD pixel row number (‘2CRV5P’ header keyword).

sector

TESS Sector number (‘SECTOR’ header keyword).

shape

Return the cube dimension shape.

time

Returns the time for all good-quality cadences.

wcs

Returns an astropy.wcs.WCS object with the World Coordinate System solution for the target pixel file.

Returns
wastropy.wcs.WCS object

WCS solution

Methods Documentation

create_threshold_mask(self, threshold=3, reference_pixel='center')

Returns an aperture mask creating using the thresholding method.

This method will identify the pixels in the TargetPixelFile which show a median flux that is brighter than threshold times the standard deviation above the overall median. The standard deviation is estimated in a robust way by multiplying the Median Absolute Deviation (MAD) with 1.4826.

If the thresholding method yields multiple contiguous regions, then only the region closest to the (col, row) coordinate specified by reference_pixel is returned. For exmaple, reference_pixel=(0, 0) will pick the region closest to the bottom left corner. By default, the region closest to the center of the mask will be returned. If reference_pixel=None then all regions will be returned.

Parameters
thresholdfloat

A value for the number of sigma by which a pixel needs to be brighter than the median flux to be included in the aperture mask.

reference_pixel: (int, int) tuple, ‘center’, or None

(col, row) pixel coordinate closest to the desired region. For example, use reference_pixel=(0,0) to select the region closest to the bottom left corner of the target pixel file. If ‘center’ (default) then the region closest to the center pixel will be selected. If None then all regions will be selected.

Returns
aperture_maskndarray

2D boolean numpy array containing True for pixels above the threshold.

cutout(self, center=None, size=5)

Cut a rectangle out of the Target Pixel File.

This methods returns a new TargetPixelFile object containing a rectangle of a given size cut out around a given center.

Parameters
center(int, int) tuple or astropy.SkyCoord

Center of the cutout. If an (int, int) tuple is passed, it will be interpreted as the (column, row) coordinates relative to the bottom-left corner of the TPF. If an astropy.SkyCoord is passed then the sky coordinate will be used instead. If None (default) then the center of the TPF will be used.

sizeint or (int, int) tuple

Number of pixels to cut out. If a single integer is passed then a square of that size will be cut. If a tuple is passed then a rectangle with dimensions (column_size, row_size) will be cut.

Returns
tpflightkurve.TargetPixelFile object

New and smaller Target Pixel File object containing only the data cut out.

estimate_centroids(self, aperture_mask='pipeline', method='moments')

Returns the flux center of an object inside aperture_mask.

Telescopes tend to smear out the light from a point-like star over multiple pixels. For this reason, it is common to estimate the position of a star by computing the geometric center of its image. Astronomers refer to this position as the centroid of the object, i.e. the term centroid is often used as a generic synonym to refer to the measured position of an object in a telescope exposure.

This function provides two methods to estimate the position of a star:

  • method='moments' will compute the “center of mass” of the light based on the 2D image moments of the pixels inside aperture_mask.

  • method='quadratic' will fit a two-dimensional, second-order polynomial to the 3x3 patch of pixels centered on the brightest pixel inside the aperture_mask, and return the peak of that polynomial. Following Vakili & Hogg 2016 (ArXiv:1610.05873, Section 3.2).

Parameters
aperture_mask‘pipeline’, ‘threshold’, ‘all’, or array-like

Which pixels contain the object to be measured, i.e. which pixels should be used in the estimation? If None or ‘all’ are passed, all pixels in the pixel file will be used. If ‘pipeline’ is passed, the mask suggested by the official pipeline will be used. If ‘threshold’ is passed, all pixels brighter than 3-sigma above the median flux will be used. Alternatively, users can pass a boolean array describing the aperture mask such that True means that the pixel will be used.

method‘moments’ or ‘quadratic’

Defines which method to use to estimate the centroids. ‘moments’ computes the centroid based on the sample moments of the data. ‘quadratic’ fits a 2D polynomial to the data and returns the coordinate of the peak of that polynomial.

Returns
columns, rowsarray, array

Arrays containing the column and row positions for the centroid for each cadence, or NaN for cadences where the estimation failed.

extract_aperture_photometry(self, aperture_mask='pipeline', centroid_method='moments')

Performs aperture photometry.

Parameters
aperture_maskarray-like, ‘pipeline’, ‘threshold’ or ‘all’

A boolean array describing the aperture such that True means that the pixel will be used. If None or ‘all’ are passed, all pixels will be used. If ‘pipeline’ is passed, the mask suggested by the official pipeline will be returned. If ‘threshold’ is passed, all pixels brighter than 3-sigma above the median flux will be used. The default behaviour is to use the TESS pipeline mask.

centroid_methodstr, ‘moments’ or ‘quadratic’

For the details on this arguments, please refer to the documentation for TargetPixelFile.estimate_centroids.

Returns
lcTessLightCurve object

Contains the summed flux within the aperture for each cadence.

get_bkg_lightcurve(self, aperture_mask=None)
get_coordinates(self, cadence='all')

Returns two 3D arrays of RA and Dec values in decimal degrees.

If cadence number is given, returns 2D arrays for that cadence. If cadence is ‘all’ returns one RA, Dec value for each pixel in every cadence. Uses the WCS solution and the POS_CORR data from TPF header.

Parameters
cadence‘all’ or int

Which cadences to return the RA Dec coordinates for.

Returns
ranumpy array, same shape as tpf.flux[cadence]

Array containing RA values for every pixel, for every cadence.

decnumpy array, same shape as tpf.flux[cadence]

Array containing Dec values for every pixel, for every cadence.

get_keyword(self, keyword, hdu=0, default=None)

Returns a header keyword value.

If the keyword is Undefined or does not exist, then return default instead.

interact(self, notebook_url='localhost:8888', max_cadences=30000, aperture_mask='pipeline', exported_filename=None, transform_func=None, ylim_func=None)

Display an interactive Jupyter Notebook widget to inspect the pixel data.

The widget will show both the lightcurve and pixel data. By default, the lightcurve shown is obtained by calling the to_lightcurve() method, unless the user supplies a custom LightCurve object. This feature requires an optional dependency, bokeh (v0.12.15 or later). This dependency can be installed using e.g. conda install bokeh.

At this time, this feature only works inside an active Jupyter Notebook, and tends to be too slow when more than ~30,000 cadences are contained in the TPF (e.g. short cadence data).

Parameters
notebook_urlstr

Location of the Jupyter notebook page (default: “localhost:8888”) When showing Bokeh applications, the Bokeh server must be explicitly configured to allow connections originating from different URLs. This parameter defaults to the standard notebook host and port. If you are running on a different location, you will need to supply this value for the application to display properly. If no protocol is supplied in the URL, e.g. if it is of the form “localhost:8888”, then “http” will be used.

max_cadencesint

Raise a RuntimeError if the number of cadences shown is larger than this value. This limit helps keep browsers from becoming unresponsive.

aperture_maskarray-like, ‘pipeline’, ‘threshold’ or ‘all’

A boolean array describing the aperture such that True means that the pixel will be used. If None or ‘all’ are passed, all pixels will be used. If ‘pipeline’ is passed, the mask suggested by the official pipeline will be returned. If ‘threshold’ is passed, all pixels brighter than 3-sigma above the median flux will be used.

exported_filename: str

An optional filename to assign to exported fits files containing the custom aperture mask generated by clicking on pixels in interact. The default adds a suffix ‘-custom-aperture-mask.fits’ to the TargetPixelFile basename.

transform_func: function

A function that transforms the lightcurve. The function takes in a LightCurve object as input and returns a LightCurve object as output. The function can be complex, such as detrending the lightcurve. In this way, the interactive selection of aperture mask can be evaluated after inspection of the transformed lightcurve. The transform_func is applied before saving a fits file. Default: None (no transform is applied).

ylim_func: function

A function that returns ylimits (low, high) given a LightCurve object. The default is to return an expanded window around the 10-90th percentile of lightcurve flux values.

Examples

To select an aperture mask for V827 Tau:

>>> import lightkurve as lk
>>> tpf = lk.search_targetpixelfile("V827 Tau", mission="K2").download()
>>> tpf.interact()

To see the full y-axis dynamic range of your lightcurve and normalize the lightcurve after each pixel selection:

>>> ylim_func = lambda lc: (0.0, lc.flux.max())
>>> transform_func = lambda lc: lc.normalize()
>>> tpf.interact(ylim_func=ylim_func, transform_func=transform_func)
interact_sky(self, notebook_url='localhost:8888', magnitude_limit=18)

Display a Jupyter Notebook widget showing Gaia DR2 positions on top of the pixels.

Parameters
notebook_urlstr

Location of the Jupyter notebook page (default: “localhost:8888”) When showing Bokeh applications, the Bokeh server must be explicitly configured to allow connections originating from different URLs. This parameter defaults to the standard notebook host and port. If you are running on a different location, you will need to supply this value for the application to display properly. If no protocol is supplied in the URL, e.g. if it is of the form “localhost:8888”, then “http” will be used.

magnitude_limitfloat

A value to limit the results in based on Gaia Gmag. Default, 18.

plot(self, ax=None, frame=0, cadenceno=None, bkg=False, aperture_mask=None, show_colorbar=True, mask_color='pink', title=None, style='lightkurve', **kwargs)

Plot the pixel data for a single frame (i.e. at a single time).

The time can be specified by frame index number (frame=0 will show the first frame) or absolute cadence number (cadenceno).

Parameters
axAxes

A matplotlib axes object to plot into. If no axes is provided, a new one will be generated.

frameint

Frame number. The default is 0, i.e. the first frame.

cadencenoint, optional

Alternatively, a cadence number can be provided. This argument has priority over frame number.

bkgbool

If True, background will be added to the pixel values.

aperture_maskndarray or str

Highlight pixels selected by aperture_mask.

show_colorbarbool

Whether or not to show the colorbar

mask_colorstr

Color to show the aperture mask

stylestr

Path or URL to a matplotlib style file, or name of one of matplotlib’s built-in stylesheets (e.g. ‘ggplot’). Lightkurve’s custom stylesheet is used by default.

kwargsdict

Keywords arguments passed to lightkurve.utils.plot_image.

Returns
axAxes

The matplotlib axes object.

show_properties(self)

Prints a description of all non-callable attributes.

Prints in order of type (ints, strings, lists, arrays, others).

to_corrector(self, method='pld')

Returns a Corrector instance to remove systematics.

Parameters
methodsstring

Currently, only “pld” is supported. This will return a PLDCorrector class instance.

Returns
correcterlightkurve.Correcter

Instance of a Corrector class, which typically provides correct() and diagnose() methods.

to_fits(self, output_fn=None, overwrite=False)

Writes the TPF to a FITS file on disk.

to_lightcurve(self, method='aperture', **kwargs)

Performs photometry on the pixel data and returns a LightCurve object.

See the docstring of aperture_photometry() for valid arguments if the method is ‘aperture’. Otherwise, see the docstring of prf_photometry() for valid arguments if the method is ‘prf’.

Parameters
method‘aperture’ or ‘prf’.

Photometry method to use.

**kwargsdict

Extra arguments to be passed to the aperture_photometry or the prf_photometry method of this class.

Returns
lcLightCurve object

Object containing the resulting lightcurve.