astrobase.fakelcs.generation module

This generates light curves of variable stars using the astrobase.lcmodels package, adds noise and observation sampling to them based on given parameters (or example light curves). See fakelcrecovery.py for functions that run a full recovery simulation.

NOTE 1: the parameters for these light curves are currently all chosen from uniform distributions, which obviously doesn’t reflect reality (some of the current parameter upper and lower limits are realistic, however). Some of these will be updated with real-life distributions as soon as I find them, especially for periods and amplitudes (along with references).

NOTE 2: To generate custom distributions, one can subclass scipy.stats.rv_continuous and override the _pdf and _cdf methods (or just the _rvs method directly to get the distributed variables if distribution’s location and scale don’t really matter). This is described here:

https://docs.scipy.org/doc/scipy/reference/generated/scipy.stats.rv_continuous.html

and doesn’t seem to be restricted to distributions described by analytic functions only. It’s probably possible to get a histogram of some complex distribution in parameter space and use a kernel-density estimator to get the PDF and CDF, e.g.:

http://scikit-learn.org/stable/modules/density.html#kernel-density.

NOTE 3: any distribution parameter below having to do with magnitudes/flux is by default in MAGNITUDES. So depth, amplitude, etc. distributions and their limits will have to be adjusted appropriately for fluxes. IT IS NOT SUFFICIENT to just set magsarefluxes = True.

FIXME: in the future, we’ll do all the amplitude, etc. distributions in differential fluxes canonically, and then take logs where appropriate if magsarefluxes = False.

FIXME: we should add RA/DEC values that are taken from GAIA if we provide a radec box for the simulation to take place in. in this way, we can parameterize blending and take it into account in the recovery as well.

FIXME: check if object coordinates end up so that two or more objects lie within a chosen blend radius. if this happens, we should check if the blender(s) are variable, and add in some fraction of their phased light curve to the blendee. if the blender(s) are not variable, add in a constant fraction of the brightness to the blendee’s light curve. the blending fraction is multiplied into the light curve of the blender(s) and the resulting flux added to the blendee’s light curve.

  • given the FWHM of the instrument, figure out the overlap
  • we need to calculate the pixel area for blendee and the sum of pixel areas covered by the blenders. this will require input kwargs for pixel size of the detector and FWHM of the star (this might need to be calculated based on the brightness of the star)

FIXME: add input from TRILEGAL produced .dat files for color and mag information. This will let us generate pulsating variables with their actual colors.

astrobase.fakelcs.generation.generate_transit_lightcurve(times, mags=None, errs=None, paramdists={'transitdepth': <scipy.stats._distn_infrastructure.rv_frozen object>, 'transitduration': <scipy.stats._distn_infrastructure.rv_frozen object>, 'transitperiod': <scipy.stats._distn_infrastructure.rv_frozen object>}, magsarefluxes=False)[source]

This generates fake planet transit light curves.

Parameters:
  • times (np.array) – This is an array of time values that will be used as the time base.
  • mags,errs (np.array) – These arrays will have the model added to them. If either is None, np.full_like(times, 0.0) will used as a substitute and the model light curve will be centered around 0.0.
  • paramdists (dict) –

    This is a dict containing parameter distributions to use for the model params, containing the following keys

    {'transitperiod', 'transitdepth', 'transitduration'}
    

    The values of these keys should all be ‘frozen’ scipy.stats distribution objects, e.g.:

    https://docs.scipy.org/doc/scipy/reference/stats.html#continuous-distributions The variability epoch will be automatically chosen from a uniform distribution between times.min() and times.max().

    The ingress duration will be automatically chosen from a uniform distribution ranging from 0.05 to 0.5 of the transitduration.

    The transitdepth will be flipped automatically as appropriate if magsarefluxes=True.

  • magsarefluxes (bool) – If the generated time series is meant to be a flux time-series, set this to True to get the correct sign of variability amplitude.
Returns:

A dict of the form below is returned:

