astrobase.astrokep module

Contains various useful tools for analyzing Kepler light curves.

astrobase.astrokep.keplerflux_to_keplermag(keplerflux, f12=174000.0)[source]

This converts the Kepler flux in electrons/sec to Kepler magnitude.

The kepler mag/flux relation is:

fkep = (10.0**(-0.4*(kepmag - 12.0)))*f12
f12 = 1.74e5 # electrons/sec
Parameters:
  • keplerflux (float or array-like) – The flux value(s) to convert to magnitudes.
  • f12 (float) – The flux value in the Kepler band corresponding to Kepler mag = 12.0.
Returns:

Magnitudes in the Kepler band corresponding to the input keplerflux flux value(s).

Return type:

np.array

astrobase.astrokep.keplermag_to_keplerflux(keplermag, f12=174000.0)[source]

This converts the Kepler mag back to Kepler flux.

Parameters:
  • keplermag (float or array-like) – The Kepler magnitude value(s) to convert to fluxes.
  • f12 (float) – The flux value in the Kepler band corresponding to Kepler mag = 12.0.
Returns:

Fluxes in the Kepler band corresponding to the input keplermag magnitude value(s).

Return type:

np.array

astrobase.astrokep.keplermag_to_sdssr(keplermag, kic_sdssg, kic_sdssr)[source]

Converts magnitude measurements in Kepler band to SDSS r band.

Parameters:
  • keplermag (float or array-like) – The Kepler magnitude value(s) to convert to fluxes.
  • kic_sdssg,kic_sdssr (float or array-like) – The SDSS g and r magnitudes of the object(s) from the Kepler Input Catalog. The .llc.fits MAST light curve file for a Kepler object contains these values in the FITS extension 0 header.
Returns:

SDSS r band magnitude(s) converted from the Kepler band magnitude.

Return type:

float or array-like

astrobase.astrokep.flux_ppm_to_magnitudes(ppm)[source]

This converts Kepler’s flux parts-per-million to magnitudes.

Mostly useful for turning PPMs reported by Kepler or TESS into millimag values to compare with ground-based surveys.

Parameters:ppm (float or array-like) – Kepler flux measurement errors or RMS values in parts-per-million.
Returns:Measurement errors or RMS values expressed in magnitudes.
Return type:float or array-like
astrobase.astrokep.read_kepler_fitslc(lcfits, headerkeys=['TIMESYS', 'BJDREFI', 'BJDREFF', 'OBJECT', 'KEPLERID', 'RA_OBJ', 'DEC_OBJ', 'EQUINOX', 'EXPOSURE', 'CDPP3_0', 'CDPP6_0', 'CDPP12_0', 'PDCVAR', 'PDCMETHD', 'CROWDSAP', 'FLFRCSAP'], datakeys=['TIME', 'TIMECORR', 'CADENCENO', 'SAP_QUALITY', 'PSF_CENTR1', 'PSF_CENTR1_ERR', 'PSF_CENTR2', 'PSF_CENTR2_ERR', 'MOM_CENTR1', 'MOM_CENTR1_ERR', 'MOM_CENTR2', 'MOM_CENTR2_ERR'], sapkeys=['SAP_FLUX', 'SAP_FLUX_ERR', 'SAP_BKG', 'SAP_BKG_ERR'], pdckeys=['PDCSAP_FLUX', 'PDCSAP_FLUX_ERR'], topkeys=['CHANNEL', 'SKYGROUP', 'MODULE', 'OUTPUT', 'QUARTER', 'SEASON', 'CAMPAIGN', 'DATA_REL', 'OBSMODE', 'PMRA', 'PMDEC', 'PMTOTAL', 'PARALLAX', 'GLON', 'GLAT', 'GMAG', 'RMAG', 'IMAG', 'ZMAG', 'D51MAG', 'JMAG', 'HMAG', 'KMAG', 'KEPMAG', 'GRCOLOR', 'JKCOLOR', 'GKCOLOR', 'TEFF', 'LOGG', 'FEH', 'EBMINUSV', 'AV', 'RADIUS', 'TMINDEX'], apkeys=['NPIXSAP', 'NPIXMISS', 'CDELT1', 'CDELT2'], appendto=None, normalize=False)[source]

