API Documentation¶
Summary¶
sntd.fitting.fit_data ([curves, snType, …]) |
The main high-level fitting function. |
sntd.simulation.createMultiplyImagedSN (…) |
Generate a multiply-imaged SN light curve set, with user-specified time delays and magnifications. |
sntd.survey_cosmo.Survey ([N, dTL, dTT, zl, …]) |
A Survey class that enables cosmology tests with assumptions of redshift distributions and time delay/lens model precision. |
sntd.curve_io.image_lc ([zpsys]) |
SNTD class that describes each image light curve of a MISN |
sntd.curve_io.MISN ([telescopename, object_name]) |
The main object for SNTD. |
sntd.curve_io.table_factory (tables[, …]) |
This function will create a new curve object using an astropy table or tables. |
sntd.ml.realizeMicro ([arand, debug, kappas, …]) |
Creates a microcaustic realization based on Wambsganss 1990 microlens code. |
sntd.ml.microcaustic_field_to_curve (field, …) |
Convolves an expanding photosphere (achromatic disc) with a microcaustic to generate a magnification curve. |
sntd.models.unresolvedMISN (curve_models[, …]) |
Wrapper class of SNCosmo.Model to provide support for unresolved lensed SN fitting. |
Fitting¶
-
sntd.fitting.
fit_data
(curves=None, snType='Ia', bands=None, models=None, params=None, bounds={}, ignore=None, constants={}, ignore_models=[], method='parallel', t0_guess=None, effect_names=[], effect_frames=[], batch_init=None, cut_time=None, force_positive_param=[], dust=None, microlensing=None, fitOrder=None, color_bands=None, color_param_ignore=[], min_points_per_band=3, identify_micro=False, min_n_bands=1, max_n_bands=None, n_cores_per_node=1, npar_cores=4, max_batch_jobs=199, max_cadence=None, fit_colors=None, fit_prior=None, par_or_batch='parallel', batch_partition=None, nbatch_jobs=None, batch_python_path=None, n_per_node=None, fast_model_selection=True, wait_for_batch=False, band_order=None, set_from_simMeta={}, guess_amplitude=True, trial_fit=False, clip_data=False, use_MLE=False, kernel='RBF', refImage='image_1', nMicroSamples=100, color_curve=None, warning_supress=True, micro_fit_bands='all', verbose=True, **kwargs)[source]¶ The main high-level fitting function.
Parameters: - curves (
MISN
) – The MISN object containing the multiple images to fit. - snType (str) – The supernova classification
- bands (
list
ofBandpass
orstr
, orBandpass
orstr
) – The band(s) to be fit - models (
list
ofModel
or str, orModel
orstr
) – The model(s) to be used for fitting to the data - params (
list
ofstr
) – The parameters to be fit for the models inside of the parameter models - bounds (
dict
) – A dictionary with parameters in params as keys and a tuple of bounds as values - ignore (
list
ofstr
) – List of parameters to ignore - constants (
dict
) – Dictionary with parameters as keys and the constant value you want to set them to as values - ignore_models (class:~list) – List of model names to ignore, usually used if you did not specify the “models” parameter and let all models for a given SN type be chosen, but you want to ignore one or more.
- method (
str
orlist
) – Needs to be ‘parallel’, ‘series’, or ‘color’, or a list containting one or more of these - t0_guess (
dict
) – Dictionary with image names (i.e. ‘image_1’,’image_2’) as keys and a guess for time of peak as values - effect_names (
list
ofstr
) – List of effect names if model contains aPropagationEffect
. - effect_frames (
list
ofstr
) – List of the frames (e.g. obs or rest) that correspond to the effects in effect_names - batch_init (
str
) – A string to be pasted into the batch python file (e.g. extra imports or filters added to sncosmo.) - cut_time (
list
) – The start and end (rest frame) phase that you want to fit in, default accept all phases. - force_positive_param (
list
) – Optional list of parameters to always make positive. - dust (
sncosmo.PropagationEffect
) – An sncosmo dust propagation effect to include in the model - microlensing (str) – If None microlensing is ignored, otherwise should be str (e.g. achromatic, chromatic)
- fitOrder (
list
) – The order you want to fit the images if using parallel method (default chooses by npoints/SNR) - color_bands (
list
) – If using multiple methods (in batch mode), the subset of bands to use for color fitting. - color_param_ignore (
list
) – If using multiple methods, parameters you may want to fit for one method but not for color method (e.g. stretch) - min_points_per_band (int) – Only accept bands to fit with this number of points fitting other criterion (e.g. minsnr)
- identify_micro (bool) – If True, function is run to attempt to identify bands where microlensing is least problematic.
- min_n_bands (int) – Checks the SN to make sure it has this number of bands (with min_points_per_band in each)
- max_n_bands (int) – The best n bands are chosen from the data.
- n_cores_per_node (int) – The number of cores to run parallelization on per node
- npar_cores (int) – The number of cores to devote to parallelization
- max_batch_jobs (int) – The maximum number of jobs allowed by your slurm task manager.
- max_cadence (int) – To clip each image of a MISN to this cadence
- fit_colors (list) – List of colors to use in color fitting (e.g. [‘bessellb-bessellv’,’bessellb-bessellr’])
- fit_prior (
MISN
or bool) – if implementing parallel method alongside others and fit_prior is True, will use output of parallel as prior for series/color. If SNTD MISN object, used as prior for series or color. - par_or_batch (str) – if providing a list of SNe, par means multiprocessing and batch means sbatch. Must supply other batch parameters if batch is chosen, so parallel is default.
- batch_partition (str) – The name of the partition for sbatch command
- nbatch_jobs (int) – number of jobs (10 jobs for 100 light curves is 10 light curves per job)
- batch_python_path (str) – path to python you want to use for batch mode (if different from current)
- n_per_node (int) – Number of SNe to fit per node (in series) in batch mode. If none, just distributes all SNe across the number of jobs you have by default.
- fast_model_selection (bool) – If you are providing a list of models and want the best fit, turning this on will make the fitter choose based on a simple minuit fit before moving to the full sntd fitting. If false, each model will be fitted with the full sntd fitting and the best will be chosen.
- wait_for_batch (bool) – if false, submits job in the background. If true, waits for job to finish (shows progress bar) and returns output.
- band_order (
list
) – If you want colors to be fit in a specific order (e.g. B-V instead of V-B depending on band order) - set_from_simMeta (
dict
) – Dictionary where keys are model parameters and values are the corresponding key in theMISN
.images.simMeta dictionary (e.g. {‘z’:’sim_redshift’} if you want to set the model redshift based on a simulated redshift in simMeta called ‘sim_redshfit’) - guess_amplitude (bool) – If True, the amplitude parameter for the model is estimated, as well as its bounds
- trial_fit (bool) – If true, a simple minuit fit is performed to locate the parameter space for nestle fits, otherwise the full parameter range in bounds is used.
- clip_data (bool) – If true, criterion like minsnr and cut_time actually will remove data from the light curve, as opposed to simply not fitting those data points.
- use_MLE (bool) – If true, uses MLE as the parameter estimator instead of the median of the nested sampling samples
- kernel (str) – The kernel to use for microlensing GPR
- refImage (str) – The name of the image you want to be the reference image (i.e. image_1,image_2, etc.)
- nMicroSamples (int) – The number of pulls from the GPR posterior you want to use for microlensing uncertainty estimation
- color_curve (
astropy.Table
) – A color curve to define the relationship between bands for parameterized light curve model. - warning_supress (bool) – Turns on or off warnings
- micro_fit_bands (str or list of str) – The band(s) to fit microlensing. All assumes achromatic, and will fit all bands together.
- verbose (bool) – Turns on/off the verbosity flag
Returns: fitted_MISN – The same MISN that was passed to fit_data, but with new fits and time delay measurements included. List if list was provided.
Return type: Examples
>>> fitCurves=sntd.fit_data(myMISN,snType='Ia', models='salt2-extended',bands=['F110W','F125W'], params=['x0','x1','t0','c'],constants={'z':1.33},bounds={'t0':(-15,15),'x1':(-2,2),'c':(0,1)}, method='parallel',microlensing=None)
- curves (
Simulation¶
-
sntd.simulation.
createMultiplyImagedSN
(sourcename, snType, redshift, z_lens=None, telescopename='telescope', objectName='object', time_delays=[10.0, 50.0], magnifications=[2.0, 1.0], sn_params={}, av_dists={}, numImages=2, cadence=5, epochs=30, clip_time=[-30, 150], bands=['F105W', 'F160W'], start_time=None, gain=200.0, skynoiseRange=(1, 1.1), timeArr=None, zpsys='ab', zp=None, hostr_v=3.1, lensr_v=3.1, microlensing_type=None, microlensing_params=[], ml_loc=[None, None], dust_model='CCM89Dust', av_host=0.3, av_lens=0, fix_luminosity=False, minsnr=0.0, scatter=True, snrFunc=None)[source]¶ Generate a multiply-imaged SN light curve set, with user-specified time delays and magnifications.
Parameters: - sourcename (
Source
or str) – The model for the spectral evolution of the source. If a string is given, it is used to retrieve aSource
from the registry. - snType (str) – The classification of the supernova
- redshift (float) – Redshift of the source
- z_lens (float) – Redshift of the lens
- telescopename (str) – The name of the telescope used for observations
- objectName (str) – The name of the simulated supernova
- time_delays (
list
offloat
) – The relative time delays for the multiple images of the supernova. Must be same length as numImages - magnifications (
list
offloat
) – The relative magnifications for the multiple images of hte supernova. Must be same length as numImages - sn_params (
dict
) – A dictionary with SN model parameters as keys, and some sort of pdf or other callable function that returns the model parameter. - av_dists (
dict
) – A dictionary with ‘host’ or ‘lens’ as keys, and some sort of pdf or other callable function that returns the relevant av parameter. - numImages (int) – The number of images to simulate
- cadence (float) – The cadence of the simulated observations (if timeArr is not defined)
- epochs (int) – The number of simulated observations (if timeArr is not defined)
- clip_time (
list
) – Rest frame phase start and end time, will clip output table to these values. - bands (
list
ofBandpass
orstr
) – The bandpass(es) used for simulated observations - start_time (float) – Start time for the leading image. If None, start will be the first value in the time array, and the peak will be this ± the relative time delay for the leading image.
- gain (float) – Gain of the telescope “obtaining” the simulated observations (if snrFunc not defined)
- skynoiseRange (
list
offloat
) – The left and right bounds of sky noise used to define observational noise (if snrFunc not defined) - timeArr (
list
offloat
) – A list of times that define the simulated observation epochs - zpsys (str or
MagSystem
) – The zero-point system used to define the photometry - zp (float or
list
offloat
) – The zero-point used to define the photometry, list if simulating multiple bandpasses. Then this list must be the same length as bands - hostr_v (float) – The r_v parameter for the host.
- lensr_v (float) – The r_v parameter for the lens.
- microlensing_type (str) – If microlensing is to be included, defines whether it is “AchromaticSplineMicrolensing” or “AchromaticMicrolensing”
- microlensing_params (
array
orlist
ofint
) – If using AchromaticSplineMicrolensing, then this params list must give three values for [nanchor, sigmadm, nspl]. If using AchromaticMicrolensing, then this must be a microcaustic defined by a 2D numpy array - ml_loc (
list
) – List containing tuple locations of SN on microlensing map (random if Nones) - dust_model (str) – The dust model to be used for simulations, see sncosmo documentation for options
- av_host (float) – The A<sub>V</sub> parameter for the simulated dust effect in the source plane
- av_lens (float) – The A<sub>V</sub> parameter for the simulated dust effect in the lens plane
- fix_luminosity (bool) – Set the luminosity of every SN to be the peak of the distribution
- minsnr (float) – A minimum SNR threshold for observations when defining uncertainty
- scatter (bool) – Boolean that decides whether Gaussian scatter is applied to simulated observations
- snrFunc (
interp1d
ordict
) – An interpolation function that defines the signal to noise ratio (SNR) as a function of magnitude in the AB system. Used to define the observations instead of telescope parameters like gain and skynoise. This can be a dictionary, so that it’s different for each filter with filters as keys and interpolation functions as values.
Returns: MISN – A MISN object containing each of the multiply-imaged SN light curves and the simulation parameters.
Return type: MISN
Examples
>>> myMISN = sntd.createMultiplyImagedSN('salt2', 'Ia', 1.33,z_lens=.53, bands=['F110W'], zp=[26.8], cadence=5., epochs=35.,skynoiseRange=(.001,.005),gain=70. , time_delays=[10., 78.], magnifications=[7,3.5], objectName='My Type Ia SN', telescopename='HST',minsnr=5.0)
- sourcename (
Cosmology¶
-
class
sntd.survey_cosmo.
Survey
(N=10, dTL=5, dTT=5, zl=0.3, zs=0.8, P=1, calc_ensemble_P=False, name='mySurvey', sys_dTL=0, **kwargs)[source]¶ A Survey class that enables cosmology tests with assumptions of redshift distributions and time delay/lens model precision.
Parameters: - N (int) – The number of discovered glSN (overruled by zl/zs below)
- dTL (float) – The percent precision on each lens model.
- dTT (float) – The percent precision on each time delay measurement
- zl (float or list) – Redshift(s) of the lens(es). If a float, assumes you want N identical SN.
- zs (float or list) – Redshift(s) of the source(s). If a float, assumes you want N identical SN.
- P (float or list) – The probability of each source/lens redshift combo, defaults to 1 (i.e. equal weight)
- calc_ensemble_P (bool) – If True, probabilities are calculated based on gaussian distributions
- name (str) – Name of your survey
-
dTc
¶ Returns the combined lens/delay uncertainty assuming statistical uncertainties
-
plot_survey_contour
(params, math_labels=None, color='#1f77b4', filled=True, confidence=[0.68, 0.95], fom=False, ax=None, alphas=[0.9, 0.3], show_legend=False, show_nestle=True, use_H0=False, **kwargs)[source]¶ Plots the contours of a nestle or gridded survey
Parameters: - params (list) – 2 parameters to plot that have been included in a survey grid or nestle survey
- math_labels (list) – list of latex labels for params (for labeling plot)
- color (str or tuple) – color for plotting in matplotlib
- filled (bool) – filled vs. outlined (if false) contours
- confidence (list or float) – confidence interval(s) to plot
- fom (bool) – if True, FOM is calculated. MAJOR WARNING: This will be enormously biased if the full contour is not being drawn, set limits accordingly
- ax (:class:'~plt.axes') – If you want the contours drawn onto existing axes
- alphas (list or float) – draws confidence intervals with this alpha, must match up with “confidence” parameter
- show_legend (bool) – If True, a legend is shown
- show_nestle (bool) – If both nestle and grid have been completed, choose nestle if true to show
- use_H0 (bool) – If true and one of the params is ‘h’, multiply by 100
-
plot_survey_gradient
(params, math_labels=None)[source]¶ Plots the gradient survey grid, where one parameter is assumed to be systematically biased
Parameters:
-
survey_fisher
(params, dx=1e-06)[source]¶ Create a fisher matrix using params, based on the overarching survey parameters.
Parameters:
-
survey_grid
(vparam_names, bounds={}, npoints=100, grad_param=None, constants={}, grad_param_bounds=None, ngrad=10, grid_param1=None, grid_param2=None, **kwargs)[source]¶ Calculate cosmological contours by varying 2 parameters in a grid.
Parameters: - vparam_names (list) – The names of parameters to vary
- bounds (dict) – Dictionary with param names as keys and list/tuple/array of bounds as values
- npoints (int) – The number of sample points
- grad_param (str) – Parameter to assume we’re to have measured wrong (see C&M 2009 Figure 4)
- constants (dict) – Constants that are not the defaults
- grad_param_bounds (dict) – Bounds for grad_param, same format as bounds
- ngrad (int) – Number of grid points to vary grad_param
- grid_param1 (iterable) – Optional choice of grid param 1 list (default uniform based on bounds)
- grid_param2 (iterable) – Optional choice of grid param 2 list (default uniform based on bounds)
Returns: - Adds to class attribute “grid” (a dictionary), with a comma-separated list of
- the vparam_names as the key and grid values as the value.
-
survey_nestle
(vparam_names, bounds, constants={}, npoints=100, **kwargs)[source]¶ Calculate cosmological contours in an MCMC-like fashion.
Parameters: - vparam_names (list) – The names of parameters to vary
- bounds (dict) – Dictionary with param names as keys and list/tuple/array of bounds as values
- npoints (int) – The number of sample points
- grad_param (str) – Parameter to assume we’re to have measured wrong (see C&M 2009 Figure 4)
- constants (dict) – Constants that are not the defaults
- grad_param_bounds (dict) – Bounds for grad_param, same format as bounds
- ngrad (int) – Number of grid points to vary grad_param
Returns: - Adds to class attribute “grid” (a dictionary), with a comma-separated list of
- the vparam_names as the key and grid values as the value.
-
class
sntd.survey_cosmo.
Fisher
(inroot='', xvar='', yvar='', fixes='', margs=[], data=[], data_is_cov=False, params=[], silent=False, name='my_fisher', cosmo_truths={'Ode0': 0.7, 'Om0': 0.3, 'h': 0.7, 'w': 0, 'w0': 0, 'wa': -1})[source]¶ Fisher class is more or less copied from Dan Coe’s arxiv paper about Fisher matrices: https://arxiv.org/pdf/0906.4123.pdf
Only some minor adjustment from me.
Parameters: - inroot (str) – base directory to read fisher matrices from if using “load” function
- xvar (str) – can optionally set an x variable as the “default for other functions”
- yvar (str) – can optionally set a y variable as the “default for other functions”
- fixes (str) – comma-separated str of parameters to fix (if fix is called)
- margs (list) – list of parameters to marginalize over (if marg is called)
- data (np.ndarray) – an input fisher matrix (NxN where N is the length of your params)
- data_is_cov (bool) – If true, assumes that “data” is actually a covariance matrix
- params (list) – Parameter names
- silent (bool) – verbosity flag for some functions
- name (str) – Name to plot in legend for figures
- cosmo_trues (dict) – True parameters to plot as dashed lines in plots
-
addpar
(param)[source]¶ Add a parameter to the list of parameters
Parameters: param (str) – The parameter to add
-
dxdyp
(xvar='', yvar='')[source]¶ Return uncertainty in two parameters and their correlation
- xvar: str
- x variable
- yvar: str
- y variable
Returns: - dx (float) – x uncertainty
- dy (float) – y uncertainty
- p (float) – rho (correlation parameter)
-
fix
(fixes=[])[source]¶ Fix parameters constant <==> Remove them from the Fisher matrix
Parameters: fixes (list) – List of parameters to fix, otherwise uses fixes attribute
-
load
(fishdir='')[source]¶ Loads existing Fisher matrix
Parameters: fishdir (str) – The directory (default inroot)
-
marg
(margs=[])[source]¶ Marginalize over variables: Remove them from the covariance matrix
Parameters: margs (list) – List of parameters to marginalize over, otherwise uses margs attribute
-
merit
(param1, param2)[source]¶ Calculates the figure of merit for a pair of parameters.
Parameters: Returns: FOM – The Figure of Merit
Return type:
-
pindex
(param)[source]¶ Index of parameter
Parameters: param (str) – Parameter you want the index of Returns: index – The index of param Return type: int
-
plot
(param1, param2, x_limits, y_limits, math_label1=None, math_label2=None, bestfit1=None, bestfit2=None, alpha=0.9, color_list=None, print_merit=True, col_order=True, show_uncertainty=False)[source]¶ Plot contours from fisher matrix. This will plot all contours from matrices that have been added together make this matrix.
Parameters: - param1 (str) – Parameter 1 to plot
- param2 (str) – Parameter 2 to plot
- xlimits (list or tuple or
ndarray
) – The x parameter limits for plotting - ylimits (list or tuple or
ndarray
) – The y parameter limits for plotting - math_label1 (str) – A latex label for axis 1
- math_label2 (str) – A latex label for axis 2
- bestfit1 (float) – The true/best fit value for parameter 1 (default self.cosmo_truths)
- bestfit2 (float) – The true/best fit value for parameter 2 (default self.cosmo_truths)
- alpha (float) – The alpha for plotting
- color_list (list) – List of colors to use for plotting
- print_merit (bool) – If True, the figure of merit it calculated and added to the legend
- show_uncertainty (bool) – If true, the final uncertainty in each parameter is printed on the plot
-
rename
(pdict1=None)[source]¶ Rename parameters given a dictionary of names & nicknames
Parameters: pdict1 (dict) – Dictionary containing old parameters as keys and new as vals
I/O¶
-
class
sntd.curve_io.
image_lc
(zpsys='AB')[source]¶ SNTD class that describes each image light curve of a MISN
-
bands
= None¶ str band names used
Type: @type
-
chisq
¶ A property that calculates the chisq based on your best fit model.
-
fits
= None¶ newDict
Contains fit information from fit_dataType: @type
-
meta
= None¶ dict
The metadata for the MISN object, intialized with an empty “info” key value pair. It’s populated by added _metachar__ characters into the header of your data file.Type: @type
-
reduced_chisq
¶ A property that calculates the reduced chisq based on your best fit model.
-
simMeta
= None¶ dict
A dictionary containing simulation metadata if this is a simulated curve objectType: @type
-
table
= None¶ astropy.table.Table
A table containing the data, used for SNCosmo functions, etc.Type: @type
-
zpsys
= None¶ str The zero-point system for this curve object
Type: @type
-
-
class
sntd.curve_io.
MISN
(telescopename='Unknown', object_name='Unknown')[source]¶ The main object for SNTD. This organizes a MISN, containing the multiple light curves in self.images, all the fits, etc.
-
add_image_lc
(myImage, key=None)[source]¶ Adds an image_lc object to the existing MISN (i.e. adds an image to a MISN)
Parameters: - myImage (
image_lc
) – The curve to add to self. - key (str) – The key you want to save this as, default is ‘image_1,image_2,etc.’
Returns: self
Return type: - myImage (
-
clip_data
(im, minsnr=-inf, mintime=-inf, maxtime=inf, peak=0, remove_bands=[], max_cadence=None, rm_NaN=True)[source]¶ Clips the data of an image based on various properties.
Parameters: - im (str) – The image to clip
- minsnr (float) – Clip based on a minimum SNR
- mintime (float) – Clip based on a minimum time (observer frame relative to peak)
- maxtime (float) – Clip based on a maximum time (observer frame relative to peak)
- peak (float) – Used in conjunction with min/max time
- remove_bands (list) – List of bands to remove from the light curve
- max_cadence (float) – Clips data so that points are spread by at least max_cadence
- rm_NaN (bool) – If True, removed NaNs from the data.
-
color_chisq
¶ A property that calculates the chisq based on your best fit color model.
-
color_reduced_chisq
¶ A property that calculates the reduced chisq based on your best fit color model.
-
color_table
(band1s, band2s, time_delays=None, referenceImage='image_1', ignore_images=[], static=False, model=None, minsnr=0.0)[source]¶ Takes the multiple images in self.images and combines the data into a single color curve using defined time delays and magnifications or best (quick) guesses.
Parameters: - band1s (str or list) – The first band(s) for color curve(s)
- band2s (str or list) – The second band(s) for color curve(s)
- time_delays (
dict
) – Dictionary with image names as keys and relative time delays as values (e.g. {‘image_1’:0,’image_2’:20}). Guessed if None. - referenceImage (str) – The image you want to be the reference (e.g. image_1, image_2, etc.)
- ignore_images (
list
) – List of images you do not want to include in the color curve. - static (bool) – Make the color curve, don’t shift the data
- model (
Model
) – If you want to use an sncosmo Model (and the guess_t0_amplitude method) to guess time delays - minsnr (float) – Cut data that don’t meet this threshold before making the color curve.
Returns: self
Return type:
-
combine_curves
(time_delays=None, magnifications=None, referenceImage='image_1', static=False, model=None, minsnr=0)[source]¶ Takes the multiple images in self.images and combines the data into a single light curve using defined time delays and magnifications or best (quick) guesses.
Parameters: - time_delays (
dict
) – Dictionary with image names as keys and relative time delays as values (e.g. {‘image_1’:0,’image_2’:20}). Guessed if None. - magnifications (
dict
) – Dictionary with image names as keys and relative magnifications as values (e.g. {‘image_1’:0,’image_2’:20}). Guessed if None. - referenceImage (str) – The image you want to be the reference (e.g. image_1, image_2, etc.)
- ignore_images (
list
) – List of images you do not want to include in the color curve. - static (bool) – Make the color curve, don’t shift the data
- model (
Model
) – If you want to use an sncosmo Model (and the guess_t0_amplitude method) to guess time delays - minsnr (float) – Cut data that don’t meet this threshold before making the color curve.
Returns: self
Return type: - time_delays (
-
meta
= None¶ dict
The metadata for the MISN object, intialized with an empty “info” key value pair. It’s populated by added _metachar__ characters into the header of your data file.Type: @type
-
object
= None¶ str Object of interest
Type: @type
-
plot_fit
(method='parallel', par_image=None)[source]¶ Makes a corner plot based on one of the fitting methods
Parameters: method (str) – parallel, series, or color Returns: figure object Return type: figure
-
plot_microlensing_fit
(show_all_samples=False)[source]¶ Shows a plot of the best-fit microlensing curve from GPR.
Parameters: show_all_samples (bool) – If True, show all GPR samples. Returns: figure object Return type: figure
-
plot_object
(bands='all', savefig=False, plot3D=False, filename='mySN', orientation='horizontal', method='separate', showModel=False, showFit=False, showMicro=False, plot_unresolved=False, **kwargs)[source]¶ - Plot the multiply-imaged SN light curves and show/save to a file.
- Each subplot shows a single-band light curve, for all images of the SN.
Parameters: - bands (str or
list
ofstr
) – ‘all’ = plot all bands; or provide a list of bands to plot - savefig (bool) – boolean to save or not save plot
- plot3D (bool) – boolean to plot in 3D with plotly
- filename (str) – if savefig is True, this is the output filename
- orientation (str) – ‘horizontal’ = all subplots are in a single row ‘vertical’ = all subplots are in a single column
- method (str) – Plots the result of separate, series, or color curve method
- showModel (bool) – If true, the underlying model before microlensing is plotted as well
- showFit (bool) – If true and it exists, the best fit model from self.images[‘image’].fits.model is plotted
- showMicro (bool) – If true and it exists, the simulated microlensing is plotted as well
- plot_unresolved (bool) – If true the individual models of an unresolved fit will be plotted
Returns: figure
Return type: ~matplotlib.pyplot.figure
-
quality_check
(min_n_bands=1, min_n_points_per_band=1, clip=False, method='parallel')[source]¶ Checks the images of a SN to make sure they pass minimum thresholds for fitting.
Parameters: - min_n_bands (int) – The minimum number of bands needed to pass
- min_n_points_per_band (int) – The minimum number of bands in a given band to pass
- clip (bool) – If True, “bad” bands are clipped in place
- method (str) – Should be parallel, series, or color. Checks all images (parallel), or the series table (series), or the color table (color)
Returns: self
Return type:
-
series_chisq
¶ A property that calculates the chisq based on your best fit series model.
-
series_reduced_chisq
¶ A property that calculates the reduced chisq based on your best fit series model.
-
table
= None¶ Table
The astropy table containing all of the data in your data fileType: @type
-
telescopename
= None¶ str Name of the telescope that the data were gathered from
Type: @type
-
-
sntd.curve_io.
read_data
(filename, **kwargs)[source]¶ - Used to read a light curve or curve object in pickle format.
- Either way, it’ll come out as a curve object.
Parameters: filename (str) – Name of the file to be read (ascii or pickle) Returns: curve Return type: curve
orMISN
-
sntd.curve_io.
write_data
(curves, filename=None, protocol=-1)[source]¶ - Used to write a MISN object to a pickle
- to be read later
Parameters: Returns: Return type:
-
sntd.curve_io.
table_factory
(tables, telescopename='Unknown', object_name='Unknown')[source]¶ This function will create a new curve object using an astropy table or tables.
Parameters: - tables (~astropy.table.Table or
list
) – Astropy table with all of your data from your data file, or a list of such tables. - telescopename (str) – Name of telescope for labeling purposes inside curve object
- object_name (str) – Name of object for labeling purposes inside curve object (e.g. SN2006jf, etc.)
Returns: curve
Return type: curve
- tables (~astropy.table.Table or
Microlensing¶
-
sntd.ml.
realizeMicro
(arand=0.25, debug=0, kappas=0.75, kappac=0.15, gamma=0.76, eps=0.6, nray=300, minmass=10, maxmass=10, power=-2.35, pixmax=5, pixminx=0, pixminy=0, pixdif=10, fracpixd=0.3, iwrite=0, verbose=False)[source]¶ Creates a microcaustic realization based on Wambsganss 1990 microlens code. All parameters are optional as they have defaults, see Wambsganss documentation for details on parameters.
-
sntd.ml.
microcaustic_field_to_curve
(field, time, zl, zs, velocity=<Quantity 10000. km / s>, M=<Quantity 1.98840987e+30 kg>, width_in_einstein_radii=10, loc='Random', plot=False, ax=None, showCurve=True, rescale=True)[source]¶ Convolves an expanding photosphere (achromatic disc) with a microcaustic to generate a magnification curve.
Parameters: - field (
numpy.ndarray
) – An opened fits file of a microcaustic, can be generated by realizeMicro - time (
numpy.array
) – Time array you want for microlensing magnification curve, explosion time is 0 - zl (float) – redshift of the lens
- zs (float) – redshift of the source
- velocity (float*
astropy.units.Unit
) – The average velocity of the expanding photosphere - M (float*
Unit
) – The mass of the deflector - width_in_einstein_radii (float) – The width of your map in units of Einstein radii
- loc (str or tuple) – Random is defualt for location of the supernova, or pixel (x,y) coordiante can be specified
- plot (bool) – If true, plots the expanding photosphere on the microcaustic
- ax (~matplotlib.pyplot.axis) – An optional axis object to plot on. If you want to show the curve, this should be a list like this: [main_ax,lower_ax]
- showCurve (bool) – If true, the microlensing curve is plotted below the microcaustic
- rescale (bool) – If true, assumes image needs to be rescaled: (x-1024)/256
Returns: - time (
numpy.array
) – The time array for the magnification curve - dmag (
numpy.array
) – The magnification curve.
- field (
Models¶
-
class
sntd.models.
unresolvedMISN
(curve_models, delays=None, magnifications=None)[source]¶ Wrapper class of SNCosmo.Model to provide support for unresolved lensed SN fitting.
-
add_effect
(effect, name, frame)[source]¶ Add a PropagationEffect to the model.
Parameters: - effect (~sncosmo.PropagationEffect) – Propagation effect.
- name (str) – Name of the effect.
- frame ({'rest', 'obs', 'free'}) –
-