snat_sim.models

The modeling module defines classes for modeling physical phenomena. This includes SNe Ia light-curves, the propagation of light through atmospheric water vapor (with and without variation in time), and the broad-band photometry of stellar catalogs.

Model Summaries

A summary of the available models is provided below:

PWVTransmissionModel

Models atmospheric transmission due to PWV at a fixed resolution

PWVModel

Model for interpolating the atmospheric water vapor at a given date and time

SNModel

An observer-frame supernova model composed of a Source and zero or more effects

ReferenceStar

Representation of modeled stellar spectra from the Goettingen Spectral Library

ReferenceCatalog

A rudimentary implementation of a reference star catalog

VariableCatalog

A reference star catalog that determines the time dependent PWV concentration along the line of sight from an underlying PWV model

Supernova models (provided by the SNModel class) are designed to closely resemble the behavior of the sncosmo package. However, unlike sncosmo.Model objects, the snat_sim.SNModel class provides support for propagation effects that vary with time. A summary of custom propagation effects provided by the snat_sim package is listed below:

StaticPWVTrans

Propagation effect for the atmospheric absorption of light due to time static PWV

SeasonalPWVTrans

Atmospheric propagation effect for a fixed PWV concentration per-season

VariablePWVTrans

Propagation effect for the atmospheric absorption of light due to time variable PWV

Simulation results are generally returned using one or more of the following object types:

ObservedCadence

The observational sampling of an astronomical light-curve

LightCurve

Abstract representation of an astronomical light-curve.

SNFitResult

Represents results from a SNModel being fit to a LightCurve

Usage Example

To ensure ease of use, supernova modeling with the snat_sim package is backward compatible with the sncosmo package and follows a similar design pattern. The SNModel class from snat_sim acts as a drop-in replacement for the Model classes built into sncosmo but provides extended functionality.

Supernova models are instantiated using a spectral template with the addition of optional observer or rest-frame effects. In the following example, atmospheric propagation effects due to precipitable water vapor are added to a Salt2 supernova model.

>>> from snat_sim import models

>>> # Create a supernova model
>>> supernova_model = models.SNModel('salt2')

>>> # Create a model for the atmosphere
>>> atm_transmission = models.StaticPWVTrans()
>>> atm_transmission.set(pwv=4)
>>> supernova_model.add_effect(effect=atm_transmission, name='Atmosphere', frame='obs')

Unlike sncosmo, the process of simulating or fitting light-curves is handled directly from the model class. We consider these two cases seperately in the following sections.

Simulating Light-Curves

To simulate a light-curve, you first need to establish the desired light-curve cadence (i.e., how the light-curve should be sampled in time):

>>> cadence = models.ObservedCadence(
...     obs_times=[-1, 0, 1],
...     bands=['sdssr', 'sdssr', 'sdssr'],
...     zp=25, zpsys='AB', skynoise=0, gain=1
... )

Light-curves can then be simulated directly from the model using the simulate_lc method. Notice that the returned type is a LightCurve object:

>>> # Here we simulate a light-curve with statistical noise
>>> light_curve_with_noise = supernova_model.simulate_lc(cadence)
>>> print(type(light_curve_with_noise))
<class 'snat_sim.models.light_curve.LightCurve'>

>>> # Here we simulate a light-curve with a fixed signal to noise ratio
>>> light_curve_fixed_snr = supernova_model.simulate_lc(cadence, fixed_snr=5)
>>> print(type(light_curve_fixed_snr))
<class 'snat_sim.models.light_curve.LightCurve'>

The LightCurve class represents astronomical light-curves and provides an easy interface for casting the data into other commonly used object types.

>>> import sncosmo
>>> from snat_sim.models import LightCurve

>>> example_data = sncosmo.load_example_data()
>>> light_curve_data = LightCurve.from_sncosmo(example_data)

>>> lc_as_dataframe = light_curve_data.to_pandas()
>>> lc_as_table = light_curve_data.to_astropy()

Fitting Light-Curves

Photometric data can be fit directly from the model using the fit_lc method. sncosmo users should note that the return signature is not the same as sncosmo.:

>>> supernova_model.set(z=.5, t0=55100.0)
>>> fit_result = supernova_model.fit_lc(
...     light_curve_data, vparam_names=['x0', 'x1', 'c'])