{'vartype': 'planet',
 'params': {'transitperiod': generated value of period,
            'transitepoch': generated value of epoch,
            'transitdepth': generated value of transit depth,
            'transitduration': generated value of transit duration,
            'ingressduration': generated value of transit ingress
                               duration},
 'times': the model times,
 'mags': the model mags,
 'errs': the model errs,
 'varperiod': the generated period of variability == 'transitperiod'
 'varamplitude': the generated amplitude of
                 variability == 'transitdepth'}

Return type:

dict

astrobase.fakelcs.generation.generate_eb_lightcurve(times, mags=None, errs=None, paramdists={'depthratio': <scipy.stats._distn_infrastructure.rv_frozen object>, 'pdepth': <scipy.stats._distn_infrastructure.rv_frozen object>, 'pduration': <scipy.stats._distn_infrastructure.rv_frozen object>, 'period': <scipy.stats._distn_infrastructure.rv_frozen object>, 'secphase': <scipy.stats._distn_infrastructure.rv_frozen object>}, magsarefluxes=False)[source]

This generates fake EB light curves.

Parameters:
  • times (np.array) – This is an array of time values that will be used as the time base.
  • mags,errs (np.array) – These arrays will have the model added to them. If either is None, np.full_like(times, 0.0) will used as a substitute and the model light curve will be centered around 0.0.
  • paramdists (dict) –

    This is a dict containing parameter distributions to use for the model params, containing the following keys

    {'period', 'pdepth', 'pduration', 'depthratio', 'secphase'}
    

    The values of these keys should all be ‘frozen’ scipy.stats distribution objects, e.g.:

    https://docs.scipy.org/doc/scipy/reference/stats.html#continuous-distributions The variability epoch will be automatically chosen from a uniform distribution between times.min() and times.max().

    The pdepth will be flipped automatically as appropriate if magsarefluxes=True.

  • magsarefluxes (bool) – If the generated time series is meant to be a flux time-series, set this to True to get the correct sign of variability amplitude.
Returns:

A dict of the form below is returned:

{'vartype': 'EB',
 'params': {'period': generated value of period,
            'epoch': generated value of epoch,
            'pdepth': generated value of priary eclipse depth,
            'pduration': generated value of prim eclipse duration,
            'depthratio': generated value of prim/sec eclipse
                          depth ratio},
 'times': the model times,
 'mags': the model mags,
 'errs': the model errs,
 'varperiod': the generated period of variability == 'period'
 'varamplitude': the generated amplitude of
                 variability == 'pdepth'}

Return type:

dict

astrobase.fakelcs.generation.generate_flare_lightcurve(times, mags=None, errs=None, paramdists={'amplitude': <scipy.stats._distn_infrastructure.rv_frozen object>, 'decayconst': <scipy.stats._distn_infrastructure.rv_frozen object>, 'nflares': [1, 5], 'risestdev': <scipy.stats._distn_infrastructure.rv_frozen object>}, magsarefluxes=False)[source]

This generates fake flare light curves.

Parameters:
  • times (np.array) – This is an array of time values that will be used as the time base.
  • mags,errs (np.array) – These arrays will have the model added to them. If either is None, np.full_like(times, 0.0) will used as a substitute and the model light curve will be centered around 0.0.
  • paramdists (dict) –

    This is a dict containing parameter distributions to use for the model params, containing the following keys

    {'amplitude', 'nflares', 'risestdev', 'decayconst'}
    

    The values of these keys should all be ‘frozen’ scipy.stats distribution objects, e.g.:

    https://docs.scipy.org/doc/scipy/reference/stats.html#continuous-distributions The flare_peak_time for each flare will be generated automatically between times.min() and times.max() using a uniform distribution.

    The amplitude will be flipped automatically as appropriate if magsarefluxes=True.

  • magsarefluxes (bool) – If the generated time series is meant to be a flux time-series, set this to True to get the correct sign of variability amplitude.
Returns:

A dict of the form below is returned:

{'vartype': 'flare',
 'params': {'amplitude': generated value of flare amplitudes,
            'nflares': generated value of number of flares,
            'risestdev': generated value of stdev of rise time,
            'decayconst': generated value of decay constant,
            'peaktime': generated value of flare peak time},
 'times': the model times,
 'mags': the model mags,
 'errs': the model errs,
 'varamplitude': the generated amplitude of
                 variability == 'amplitude'}