This extracts the light curve from a single Kepler or K2 LC FITS file.

This works on the light curves available at MAST:

  • kplr{kepid}-{somedatething}_llc.fits files from the Kepler mission
  • ktwo{epicid}-c{campaign}_llc.fits files from the K2 mission
Parameters:
  • lcfits (str) – The filename of a MAST Kepler/K2 light curve FITS file.
  • headerkeys (list) – A list of FITS header keys that will be extracted from the FITS light curve file. These describe the observations. The default value for this is given in LCHEADERKEYS above.
  • datakeys (list) – A list of FITS column names that correspond to the auxiliary measurements in the light curve. The default is LCDATAKEYS above.
  • sapkeys (list) – A list of FITS column names that correspond to the SAP flux measurements in the light curve. The default is LCSAPKEYS above.
  • pdckeys (list) – A list of FITS column names that correspond to the PDC flux measurements in the light curve. The default is LCPDCKEYS above.
  • topkeys (list) – A list of FITS header keys that describe the object in the light curve. The default is LCTOPKEYS above.
  • apkeys (list) – A list of FITS header keys that describe the flux measurement apertures used by the Kepler/K2 pipeline. The default is LCAPERTUREKEYS above.
  • appendto (lcdict or None) – If appendto is an lcdict, will append measurements of this lcdict to that lcdict. This is used for consolidating light curves for the same object across different files (quarters). The appending does not care about the time order. To consolidate light curves in time order, use consolidate_kepler_fitslc below.
  • normalize (bool) – If True, then each component light curve’s SAP_FLUX and PDCSAP_FLUX measurements will be normalized to 1.0 by dividing out the median flux for the component light curve.
Returns:

Returns an lcdict (this is useable by most astrobase functions for LC processing).

Return type:

lcdict

astrobase.astrokep.consolidate_kepler_fitslc(keplerid, lcfitsdir, normalize=True, headerkeys=['TIMESYS', 'BJDREFI', 'BJDREFF', 'OBJECT', 'KEPLERID', 'RA_OBJ', 'DEC_OBJ', 'EQUINOX', 'EXPOSURE', 'CDPP3_0', 'CDPP6_0', 'CDPP12_0', 'PDCVAR', 'PDCMETHD', 'CROWDSAP', 'FLFRCSAP'], datakeys=['TIME', 'TIMECORR', 'CADENCENO', 'SAP_QUALITY', 'PSF_CENTR1', 'PSF_CENTR1_ERR', 'PSF_CENTR2', 'PSF_CENTR2_ERR', 'MOM_CENTR1', 'MOM_CENTR1_ERR', 'MOM_CENTR2', 'MOM_CENTR2_ERR'], sapkeys=['SAP_FLUX', 'SAP_FLUX_ERR', 'SAP_BKG', 'SAP_BKG_ERR'], pdckeys=['PDCSAP_FLUX', 'PDCSAP_FLUX_ERR'], topkeys=['CHANNEL', 'SKYGROUP', 'MODULE', 'OUTPUT', 'QUARTER', 'SEASON', 'CAMPAIGN', 'DATA_REL', 'OBSMODE', 'PMRA', 'PMDEC', 'PMTOTAL', 'PARALLAX', 'GLON', 'GLAT', 'GMAG', 'RMAG', 'IMAG', 'ZMAG', 'D51MAG', 'JMAG', 'HMAG', 'KMAG', 'KEPMAG', 'GRCOLOR', 'JKCOLOR', 'GKCOLOR', 'TEFF', 'LOGG', 'FEH', 'EBMINUSV', 'AV', 'RADIUS', 'TMINDEX'], apkeys=['NPIXSAP', 'NPIXMISS', 'CDELT1', 'CDELT2'])[source]

This gets all Kepler/K2 light curves for the given keplerid in lcfitsdir.

Searches recursively in lcfitsdir for all of the files belonging to the specified keplerid. Sorts the light curves by time. Returns an lcdict. This is meant to be used to consolidate light curves for a single object across Kepler quarters.