>>> print(type(fit_result))
<class 'snat_sim.models.supernova.SNFitResult'>

The SNFitResult object is similar to the Result class from sncosmo but is not backward compatible. SNFitResult instances provide access to fit results as floats or pandas objects (e.g., pandas.Series or pandas.DataFrame) depending on the value.

>>> print(fit_result)
     success: True
     message: Minimization exited successfully.
       ncall: 77
        nfit: 1
       chisq: 35.537
        ndof: 37
 param_names: ['z', 't0', 'x0', 'x1', 'c', 'Atmospherepwv']
  parameters: [5.000e-01 5.510e+04 1.194e-05 4.257e-01 2.483e-01 4.000e+00]
vparam_names: ['x0', 'x1', 'c']
      errors: [3.830e-07 3.173e-01 2.931e-02]
  covariance:
    [[ 1.467e-13 -6.472e-08 -7.966e-09]
     [-6.472e-08  1.007e-01  1.174e-03]
     [-7.966e-09  1.174e-03  8.590e-04]]

Fit result objects are also capable of calculating the variance in the distance modulus given the alpha/beta standardization parameters:

>>> from snat_sim import constants as const

>>> # The covariance matrix used when determining error values
>>> fit_result.salt_covariance_linear()
          mB        x1         c
mB  0.001213  0.005886  0.000725
x1  0.005886  0.100684  0.001174
c   0.000725  0.001174  0.000859

>>> # Here we use alpha and beta parameters from Betoule et al. 2014
>>> mu_variance = fit_result.mu_variance_linear(
...     alpha=const.betoule_alpha, beta=const.betoule_beta)

>>> print(f'{mu_variance: .5f}')
 0.00762

Calibrating Light-Curves

Your use case may involve calibrating simulated light-curves relative to a stellar reference catalog. The spectrum for individual spectral types can be retreived using the ReferenceStar class:

>>> from snat_sim.models import reference_star

>>> g2_star = reference_star.ReferenceStar('G2')
>>> print(g2_star.to_pandas())
3000.000     4.960049e+17
3000.006     4.659192e+17
3000.012     4.304657e+17
3000.018     3.751426e+17
3000.024     2.847191e+17
                 ...
11999.920    1.366567e+18
11999.940    1.366673e+18
11999.960    1.366418e+18
11999.980    1.365863e+18
12000.000    1.365315e+18
Length: 933333, dtype: float32

A ReferenceCatalog is used to represent a collection of stars with different stellar types. Catalog instances can be used to calibrate supernoca light-curves.

>>> reference_catalog = reference_star.ReferenceCatalog('G2', 'M5')
>>> print(reference_catalog.calibrate_lc(light_curve_data, pwv=4))

Module Docs

Reference Stars API

Models for the spectro-photometric flux of stars and their combination into reference catalogs.

class snat_sim.models.reference_star.ReferenceStar(spectral_type)[source]

Representation of modeled stellar spectra from the Goettingen Spectral Library

__init__(spectral_type)[source]

Load a spectrum for the given spectral type

Flux values are returned in phot/cm2/s/angstrom and are index by wavelength values in Angstroms.

Parameters:

spectral_type (str) – Spectral type (e.g., G2)

to_pandas()[source]

Return the spectral template as a pandas.Series object

Return type:

Series

static get_available_types()[source]

Return the spectral types available on disk

Return type:

List[str]

Returns:

A list fo spectral types

flux_scaling_dataframe()[source]

Retrieve pre-tabulated scale factors uses when calibrating flux values for PWV absorption

Return type:

DataFrame

Returns:

A DataFrame indexed by PWV with columns for flux

flux(band, pwv)[source]

Return the reference star flux values

Parameters:
  • band (str) – Band to get flux for

  • pwv (Union[float, float16, float32, float64, float128, int, int8, int16, int32, int64, uint8, uint16, uint32, uint64, longlong, ulonglong, array]) – PWV concentration along the line of sight

Return type:

ndarray

Returns:

The normalized flux at the given PWV value(s)

norm_flux(band, pwv)[source]

Return the normalized reference star flux values

Parameters:
  • band (str) – Band to get flux for

  • pwv (Union[float, float16, float32, float64, float128, int, int8, int16, int32, int64, uint8, uint16, uint32, uint64, longlong, ulonglong, array]) – PWV concentration along the line of sight