Return type:

dict

astrobase.fakelcs.generation.generate_sinusoidal_lightcurve(times, mags=None, errs=None, paramdists={'amplitude': <scipy.stats._distn_infrastructure.rv_frozen object>, 'fourierorder': [2, 10], 'period': <scipy.stats._distn_infrastructure.rv_frozen object>, 'phioffset': 0.0}, magsarefluxes=False)[source]

This generates fake sinusoidal light curves.

This can be used for a variety of sinusoidal variables, e.g. RRab, RRc, Cepheids, Miras, etc. The functions that generate these model LCs below implement the following table:

## FOURIER PARAMS FOR SINUSOIDAL VARIABLES
#
# type        fourier           period [days]
#             order    dist     limits         dist

# RRab        8 to 10  uniform  0.45--0.80     uniform
# RRc         3 to 6   uniform  0.10--0.40     uniform
# HADS        7 to 9   uniform  0.04--0.10     uniform
# rotator     2 to 5   uniform  0.80--120.0    uniform
# LPV         2 to 5   uniform  250--500.0     uniform

FIXME: for better model LCs, figure out how scipy.signal.butter works and low-pass filter using scipy.signal.filtfilt.

Parameters:
  • times (np.array) – This is an array of time values that will be used as the time base.
  • mags,errs (np.array) – These arrays will have the model added to them. If either is None, np.full_like(times, 0.0) will used as a substitute and the model light curve will be centered around 0.0.
  • paramdists (dict) –

    This is a dict containing parameter distributions to use for the model params, containing the following keys

    {'period', 'fourierorder', 'amplitude', 'phioffset'}
    

    The values of these keys should all be ‘frozen’ scipy.stats distribution objects, e.g.:

    https://docs.scipy.org/doc/scipy/reference/stats.html#continuous-distributions The variability epoch will be automatically chosen from a uniform distribution between times.min() and times.max().

    The amplitude will be flipped automatically as appropriate if magsarefluxes=True.

  • magsarefluxes (bool) – If the generated time series is meant to be a flux time-series, set this to True to get the correct sign of variability amplitude.
Returns:

A dict of the form below is returned:

{'vartype': 'sinusoidal',
 'params': {'period': generated value of period,
            'epoch': generated value of epoch,
            'amplitude': generated value of amplitude,
            'fourierorder': generated value of fourier order,
            'fourieramps': generated values of fourier amplitudes,
            'fourierphases': generated values of fourier phases},
 'times': the model times,
 'mags': the model mags,
 'errs': the model errs,
 'varperiod': the generated period of variability == 'period'
 'varamplitude': the generated amplitude of
                 variability == 'amplitude'}

Return type:

dict

astrobase.fakelcs.generation.generate_rrab_lightcurve(times, mags=None, errs=None, paramdists={'amplitude': <scipy.stats._distn_infrastructure.rv_frozen object>, 'fourierorder': [8, 11], 'period': <scipy.stats._distn_infrastructure.rv_frozen object>, 'phioffset': 3.141592653589793}, magsarefluxes=False)[source]

This generates fake RRab light curves.

Parameters:
  • times (np.array) – This is an array of time values that will be used as the time base.
  • mags,errs (np.array) – These arrays will have the model added to them. If either is None, np.full_like(times, 0.0) will used as a substitute and the model light curve will be centered around 0.0.
  • paramdists (dict) –

    This is a dict containing parameter distributions to use for the model params, containing the following keys

    {'period', 'fourierorder', 'amplitude'}
    

    The values of these keys should all be ‘frozen’ scipy.stats distribution objects, e.g.:

    https://docs.scipy.org/doc/scipy/reference/stats.html#continuous-distributions The variability epoch will be automatically chosen from a uniform distribution between times.min() and times.max().

    The amplitude will be flipped automatically as appropriate if magsarefluxes=True.

  • magsarefluxes (bool) – If the generated time series is meant to be a flux time-series, set this to True to get the correct sign of variability amplitude.