NOTE: keplerid is an integer (without the leading zeros). This is usually the KIC ID.

NOTE: if light curve time arrays contain nans, these and their associated measurements will be sorted to the end of the final combined arrays.

Parameters:
  • keplerid (int) – The Kepler ID of the object to consolidate LCs for, as an integer without any leading zeros. This is usually the KIC or EPIC ID.
  • lcfitsdir (str) – The directory to look in for LCs of the specified object.
  • normalize (bool) – If True, then each component light curve’s SAP_FLUX and PDCSAP_FLUX measurements will be normalized to 1.0 by dividing out the median flux for the component light curve.
  • headerkeys (list) – A list of FITS header keys that will be extracted from the FITS light curve file. These describe the observations. The default value for this is given in LCHEADERKEYS above.
  • datakeys (list) – A list of FITS column names that correspond to the auxiliary measurements in the light curve. The default is LCDATAKEYS above.
  • sapkeys (list) – A list of FITS column names that correspond to the SAP flux measurements in the light curve. The default is LCSAPKEYS above.
  • pdckeys (list) – A list of FITS column names that correspond to the PDC flux measurements in the light curve. The default is LCPDCKEYS above.
  • topkeys (list) – A list of FITS header keys that describe the object in the light curve. The default is LCTOPKEYS above.
  • apkeys (list) – A list of FITS header keys that describe the flux measurement apertures used by the Kepler/K2 pipeline. The default is LCAPERTUREKEYS above.
Returns:

Returns an lcdict (this is useable by most astrobase functions for LC processing).

Return type:

lcdict

astrobase.astrokep.read_k2sff_lightcurve(lcfits)[source]

This reads a K2 SFF (Vandenberg+ 2014) light curve into an lcdict.

Use this with the light curves from the K2 SFF project at MAST.

Parameters:lcfits (str) – The filename of the FITS light curve file downloaded from MAST.
Returns:Returns an lcdict (this is useable by most astrobase functions for LC processing).
Return type:lcdict
astrobase.astrokep.kepler_lcdict_to_pkl(lcdict, outfile=None)[source]

This writes the lcdict to a Python pickle.

Parameters:
  • lcdict (lcdict) – This is the input lcdict to write to a pickle.
  • outfile (str or None) – If this is None, the object’s Kepler ID/EPIC ID will determined from the lcdict and used to form the filename of the output pickle file. If this is a str, the provided filename will be used.
Returns:

The absolute path to the written pickle file.

Return type:

str

astrobase.astrokep.read_kepler_pklc(picklefile)[source]

This turns the pickled lightcurve file back into an lcdict.

Parameters:picklefile (str) – The path to a previously written Kepler LC picklefile generated by kepler_lcdict_to_pkl above.
Returns:Returns an lcdict (this is useable by most astrobase functions for LC processing).
Return type:lcdict
astrobase.astrokep.stitch_kepler_lcdict(lcdict)[source]

This stitches Kepler light curves together across quarters.

FIXME: implement this.

Parameters:lcdict (lcdict) – An lcdict produced by consolidate_kepler_fitslc. The flux measurements between quarters will be stitched together.
Returns:Returns an lcdict (this is useable by most astrobase functions for LC processing). The flux measurements will have been shifted to form a seamless light curve across quarters suitable for long-term variability investigation.
Return type:lcdict
astrobase.astrokep.filter_kepler_lcdict(lcdict, filterflags=True, nanfilter='sap, pdc', timestoignore=None)[source]

This filters the Kepler lcdict, removing nans and bad observations.

By default, this function removes points in the Kepler LC that have ANY quality flags set.

Parameters:
  • lcdict (lcdict) – An lcdict produced by consolidate_kepler_fitslc or read_kepler_fitslc.
  • filterflags (bool) – If True, will remove any measurements that have non-zero quality flags present. This usually indicates an issue with the instrument or spacecraft.
  • nanfilter ({'sap','pdc','sap,pdc'}) – Indicates the flux measurement type(s) to apply the filtering to.
  • timestoignore (list of tuples or None) –

    This is of the form:

    [(time1_start, time1_end), (time2_start, time2_end), ...]
    

    and indicates the start and end times to mask out of the final lcdict. Use this to remove anything that wasn’t caught by the quality flags.