Return type:

ndarray

Returns:

The normalized flux at the given PWV value(s)

class snat_sim.models.reference_star.ReferenceCatalog(*spectral_types)[source]

A rudimentary implementation of a reference star catalog

__init__(*spectral_types)[source]

Create a reference star catalog composed of the given spectral types

Parameters:

*spectral_types – Spectral types for the catalog (e.g., ‘G2’, ‘M5’, ‘K2’)

average_norm_flux(band, pwv)[source]

Return the average normalized reference star flux

Parameters:
  • band (str) – Band to get flux for

  • pwv (Union[float, float16, float32, float64, float128, int, int8, int16, int32, int64, uint8, uint16, uint32, uint64, longlong, ulonglong, Collection, ndarray]) – PWV concentration along the line of sight

Return type:

ndarray

Returns:

The normalized flux at the given PWV value(s)

calibrate_lc(light_curve, pwv)[source]

Divide normalized reference flux from a light-curve

Recalibrate flux values using the average change in flux of a collection of reference stars.

Parameters:
  • light_curve (LightCurve) – Light curve to calibrate

  • pwv (Union[float, float16, float32, float64, float128, int, int8, int16, int32, int64, uint8, uint16, uint32, uint64, longlong, ulonglong, ndarray]) – PWV concentration along the line of sight

Return type:

LightCurve

Returns:

A modified copy of the light_curve argument

class snat_sim.models.reference_star.VariableCatalog(pwv_model, *spectral_types)[source]

A reference star catalog that determines the time dependent PWV concentration along the line of sight from an underlying PWV model

__init__(pwv_model, *spectral_types)[source]

Create a reference star catalog composed of the given spectral types and a PWV model

Parameters:
  • *spectral_types – Spectral types for the catalog (e.g., ‘G2’, ‘M5’, ‘K2’)

  • pwv_model (PWVModel) – The PWV model to determine the zenith PWV concentration from

average_norm_flux(band, time, ra, dec, lat=-30.244573, lon=-70.7499537, alt=1024, time_format='mjd')[source]

Return the average normalized reference star flux

Parameters:
  • band (str) – Band to get flux for

  • time (Union[float, ndarray, Collection]) – Time at which the target is observed

  • ra (float) – Right Ascension of the target (Deg)

  • dec (float) – Declination of the target (Deg)

  • lat (float) – Latitude of the observer (Deg)

  • lon (float) – Longitude of the observer (Deg)

  • alt (float) – Altitude of the observer (m)

  • time_format (str) – Astropy supported format of the time value (Default: ‘mjd’)

Return type:

ndarray

Returns:

The normalized flux at the given PWV value(s)

calibrate_lc(light_curve, ra, dec, lat=-30.244573, lon=-70.7499537, alt=1024, time_format='mjd')[source]

Divide normalized reference flux from a light-curve

Recalibrate flux values using the average change in flux of a collection of reference stars.

Parameters:
  • light_curve (LightCurve) – Yhr light-curve to calibrate

  • ra (float) – Right Ascension of the target (Deg)

  • dec (float) – Declination of the target (Deg)

  • lat (float) – Latitude of the observer (Deg)

  • lon (float) – Longitude of the observer (Deg)

  • alt (float) – Altitude of the observer (m)

  • time_format (str) – Astropy supported format of the time value (Default: ‘mjd’)

Return type:

LightCurve

Returns:

A modified copy of the light_curve argument

PWV Effects API

Modeling functionality for Precipitable Water Vapor (PWV) and related observational effects.

class snat_sim.models.pwv.PWVModel(pwv_series)[source]

Model for interpolating the atmospheric water vapor at a given date and time

__init__(pwv_series)[source]

Build a model for time variable PWV by interpolating from a given PWV time series

Parameters:

pwv_series (Series) – PWV values with a datetime index

static from_suominet_receiver(receiver, year, supp_years=None)[source]

Construct a PWVModel instance using data from a SuomiNet receiver

Parameters:
  • receiver (GPSReceiver) – GPS receiver to access data from

  • year (int) – Year to use data from when building the model

  • supp_years (Optional[Collection[int]]) – Years to supplement data with when missing from year

Return type:

PWVModel

Returns:

A PWV model that interpolates from data taken by the given receiver

