Scripting

FRBGui ships with the following modules that can be used for scripting outside of the graphical interface:

  • driftrate.py : Utilities for manipulating FRB waterfalls and performing measurements on them.

  • driftlaw.py : Utilities for processing FRBGui measurements for the purpose of studying spectro-temporal relationships.

These can be accessed with

1import driftrate
2import driftlaw

This page is currently under construction and may be incomplete.

frbgui module

The frbgui function itself has multiple options that can be used to configure the startup of FRBGui

As stated in the Overview, FRBGui can be called from a script like so:

1from frbgui import frbgui
frbgui.frbgui(filefilter='*.npz', datadir='', maskfile=None, dmrange=None, numtrials=10, dmstep=0.1, winwidth=1700, winheight=850, winpos=(0, 0))

Start FRBGui window. Every option has a default so all parameters are optional.

Parameters:
  • filefilter (str) – A glob pattern to filter for by default. e.g. ‘*.npz’

  • datadir (str) – path to the directory to start frbgui in

  • maskfile (str) – if masks were exported previously, they can be loaded by default by adding a path to the mask file here

  • dmrange (tuple) – a list or tuple of the default DM range to use (e.g. [555, 575])

  • numtrials (int) – number of measurements to perform over dmrange

  • dmstep (float) – size of DM step to use. Will overwrite numtrials

  • windiwdth (int) – width in pixels of FRBGui window

  • winheight (int) – height in pixels of FRBGui window

  • winpos (tuple) – (x,y) coordinates on the screen of where to start FRBGui

driftrate module

This module contains Utilities for manipulating FRB waterfalls (such as dedispersion and subsampling) and performing measurements on them.

driftrate.autocorr2d(data, maskcorrpeak=True)

Returns a 2D autocorrelation computed via an intermediate FFT

Parameters:
  • data (np.ndarray) – 2D array of size MxN

  • maskcorrpeak (bool, optional) – Set to False to keep the zero lag peak

Returns:

2D array of size 2M-1 x 2N-1. The autocorrelation of data.

Return type:

np.ndarray

driftrate.cropwfall(wfall, twidth=150, pkidx=None)

Crops a waterfall and attempts to center the burst in the array.

Parameters:
  • wfall (np.ndarray) – the array to be cropped of size M x N

  • twidth (int) – the number of channels on either side of the burst. Total width is therefore 2*twidth

  • pkidx (int, optional) – the index of the waterfall to center. If None will use np.argmax of the timeseries as the index to center.

Returns:

cropped waterfall of size M x (2*twidth)

Return type:

np.ndarray

driftrate.dedisperse(intensity, DM, nu_low, df_mhz, dt_ms, cshift=0, a_dm=4149377.59336)

Incoherently dedisperse a dynamic spectrum to the specified DM.

Computes the time delay at each frequency and correspondingly rolls the data

\[\Delta t = - a\Big(\frac{1}{\nu^2_i} - \frac{1}{\nu^2_\text{high}}\Big)\text{DM}\]

By default uses the pulsar community dispersion constant of

\(a = 4.14937759336e6 \quad\text{MHz}^2 \text{cm}^3 \text{pc}^-1 \text{ms}\).

Parameters:
  • intensity (np.ndarray) – the dynamic spectrum (or waterfall) to dedisperse. intensity[0] should correspond to the lowest frequency.

  • DM (float) – the dispersion measure in pc/cm^3

  • nu_low (float) – the frequency at the bottom of the band

  • df_mhz (float) – the frequency resolution of the channels in MHZ

  • dt_ms (float) – the time resolution of channels in ms

  • cshift (int, optional) – additional horizontal shift to add after dedispersion in number of channels

  • a_dm (float) – the dispersion constant to use when dedispersing. See units above.

Returns:

the dedispersed 2d intensity array

Return type:

np.ndarray

driftrate.exportresults(results)

Creates a dataframe of measurement results.

See driftrate.columns.

Parameters:

results (list) – the results list

driftrate.findCenter(burstwindow)

Find the center index with uncertainties of the peak of the waterfall. This is done by time averaging the waterfall to obtain a spectrum and then finding the average channel weighted by the square of the intensity.

This is calculated as

\[\bar{\nu_i} = \frac{\sum_{\nu_i} \nu_i I^2(\nu_i)}{\sum_{\nu_i} I^2(\nu_i)}\]

The full uncertainty (after taking into account units and standard error propagation) is given by

\[\delta \bar{\nu} = \sqrt{\nu_{\text{res}}^2/12 + 4\sigma_I^2\Big(\nu_{\text{res}}\frac{\sum_{\nu} (\nu - \bar{\nu})I(\nu) }{\sum_{\nu} I^2(\nu)}\Big)^2}\]
Parameters:

burstwindow (np.ndarray) – 2d array of the waterfall

Returns:

(meanfreqi, errorsum) - tuple of the mean peak frequency index and summation term of the uncertainty.

You may then calculate the mean frequency in MHz with

center_f = meanfreqi*fres_MHz + lowest_freq

The uncertainty on center_f is then

center_f_err = np.sqrt(fres_MHz**2/12 + 4*wfallsigma**2*(errorsum*fres_MHz)**2)

wfallsigma (\(\sigma_I\)) is the standard deviation of the intensity noise. This can be sampled from the noise in the waterfall that doesn’t include the burst. For example:

