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.
A summary of the available models is provided below:
|
Models atmospheric transmission due to PWV at a fixed resolution |
|
Model for interpolating the atmospheric water vapor at a given date and time |
|
An observer-frame supernova model composed of a Source and zero or more effects |
|
Representation of modeled stellar spectra from the Goettingen Spectral Library |
|
A rudimentary implementation of a reference star catalog |
|
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:
|
Propagation effect for the atmospheric absorption of light due to time static PWV |
|
Atmospheric propagation effect for a fixed PWV concentration per-season |
|
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:
|
The observational sampling of an astronomical light-curve |
|
Abstract representation of an astronomical light-curve. |
|
Represents results from a |
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.
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()
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
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))
Models for the spectro-photometric flux of stars and their combination into reference catalogs.
Representation of modeled stellar spectra from the Goettingen Spectral Library
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.
spectral_type (str
) – Spectral type (e.g., G2)
Return the spectral types available on disk
List
[str
]
A list fo spectral types
Retrieve pre-tabulated scale factors uses when calibrating flux values for PWV absorption
DataFrame
A DataFrame indexed by PWV with columns for flux
Return the reference star flux values
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
ndarray
The normalized flux at the given PWV value(s)
Return the normalized reference star flux values
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
ndarray
The normalized flux at the given PWV value(s)
A rudimentary implementation of a reference star catalog
Create a reference star catalog composed of the given spectral types
*spectral_types – Spectral types for the catalog (e.g., ‘G2’, ‘M5’, ‘K2’)
Return the average normalized reference star flux
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
ndarray
The normalized flux at the given PWV value(s)
Divide normalized reference flux from a light-curve
Recalibrate flux values using the average change in flux of a collection of reference stars.
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
A modified copy of the light_curve
argument
A reference star catalog that determines the time dependent PWV concentration along the line of sight from an underlying PWV model
Create a reference star catalog composed of the given spectral types and a PWV model
*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
Return the average normalized reference star flux
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’)
ndarray
The normalized flux at the given PWV value(s)
Divide normalized reference flux from a light-curve
Recalibrate flux values using the average change in flux of a collection of reference stars.
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’)
A modified copy of the light_curve
argument
Modeling functionality for Precipitable Water Vapor (PWV) and related observational effects.
Model for interpolating the atmospheric water vapor at a given date and time
Build a model for time variable PWV by interpolating from a given PWV time series
pwv_series (Series
) – PWV values with a datetime index
Construct a PWVModel
instance using data from a SuomiNet receiver
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
A PWV model that interpolates from data taken by the given receiver
Calculate the airmass through which a target is observed
Default latitude, longitude, and altitude are set to the Rubin Observatory.
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
Airmass in units of Sec(z)
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.
time – The time to interpolate PWV for
time_format – Astropy supported format of the time value (Default: ‘mjd’)
The PWV at zenith for the given time(s)
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
.
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’)
Union
[float
, array
]
The PWV concentration along the line of sight to the target
Models atmospheric transmission due to PWV at a fixed resolution
Instantiate a PWV transmission model at the given resolution
Transmission values are determined using the v1_transmission
model
from the pwv_kpno
package.
resolution (Optional
[float
]) – Resolution to bin the atmospheric model to
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.
pwv – Line of sight PWV to interpolate for
wave – Wavelengths to evaluate transmission (Defaults to samp_wave
attribute)
The interpolated transmission at the given wavelengths / resolution
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.
pwv: Atmospheric concentration of PWV along line of sight in mm
Time independent atmospheric transmission effect due to PWV
Setting the transmission_res
argument to None
results in the
highest available transmission model available.
transmission_res (float) – Reduce the underlying transmission model by binning to the given resolution
Resolution used when binning the underlying atmospheric transmission model
float
Base class for building time-variable PWV propagation effects
Time variable atmospheric transmission effect due to PWV
Setting the transmission_res
argument to None
results in the
highest available transmission model available.
pwv: Atmospheric concentration of PWV along line of sight in mm
transmission_res (float
) – Reduce the underlying transmission model by binning to the given resolution
The PWV concentration used by the propagation effect at a given time
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
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
)
An array of PWV values in units of mm
Propagate the flux through the atmosphere
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
ndarray
An array of flux values after suffering from PWV absorption
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).
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)
Time variable atmospheric transmission due to PWV
Setting the transmission_res
argument to None
results in the
highest available transmission model available.
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
The PWV concentration used by the propagation effect at a given time
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
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
)
An array of PWV values in units of mm
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.
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
Time variable atmospheric transmission due to PWV that changes per season
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
The PWV concentration used by the propagation effect at a given time
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
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
)
An array of PWV values in units of mm
Modeling functionality for the spectro-photometric time-sampling of supernovae.
An observer-frame supernova model composed of a Source and zero or more effects
Create an SNModel` instance from a sncosmo.Model
instance
model (Model
) – The sncosmo model to build from
An SNModel
object
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.
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
The simulated light-curve
Fit model parameters to photometric data via a chi-squared minimization
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.
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.
The fit result and a copy of the model set to the fitted parameters
Represents results from a SNModel
being fit to a LightCurve
The names of the model parameters
List
[str
]
The model parameters
Series
List of parameter names varied in the fit
List
[str
]
The values of the varied parameters
Series
The covariance matrix
Optional
[DataFrame
]
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.
x0_truth (Optional
[float
]) – Optionally assert an alternative x0 value
DataFrame
The covariance matrix asd a pandas DataFrame
Calculate the variance in distance modulus
Determined using the covariance matrix of apparent magnitude and
salt2 parameters (See the salt_covariance_linear
method).
alpha (float
) – The stretch correction factor
beta (float
) – The color excess correction factor
float
The variance in mu
Base class for propagation effects that vary with time
Similar to sncosmo.PropagationEffect
class, but the propagate
method accepts a time
argument.
Propagate the flux through the atmosphere
wave (ndarray
) – An array of wavelength values
flux (ndarray
) – An array of flux values
time (ndarray
) – Array of time values
ndarray
An array of flux values after suffering propagation effects
Abstract representations of astronomical time-series data
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.
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
Return the observational cadence as an astropy.Table
The returned table of observations is formatted for use with with
the sncosmo
package.
Table
An astropy table representing the observational cadence in sncosmo
format
Abstract representation of an astronomical light-curve.
An astronomical light-curve
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
Create a LightCurve
instance from an astropy table in the SNCosmo format
data (Table
) – A table in the sncosmo format
A LightCurve
instance
Return the light-curve data as an astropy Table
formatted for use with sncosmo
Table
An Table
instance formatted for use with sncosmo
Return the light-curve data as a pandas DataFrame
DataFrame
A DataFrame
instance with the light-curve data