static calc_airmass(time, ra, dec, lat=-30.244573, lon=-70.7499537, alt=1024, time_format='mjd', raise_below_horizon=True)[source]

Calculate the airmass through which a target is observed

Default latitude, longitude, and altitude are set to the Rubin Observatory.

Parameters:
  • time – Time at which the target is observed

  • ra – Right Ascension of the target (Deg)

  • dec – Declination of the target (Deg)

  • lat – Latitude of the observer (Deg)

  • lon – Longitude of the observer (Deg)

  • alt – Altitude of the observer (m)

  • time_format – Astropy supported format of the time value (Default: ‘mjd’)

  • raise_below_horizon – If true, raise a ValueError for an airmasses less than 1

Returns:

Airmass in units of Sec(z)

pwv_zenith(time, time_format='mjd')[source]

Interpolate the PWV at zenith as a function of time

The time_format argument can be set to None when passing datetime objects for time instead of numerical values.

Parameters:
  • time – The time to interpolate PWV for

  • time_format – Astropy supported format of the time value (Default: ‘mjd’)

Returns:

The PWV at zenith for the given time(s)

pwv_los(time, ra, dec, lat=-30.244573, lon=-70.7499537, alt=1024, time_format='mjd')[source]

Interpolate the PWV along the line of sight for the given time

The time_format argument can be set to None when passing datetime objects instead of numerical values for time.

Parameters:
  • time – Time at which the target is observed

  • ra – Right Ascension of the target (Deg)

  • dec – Declination of the target (Deg)

  • lat – Latitude of the observer (Deg)

  • lon – Longitude of the observer (Deg)

  • alt – Altitude of the observer (m)

  • time_format – Astropy supported format of the time value (Default: ‘mjd’)

Return type:

Union[float, array]

Returns:

The PWV concentration along the line of sight to the target

seasonal_averages()[source]

Calculate the average PWV in each season

Assumes seasons based on equinox and solstice dates in the year 2020.

Return type:

Dict[str, float]

Returns:

A dictionary with the average PWV in each season (in mm)

class snat_sim.models.pwv.PWVTransmissionModel(resolution=None)[source]

Models atmospheric transmission due to PWV at a fixed resolution

__init__(resolution=None)[source]

Instantiate a PWV transmission model at the given resolution

Transmission values are determined using the v1_transmission model from the pwv_kpno package.

Parameters:

resolution (Optional[float]) – Resolution to bin the atmospheric model to

calc_transmission(pwv, wave=None)[source]

Evaluate the transmission model at the given wavelengths

Returns a Series object if pwv is a scalar, and a DataFrame object if pwv is an array. Wavelengths are expected in angstroms.

Parameters:
  • pwv – Line of sight PWV to interpolate for

  • wave – Wavelengths to evaluate transmission (Defaults to samp_wave attribute)

Returns:

The interpolated transmission at the given wavelengths / resolution

class snat_sim.models.pwv.StaticPWVTrans(transmission_res=5)[source]

Propagation effect for the atmospheric absorption of light due to time static PWV

The value of the pwv parameter reflects the total PWV concentration along the line of sight. It is NOT scaled using the airmass through which an object is observed.

Effect Parameters:

pwv: Atmospheric concentration of PWV along line of sight in mm

__init__(transmission_res=5)[source]

Time independent atmospheric transmission effect due to PWV

Setting the transmission_res argument to None results in the highest available transmission model available.

Parameters:

transmission_res (float) – Reduce the underlying transmission model by binning to the given resolution

property transmission_res: float

Resolution used when binning the underlying atmospheric transmission model

Return type:

float

propagate(wave, flux, *args)[source]

Propagate the flux through the atmosphere

Parameters:
  • wave (ndarray) – A 1D array of wavelength values

  • flux (ndarray) – An array of flux values

Return type:

ndarray

Returns:

An array of flux values after suffering from PWV absorption

class snat_sim.models.pwv.AbstractVariablePWVEffect(transmission_res=5)[source]

Base class for building time-variable PWV propagation effects

__init__(transmission_res=5)[source]

Time variable atmospheric transmission effect due to PWV

Setting the transmission_res argument to None results in the highest available transmission model available.

Effect Parameters:

pwv: Atmospheric concentration of PWV along line of sight in mm

Parameters:

transmission_res (float) – Reduce the underlying transmission model by binning to the given resolution