Returns:

Returns an lcdict (this is useable by most astrobase functions for LC processing). The lcdict is filtered IN PLACE!

Return type:

lcdict

astrobase.astrokep.epd_kepler_lightcurve(lcdict, xccol='mom_centr1', yccol='mom_centr2', timestoignore=None, filterflags=True, writetodict=True, epdsmooth=5)[source]

This runs EPD on the Kepler light curve.

Following Huang et al. 2015, we fit the following EPD function to a smoothed light curve, and then subtract it to obtain EPD corrected magnitudes:

f = c0 +
    c1*sin(2*pi*x) + c2*cos(2*pi*x) + c3*sin(2*pi*y) + c4*cos(2*pi*y) +
    c5*sin(4*pi*x) + c6*cos(4*pi*x) + c7*sin(4*pi*y) + c8*cos(4*pi*y) +
    c9*bgv + c10*bge

By default, this function removes points in the Kepler LC that have ANY quality flags set.

Parameters:
  • lcdict (lcdict) – An lcdict produced by consolidate_kepler_fitslc or read_kepler_fitslc.
  • xcol,ycol (str) – Indicates the x and y coordinate column names to use from the Kepler LC in the EPD fit.
  • timestoignore (list of tuples) –

    This is of the form:

    [(time1_start, time1_end), (time2_start, time2_end), ...]
    

    and indicates the start and end times to mask out of the final lcdict. Use this to remove anything that wasn’t caught by the quality flags.

  • filterflags (bool) – If True, will remove any measurements that have non-zero quality flags present. This usually indicates an issue with the instrument or spacecraft.
  • writetodict (bool) –

    If writetodict is True, adds the following columns to the lcdict:

    epd_time = time array
    epd_sapflux = uncorrected flux before EPD
    epd_epdsapflux = corrected flux after EPD
    epd_epdsapcorr = EPD flux corrections
    epd_bkg = background array
    epd_bkg_err = background errors array
    epd_xcc = xcoord array
    epd_ycc = ycoord array
    epd_quality = quality flag array
    

    and updates the ‘columns’ list in the lcdict as well.

  • epdsmooth (int) – Sets the number of light curve points to smooth over when generating the EPD fit function.
Returns:

Returns a tuple of the form: (times, epdfluxes, fitcoeffs, epdfit)

Return type:

tuple

astrobase.astrokep.rfepd_kepler_lightcurve(lcdict, xccol='mom_centr1', yccol='mom_centr2', timestoignore=None, filterflags=True, writetodict=True, epdsmooth=23, decorr='xcc, ycc', nrftrees=200)[source]

This uses a RandomForestRegressor to fit and decorrelate Kepler light curves.

Fits the X and Y positions, the background, and background error.

By default, this function removes points in the Kepler LC that have ANY quality flags set.

Parameters:
  • lcdict (lcdict) – An lcdict produced by consolidate_kepler_fitslc or read_kepler_fitslc.
  • xcol,ycol (str) – Indicates the x and y coordinate column names to use from the Kepler LC in the EPD fit.
  • timestoignore (list of tuples) –

    This is of the form:

    [(time1_start, time1_end), (time2_start, time2_end), ...]
    

    and indicates the start and end times to mask out of the final lcdict. Use this to remove anything that wasn’t caught by the quality flags.

  • filterflags (bool) – If True, will remove any measurements that have non-zero quality flags present. This usually indicates an issue with the instrument or spacecraft.
  • writetodict (bool) –

    If writetodict is True, adds the following columns to the lcdict:

    rfepd_time = time array
    rfepd_sapflux = uncorrected flux before EPD
    rfepd_epdsapflux = corrected flux after EPD
    rfepd_epdsapcorr = EPD flux corrections
    rfepd_bkg = background array
    rfepd_bkg_err = background errors array
    rfepd_xcc = xcoord array
    rfepd_ycc = ycoord array
    rfepd_quality = quality flag array
    

    and updates the ‘columns’ list in the lcdict as well.

  • epdsmooth (int) – Sets the number of light curve points to smooth over when generating the EPD fit function.
  • decorr ({'xcc,ycc','bgv,bge','xcc,ycc,bgv,bge'}) – Indicates whether to use the x,y coords alone; background value and error alone; or x,y coords and background value, error in combination as the features to training the RandomForestRegressor on and perform the fit.
  • nrftrees (int) – The number of trees to use in the RandomForestRegressor.