Returns:

A dict of the form below is returned:

{'vartype': 'RRab',
 'params': {'period': generated value of period,
            'epoch': generated value of epoch,
            'amplitude': generated value of amplitude,
            'fourierorder': generated value of fourier order,
            'fourieramps': generated values of fourier amplitudes,
            'fourierphases': generated values of fourier phases},
 'times': the model times,
 'mags': the model mags,
 'errs': the model errs,
 'varperiod': the generated period of variability == 'period'
 'varamplitude': the generated amplitude of
                 variability == 'amplitude'}

Return type:

dict

astrobase.fakelcs.generation.generate_rrc_lightcurve(times, mags=None, errs=None, paramdists={'amplitude': <scipy.stats._distn_infrastructure.rv_frozen object>, 'fourierorder': [2, 3], 'period': <scipy.stats._distn_infrastructure.rv_frozen object>, 'phioffset': 4.71238898038469}, magsarefluxes=False)[source]

This generates fake RRc light curves.

Parameters:
  • times (np.array) – This is an array of time values that will be used as the time base.
  • mags,errs (np.array) – These arrays will have the model added to them. If either is None, np.full_like(times, 0.0) will used as a substitute and the model light curve will be centered around 0.0.
  • paramdists (dict) –

    This is a dict containing parameter distributions to use for the model params, containing the following keys

    {'period', 'fourierorder', 'amplitude'}
    

    The values of these keys should all be ‘frozen’ scipy.stats distribution objects, e.g.:

    https://docs.scipy.org/doc/scipy/reference/stats.html#continuous-distributions The variability epoch will be automatically chosen from a uniform distribution between times.min() and times.max().

    The amplitude will be flipped automatically as appropriate if magsarefluxes=True.

  • magsarefluxes (bool) – If the generated time series is meant to be a flux time-series, set this to True to get the correct sign of variability amplitude.
Returns:

A dict of the form below is returned:

{'vartype': 'RRc',
 'params': {'period': generated value of period,
            'epoch': generated value of epoch,
            'amplitude': generated value of amplitude,
            'fourierorder': generated value of fourier order,
            'fourieramps': generated values of fourier amplitudes,
            'fourierphases': generated values of fourier phases},
 'times': the model times,
 'mags': the model mags,
 'errs': the model errs,
 'varperiod': the generated period of variability == 'period'
 'varamplitude': the generated amplitude of
                 variability == 'amplitude'}

Return type:

dict

astrobase.fakelcs.generation.generate_hads_lightcurve(times, mags=None, errs=None, paramdists={'amplitude': <scipy.stats._distn_infrastructure.rv_frozen object>, 'fourierorder': [5, 10], 'period': <scipy.stats._distn_infrastructure.rv_frozen object>, 'phioffset': 3.141592653589793}, magsarefluxes=False)[source]

This generates fake HADS light curves.

Parameters:
  • times (np.array) – This is an array of time values that will be used as the time base.
  • mags,errs (np.array) – These arrays will have the model added to them. If either is None, np.full_like(times, 0.0) will used as a substitute and the model light curve will be centered around 0.0.
  • paramdists (dict) –

    This is a dict containing parameter distributions to use for the model params, containing the following keys

    {'period', 'fourierorder', 'amplitude'}
    

    The values of these keys should all be ‘frozen’ scipy.stats distribution objects, e.g.:

    https://docs.scipy.org/doc/scipy/reference/stats.html#continuous-distributions The variability epoch will be automatically chosen from a uniform distribution between times.min() and times.max().

    The amplitude will be flipped automatically as appropriate if magsarefluxes=True.

  • magsarefluxes (bool) – If the generated time series is meant to be a flux time-series, set this to True to get the correct sign of variability amplitude.
Returns:

A dict of the form below is returned:

{'vartype': 'HADS',
 'params': {'period': generated value of period,
            'epoch': generated value of epoch,
            'amplitude': generated value of amplitude,
            'fourierorder': generated value of fourier order,
            'fourieramps': generated values of fourier amplitudes,
            'fourierphases': generated values of fourier phases},
 'times': the model times,
 'mags': the model mags,
 'errs': the model errs,
 'varperiod': the generated period of variability == 'period'
 'varamplitude': the generated amplitude of
                 variability == 'amplitude'}

Return type:

dict

astrobase.fakelcs.generation.generate_rotator_lightcurve(times, mags=None, errs=None, paramdists={'amplitude': <scipy.stats._distn_infrastructure.rv_frozen object>, 'fourierorder': [2, 3], 'period': <scipy.stats._distn_infrastructure.rv_frozen object>, 'phioffset': 4.71238898038469}, magsarefluxes=False)[source]

This generates fake rotator light curves.

Parameters:
  • times (np.array) – This is an array of time values that will be used as the time base.
  • mags,errs (np.array) – These arrays will have the model added to them. If either is None, np.full_like(times, 0.0) will used as a substitute and the model light curve will be centered around 0.0.
  • paramdists (dict) –

    This is a dict containing parameter distributions to use for the model params, containing the following keys

    {'period', 'fourierorder', 'amplitude'}
    

    The values of these keys should all be ‘frozen’ scipy.stats distribution objects, e.g.:

    https://docs.scipy.org/doc/scipy/reference/stats.html#continuous-distributions The variability epoch will be automatically chosen from a uniform distribution between times.min() and times.max().

    The amplitude will be flipped automatically as appropriate if magsarefluxes=True.

  • magsarefluxes (bool) – If the generated time series is meant to be a flux time-series, set this to True to get the correct sign of variability amplitude.
Returns:

A dict of the form below is returned:

{'vartype': 'rotator',
 'params': {'period': generated value of period,
            'epoch': generated value of epoch,
            'amplitude': generated value of amplitude,
            'fourierorder': generated value of fourier order,
            'fourieramps': generated values of fourier amplitudes,
            'fourierphases': generated values of fourier phases},
 'times': the model times,
 'mags': the model mags,
 'errs': the model errs,
 'varperiod': the generated period of variability == 'period'
 'varamplitude': the generated amplitude of
                 variability == 'amplitude'}

Return type:

dict

astrobase.fakelcs.generation.generate_lpv_lightcurve(times, mags=None, errs=None, paramdists={'amplitude': <scipy.stats._distn_infrastructure.rv_frozen object>, 'fourierorder': [2, 3], 'period': <scipy.stats._distn_infrastructure.rv_frozen object>, 'phioffset': 4.71238898038469}, magsarefluxes=False)[source]

This generates fake long-period-variable (LPV) light curves.

Parameters:
  • times (np.array) – This is an array of time values that will be used as the time base.
  • mags,errs (np.array) – These arrays will have the model added to them. If either is None, np.full_like(times, 0.0) will used as a substitute and the model light curve will be centered around 0.0.
  • paramdists (dict) –

    This is a dict containing parameter distributions to use for the model params, containing the following keys

    {'period', 'fourierorder', 'amplitude'}
    

    The values of these keys should all be ‘frozen’ scipy.stats distribution objects, e.g.:

    https://docs.scipy.org/doc/scipy/reference/stats.html#continuous-distributions The variability epoch will be automatically chosen from a uniform distribution between times.min() and times.max().

    The amplitude will be flipped automatically as appropriate if magsarefluxes=True.

  • magsarefluxes (bool) – If the generated time series is meant to be a flux time-series, set this to True to get the correct sign of variability amplitude.
Returns:

A dict of the form below is returned:

{'vartype': 'LPV',
 'params': {'period': generated value of period,
            'epoch': generated value of epoch,
            'amplitude': generated value of amplitude,
            'fourierorder': generated value of fourier order,
            'fourieramps': generated values of fourier amplitudes,
            'fourierphases': generated values of fourier phases},
 'times': the model times,
 'mags': the model mags,
 'errs': the model errs,
 'varperiod': the generated period of variability == 'period'
 'varamplitude': the generated amplitude of
                 variability == 'amplitude'}