abstract assumed_pwv(time)[source]

The PWV concentration used by the propagation effect at a given time

Parameters:

time (TypeVar(FloatColl, Union[float, float16, float32, float64, float128, int, int8, int16, int32, int64, uint8, uint16, uint32, uint64, longlong, ulonglong], Collection[Union[float, float16, float32, float64, float128, int, int8, int16, int32, int64, uint8, uint16, uint32, uint64, longlong, ulonglong]], ndarray)) – Time to get the PWV concentration for

Return type:

TypeVar(FloatColl, Union[float, float16, float32, float64, float128, int, int8, int16, int32, int64, uint8, uint16, uint32, uint64, longlong, ulonglong], Collection[Union[float, float16, float32, float64, float128, int, int8, int16, int32, int64, uint8, uint16, uint32, uint64, longlong, ulonglong]], ndarray)

Returns:

An array of PWV values in units of mm

propagate(wave, flux, time)[source]

Propagate the flux through the atmosphere

Parameters:
  • wave (ndarray) – An array of wavelength values

  • flux (ndarray) – An array of flux values

  • time (Union[float, ndarray]) – Array of time values to determine PWV for

Return type:

ndarray

Returns:

An array of flux values after suffering from PWV absorption

class snat_sim.models.pwv.VariablePWVTrans(pwv_model, time_format='mjd', transmission_res=5.0)[source]

Propagation effect for the atmospheric absorption of light due to time variable PWV

The value of the pwv parameter reflects the total PWV concentration at zenith. This value is scaled using the airmass as calculated using the position of the observer (lat, lon, and alt parameters) and the object (ra and dec parameters).

Effect Parameters:

ra: Target Right Ascension in degrees dec: Target Declination in degrees lat: Observer latitude in degrees (defaults to location of VRO) lon: Observer longitude in degrees (defaults to location of VRO) alt: Observer altitude in meters (defaults to height of VRO)

__init__(pwv_model, time_format='mjd', transmission_res=5.0)[source]

Time variable atmospheric transmission due to PWV

Setting the transmission_res argument to None results in the highest available transmission model available.

Parameters:
  • pwv_model (PWVModel) – Returns PWV at zenith for a given time value and time format

  • time_format (str) – Astropy recognized time format used by the pwv_interpolator

  • transmission_res (float) – Reduce the underlying transmission model by binning to the given resolution

assumed_pwv(time)[source]

The PWV concentration used by the propagation effect at a given time

Parameters:

time (TypeVar(FloatColl, Union[float, float16, float32, float64, float128, int, int8, int16, int32, int64, uint8, uint16, uint32, uint64, longlong, ulonglong], Collection[Union[float, float16, float32, float64, float128, int, int8, int16, int32, int64, uint8, uint16, uint32, uint64, longlong, ulonglong]], ndarray)) – Time to get the PWV concentration for

Return type:

TypeVar(FloatColl, Union[float, float16, float32, float64, float128, int, int8, int16, int32, int64, uint8, uint16, uint32, uint64, longlong, ulonglong], Collection[Union[float, float16, float32, float64, float128, int, int8, int16, int32, int64, uint8, uint16, uint32, uint64, longlong, ulonglong]], ndarray)

Returns:

An array of PWV values in units of mm

class snat_sim.models.pwv.SeasonalPWVTrans(time_format='mjd', transmission_res=5.0)[source]

Atmospheric propagation effect for a fixed PWV concentration per-season

The value of the winter, summer, spring and fall parameters reflect the total PWV concentration at zenith during those time periods. These values are scaled using the airmass calculated using the position of the observer (lat, lon, and alt parameters) and the object (ra and dec parameters).

Season names are defined for an individual in the southern hemisphere.

Effect Parameters:

ra: Target Right Ascension in degrees dec: Target Declination in degrees lat: Observer latitude in degrees (defaults to location of VRO) lon: Observer longitude in degrees (defaults to location of VRO) alt: Observer altitude in meters (defaults to height of VRO) winter: PWV concentration in the winter spring: PWV concentration in the spring summer: PWV concentration in the summer fall: PWV concentration in the fall

__init__(time_format='mjd', transmission_res=5.0)[source]

Time variable atmospheric transmission due to PWV that changes per season