Returns:

Returns a tuple of the form: (times, corrected_fluxes, flux_corrections)

Return type:

tuple

astrobase.astrokep.detrend_centroid(lcd, detrend='legendre', sigclip=None, mingap=0.5)[source]

Detrends the x and y coordinate centroids for a Kepler light curve.

Given an lcdict for a single quarter of Kepler data, returned by read_kepler_fitslc, this function returns this same dictionary, appending detrended centroid_x and centroid_y values.

Here “detrended” means “finite, SAP quality flag set to 0, sigma clipped, timegroups selected based on mingap day gaps, then fit vs time by a legendre polynomial of lowish degree”.

Parameters:
  • lcd (lcdict) – An lcdict generated by the read_kepler_fitslc function.
  • detrend ({'legendre'}) – Method by which to detrend the LC. ‘legendre’ is the only thing implemented at the moment.
  • sigclip (None or float or int or sequence of floats/ints) – Determines the type and amount of sigma-clipping done on the light curve to remove outliers. If None, no sigma-clipping is performed. If a two element sequence of floats/ints, the first element corresponds to the fainter sigma-clip limit, and the second element corresponds to the brighter sigma-clip limit.
  • mingap (float) – Number of days by which to define “timegroups” (for individual fitting each of timegroup, and to eliminate “burn-in” of Kepler spacecraft. For long cadence data, 0.5 days is typical.
Returns:

This is of the form (lcd, errflag), where:

lcd : an lcdict with the new key lcd[‘centroids’], containing the detrended times, (centroid_x, centroid_y) values, and their errors.

errflag : boolean error flag, could be raised at various points.

Return type:

tuple

astrobase.astrokep.get_centroid_offsets(lcd, t_ing_egr, oot_buffer_time=0.1, sample_factor=3)[source]

After running detrend_centroid, this gets positions of centroids during transits, and outside of transits.

These positions can then be used in a false positive analysis.

This routine requires knowing the ingress and egress times for every transit of interest within the quarter this routine is being called for. There is currently no astrobase routine that automates this for periodic transits (it must be done in a calling routine).

To get out of transit centroids, this routine takes points outside of the “buffer” set by oot_buffer_time, sampling 3x as many points on either side of the transit as are in the transit (or however many are specified by sample_factor).

Parameters:
  • lcd (lcdict) – An lcdict generated by the read_kepler_fitslc function. We assume that the detrend_centroid function has been run on this lcdict.
  • t_ing_egr (list of tuples) –

    This is of the form:

    [(ingress time of i^th transit, egress time of i^th transit)]
    

    for i the transit number index in this quarter (starts at zero at the beginning of every quarter). Assumes units of BJD.

  • oot_buffer_time (float) – Number of days away from ingress and egress times to begin sampling “out of transit” centroid points. The number of out of transit points to take per transit is 3x the number of points in transit.
  • sample_factor (float) – The size of out of transit window from which to sample.
Returns:

This is a dictionary keyed by transit number (i.e., the same index as t_ing_egr), where each key contains the following value:

{'ctd_x_in_tra':ctd_x_in_tra,
 'ctd_y_in_tra':ctd_y_in_tra,
 'ctd_x_oot':ctd_x_oot,
 'ctd_y_oot':ctd_y_oot,
 'npts_in_tra':len(ctd_x_in_tra),
 'npts_oot':len(ctd_x_oot),
 'in_tra_times':in_tra_times,
 'oot_times':oot_times}

Return type:

dict