Return type:

dict

astrobase.fakelcs.generation.generate_cepheid_lightcurve(times, mags=None, errs=None, paramdists={'amplitude': <scipy.stats._distn_infrastructure.rv_frozen object>, 'fourierorder': [8, 11], 'period': <scipy.stats._distn_infrastructure.rv_frozen object>, 'phioffset': 3.141592653589793}, magsarefluxes=False)[source]

This generates fake Cepheid light curves.

Parameters:
  • times (np.array) – This is an array of time values that will be used as the time base.
  • mags,errs (np.array) – These arrays will have the model added to them. If either is None, np.full_like(times, 0.0) will used as a substitute and the model light curve will be centered around 0.0.
  • paramdists (dict) –

    This is a dict containing parameter distributions to use for the model params, containing the following keys

    {'period', 'fourierorder', 'amplitude'}
    

    The values of these keys should all be ‘frozen’ scipy.stats distribution objects, e.g.:

    https://docs.scipy.org/doc/scipy/reference/stats.html#continuous-distributions The variability epoch will be automatically chosen from a uniform distribution between times.min() and times.max().

    The amplitude will be flipped automatically as appropriate if magsarefluxes=True.

  • magsarefluxes (bool) – If the generated time series is meant to be a flux time-series, set this to True to get the correct sign of variability amplitude.
Returns:

A dict of the form below is returned:

{'vartype': 'cepheid',
 'params': {'period': generated value of period,
            'epoch': generated value of epoch,
            'amplitude': generated value of amplitude,
            'fourierorder': generated value of fourier order,
            'fourieramps': generated values of fourier amplitudes,
            'fourierphases': generated values of fourier phases},
 'times': the model times,
 'mags': the model mags,
 'errs': the model errs,
 'varperiod': the generated period of variability == 'period'
 'varamplitude': the generated amplitude of
                 variability == 'amplitude'}

Return type:

dict

astrobase.fakelcs.generation.make_fakelc(lcfile, outdir, magrms=None, randomizemags=True, randomizecoords=False, lcformat='hat-sql', lcformatdir=None, timecols=None, magcols=None, errcols=None)[source]

This preprocesses an input real LC and sets it up to be a fake LC.

Parameters:
  • lcfile (str) – This is an input light curve file that will be used to copy over the time-base. This will be used to generate the time-base for fake light curves to provide a realistic simulation of the observing window function.
  • outdir (str) – The output directory where the the fake light curve will be written.
  • magrms (dict) – This is a dict containing the SDSS r mag-RMS (SDSS rmag-MAD preferably) relation based on all light curves that the input lcfile is from. This will be used to generate the median mag and noise corresponding to the magnitude chosen for this fake LC.
  • randomizemags (bool) – If this is True, then a random mag between the first and last magbin in magrms will be chosen as the median mag for this light curve. This choice will be weighted by the mag bin probability obtained from the magrms kwarg. Otherwise, the median mag will be taken from the input lcfile’s lcdict[‘objectinfo’][‘sdssr’] key or a transformed SDSS r mag generated from the input lcfile’s lcdict[‘objectinfo’][‘jmag’], [‘hmag’], and [‘kmag’] keys. The magrms relation for each magcol will be used to generate Gaussian noise at the correct level for the magbin this light curve’s median mag falls into.
  • randomizecoords (bool) – If this is True, will randomize the RA, DEC of the output fake object and not copy over the RA/DEC from the real input object.
  • lcformat (str) – This is the formatkey associated with your input real light curve format, which you previously passed in to the lcproc.register_lcformat function. This will be used to look up how to find and read the light curve specified in lcfile.
  • lcformatdir (str or None) – If this is provided, gives the path to a directory when you’ve stored your lcformat description JSONs, other than the usual directories lcproc knows to search for them in. Use this along with lcformat to specify an LC format JSON file that’s not currently registered with lcproc.
  • timecols (list of str or None) – The timecol keys to use from the input lcdict in generating the fake light curve. Fake LCs will be generated for each each timecol/magcol/errcol combination in the input light curve.
  • magcols (list of str or None) – The magcol keys to use from the input lcdict in generating the fake light curve. Fake LCs will be generated for each each timecol/magcol/errcol combination in the input light curve.
  • errcols (list of str or None) – The errcol keys to use from the input lcdict in generating the fake light curve. Fake LCs will be generated for each each timecol/magcol/errcol combination in the input light curve.