Parameters:
  • time_format (str) – Astropy recognized time format used by the pwv_interpolator

  • transmission_res (float) – Reduce the underlying transmission model by binning to the given resolution

assumed_pwv(time)[source]

The PWV concentration used by the propagation effect at a given time

Parameters:

time (TypeVar(FloatColl, Union[float, float16, float32, float64, float128, int, int8, int16, int32, int64, uint8, uint16, uint32, uint64, longlong, ulonglong], Collection[Union[float, float16, float32, float64, float128, int, int8, int16, int32, int64, uint8, uint16, uint32, uint64, longlong, ulonglong]], ndarray)) – Time to get the PWV concentration for

Return type:

TypeVar(FloatColl, Union[float, float16, float32, float64, float128, int, int8, int16, int32, int64, uint8, uint16, uint32, uint64, longlong, ulonglong], Collection[Union[float, float16, float32, float64, float128, int, int8, int16, int32, int64, uint8, uint16, uint32, uint64, longlong, ulonglong]], ndarray)

Returns:

An array of PWV values in units of mm

static from_pwv_model(pwv_model)[source]

Create a new instance using the averaged per-season PWV from a PWV model

Supernova API

Modeling functionality for the spectro-photometric time-sampling of supernovae.

class snat_sim.models.supernova.SNModel(source, effects=None, effect_names=None, effect_frames=None)[source]

An observer-frame supernova model composed of a Source and zero or more effects

static from_sncosmo(model)[source]