wfallsigma = np.std( burstwindow[:, 0:burstwindow.shape[1]//20] )

Return type:

tuple

driftrate.fitdatagaussiannlsq(data, extents, p0=[], sigma=0, bounds=(-inf, inf))

Fit 2d gaussian to data with the non-linear least squares algorithm implemented by scipy.optimize.curve_fit.

Uses the physical units as the data’s coordinates.

Parameters:
  • data (np.ndarray) – 2D array of the data to fit on

  • extents (tuple) – Tuple of the initial time, final time, initial frequency, and final frequency (ti, tf, nu_i, nu_f). See getExtents().

  • p0 (list, optional) – initial guess to used, input as a list of floats corresponding to [amplitude, x0, y0, sigmax, sigmay, theta]. See twoD_Gaussian().

  • sigma (float, optional) – uncertainty to use during fitting. See Scipy Documentation for more information

  • bounds (tuple, optional) –

    lower and upper bounds on parameters. See Scipy Documentation for more information

Returns:

tuple of (popt, pcov) where popt is the list of gaussian parameters found and pcov is the covariance matrix.

Return type:

tuple

driftrate.fitgaussiannlsq(data, p0=[], sigma=0, bounds=(-inf, inf))

Deprecated function. See fitdatagaussianlsq()

driftrate.getDataCoords(extents, shape)
driftrate.getExtents(wfall, df: float = 1.0, dt: float = 1.0, lowest_freq: float = 1.0, lowest_time: float = 0.0)

Given a waterfall’s time and frequency resolutions, return the extents of the axes as well as the axes of its autocorrelation.

Convenience function for plt.imshow’s extent keyword

Parameters:
  • wfall (np.ndarray) – 2D array of the waterfall

  • df (float) – frequency resolution

  • dt (float) – time resolution

  • lowest_freq (float) – lowest frequency in the band

Returns:

(extents, corrextents) where extents is a list of the waterfall extents and corrextents is a list of the corresponding autocorrelation extents

Return type:

tuple

driftrate.getSubbursts(wfall_cr, df, dt, lowest_freq, regions)

Split a waterfall into regions.

Specify regions with a dictionary. For example

regions = {"a": [0, 3], "b": [3, 5.5]}

where the arrays are the timestamps in milliseconds.

Regions should be named as “a”, “b”, “c”, etc. See driftrate.subburst_suffixes.

Meant for use by frbgui.

Parameters:
  • wfall_cr (np.ndarray) – waterfall to split

  • df (float) – frequency resolution

  • dt (float) – time resolution

  • lowest_freq (float) – lowest frequency

  • regions (dict) – Dictionary of {“a”: [start_ms, end_ms], etc..}

Returns:

A dictionary of cropped waterfalls. Follows {“a”: np.ndarray}

Return type:

dict

driftrate.makeDataFitmap(popt, corr, extents)

Make a fitmap using physical coordinates

Parameters:
  • popt (list) – List of 2D gaussian model parameters. [amplitude, xo, yo, sigmax, sigmay, theta)]

  • corr (np.ndarray) – the 2d autocorrelation. Can simply pass a 2d array of the same size

  • extents (list) – Return value of getExtents()

Returns:

a 2d array of the 2d gaussian specified by popt

Return type:

np.ndarray

driftrate.makeFitmap(popt, corr)

Deprecated function. See makeDataFitmap()

driftrate.moments(data)

Returns (height, x, y, width_x, width_y) initial gaussian parameters of a 2D distribution by calculating its moments

Parameters:

data (np.ndarray) – the data to compute moments of

Returns:

(height, x, y, sigma_x, sigma_y)

Return type:

list

driftrate.plotResults(resultsfile, datafiles=[], masks=None, figsize=(14, 16), nrows=6, ncols=4, clip=1, snrLines=False, show=False)

Given a results CSV produced by FRBGui will plot fit results by burst at each DM in a stampcard and save a PDF with the same name.

Parameters:
  • resultsfile (str) – Filename of CSV file

  • datafiles (list[str]) – list of burst filenames in FRBGui’s .npz Burst Format

  • masks (str) – filename to FRBGui maskfile

  • figsize (tuple(float)) – (width, height) of figure

  • nrows (int) – number of rows in stampcard

  • ncols (int) – number of cols in stampcard

  • clip (int) – factor to clip the color scale of the autocorrelation figure for improving display SNR

  • snrLines (bool) – If true plots 1 sigma lines around the autocorrelation

  • show (bool) – Show the figure in a window if true. Displays a pdf otherwise.

Returns:

True if completed. Saves a PDF file when completed.

Return type:

bool

driftrate.plotStampcard(loadfunc, fileglob='*.npy', figsize=(14, 16), nrows=6, ncols=4, twidth=150)

Plot bursts and their autocorrelations in a stampcard

Optionally do custom processing to find the slope and plot the resulting fit as well.

Parameters:
  • loadfunc (function) – a function that accepts “filename” to a waterfall as an argument and loads the waterfall and returns (subfall, pkidx, wfall) where wfall is the loaded waterfall, subfall is the subbanded waterfall that you want to see, and pkidx is the index in the timeseries of the data with peak intensity. Implement this yourself to load your data files

  • fileglob (str) – a glob that matches on your data files (e.g. “*.npy”)

  • figsize (tuple[float]) – Tuple of (width, height) for the resulting figure

  • nrows (int) – number of rows in the stampcard

  • ncols (int) – number of columns in the stampcard

  • twidth (int) – number of channels on either side of the burst to plot

Returns:

shows a figure. Useful in jupyterlab for example.

Return type:

None

driftrate.processBurst(burstwindow, fres_MHz, tres_ms, lowest_freq, burstkey=1, p0=[], popt_custom=[], bounds=(-inf, inf), nclip=None, clip=None, plot=False, corrsigma=None, wfallsigma=None, maskcorrpeak=True, verbose=True, lowest_time=0)

Given a waterfall of a burst, will use the 2d autocorrelation+gaussian fitting method to perform spectro-temporal measurements of the burst

Can optionally plot the measurement.

Parameters:
  • burstwindow (np.ndarray) – 2d array of the burst waterfall

  • fres_MHz (float) – frequency resolution in MHz

  • tres_ms (float) – time resolution in ms

  • lowest_freq (float) – lowest frequency in MHz

  • burstkey (int, optional) – Burst number. Used in plot title

  • p0 (list, optional) – Initial 2d gaussian guess to use. [amplitude, x0, y0, sigmax, sigmay, theta]

  • popt_custom (list, optional) – Override the fit and plot your own 2D gaussian by placing your popt list here

  • bounds (tuple, optional) – parameter bounds. See fitdatagaussianlsq()

  • nclip (float, optional) – minimum clip value of autocorrelation. Applied before fitting.

  • clip (float, optional) – maximum clip value of autocorrelation. Applied before fitting.

  • plot (bool, optional) – If true will display a diagnostic plot of the fit

  • corrsigma (float, optional) – Standard deviation of noise of autocorrelation. Used when fitting.

  • wfallsigma (float, optional) – Standard deviation of noise of waterfall. Used when fitting.

  • verbose (bool, optional) – Set to False to limit console output

  • maskcorrpeak (bool, optional) – Set to False to keep the zero lag correlation peak

  • lowest_time (float, optional) – starting time (ms) of waterfall. Default 0.

Returns:

(slope, slope_error, popt, perr, theta, red_chisq, center_f, center_f_err, fitmap)

Return type:

tuple

driftrate.processDMRange(burstname, wfall, burstdm, dmrange, fres_MHz, tres_ms, lowest_freq, p0=[], corrsigma=None, wfallsigma=None, tqdmout=None, progress_cb=None, lowest_time=0)

Process burst measurements over a range of DMs.

Dedisperses to each DM in dmrange and calls processBurst() on the resulting waterfall. Returns a list of all measurements as well as a dataframe. Columns of the dataframe are given by driftrate.columns

Parameters:
  • burstname (str) – name of burst to use in results dataframe

  • wfall (np.ndarray) – 2d array of burst waterfall

  • burstdm (float) – the DM of the burst in wfall

  • dmrange (list[float]) – a list of the DMs you would measurements at

  • fres_MHz (float) – frequency resolution in MHz

  • tres_ms (float) – time resolution in milliseconds

  • lowest_freq (float) – lowest frequency of the waterfall

  • p0 (list, optional) – optionally provide an initial guess for the 2d gaussian fit with list of [amplitude, x0, y0, sigmax, sigmay, theta]

  • corrsigma (float, optional) – standard deviation of the burst’s autocorrelation

  • wfallsigma (float, optional) – standard deviation of the burst waterfall

  • tqdmout (object, optional) – output for TQDM’s progress bar

  • progress_cb (function, optional) – callback to run after each DM is processed. Called with ( #of DMs processed) / len(dmrange) and a string of “{#}/{len(dmrange)}”

  • lowest_time (float, optional) – the starting time of the waterfall, passed to processBurst().

Returns:

(list of measurement results, dataframe of measurement results)

Return type:

tuple

driftrate.readRegions(resultsdf)

Parse regions out of a results csv produced by FRBGui

Parameters:

resultsdf (pd.DataFrame) – dataframe of results. Can load with pd.read_csv(filename).set_index('name')

driftrate.scilabel(num, err)

Utility for pretty scientific notation labels

e.g. (6.0 +/- 0.3) $times 10^2$ MHz/ms

Parameters:
  • num (float) – the number to label

  • err (float) – the uncertainty of num

Returns:

the formatted label string

Return type:

str

driftrate.structureParameter(wfall, dt, tstart, tend)

wip. Compute the structure parameter of a burst.

See eq. 1 in gajjar et al. 2018

Parameters:
  • wfall (np.ndarray) – burst waterfall

  • dt (float) – time resolution

  • tstart (int) – time channel to start at

  • tend (int) – time channel to end at

Returns:

the structure parameter

Return type:

float

driftrate.subband(wfall, nsub)

Downsample frequency channels of a waterfall.

See subsample() for more general method.

Parameters:
  • wfall (np.ndarray) – 2d array

  • nsub (int) – number of channels. nsub should evenly divide wfall.shape[0]

Returns:

2d subbanded array.

Return type:

np.ndarray

driftrate.subsample(m, nfreq, ntime)

Subsample a waterfall in time and frequency

Parameters:
  • m (np.ndarray) – 2d array to subsample.

  • nfreq (int) – the number of frequency channels desired. Should evenly divide m.shape[0]

  • ntime (int) – the number of time channels desired. Should evenly divide m.shape[1]

Returns:

the subsampled array

Return type:

np.ndarray

driftrate.subtractbg(wfall, tleft: int = 0, tright: int = 1)

Subtract a background sample from a waterfall

This will compute the mean along the time axis to produce an array of noise as a function of frequency and subtract that array from the entire waterfall.

Avoid sampling the burst if possible to improve accuracy of measurements.

Parameters:
  • wfall (np.ndarray) – 2d array of burst waterfall

  • tleft (int) – time channel to start sample from

  • tright (int) – time channel to end sample from

Returns:

the waterfall subtracted by the background sample

Return type:

np.ndarray

driftrate.twoD_Gaussian(point, amplitude, xo, yo, sigma_x, sigma_y, theta)

2D rotatable Gaussian Model function. Used as the model function in scipy.optimize.curve_fit

Parameters:
  • point (tuple) – (y, x) point to evaluate the model at

  • amplitude (float) – amplitude of the 2D gaussian

  • xo (float) – central x position of the gaussian

  • yo (float) – central y position of the gaussian

  • sigma_x (float) – standard deviation in x-direction

  • sigma_y (float) – standard deviation in y-direction

  • theta (float) – angle of gaussian

Returns:

1d array of raveled 2d gaussian data

Return type:

list

driftrate.updatenpz(npz, field, val)

Update a field in a npz file and resave it.

Parameters:
  • npz (str) – filename of .npz file

  • field (str) – the field you would like to update

  • val (Any) – the value to update the field to

Returns:

None

driftlaw module

driftlaw.atanmodel(B, x)
driftlaw.bakeMeasurements(sources, names, exclusions, targetDMs, logging=True, tagColors=['r', 'r', 'b', 'g', 'y', 'c', 'tomato', 'c'], showDMtraces=False, exclude_set=None)

Process completed measurements performed over a range of DMs for the purposes of analysing them in the context of the triggered relativistic dynamical model (TRDM).

Performs fits of the sub-burst slope law (i.e. slope = A/duration, A is a constant) by grouping measurements by DM. Provides detailed logging information.

The logging information of this function should be carefully reviewed to ensure the science and questions you are investigating are not being affected. Reviewing the information in this function should provide insight into how your bursts as a cohort are behaving as the DM changes.

The outputs (the “baked measurements”) of this function are used to produce figures of spectro-temporal properties of the measurements and explore their relationships.

Specifically:

  1. Filter measurements based on uncertainties and/or unphysical values. Measurements with uncertainty greater than 40% of the measurement value are discarded. Measurements with uncertainties greater than 10^8 are discarded. Measurements where the sub-burst slope is positive are discarded due to assumed over-dedispersion. Measurements with no valid fit (negative gaussian amplitude) are excluded.

  2. Color data points based on spreadsheets.

  3. Finds range of fits by source.

  4. Tags data that is at the target DM and splits up measurements by DM. This can be then used to plot burst measurements at a “representative” DM.

  5. Compute ranges of measurements and use these as conservative uncertainty bars on plots of bursts at their representative DM.

Example usage: .. code-block:: python

bakeddata, fitdata, sources, extradata = driftlaw.bakeMeasurements(sources, names, exclusions, targetDMs, logging=False, tagColors=tagColors)

See here for a (nontrivial) usage example.

Parameters:
  • sources (list[pd.DataFrame]) – list of DataFrames of frbgui results spreadsheets. Spreadsheets should be complete, that is, they are the output of computeModelDetails().

  • names (list[str]) – list of names for each dataset listed in sources. used in figures. One name per data source.

  • exclusions (list[list[str]]) – list of bursts to exclude by name. One list per data source.

  • targetDMs (list[float]) – Representative DM to display each burst at. One per data source. After running this function once it is a good idea to run getOptimalDMs() and use that list as your representative DMs. getOptimalDMs() requires the baked data so cannot be run first.

  • logging (bool, optional) – Set to False to disable logging. Do this after you’ve reviewed the logging output.

  • tagColors (list[str], optional) – list of matplotlib colors for plotting purposes. One per data source

  • showDMtraces (bool, optional) – if True will display plots of Slope vs. duration for each burst over the DM range. Useful for looking at the change in measurements as DM varies.

Returns:

(bakeddata, fitdata, sources, extradata).

  • bakkeddata is an dictionary containing “frames”, “labels”, and “colors”. Used in plotSlopeVsDuration() to produce figures of measurements at a representative DM with uncertainties estimated from the DM range used.

  • fitdata is a DataFrame of the fit results (i.e. fitting slope = A/duration) at each DM and is used by getOptimalDMs(). It contains the following columns and information:

'name': # source name
'DM': # DM fit was performed at
'param': # The fit parameter to the fit A/duration
'err': # the error on A ,
'red_chisq': # the reduced chisquared,
'numbursts': # the number of bursts used, after measurement exclusions
  • sources is the updated list of dataframes after measurement exclusions.

  • extradata Some extra information regarding the sub-burst slope vs frequency relationship.

Return type:

tuple

driftlaw.cleanAngle(row)
driftlaw.computeModelDetails(frame, channelSpaceDuration=False)

Takes a dataframe and computes columns related to the triggered relativistic dynamical model from the 2d gaussian parameters found in the measurement stage.

Importantly, computes the sub-burst slope, duration, bandwidth, and uncertainties associated with these.

Can be run even if the columns this function computes already exist and can serve to recalculate or update computed values.

Parameters:

frame (pd.DataFrame) – the results dataframe from FRBGui or driftrate.processDMRange()

Returns:

the dataframe with new computed columns

Return type:

pd.DataFrame

driftlaw.fitodr(frame, beta0=[1000], errorfunc=<function log_error>, log=True)
driftlaw.fitreciprocal(x, data, sigma=1)
driftlaw.fitreciprocal_log(x, data, sigma=1, loglog=False)
driftlaw.getOptimalDMs(fitdata, log=False)

Using the fitdata found in bakeMeasurements(), identify which DM results in the lowest reduced chi-squared when fitting the sub-burst slope law while using all the bursts.

Thus this method hypothesizes that the sub-burst slope law must hold and uses that to define an “optimal” DM.

Because some measurements are excluded, some DMs (usually on the higher end of the range) may lead to entirely excluding a burst from analysis. Here we choose to interpret this as the DM being too high and unphysical.

That is, in addition to having the best fit to the sub-burst slope law, the “optimal” DM is also defined to be the DM where all bursts in our sample have physically valid (under our assumptions) measurements.

Parameters:
  • fitdata (pd.Dataframe) – Dataframe output from bakeMeasurements()

  • log (bool, optional) – if True, ouptut tables of the fitdata with their reduced chi-squareds and number of remaining bursts by DM. Can be used to manually pick an “optimal” DM.

Returns:

list of optimal DMs by source. Use this as the input to py:meth:bakeMeasurements to produce publication figures of your measurements at a representative DM that is the optimal DM as defined here.

Return type:

list

driftlaw.limitedDMrangeerror(frame)
driftlaw.limitedDMrangeerror_odr(frame)

The range errors are asymmetric. Take the largest error

driftlaw.limitedDMrangeerror_recip(frame)
driftlaw.limitedDMsloperanges(fitdf, source, threshold=0)

Like sloperanges but only for DMs where all the bursts in a sample have valid measurements

driftlaw.log_error(frame)

see modelerror()

driftlaw.log_log(x, k, b)
driftlaw.modelerror(frame)
driftlaw.modelerror_recip(frame)
driftlaw.offset_atanmodel(B, x, zero_ddm_fit=6.554)
driftlaw.plotSlopeVsDuration(frames=[], labels=[], title=None, logscale=True, annotatei=0, markers=['o', '^', 'v', 'd', 'p'], hidefit=[], hidefitlabel=False, fitlines=['r-', 'b--', 'g-.'], fitextents=None, figsize=(17, 9), errorfunc=<function modelerror>, fiterrorfunc=<function rangelog_error>, dmtrace=False, ax=None)

Plot the normalized sub-burst slope vs. sub-burst duration for a list of baked measurements.

Finds a fit of the form

\[\frac{1}{\nu} \Big|\frac{d\nu}{dt}\Big| = \frac{A}{t}\]

Bake your spectro-temporal measurements with bakeMeasurements() first and use the outputs as the frames and labels parameters needed here.

Allows you to define how the uncertainties are calculated.

The following is an example of how to call this function using bakeddata. In this example the uncertainties are estimated from the range of measurements obtained over the range of DMs that include all bursts after exclusions are performed.

ax, _ = driftlaw.plotSlopeVsDuration(bakeddata['frames'], bakeddata['labels'], title="My Favorite Repeating FRB", logscale=True, errorfunc=driftlaw.limitedDMrangeerror, fiterrorfunc=driftlaw.log_error)

See here for a (nontrivial) usage example.

Parameters:
  • frames (list[pd.DataFrame]) – list of dataframes of measurements. Produce this with bakeMeasurements()

  • labels (list[str]) – Labels to use on plot. One per source dataframe.

  • title (str, optional) – title of plot

  • logscale (bool, optional) – if True use a logscale

  • annotatei (int or list[int], optional) – index or list of indices of frames. Indices which frames’ bursts you would like to see annotated. Useful for identifying outliers for further review.

  • markers (str or list[str], optional) – matplotlib plot markers

  • hidefit (int or list[int], optional) – index or indices of fits to hide. Starts from 0.

  • hidefitlabel (bool, optional) – If True will not display a label in the legend of fits

  • fitlines (str or list[str], optional) – matplotlib linestyles used for fits

  • fitextents (tuple, optional) – (min duration, max duration) extents over which to display the fits

  • figsize (tuple, optional) – matplotlib figure size (width, height)

  • errorfunc (function, optional) – the function you want to evaluate your measurement uncertainties. Used as part of pd.scatter. If unsure, use driftlaw.limitedDMrangeerror() or the default value.

  • fiterrorfunc (function, optional) – function for evaluating fit uncertainties. If unsure use driftlaw.log_error() or the default value.

  • dmtrace (bool, optional) – if True create an additional plot of norm. slope vs. duration for each frame showing the trace of the point over the range of DMs used to measure.

  • ax (matplotlib.axes.Axes, optional) – the axis to plot on.

Returns:

(ax, fits). ax is the figure axes and can be used to make additional modifications to the figure. fits is a list of [label, param, err, red_chisq, residuals, len(frame)].

Return type:

tuple

driftlaw.rangeerror(frame)

These ranges are not errors in the statistical sense. they are the min/max values, which should be larger than the real errors. So this is extremely conservative while also being easier to compute.

The strange shape of the returned value is due to a quirk in the way pandas handles asymmetric errors.

driftlaw.rangeerror_odr(frame)

The range errors are asymmetric. Take the largest error

driftlaw.rangelog_error(frame)

The range errors are asymmetric. Average the error

driftlaw.reciprocal(x, a)
driftlaw.reciprocal_log(x, b)
driftlaw.reciprocal_odr(B, x)
driftlaw.reciprocal_odr_log(B, x)
driftlaw.sloperanges(source)

Given all burst and model data at different trial DMs, computes the range of slopes durations across the range of trial DMs

Used as part of estimating uncertainties in bakeMeasurements().