Returns:

A tuple of the following form is returned:

(fakelc_fpath,
 fakelc_lcdict['columns'],
 fakelc_lcdict['objectinfo'],
 fakelc_lcdict['moments'])

Return type:

tuple

astrobase.fakelcs.generation.collection_worker(task)[source]

This wraps process_fakelc for make_fakelc_collection below.

Parameters:task (tuple) –

This is of the form:

task[0] = lcfile
task[1] = outdir
task[2] = magrms
task[3] = dict with keys: {'lcformat', 'timecols', 'magcols',
                           'errcols', 'randomizeinfo'}
Returns:This returns a tuple of the form:
(fakelc_fpath,
 fakelc_lcdict['columns'],
 fakelc_lcdict['objectinfo'],
 fakelc_lcdict['moments'])
Return type:tuple
astrobase.fakelcs.generation.make_fakelc_collection(lclist, simbasedir, magrmsfrom, magrms_interpolate='quadratic', magrms_fillvalue='extrapolate', maxlcs=25000, maxvars=2000, randomizemags=True, randomizecoords=False, vartypes=('EB', 'RRab', 'RRc', 'cepheid', 'rotator', 'flare', 'HADS', 'planet', 'LPV'), lcformat='hat-sql', lcformatdir=None, timecols=None, magcols=None, errcols=None)[source]

This prepares light curves for the recovery sim.

Collects light curves from lclist using a uniform sampling among them. Copies them to the simbasedir, zeroes out their mags and errs but keeps their time bases, also keeps their RMS and median mags for later use. Calculates the mag-rms relation for the entire collection and writes that to the simbasedir as well.

The purpose of this function is to copy over the time base and mag-rms relation of an existing light curve collection to use it as the basis for a variability recovery simulation.

This returns a pickle written to the simbasedir that contains all the information for the chosen ensemble of fake light curves and writes all generated light curves to the simbasedir/lightcurves directory. Run the add_variability_to_fakelc_collection function after this function to add variability of the specified type to these generated light curves.

Parameters:
  • lclist (list of str) – This is a list of existing project light curves. This can be generated from astrobase.lcproc.catalogs.make_lclist() or similar.
  • simbasedir (str) – This is the directory to where the fake light curves and their information will be copied to.
  • magrmsfrom (str or dict) –

    This is used to generate magnitudes and RMSes for the objects in the output collection of fake light curves. This arg is either a string pointing to an existing pickle file that must contain a dict or a dict variable that MUST have the following key-vals at a minimum:

    {'<magcol1_name>': {
          'binned_sdssr_median': array of median mags for each magbin
          'binned_lcmad_median': array of LC MAD values per magbin
     },
     '<magcol2_name>': {
          'binned_sdssr_median': array of median mags for each magbin
          'binned_lcmad_median': array of LC MAD values per magbin
     },
     .
     .
     ...}
    

    where magcol1_name, etc. are the same as the magcols listed in the magcols kwarg (or the default magcols for the specified lcformat). Examples of the magrmsfrom dict (or pickle) required can be generated by the astrobase.lcproc.varthreshold.variability_threshold() function.

  • magrms_interpolate,magrms_fillvalue (str) –

    These are arguments that will be passed directly to the scipy.interpolate.interp1d function to generate interpolating functions for the mag-RMS relation. See:

    https://docs.scipy.org/doc/scipy/reference/generated/scipy.interpolate.interp1d.html

    for details.

  • maxlcs (int) – This is the total number of light curves to choose from lclist and generate as fake LCs.
  • maxvars (int) – This is the total number of fake light curves that will be marked as variable.
  • vartypes (list of str) – This is a list of variable types to put into the collection. The vartypes for each fake variable star will be chosen uniformly from this list.
  • lcformat (str) – This is the formatkey associated with your input real light curves’ format, which you previously passed in to the lcproc.register_lcformat function. This will be used to look up how to find and read the light curves specified in lclist.
  • lcformatdir (str or None) – If this is provided, gives the path to a directory when you’ve stored your lcformat description JSONs, other than the usual directories lcproc knows to search for them in. Use this along with lcformat to specify an LC format JSON file that’s not currently registered with lcproc.
  • timecols (list of str or None) – The timecol keys to use from the input lcdict in generating the fake light curve. Fake LCs will be generated for each each timecol/magcol/errcol combination in the input light curves.
  • magcols (list of str or None) – The magcol keys to use from the input lcdict in generating the fake light curve. Fake LCs will be generated for each each timecol/magcol/errcol combination in the input light curves.
  • errcols (list of str or None) – The errcol keys to use from the input lcdict in generating the fake light curve. Fake LCs will be generated for each each timecol/magcol/errcol combination in the input light curves.