Create an SNModel` instance from a sncosmo.Model instance

Parameters:

model (Model) – The sncosmo model to build from

Return type:

SNModel

Returns:

An SNModel object

simulate_lc(cadence, scatter=True, fixed_snr=None)[source]

Simulate a SN light-curve

If scatter is True, then simulated flux values include an added random component drawn from a normal distribution with a standard deviation equal to the error of the observation.

Parameters:
  • cadence (ObservedCadence) – Observational cadence to evaluate the light-curve with

  • scatter (bool) – Whether to add random noise to the flux values

  • fixed_snr (Optional[float]) – Optionally simulate the light-curve using a fixed signal to noise ratio

Return type:

LightCurve

Returns:

The simulated light-curve

fit_lc(data=None, vparam_names=(), bounds=None, method='minuit', guess_amplitude=True, guess_t0=True, guess_z=True, minsnr=5.0, modelcov=False, maxcall=10000, phase_range=None, wave_range=None)[source]

Fit model parameters to photometric data via a chi-squared minimization

Fitting behavior:
  • Parameters of the parent instance are not modified during the fit

  • If modelcov is enabled, the fit is performed multiple times until convergence.

  • The t0 parameter has a default fitting boundary such that the latest phase of the model lines up with the earliest data point and the earliest phase of the model lines up with the latest data point.

Parameters:
  • data – Table of photometric data.

  • vparam_names – Model parameters to vary in the fit.

  • bounds – Bounded range for each parameter. Keys should be parameter names, values are tuples.

  • guess_amplitude – Whether to guess the amplitude from the data.

  • guess_t0 – Whether to guess t0. Only has an effect when fitting t0.

  • guess_z – Whether to guess z (redshift). Only has an effect when fitting redshift.

  • minsnr – When guessing amplitude and t0, only use data with signal-to-noise ratio greater than this value.

  • method – Minimization method to use. Either “minuit” or “emcee”

  • modelcov – Include model covariance when calculating chisq. Default is False.

  • maxcall – Maximum number of chi-square iterations to evaluate when fitting.

  • phase_range – If given, discard data outside this range of phases.

  • wave_range – If given, discard data with bandpass effective wavelengths outside this range.

Returns:

The fit result and a copy of the model set to the fitted parameters

class snat_sim.models.supernova.SNFitResult(*args, **kwargs)[source]

Represents results from a SNModel being fit to a LightCurve

property param_names: List[str]

The names of the model parameters

Return type:

List[str]

property parameters: Series

The model parameters

Return type:

Series

property vparam_names: List[str]

List of parameter names varied in the fit

Return type:

List[str]

property vparameters: Series

The values of the varied parameters

Return type:

Series

property covariance: Optional[DataFrame]

The covariance matrix

Return type:

Optional[DataFrame]

salt_covariance_linear(x0_truth=None)[source]

The covariance matrix of apparent magnitude and salt2 parameters

Will raise an error if the x0, x1 and c parameters are not varied in the fit.

Parameters:

x0_truth (Optional[float]) – Optionally assert an alternative x0 value

Return type:

DataFrame

Returns:

The covariance matrix asd a pandas DataFrame

mu_variance_linear(alpha, beta)[source]

Calculate the variance in distance modulus

Determined using the covariance matrix of apparent magnitude and salt2 parameters (See the salt_covariance_linear method).

Parameters:
  • alpha (float) – The stretch correction factor

  • beta (float) – The color excess correction factor

Return type:

float

Returns:

The variance in mu

class snat_sim.models.supernova.VariablePropagationEffect[source]

Base class for propagation effects that vary with time

Similar to sncosmo.PropagationEffect class, but the propagate method accepts a time argument.

abstract propagate(wave, flux, time)[source]

Propagate the flux through the atmosphere

Parameters:
  • wave (ndarray) – An array of wavelength values

  • flux (ndarray) – An array of flux values

  • time (ndarray) – Array of time values

Return type:

ndarray

Returns:

An array of flux values after suffering propagation effects

Light-Curve API

Abstract representations of astronomical time-series data

class snat_sim.models.light_curve.ObservedCadence(obs_times, bands, skynoise=<property object>, zp=<property object>, zpsys=<property object>, gain=<property object>)[source]

The observational sampling of an astronomical light-curve

The zero-point, zero point system, and gain arguments can be a collection of values (one per obs_time value), or a single value to apply at all observation times.

Parameters:
  • obs_times (Collection[float]) – Array of observation times for the light-curve

  • bands (Collection[str]) – Array of bands for each observation

  • zp (TypeVar(FloatColl, Union[float, float16, float32, float64, float128, int, int8, int16, int32, int64, uint8, uint16, uint32, uint64, longlong, ulonglong], Collection[Union[float, float16, float32, float64, float128, int, int8, int16, int32, int64, uint8, uint16, uint32, uint64, longlong, ulonglong]], ndarray)) – The zero-point or an array of zero-points for each observation

  • zpsys (Union[str, Collection[str]]) – The zero-point system or an array of zero-point systems

  • gain (TypeVar(FloatColl, Union[float, float16, float32, float64, float128, int, int8, int16, int32, int64, uint8, uint16, uint32, uint64, longlong, ulonglong], Collection[Union[float, float16, float32, float64, float128, int, int8, int16, int32, int64, uint8, uint16, uint32, uint64, longlong, ulonglong]], ndarray)) – The simulated gain or an array of gain values

to_sncosmo()[source]

Return the observational cadence as an astropy.Table

The returned table of observations is formatted for use with with the sncosmo package.

Return type:

Table

Returns:

An astropy table representing the observational cadence in sncosmo format

__init__(obs_times, bands, skynoise=<property object>, zp=<property object>, zpsys=<property object>, gain=<property object>)
class snat_sim.models.light_curve.LightCurve(time, band, flux, fluxerr, zp, zpsys, phot_flag=None)[source]

Abstract representation of an astronomical light-curve.

__init__(time, band, flux, fluxerr, zp, zpsys, phot_flag=None)[source]

An astronomical light-curve

Parameters:
  • time (Collection[float]) – Time values for each observation in a numerical format (e.g. JD or MJD)

  • band (Collection[str]) – The band each observation was performed in

  • flux (Collection[float]) – The flux of each observation

  • fluxerr (Collection[float]) – The error in each flux value

  • zp (Collection[float]) – The zero-point of each observation

  • zpsys (Collection[str]) – The zero-point system of each observation

  • phot_flag (Optional[Collection[int]]) – Optional flag for each photometric observation

static from_sncosmo(data)[source]

Create a LightCurve instance from an astropy table in the SNCosmo format

Parameters:

data (Table) – A table in the sncosmo format

Return type:

LightCurve

Returns:

A LightCurve instance

to_astropy()[source]

Return the light-curve data as an astropy Table formatted for use with sncosmo

Return type:

Table

Returns:

An Table instance formatted for use with sncosmo

to_pandas()[source]

Return the light-curve data as a pandas DataFrame

Return type:

DataFrame

Returns:

A DataFrame instance with the light-curve data

copy()[source]

Return a copy of the instance

Return type:

LightCurve