Returns:

Returns the string file name of a pickle containing all of the information for the fake LC collection that has been generated.

Return type:

str

astrobase.fakelcs.generation.add_fakelc_variability(fakelcfile, vartype, override_paramdists=None, magsarefluxes=False, overwrite=False)[source]

This adds variability of the specified type to the fake LC.

The procedure is (for each magcol):

  • read the fakelcfile, get the stored moments and vartype info
  • add the periodic variability specified in vartype and varparamdists. if vartype == None, then do nothing in this step. If override_vartype is not None, override stored vartype with specified vartype. If override_varparamdists provided, override with specified varparamdists. NOTE: the varparamdists must make sense for the vartype, otherwise, weird stuff will happen.
  • add the median mag level stored in fakelcfile to the time series
  • add Gaussian noise to the light curve as specified in fakelcfile
  • add a varinfo key and dict to the lcdict with varperiod, varepoch, varparams
  • write back to fake LC pickle
  • return the varinfo dict to the caller
Parameters:
  • fakelcfile (str) – The name of the fake LC file to process.
  • vartype (str) – The type of variability to add to this fake LC file.
  • override_paramdists (dict) – A parameter distribution dict as in the generate_XX_lightcurve functions above. If provided, will override the distribution stored in the input fake LC file itself.
  • magsarefluxes (bool) – Sets if the variability amplitude is in fluxes and not magnitudes.
  • overwite (bool) – This overwrites the input fake LC file with a new variable LC even if it’s been processed before.
Returns:

A dict of the following form is returned:

{'objectid':lcdict['objectid'],
 'lcfname':fakelcfile,
 'actual_vartype':vartype,
 'actual_varparams':lcdict['actual_varparams']}

Return type:

dict

astrobase.fakelcs.generation.add_variability_to_fakelc_collection(simbasedir, override_paramdists=None, overwrite_existingvar=False)[source]

This adds variability and noise to all fake LCs in simbasedir.

If an object is marked as variable in the fakelcs-info.pkl file in simbasedir, a variable signal will be added to its light curve based on its selected type, default period and amplitude distribution, the appropriate params, etc. the epochs for each variable object will be chosen uniformly from its time-range (and may not necessarily fall on a actual observed time). Nonvariable objects will only have noise added as determined by their params, but no variable signal will be added.

Parameters:
  • simbasedir (str) – The directory containing the fake LCs to process.
  • override_paramdists (dict) –

    This can be used to override the stored variable parameters in each fake LC. It should be a dict of the following form:

    {'<vartype1>': {'<param1>: a scipy.stats distribution function or
                               the np.random.randint function,
                    .
                    .
                    .
                    '<paramN>: a scipy.stats distribution function
                               or the np.random.randint function}
    

    for any vartype in VARTYPE_LCGEN_MAP. These are used to override the default parameter distributions for each variable type.

  • overwrite_existingvar (bool) – If this is True, then will overwrite any existing variability in the input fake LCs in simbasedir.
Returns:

This returns a dict containing the fake LC filenames as keys and variability info for each as values.

Return type:

dict