skyllh.analyses.i3.publicdata_ps package
Submodules
skyllh.analyses.i3.publicdata_ps.aeff module
- skyllh.analyses.i3.publicdata_ps.aeff.load_effective_area_array(pathfilenames)
Loads the (nbins_decnu, nbins_log10enu)-shaped 2D effective area array from the given data file.
- Parameters:
pathfilename (str | list of str) – The file name of the data file.
- Returns:
aeff_decnu_log10enu ((nbins_decnu, nbins_log10enu)-shaped 2D ndarray) – The ndarray holding the effective area for each (dec_nu,log10(E_nu/GeV)) bin.
decnu_binedges_lower ((nbins_decnu,)-shaped ndarray) – The ndarray holding the lower bin edges of the dec_nu axis.
decnu_binedges_upper ((nbins_decnu,)-shaped ndarray) – The ndarray holding the upper bin edges of the dec_nu axis.
log10_enu_binedges_lower ((nbins_log10enu,)-shaped ndarray) – The ndarray holding the lower bin edges of the log10(E_nu/GeV) axis.
log10_enu_binedges_upper ((nbins_log10enu,)-shaped ndarray) – The ndarray holding the upper bin edges of the log10(E_nu/GeV) axis.
- class skyllh.analyses.i3.publicdata_ps.aeff.PDAeff(pathfilenames, src_dec=None, min_log10enu=None, max_log10enu=None, **kwargs)
Bases:
object
This class provides a representation of the effective area provided by the public data.
Creates an effective area instance by loading the effective area data from the given file.
- Parameters:
pathfilenames (str | list of str) – The path file names of the effective area data file(s) which should be used for this public data effective area instance.
src_dec (float | None) – The source declination in radians for which detection probabilities should get pre-calculated using the
get_detection_prob_for_decnu
method.min_log10enu (float | None) – The minimum log10(E_nu/GeV) value that should be used for calculating the detection probability. If None, the lowest available neutrino energy bin edge of the effective area is used.
max_log10enu (float | None) – The maximum log10(E_nu/GeV) value that should be used for calculating the detection probability. If None, the highest available neutrino energy bin edge of the effective area is used.
- property decnu_binedges
(read-only) The bin edges of the neutrino declination axis in radians.
- property sin_decnu_binedges
(read-only) The sin of the bin edges of the neutrino declination in radians.
- property decnu_bincenters
(read-only) The bin center values of the neutrino declination axis in radians.
- property n_decnu_bins
(read-only) The number of bins of the neutrino declination axis.
- property log10_enu_binedges
(read-only) The bin edges of the log10(E_nu/GeV) neutrino energy axis.
- property log10_enu_binedges_lower
(read-only) The lower binedges of the log10(E_nu/GeV) neutrino energy axis.
- property log10_enu_binedges_upper
(read-only) The upper binedges of the log10(E_nu/GeV) neutrino energy axis.
- property log10_enu_bincenters
(read-only) The bin center values of the log10(E_nu/GeV) neutrino energy axis.
- property n_log10_enu_bins
(read-only) The number of bins of the log10 neutrino energy axis.
- property aeff_decnu_log10enu
(read-only) The effective area in cm^2 as (n_decnu,n_log10enu)-shaped 2D numpy ndarray.
- create_sin_decnu_log10_enu_spline()
DEPRECATED! Creates a FctSpline2D object representing a 2D spline of the effective area in sin(dec_nu)-log10(E_nu/GeV)-space.
- Returns:
spl (FctSpline2D instance) – The FctSpline2D instance representing a spline in the sin(dec_nu)-log10(E_nu/GeV)-space.
- get_aeff_for_decnu(decnu)
Retrieves the effective area as function of log10_enu.
- Parameters:
decnu (float) – The true neutrino declination.
- Returns:
aeff ((n,)-shaped numpy ndarray) – The effective area in cm^2 for the given true neutrino declination as a function of log10 true neutrino energy.
- get_detection_prob_for_decnu(decnu, enu_min, enu_max, enu_range_min, enu_range_max)
Calculates the detection probability for given true neutrino energy ranges for a given neutrino declination.
- Parameters:
decnu (float) – The neutrino declination in radians.
enu_min (float | ndarray of float) – The minimum energy in GeV.
enu_max (float | ndarray of float) – The maximum energy in GeV.
enu_range_min (float) – The minimum energy in GeV of the entire energy range.
enu_range_max (float) – The maximum energy in GeV of the entire energy range.
- Returns:
det_prob (ndarray of float) – The neutrino energy detection probabilities for the given true enegry ranges.
skyllh.analyses.i3.publicdata_ps.backgroundpdf module
- class skyllh.analyses.i3.publicdata_ps.backgroundpdf.PDBackgroundI3EnergyPDF(data_logE, data_sinDec, data_mcweight, data_physicsweight, logE_binning, sinDec_binning, smoothing_filter, kde_smoothing=False, **kwargs)
Bases:
EnergyPDF
,IsBackgroundPDF
,UsesBinning
This is the base class for an IceCube specific energy background PDF for the public data.
IceCube energy PDFs depend solely on the energy and the zenith angle, and hence, on the declination of the event.
The IceCube energy PDF is modeled as a 1d histogram in energy, but for different sin(declination) bins, hence, stored as a 2d histogram.
Creates a new IceCube energy PDF object for the public data.
- Parameters:
data_logE (instance of ndarray) – The 1d ndarray holding the log10(E) values of the events.
data_sinDec (instance of ndarray) – The 1d ndarray holding the sin(dec) values of the events.
data_mcweight (instance of ndarray) – The 1d ndarray holding the monte-carlo weights of the events. The final data weight will be the product of data_mcweight and data_physicsweight.
data_physicsweight (instance of ndarray) – The 1d ndarray holding the physics weights of the events. The final data weight will be the product of data_mcweight and data_physicsweight.
logE_binning (instance of BinningDefinition) – The binning definition for the log10(E) axis.
sinDec_binning (instance of BinningDefinition) – The binning definition for the sin(declination) axis.
smoothing_filter (instance of SmoothingFilter | None) – The smoothing filter to use for smoothing the energy histogram. If None, no smoothing will be applied.
kde_smoothing (bool) – Apply a kde smoothing to the energy pdf for each bin in sin(dec). This is useful for signal injections, because it ensures that the background is not zero when injecting high energy events. Default: False.
- property hist_smoothing_method
The instance of HistSmoothingMethod defining the smoothing filter of the energy PDF histogram.
- property hist
(read-only) The 2D logE-sinDec histogram array.
- property hist_mask_mc_covered
(read-only) The boolean ndarray holding the mask of the 2D histogram bins for which there is monte-carlo coverage.
- property hist_mask_mc_covered_zero_physics
(read-only) The boolean ndarray holding the mask of the 2D histogram bins for which there is monte-carlo coverage but zero physics contribution.
- property hist_mask_mc_covered_with_physics
(read-only) The boolean ndarray holding the mask of the 2D histogram bins for which there is monte-carlo coverage and has physics contribution.
- initialize_for_new_trial(tdm, tl=None, **kwargs)
Pre-compute the probability densitiy values of the trial data, which has to be done only once for a particular trial data.
- assert_is_valid_for_trial_data(tdm, tl=None)
Checks if this energy PDF covers the entire value range of the trail data events.
- Parameters:
tdm (instance of TrialDataManager) –
The instance of TrialDataManager holding the trial data events. The following data fields need to exist:
- log_energyfloat
The base-10 logarithm of the reconstructed energy value.
- sin_decfloat
The sine of the declination value of the event.
tl (instance of TimeLord | None) – The optional instance of TimeLord to measure timing information.
- Raises:
ValueError – If parts of the trial data is outside the value range of this PDF.
- get_pd(tdm, params_recarray=None, tl=None)
Calculates the energy probability density (in 1/log10(E/GeV)) of each trial data event.
- Parameters:
tdm (instance of TrialDataManager) –
The TrialDataManager instance holding the data events for which the probability should be calculated for. The following data fields must exist:
- log_energyfloat
The base-10 logarithm of the energy value of the event.
- sin_decfloat
The sin(declination) value of the event.
params_recarray (None) – Unused interface parameter.
tl (instance of TimeLord | None) – The optional TimeLord instance that should be used to measure timing information.
- Returns:
pd (instance of ndarray) – The (N_selected_events,)-shaped numpy ndarray holding the energy probability density value for each trial data event.
grads (dict) – The dictionary holding the gradients of the probability density w.r.t. each global fit parameter. By definition this PDF does not depend on any fit parameter, hence, this dictionary is empty.
- class skyllh.analyses.i3.publicdata_ps.backgroundpdf.PDDataBackgroundI3EnergyPDF(data_exp, logE_binning, sinDec_binning, smoothing_filter=None, kde_smoothing=False, **kwargs)
Bases:
PDBackgroundI3EnergyPDF
This is the IceCube energy background PDF, which gets constructed from the experimental data of the public data.
Constructs a new IceCube energy background PDF from experimental data.
- Parameters:
data_exp (instance of DataFieldRecordArray) –
The array holding the experimental data. The following data fields must exist:
- log_energyfloat
The base-10 logarithm of the reconstructed energy value of the data event.
- sin_decfloat
The sine of the reconstructed declination of the data event.
logE_binning (instance of BinningDefinition) – The binning definition for the binning in log10(E).
sinDec_binning (instance of BinningDefinition) – The binning definition for the sin(declination).
smoothing_filter (instance of SmoothingFilter | None) – The smoothing filter to use for smoothing the energy histogram. If None, no smoothing will be applied.
- class skyllh.analyses.i3.publicdata_ps.backgroundpdf.PDMCBackgroundI3EnergyPDF(pdf_log10emu_sindecmu, log10emu_binning, sindecmu_binning, **kwargs)
Bases:
EnergyPDF
,IsBackgroundPDF
,UsesBinning
This class provides a background energy PDF constructed from the public data and a monte-carlo background flux model.
Constructs a new background energy PDF with the given PDF data and binning.
- Parameters:
pdf_log10emu_sindecmu (instance of numpy ndarray) – The (n_log10emu, n_sindecmu)-shaped 2D numpy ndarray holding the PDF values in unit 1/log10(E_mu/GeV). A copy of this data will be created and held within this class instance.
log10emu_binning (instance of BinningDefinition) – The binning definition for the binning in log10(E_mu/GeV).
sindecmu_binning (instance of BinningDefinition) – The binning definition for the binning in sin(dec_mu).
- assert_is_valid_for_trial_data(tdm, tl=None)
Checks if this energy PDF covers the entire value range of the trail data events.
- Parameters:
tdm (instance of TrialDataManager) –
The instance of TrialDataManager holding the trial data events. The following data fields need to exist:
- log_energyfloat
The base-10 logarithm of the reconstructed energy value.
- sin_decfloat
The sine of the declination value of the event.
tl (instance of TimeLord | None) – The optional instance of TimeLord to measure timing information.
- Raises:
ValueError – If parts of the trial data is outside the value range of this PDF.
- get_pd(tdm, params_recarray=None, tl=None)
Gets the probability density for the given trial data events.
- Parameters:
tdm (instance of TrialDataManager) –
The instance of TrialDataManager holding the trial data events. The following data fields need to exist:
- log_energyfloat
The base-10 logarithm of the reconstructed energy value.
- sin_decfloat
The sine of the declination value of the event.
params_recarray (None) – Unused interface argument.
tl (instance of TimeLord | None) – The optional instance of TimeLord that should be used to measure timing information.
- Returns:
pd (instance of ndarray) – The (N_selected_events,)-shaped numpy ndarray holding the probability density value for each event.
grads (dict) – The dictionary holding the gradients of the probability density w.r.t. each global fit parameter. By definition this PDF does not depend on any fit parameter, hence, this dictionary is empty.
skyllh.analyses.i3.publicdata_ps.bkg_flux module
- skyllh.analyses.i3.publicdata_ps.bkg_flux.get_dOmega(dec_min, dec_max)
Calculates the solid angle given two declination angles.
- Parameters:
dec_min (float | array of float) – The smaller declination angle.
dec_max (float | array of float) – The larger declination angle.
- Returns:
solidangle (float | array of float) – The solid angle corresponding to the two given declination angles.
- skyllh.analyses.i3.publicdata_ps.bkg_flux.southpole_zen2dec(zen)
Converts zenith angles at the South Pole to declination angles.
- Parameters:
zen ((n,)-shaped 1d numpy ndarray) – The numpy ndarray holding the zenith angle values in radians.
- Returns:
dec ((n,)-shaped 1d numpy ndarray) – The numpy ndarray holding the declination angle values in radians.
- skyllh.analyses.i3.publicdata_ps.bkg_flux.get_flux_atmo_decnu_log10enu(flux_pathfilename, log10_enu_max=9)
Constructs the atmospheric flux map function f_atmo(log10(E_nu/GeV),dec_nu) in unit 1/(GeV cm^2 sr s).
- Parameters:
flux_pathfilename (str) – The pathfilename of the file containing the MCEq fluxes.
log10_enu_max (float) – The log10(E/GeV) value of the maximum neutrino energy to be considered.
- Returns:
flux_atmo ((n_dec, n_e_grid)-shaped 2D numpy ndarray) – The numpy ndarray holding the the atmospheric neutrino flux function in unit 1/(GeV cm^2 sr s).
decnu_binedges ((n_decnu+1,)-shaped 1D numpy ndarray) – The numpy ndarray holding the dec_nu bin edges.
log10_enu_binedges ((n_enu+1,)-shaped 1D numpy ndarray) – The numpy ndarray holding the neutrino energy bin edges in log10.
- skyllh.analyses.i3.publicdata_ps.bkg_flux.get_flux_astro_decnu_log10enu(decnu_binedges, log10_enu_binedges)
Constructs the astrophysical neutrino flux function f_astro(log10(E_nu/GeV),dec_nu) in unit 1/(GeV cm^2 sr s).
It uses the best fit from the IceCube publication [1].
- Parameters:
decnu_binedges ((n_decnu+1,)-shaped 1D numpy ndarray) – The numpy ndarray holding the dec_nu bin edges.
log10_enu_binedges ((n_enu+1,)-shaped 1D numpy ndarray) – The numpy ndarray holding the log10 values of the neutrino energy bin edges in GeV.
- Returns:
f_astro ((n_decnu, n_log10enu)-shaped 2D numpy ndarray) – The numpy ndarray holding the astrophysical flux values in unit 1/(GeV cm^2 sr s).
References
- skyllh.analyses.i3.publicdata_ps.bkg_flux.convert_flux_bkg_to_pdf_bkg(f_bkg, decnu_binedges, log10_enu_binedges)
Converts the given background flux function f_bkg into a background flux PDF in unit 1/(log10(E/GeV) rad).
- Parameters:
f_bkg ((n_decnu, n_enu)-shaped 2D numpy ndarray) – The numpy ndarray holding the background flux values in unit 1/(GeV cm^2 s sr).
decnu_binedges ((n_decnu+1,)-shaped 1D numpy ndarray) – The numpy ndarray holding the dec_nu bin edges in radians.
log10_enu_binedges ((n_enu+1,)-shaped 1D numpy ndarray) – The numpy ndarray holding the log10 values of the neutrino energy bin edges in GeV.
- Returns:
p_bkg ((n_decnu, n_enu)-shaped 2D numpy ndarray) – The numpy ndarray holding the background flux pdf values.
- skyllh.analyses.i3.publicdata_ps.bkg_flux.get_pd_atmo_decnu_Enu(flux_pathfilename, log10_true_e_max=9)
Constructs the atmospheric neutrino PDF p_atmo(E_nu,dec_nu) in unit 1/(GeV rad).
- Parameters:
flux_pathfilename (str) – The pathfilename of the file containing the MCEq flux.
log10_true_e_max (float) – The log10(E/GeV) value of the maximum true energy to be considered.
- Returns:
pd_atmo ((n_dec, n_e_grid)-shaped 2D numpy ndarray) – The numpy ndarray holding the the atmospheric neutrino PDF in unit 1/(GeV rad).
decnu_binedges ((n_decnu+1,)-shaped 1D numpy ndarray) – The numpy ndarray holding the dec_nu bin edges.
log10_e_grid_edges ((n_e_grid+1,)-shaped 1D numpy ndarray) – The numpy ndarray holding the energy bin edges in log10.
- skyllh.analyses.i3.publicdata_ps.bkg_flux.get_pd_atmo_E_nu_sin_dec_nu(flux_pathfilename)
Constructs the atmospheric energy PDF p_atmo(E_nu|sin(dec_nu)) in unit 1/GeV.
- Parameters:
flux_pathfilename (str) – The pathfilename of the file containing the MCEq flux.
- Returns:
pd_atmo ((n_sin_dec, n_e_grid)-shaped 2D numpy ndarray) – The numpy ndarray holding the the atmospheric energy PDF in unit 1/GeV.
sin_dec_binedges (numpy ndarray) – The (n_sin_dec+1,)-shaped 1D numpy ndarray holding the sin(dec) bin edges.
log10_e_grid_edges (numpy ndarray) – The (n_e_grid+1,)-shaped 1D numpy ndarray holding the energy bin edges in log10.
- skyllh.analyses.i3.publicdata_ps.bkg_flux.get_pd_astro_E_nu_sin_dec_nu(sin_dec_binedges, log10_e_grid_edges)
Constructs the astrophysical energy PDF p_astro(E_nu|sin(dec_nu)) in unit 1/GeV. It uses the best fit from the IceCube publication [1].
- Parameters:
sin_dec_binedges ((n_sin_dec+1,)-shaped 1D numpy ndarray) – The numpy ndarray holding the sin(dec) bin edges.
log10_e_grid_edges ((n_e_grid+1,)-shaped 1D numpy ndarray) – The numpy ndarray holding the log10 values of the energy bin edges in GeV of the energy grid.
- Returns:
pd_astro ((n_sin_dec, n_e_grid)-shaped 2D numpy ndarray) – The numpy ndarray holding the energy probability density values p(E_nu|sin_dec_nu) in unit 1/GeV.
References
- skyllh.analyses.i3.publicdata_ps.bkg_flux.get_pd_bkg_E_nu_sin_dec_nu(pd_atmo, pd_astro, log10_e_grid_edges)
Constructs the total background flux probability density p_bkg(E_nu|sin(dec_nu)) in unit 1/GeV.
- Parameters:
pd_atmo ((n_sin_dec, n_e_grid)-shaped 2D numpy ndarray) – The numpy ndarray holding the probability density values p(E_nu|sin(dec_nu)) in 1/GeV of the atmospheric flux.
pd_astro ((n_sin_dec, n_e_grid)-shaped 2D numpy ndarray) – The numpy ndarray holding the probability density values p(E_nu|sin(dec_nu)) in 1/GeV of the astrophysical flux.
log10_e_grid_edges ((n_e_grid+1,)-shaped numpy ndarray) – The numpy ndarray holding the log10 values of the energy grid bin edges in GeV.
- Returns:
pd_bkg ((n_sin_dec, n_e_grid)-shaped 2D numpy ndarray) – The numpy ndarray holding total background probability density values p_bkg(E_nu|sin(dec_nu)) in unit 1/GeV.
skyllh.analyses.i3.publicdata_ps.detsigyield module
- class skyllh.analyses.i3.publicdata_ps.detsigyield.PDSingleParamFluxPointLikeSourceI3DetSigYieldBuilder(param_grid, spline_order_sinDec=2, spline_order_param=2, ncpu=None, **kwargs)
Bases:
SingleParamFluxPointLikeSourceI3DetSigYieldBuilder
This detector signal yield builder class constructs a detector signal yield instance for a variable flux model of a single parameter, assuming a point-like source. It constructs a two-dimensional spline function in sin(dec) and the parameter, using a
scipy.interpolate.RectBivariateSpline
. Hence, the detector signal yield can vary with the declination and the parameter of the flux model.It is tailored to the IceCube detector at the South Pole, where the effective area depends solely on the zenith angle, and hence on the declination, of the source.
It takes the effective area for the detector signal yield from the auxilary detector effective area data file given by the public data.
Creates a new IceCube detector signal yield builder instance for a flux model with a single parameter. It requires the effective area from the public data, and a parameter grid to compute the parameter dependency of the detector signal yield.
- Parameters:
param_grid (instance of ParameterGrid) – The instance of ParameterGrid which defines the grid of parameter values.
spline_order_sinDec (int) – The order of the spline function for the logarithmic values of the detector signal yield along the sin(dec) axis. The default is 2.
spline_order_param (int) – The order of the spline function for the logarithmic values of the detector signal yield along the parameter axis. The default is 2.
ncpu (int | None) – The number of CPUs to utilize. If set to
None
, global setting will take place.
- assert_types_of_construct_detsigyield_arguments(shgs, **kwargs)
Checks the correct types of the arguments for the
construct_detsigyield
method.
- construct_detsigyield(dataset, data, shg, ppbar=None)
Constructs a detector signal yield 2-dimensional log spline function for the given flux model with varying parameter values.
- Parameters:
dataset (instance of Dataset) – The Dataset instance holding the sin(dec) binning definition.
data (instance of DatasetData) – The instance of DatasetData holding the monte-carlo event data. This implementation loads the effective area from the provided public data and hence does not need monte-carlo data.
shg (instance of SourceHypoGroup) – The instance of SourceHypoGroup (i.e. sources and flux model) for which the detector signal yield should get constructed.
ppbar (ProgressBar instance | None) – The instance of ProgressBar of the optional parent progress bar.
- Returns:
detsigyield (instance of SingleParamFluxPointLikeSourceI3DetSigYield) – The DetSigYield instance for a point-like source with a flux model of a single parameter.
skyllh.analyses.i3.publicdata_ps.mcbkg_ps module
The mcbkg_ps analysis is a multi-dataset time-integrated single source analysis with a two-component likelihood function using a spacial and an energy event PDF. It initializes the background energy pdf using auxiliary fluxes and pdfs, which are generated by running scripts/mceq_atm_bkg.py script.
- skyllh.analyses.i3.publicdata_ps.mcbkg_ps.TXS_location()
- skyllh.analyses.i3.publicdata_ps.mcbkg_ps.create_analysis(datasets, source, refplflux_Phi0=1, refplflux_E0=1000.0, refplflux_gamma=2, ns_seed=100, ns_min=0, ns_max=1000.0, gamma_seed=3, gamma_min=1, gamma_max=5, minimizer_impl='LBFGS', cut_sindec=None, spl_smooth=None, cap_ratio=False, compress_data=False, keep_data_fields=None, evt_sel_delta_angle_deg=10, efficiency_mode=None, tl=None, ppbar=None, logger_name=None)
Creates the Analysis instance for this particular analysis.
- Parameters:
datasets (list of Dataset instances) – The list of Dataset instances, which should be used in the analysis.
source (PointLikeSource instance) – The PointLikeSource instance defining the point source position.
refplflux_Phi0 (float) – The flux normalization to use for the reference power law flux model.
refplflux_E0 (float) – The reference energy to use for the reference power law flux model.
refplflux_gamma (float) – The spectral index to use for the reference power law flux model.
ns_seed (float) – Value to seed the minimizer with for the ns fit.
ns_min (float) – Lower bound for ns fit.
ns_max (float) – Upper bound for ns fit.
gamma_seed (float | None) – Value to seed the minimizer with for the gamma fit. If set to None, the refplflux_gamma value will be set as gamma_seed.
gamma_min (float) – Lower bound for gamma fit.
gamma_max (float) – Upper bound for gamma fit.
minimizer_impl (str) – Minimizer implementation to be used. Supported options are
"LBFGS"
(L-BFG-S minimizer used from thescipy.optimize
module), or"minuit"
(Minuit minimizer used by theiminuit
module). Default: “LBFGS”.cut_sindec (list of float | None) – sin(dec) values at which the energy cut in the southern sky should start. If None, np.sin(np.radians([-2, 0, -3, 0, 0])) is used.
spl_smooth (list of float) – Smoothing parameters for the 1D spline for the energy cut. If None, [0., 0.005, 0.05, 0.2, 0.3] is used.
compress_data (bool) – Flag if the data should get converted from float64 into float32.
keep_data_fields (list of str | None) – List of additional data field names that should get kept when loading the data.
evt_sel_delta_angle_deg (float) – The delta angle in degrees for the event selection optimization methods.
efficiency_mode (str | None) –
The efficiency mode the data should get loaded with. Possible values are:
'memory'
:The data will be load in a memory efficient way. This will require more time, because all data records of a file will be loaded sequentially.
'time'
:The data will be loaded in a time efficient way. This will require more memory, because each data file gets loaded in memory at once.
The default value is
'time'
. If set toNone
, the default value will be used.tl (TimeLord instance | None) – The TimeLord instance to use to time the creation of the analysis.
ppbar (ProgressBar instance | None) – The instance of ProgressBar for the optional parent progress bar.
logger_name (str | None) – The name of the logger to be used. If set to
None
,__name__
will be used.
- Returns:
ana (instance of SingleSourceMultiDatasetLLHRatioAnalysis) – The Analysis instance for this analysis.
skyllh.analyses.i3.publicdata_ps.pdfratio module
- class skyllh.analyses.i3.publicdata_ps.pdfratio.PDSigSetOverBkgPDFRatio(sig_pdf_set, bkg_pdf, cap_ratio=False, **kwargs)
Bases:
SigSetOverBkgPDFRatio
Creates a PDFRatio instance for the public data. It takes a signal PDF set for different discrete gamma values.
- Parameters:
sig_pdf_set (instance of PDSignalEnergyPDFSet) – The PDSignalEnergyPDFSet instance holding the set of signal energy PDFs.
bkg_pdf (instance of PDDataBackgroundI3EnergyPDF) – The PDDataBackgroundI3EnergyPDF instance holding the background energy PDF.
cap_ratio (bool) – Switch whether the S/B PDF ratio should get capped where no background is available. Default is False.
- property cap_ratio
Boolean switch whether to cap the ratio where no background information is available (True) or use the smallest possible floating point number greater than zero as background pdf value (False).
- get_ratio(tdm, src_params_recarray, tl=None)
Calculates the PDF ratio values for all events and sources.
- Parameters:
tdm (instance of TrialDataManager) – The instance of TrialDataManager holding the trial data events for which the PDF ratio values should get calculated.
src_params_recarray (instance of numpy structured ndarray | None) – The (N_sources,)-shaped numpy structured ndarray holding the parameter names and values of the sources. See the documentation of the
skyllh.core.parameters.ParameterModelMapper.create_src_params_recarray()
method for more information.tl (instance of TimeLord | None) – The optional instance of TimeLord that should be used to measure timing information.
- Returns:
ratios (instance of ndarray) – The (N_values,)-shaped 1d numpy ndarray of float holding the PDF ratio value for each trial event and source.
- get_gradient(tdm, src_params_recarray, fitparam_id, tl=None)
Retrieves the PDF ratio gradient for the global fit parameter
fitparam_id
for each trial data event and source, given the given set of parameterssrc_params_recarray
for each source.- Parameters:
tdm (instance of TrialDataManager) – The instance of TrialDataManager holding the trial data events for which the PDF ratio gradient values should get calculated.
src_params_recarray (instance of numpy structured ndarray | None) – The (N_sources,)-shaped numpy structured ndarray holding the parameter names and values of the sources. See the documentation of the
skyllh.core.parameters.ParameterModelMapper.create_src_params_recarray()
method for more information.fitparam_id (int) – The name of the fit parameter for which the gradient should get calculated.
tl (instance of TimeLord | None) – The optional TimeLord instance that should be used to measure timing information.
- Returns:
grad (instance of ndarray) – The (N_values,)-shaped numpy ndarray holding the gradient values for all sources and trial events w.r.t. the given global fit parameter.
skyllh.analyses.i3.publicdata_ps.signal_generator module
- class skyllh.analyses.i3.publicdata_ps.signal_generator.PDDatasetSignalGenerator(shg_mgr, ds, ds_idx, energy_cut_spline=None, cut_sindec=None, **kwargs)
Bases:
SignalGenerator
This class implements a signal generator for a single public data dataset.
Creates a new instance of the signal generator for generating signal events from a specific public data dataset.
- Parameters:
shg_mgr (instance of SourceHypoGroupManager) – The instance of SourceHypoGroupManager defining the source hypothesis groups.
ds (instance of Dataset) – The instance of Dataset for which signal events should get generated.
ds_idx (int) – The index of the dataset.
energy_cut_spline (scipy.interpolate.UnivariateSpline) – A spline of E(sin_dec) that defines the declination dependent energy cut in the IceCube southern sky.
cut_sindec (float) – The sine of the declination to start applying the energy cut. The cut will be applied from this declination down.
- change_shg_mgr(shg_mgr)
Changes the source hypothesis group manager. This will recreate the internal source dependent data structures.
- static create_energy_filter_mask(events, spline, cut_sindec, logger)
Creates a mask for cutting all events below
cut_sindec
that have an energy smaller than the energy spline at their declination.- Parameters:
events (instance of DataFieldRecordArray) – The instance of DataFieldRecordArray holding the generated signal events.
spline (instance of scipy.interpolate.UnivariateSpline) – A spline of E(sin_dec) that defines the declination dependent energy cut in the IceCube southern sky.
cut_sindec (float) – The sine of the declination to start applying the energy cut. The cut will be applied from this declination down.
logger (instance of logging.Logger) – The Logger instance.
- Returns:
filter_mask (instance of numpy ndarray) – The (len(events),)-shaped numpy ndarray with the mask of the events to cut.
- generate_signal_events_for_source(rss, src_idx, n_events)
Generates
n_events
signal events for the given source location and flux model.- Parameters:
rss (instance of RandomStateService) – The instance of RandomStateService providing the random number generator state.
src_idx (int) – The index of the source.
n_events (int) – Number of signal events to be generated.
- Returns:
events (instance of DataFieldRecordArray) – The numpy record array holding the event data. It contains the following data fields:
’isvalid’
’log_true_energy’
’log_energy’
’dec’
’ra’
’ang_err’
- generate_signal_events(rss, mean, poisson=True, src_detsigyield_weights_service=None, **kwargs)
Generates
mean
number of signal events.- Parameters:
rss (instance of RandomStateService) – The instance of RandomStateService providing the random number generator state.
mean (int | float) – The mean number of signal events. If the
poisson
argument is set to True, the actual number of generated signal events will be drawn from a Poisson distribution with this given mean value of signal events.poisson (bool) – If set to True, the actual number of generated signal events will be drawn from a Poisson distribution with the given mean value of signal events. If set to False, the argument
mean
specifies the actual number of generated signal events.src_detsigyield_weights_service (instance of SrcDetSigYieldWeightsService) – The instance of SrcDetSigYieldWeightsService providing the weighting of the sources within the detector.
- Returns:
n_signal (int) – The number of generated signal events.
signal_events_dict (dict of DataFieldRecordArray) – The dictionary holding the DataFieldRecordArray instances with the generated signal events. Each key of this dictionary represents the dataset index for which the signal events have been generated.
- class skyllh.analyses.i3.publicdata_ps.signal_generator.TimeDependentPDDatasetSignalGenerator(shg_mgr, ds, ds_idx, livetime, time_flux_profile, energy_cut_spline=None, cut_sindec=None, **kwargs)
Bases:
PDDatasetSignalGenerator
This time dependent signal generator for a public PS dataset generates events using the
PDDatasetSignalGenerator
class. It then draws times for each event and adds them to the event array.- Parameters:
shg_mgr (instance of SourceHypoGroupManager) – The instance of SourceHypoGroupManager that defines the list of source hypothesis groups, i.e. the list of sources.
ds (instance of Dataset) – The instance of Dataset for which signal events should get generated.
ds_idx (int) – The index of the dataset.
livetime (instance of Livetime) – The instance of Livetime providing the live-time information of the dataset.
time_flux_profile (instance of TimeFluxProfile) –
The instance of TimeFluxProfile providing the time profile of the source(s).
Note:
At this time the some time profile will be used for all sources!
energy_cut_spline (scipy.interpolate.UnivariateSpline) – A spline of E(sin_dec) that defines the declination dependent energy cut in the IceCube southern sky.
cut_sindec (float) – The sine of the declination to start applying the energy cut. The cut will be applied from this declination down.
- property livetime
The instance of Livetime providing the live-time information of the dataset.
- generate_signal_events(rss, mean, poisson=True, src_detsigyield_weights_service=None, **kwargs)
Generates
mean
number of signal events with times.- Parameters:
rss (instance of RandomStateService) – The instance of RandomStateService providing the random number generator state.
mean (int | float) – The mean number of signal events. If the
poisson
argument is set to True, the actual number of generated signal events will be drawn from a Poisson distribution with this given mean value of signal events.poisson (bool) – If set to True, the actual number of generated signal events will be drawn from a Poisson distribution with the given mean value of signal events. If set to False, the argument
mean
specifies the actual number of generated signal events.src_detsigyield_weights_service (instance of SrcDetSigYieldWeightsService) – The instance of SrcDetSigYieldWeightsService providing the weighting of the sources within the detector.
- Returns:
n_signal (int) – The number of generated signal events.
signal_events_dict (dict of DataFieldRecordArray) – The dictionary holding the DataFieldRecordArray instances with the generated signal events. Each key of this dictionary represents the dataset index for which the signal events have been generated.
skyllh.analyses.i3.publicdata_ps.signalpdf module
- class skyllh.analyses.i3.publicdata_ps.signalpdf.PDSignalEnergyPDF(f_e_spl, **kwargs)
Bases:
PDF
,IsSignalPDF
This class provides a signal energy PDF for a spectrial index value.
Creates a new signal energy PDF instance for a particular spectral index value.
- Parameters:
f_e_spl (instance of FctSpline1D) – The instance of FctSpline1D representing the spline of the energy PDF.
- assert_is_valid_for_trial_data(tdm, tl=None)
This method is supposed to check if this PDF is valid for all the given trial data. This means, it needs to check if there is a PDF value for each trial data event that will be used in the likelihood evaluation. This is just a seatbelt. The method must raise a
ValueError
if the PDF is not valid for the given trial data.- Parameters:
tdm (instance of TrialDataManager) – The instance of TrialDataManager holding the trial data events.
tl (instance of TimeLord | None) – The optional instance of TimeLord for measuring timing information.
- Raises:
ValueError – If some of the trial data is outside the PDF’s value space.
- get_pd_by_log10_reco_e(log10_reco_e, tl=None)
Calculates the probability density for the given log10(E_reco/GeV) values using the spline representation of the PDF.
- Parameters:
log10_reco_e (instance of ndarray) – The (n_log10_reco_e,)-shaped numpy ndarray holding the log10(E_reco/GeV) values for which the energy PDF should get evaluated.
tl (instance of TimeLord | None) – The optional instance of TimeLord that should be used to measure timing information.
- Returns:
pd (instance of numpy ndarray) – The (n_log10_reco_e,)-shaped numpy ndarray with the probability density for each energy value.
- get_pd(tdm, params_recarray=None, tl=None)
Calculates the probability density for all given trial data events and sources.
- Parameters:
tdm (instance of TrialDataManager) –
The instance of TrialDataManager holding the trial data events for which the probability density should be looked up. The following data fields must be present:
- log_energyfloat
The base-10 logarithm of the reconstructed energy.
params_recarray (None) – Unused interface argument.
tl (instance of TimeLord | None) – The optional instance of TimeLord that should be used to measure timing information.
- Returns:
pd (instance of ndarray) – The (N_values,)-shaped numpy ndarray holding the probability density for each trial data event and source.
grads (dict) – The dictionary holding the gradient values for each global fit parameter. By definition this PDF does not depend on any fit parameters, hence, this is an empty dictionary.
- class skyllh.analyses.i3.publicdata_ps.signalpdf.PDSignalEnergyPDFSet(ds, src_dec, fluxmodel, param_grid_set, ncpu=None, ppbar=None, **kwargs)
Bases:
PDFSet
,IsSignalPDF
,PDF
,IsParallelizable
This class provides a signal energy PDF set using the public data. It creates a set of PDSignalEnergyPDF instances, one for each spectral index value on a grid.
Creates a new PDSignalEnergyPDFSet instance for the public data.
- Parameters:
ds (instance of Dataset) – The instance of Dataset that defines the dataset of the public data.
src_dec (float) – The declination of the source in radians.
fluxmodel (instance of FactorizedFluxModel) – The instance of FactorizedFluxModel that defines the source’s flux model.
param_grid_set (instance of ParameterGrid | instance of ParameterGridSet) – The parameter grid set defining the grids of the parameters this energy PDF set depends on.
ncpu (int | None) – The number of CPUs to utilize. Global setting will take place if not specified, i.e. set to None.
ppbar (instance of ProgressBar | None) – The instance of ProgressBar for the optional parent progress bar.
skyllh.analyses.i3.publicdata_ps.smearing_matrix module
- skyllh.analyses.i3.publicdata_ps.smearing_matrix.load_smearing_histogram(pathfilenames)
Loads the 5D smearing histogram from the given data file.
- Parameters:
pathfilenames (str | list of str) – The file name of the data file.
- Returns:
histogram (5d ndarray) – The 5d histogram array holding the probability values of the smearing matrix. The axes are (true_e, true_dec, reco_e, psi, ang_err).
true_e_bin_edges (1d ndarray) – The ndarray holding the bin edges of the true energy axis.
true_dec_bin_edges (1d ndarray) – The ndarray holding the bin edges of the true declination axis in radians.
reco_e_lower_edges (3d ndarray) – The 3d ndarray holding the lower bin edges of the reco energy axis. For each pair of true_e and true_dec different reco energy bin edges are provided. The shape is (n_true_e, n_true_dec, n_reco_e).
reco_e_upper_edges (3d ndarray) – The 3d ndarray holding the upper bin edges of the reco energy axis. For each pair of true_e and true_dec different reco energy bin edges are provided. The shape is (n_true_e, n_true_dec, n_reco_e).
psi_lower_edges (4d ndarray) – The 4d ndarray holding the lower bin edges of the psi axis in radians. The shape is (n_true_e, n_true_dec, n_reco_e, n_psi).
psi_upper_edges (4d ndarray) – The 4d ndarray holding the upper bin edges of the psi axis in radians. The shape is (n_true_e, n_true_dec, n_reco_e, n_psi).
ang_err_lower_edges (5d ndarray) – The 5d ndarray holding the lower bin edges of the angular error axis in radians. The shape is (n_true_e, n_true_dec, n_reco_e, n_psi, n_ang_err).
ang_err_upper_edges (5d ndarray) – The 5d ndarray holding the upper bin edges of the angular error axis in radians. The shape is (n_true_e, n_true_dec, n_reco_e, n_psi, n_ang_err).
- class skyllh.analyses.i3.publicdata_ps.smearing_matrix.PDSmearingMatrix(pathfilenames, **kwargs)
Bases:
object
This class is a helper class for dealing with the smearing matrix provided by the public data.
Creates a smearing matrix instance by loading the smearing matrix from the given file.
- property n_log10_true_e_bins
(read-only) The number of log10 true energy bins.
- property true_e_bin_edges
(read-only) The (n_true_e+1,)-shaped 1D numpy ndarray holding the bin edges of the true energy.
Depricated! Use log10_true_enu_binedges instead!
- property true_e_bin_centers
(read-only) The (n_true_e,)-shaped 1D numpy ndarray holding the bin center values of the true energy.
- property log10_true_enu_binedges
(read-only) The (n_log10_true_enu+1,)-shaped 1D numpy ndarray holding the bin edges of the log10 true neutrino energy.
- property n_true_dec_bins
(read-only) The number of true declination bins.
- property true_dec_bin_edges
(read-only) The (n_true_dec+1,)-shaped 1D numpy ndarray holding the bin edges of the true declination.
- property true_dec_bin_centers
(read-only) The (n_true_dec,)-shaped 1D ndarray holding the bin center values of the true declination.
- property log10_reco_e_binedges_lower
(read-only) The upper bin edges of the log10 reco energy axes.
- property log10_reco_e_binedges_upper
(read-only) The upper bin edges of the log10 reco energy axes.
- property min_log10_reco_e
(read-only) The minimum value of the reconstructed energy axis.
- property max_log10_reco_e
(read-only) The maximum value of the reconstructed energy axis.
- property min_log10_psi
(read-only) The minimum log10 value of the psi axis.
- property max_log10_psi
(read-only) The maximum log10 value of the psi axis.
- property pdf
(read-only) The probability-density-function P(E_reco,psi,ang_err|E_nu,dec_nu), which, by definition, is the histogram property divided by the 3D bin volumes for E_reco, psi, and ang_err.
- get_true_dec_idx(true_dec)
Returns the true declination index for the given true declination value.
- Parameters:
true_dec (float) – The true declination value in radians.
- Returns:
true_dec_idx (int) – The index of the declination bin for the given declination value.
- get_log10_true_e_idx(log10_true_e)
Returns the bin index for the given true log10 energy value.
- Parameters:
log10_true_e (float) – The log10 value of the true energy.
- Returns:
log10_true_e_idx (int) – The index of the true log10 energy bin for the given log10 true energy value.
- get_reco_e_idx(true_e_idx, true_dec_idx, reco_e)
Returns the bin index for the given reco energy value given the given true energy and true declination bin indices.
- Parameters:
true_e_idx (int) – The index of the true energy bin.
true_dec_idx (int) – The index of the true declination bin.
reco_e (float) – The reco energy value for which the bin index should get returned.
- Returns:
reco_e_idx (int | None) – The index of the reco energy bin the given reco energy value falls into. It returns None if the value is out of range.
- get_psi_idx(true_e_idx, true_dec_idx, reco_e_idx, psi)
Returns the bin index for the given psi value given the true energy, true declination and reco energy bin indices.
- Parameters:
true_e_idx (int) – The index of the true energy bin.
true_dec_idx (int) – The index of the true declination bin.
reco_e_idx (int) – The index of the reco energy bin.
psi (float) – The psi value in radians for which the bin index should get returned.
- Returns:
psi_idx (int | None) – The index of the psi bin the given psi value falls into. It returns None if the value is out of range.
- get_ang_err_idx(true_e_idx, true_dec_idx, reco_e_idx, psi_idx, ang_err)
Returns the bin index for the given angular error value given the true energy, true declination, reco energy, and psi bin indices.
- Parameters:
true_e_idx (int) – The index of the true energy bin.
true_dec_idx (int) – The index of the true declination bin.
reco_e_idx (int) – The index of the reco energy bin.
psi_idx (int) – The index of the psi bin.
ang_err (float) – The angular error value in radians for which the bin index should get returned.
- Returns:
ang_err_idx (int | None) – The index of the angular error bin the given angular error value falls into. It returns None if the value is out of range.
- get_true_log_e_range_with_valid_log_e_pdfs(dec_idx)
Determines the true log energy range for which log_e PDFs are available for the given declination bin.
- Parameters:
dec_idx (int) – The declination bin index.
- Returns:
min_log_true_e (float) – The minimum true log energy value.
max_log_true_e (float) – The maximum true log energy value.
- get_log_e_pdf(log_true_e_idx, dec_idx)
Retrieves the log_e PDF from the given true energy bin index and source bin index. Returns (None, None, None, None) if any of the bin indices are less then zero, or if the sum of all pdf bins is zero.
- Parameters:
log_true_e_idx (int) – The index of the true energy bin.
dec_idx (int) – The index of the declination bin.
- Returns:
pdf (1d ndarray) – The log_e pdf values.
lower_bin_edges (1d ndarray) – The lower bin edges of the energy pdf histogram.
upper_bin_edges (1d ndarray) – The upper bin edges of the energy pdf histogram.
bin_widths (1d ndarray) – The bin widths of the energy pdf histogram.
- get_psi_pdf(log_true_e_idx, dec_idx, log_e_idx)
Retrieves the psi PDF from the given true energy bin index, the source bin index, and the log_e bin index. Returns (None, None, None, None) if any of the bin indices are less then zero, or if the sum of all pdf bins is zero.
- Parameters:
log_true_e_idx (int) – The index of the true energy bin.
dec_idx (int) – The index of the declination bin.
log_e_idx (int) – The index of the log_e bin.
- Returns:
pdf (1d ndarray) – The psi pdf values.
lower_bin_edges (1d ndarray) – The lower bin edges of the psi pdf histogram.
upper_bin_edges (1d ndarray) – The upper bin edges of the psi pdf histogram.
bin_widths (1d ndarray) – The bin widths of the psi pdf histogram.
- get_ang_err_pdf(log_true_e_idx, dec_idx, log_e_idx, psi_idx)
Retrieves the angular error PDF from the given true energy bin index, the source bin index, the log_e bin index, and the psi bin index. Returns (None, None, None, None) if any of the bin indices are less then zero, or if the sum of all pdf bins is zero.
- Parameters:
log_true_e_idx (int) – The index of the true energy bin.
dec_idx (int) – The index of the declination bin.
log_e_idx (int) – The index of the log_e bin.
psi_idx (int) – The index of the psi bin.
- Returns:
pdf (1d ndarray) – The ang_err pdf values.
lower_bin_edges (1d ndarray) – The lower bin edges of the ang_err pdf histogram.
upper_bin_edges (1d ndarray) – The upper bin edges of the ang_err pdf histogram.
bin_widths (1d ndarray) – The bin widths of the ang_err pdf histogram.
- sample_log_e(rss, dec_idx, log_true_e_idxs)
Samples log energy values for the given source declination and true energy bins.
- Parameters:
rss (instance of RandomStateService) – The RandomStateService which should be used for drawing random numbers from.
dec_idx (int) – The index of the source declination bin.
log_true_e_idxs (1d ndarray of int) – The bin indices of the true energy bins.
- Returns:
log_e_idx (1d ndarray of int) – The bin indices of the log_e pdf corresponding to the sampled log_e values.
log_e (1d ndarray of float) – The sampled log_e values.
- sample_psi(rss, dec_idx, log_true_e_idxs, log_e_idxs)
Samples psi values for the given source declination, true energy bins, and log_e bins.
- Parameters:
rss (instance of RandomStateService) – The RandomStateService which should be used for drawing random numbers from.
dec_idx (int) – The index of the source declination bin.
log_true_e_idxs (1d ndarray of int) – The bin indices of the true energy bins.
log_e_idxs (1d ndarray of int) – The bin indices of the log_e bins.
- Returns:
psi_idx (1d ndarray of int) – The bin indices of the psi pdf corresponding to the sampled psi values.
psi (1d ndarray of float) – The sampled psi values in radians.
- sample_ang_err(rss, dec_idx, log_true_e_idxs, log_e_idxs, psi_idxs)
Samples ang_err values for the given source declination, true energy bins, log_e bins, and psi bins.
- Parameters:
rss (instance of RandomStateService) – The RandomStateService which should be used for drawing random numbers from.
dec_idx (int) – The index of the source declination bin.
log_true_e_idxs (1d ndarray of int) – The bin indices of the true energy bins.
log_e_idxs (1d ndarray of int) – The bin indices of the log_e bins.
psi_idxs (1d ndarray of int) – The bin indices of the psi bins.
- Returns:
ang_err_idx (1d ndarray of int) – The bin indices of the angular error pdf corresponding to the sampled angular error values.
ang_err (1d ndarray of float) – The sampled angular error values in radians.
skyllh.analyses.i3.publicdata_ps.time_dependent_ps module
Setup the time-dependent analysis. For now this only works on a single dataset.
- skyllh.analyses.i3.publicdata_ps.time_dependent_ps.create_signal_time_pdf(cfg, grl, gauss=None, box=None)
Creates the signal time PDF, either a gaussian or a box shaped PDF.
- Parameters:
cfg (instance of Config) – The instance of Config holding the local configuration.
grl (instance of numpy structured ndarray) – The structured numpy ndarray holding the good-run-list data.
gauss (dict | None) – None or dictionary with {“mu”: float, “sigma”: float}.
box (dict | None) – None or dictionary with {“start”: float, “stop”: float}.
- Returns:
pdf (instance of PDF) – The created time PDF instance.
- skyllh.analyses.i3.publicdata_ps.time_dependent_ps.change_signal_time_pdf_of_llhratio_function(ana, gauss=None, box=None)
Changes the signal time PDF of the log-likelihood ratio function.
- Parameters:
ana (instance of SingleSourceMultiDatasetLLHRatioAnalysis) – The analysis instance.
gauss (dict | None) – None or dictionary with {“mu”: float, “sigma”: float}.
box (dict | None) – None or dictionary with {“start”: float, “stop”: float}.
- skyllh.analyses.i3.publicdata_ps.time_dependent_ps.get_energy_spatial_signal_over_background(ana, fitparam_values, tl=None)
Returns the signal over background ratio for (spatial_signal * energy_signal) / (spatial_background * energy_background).
- Parameters:
ana (instance of SingleSourceMultiDatasetLLHRatioAnalysis) – The analysis instance.
fitparam_values (instance of ndarray) – The (N_fitparams,)-shaped numpy ndarray holding the values of the global fit parameters, e.g. ns and gamma.
tl (instance of TimeLord | None) – The optional instance of TimeLord for measuring timing behavior.
- Returns:
ratio (1d ndarray) – Product of spatial and energy signal over background pdfs.
- skyllh.analyses.i3.publicdata_ps.time_dependent_ps.change_fluxmodel_gamma(ana, gamma)
Sets the given gamma value to the flux model of the single source.
- Parameters:
ana (instance of SingleSourceMultiDatasetLLHRatioAnalysis) – The analysis that should be used.
gamma (float) – Spectral index for the flux model.
- skyllh.analyses.i3.publicdata_ps.time_dependent_ps.change_time_flux_profile_params(ana, params)
Changes the parameters of the source’s time flux profile.
- Parameters:
ana (instance of SingleSourceMultiDatasetLLHRatioAnalysis) – The analysis that should be used.
params (dict) – The dictionary with the parameter names and values to be set.
- skyllh.analyses.i3.publicdata_ps.time_dependent_ps.calculate_TS(ana, em_results, rss)
Calculate the best TS value from the expectation maximization gamma scan results.
- Parameters:
ana (instance of SingleSourceMultiDatasetLLHRatioAnalysis) – The analysis that should be used.
em_results (instance of structured ndarray) – The numpy structured ndarray holding the EM results (from the gamma scan).
rss (instance of RandomStateService) – The instance of RandomStateService that should be used to generate random numbers from.
- Returns:
max_TS (float) – The maximal TS value of all maximized time hypotheses.
best_em_result (instance of numpy structured ndarray) – The row of
em_results
that corresponds to the best fit.best_fitparam_values (instance of numpy ndarray) – The instance of numpy ndarray holding the fit parameter values of the overall best fit result.
- skyllh.analyses.i3.publicdata_ps.time_dependent_ps.run_gamma_scan_for_single_flare(ana, remove_time=None, gamma_min=1, gamma_max=5, n_gamma=51, ppbar=None)
Runs
em_fit
for different gamma values in the signal energy PDF.- Parameters:
ana (instance of SingleSourceMultiDatasetLLHRatioAnalysis) – The analysis that should be used.
remove_time (float) – Time information of event that should be removed.
gamma_min (float) – Lower bound for gamma scan.
gamma_max (float) – Upper bound for gamma scan.
n_gamma (int) – Number of steps for gamma scan.
ppbar (instance of ProgressBar | None) – The optional parent instance of ProgressBar.
- Returns:
em_results (instance of numpy structured ndarray) – The numpy structured ndarray with fields
- gammafloat
The spectral index value.
- mufloat
The determined mean value of the gauss curve.
- sigmafloat
The determined standard deviation of the gauss curve.
- ns_emfloat
The scaling factor of the flare.
- skyllh.analyses.i3.publicdata_ps.time_dependent_ps.unblind_single_flare(ana, remove_time=None)
Run EM for a single flare on unblinded data. Similar to the original analysis, remove the alert event.
- Parameters:
ana (instance of SingleSourceMultiDatasetLLHRatioAnalysis) – The analysis that should be used.
remove_time (float) – Time of the event that should be removed. In the case of the TXS analysis:
remove_time=TXS_0506_PLUS056_ALERT_TIME
.
- Returns:
max_TS (float) – The maximal TS value of all maximized time hypotheses.
best_em_result (instance of numpy structured ndarray) – The EM result from the gamma scan corresponding to the best fit.
best_fitparam_values (instance of numpy ndarray) – The instance of numpy ndarray holding the fit parameter values of the overall best fit result.
- skyllh.analyses.i3.publicdata_ps.time_dependent_ps.do_trial_with_em(ana, rss, mean_n_sig=0, gamma_src=2, gamma_min=1, gamma_max=5, n_gamma=21, gauss=None, box=None, tl=None, ppbar=None)
Performs a trial using the expectation maximization algorithm. It runs a gamma scan and does the EM for each gamma value.
- Parameters:
ana (instance of SingleSourceMultiDatasetLLHRatioAnalysis) – The analysis instance that should be used to perform the trial.
rss (instance of RandomStateService) – The instance of RandomStateService that should be used to generate random numbers.
mean_n_sig (float) – The mean number of signal events that should be generated.
gamma_src (float) – The spectral index of the source.
gamma_min (float) – Lower bound of the gamma scan.
gamma_max (float) – Upper bound of the gamma scan.
n_gamma (int) – Number of steps of the gamma scan.
gauss (dict | None) – Properties of the Gaussian time PDF. None or dictionary with {“mu”: float, “sigma”: float}.
box (dict | None) – Properties of the box time PDF. None or dictionary with {“start”: float, “stop”: float}.
tl (instance of TimeLord | None) – The optional instance of TimeLord to measure timing information.
ppbar (instance of ProgressBar | None) – The optional parent instance of ProgressBar.
- Returns:
trial (instance of structured ndarray) – The numpy structured ndarray of length 1 with the following fields:
- seednumpy.int64
The seed value used to generate the trial.
- mean_n_signumpy.float64
The mean number of signal events of the trial.
- n_signumpy.int64
The actual number of signal events in the trial.
- gamma_srcnumpy.float64
The spectral index of the source.
- mu_signumpy.float64
The mean value of the Gaussian time PDF of the source.
- sigma_signumpy.float64
The sigma value of the Gaussian time PDF of the source.
- start_signumpy.float64
The start time of the box time PDF of the source.
- stop_signumpy.float64
The stop time of the box time PDF of the source.
- tsnumpy.float64
The test-statistic value of the trial.
- ns_fitnumpy.float64
The fitted number of signal events.
- ns_emnumpy.float64
The scaling factor of the flare.
- gamma_fitnumpy.float64
The fitted spectral index of the trial.
- gamma_emnumpy.float64
The spectral index of the best EM trial.
- mu_fitnumpy.float64
The fitted mean value of the Gaussian time PDF.
- sigma_fitnumpy.float64
The fitted sigma value of the Gaussian time PDF.
- skyllh.analyses.i3.publicdata_ps.time_dependent_ps.do_trials_with_em(ana, n=1000, ncpu=None, seed=1, mean_n_sig=0, gamma_src=2, gamma_min=1, gamma_max=5, n_gamma=21, gauss=None, box=None, tl=None, ppbar=None)
Performs
n_trials
trials using the expectation maximization algorithm. For each trial it runs a gamma scan and does the EM for each gamma value.- Parameters:
ana (instance of SingleSourceMultiDatasetLLHRatioAnalysis) – The analysis instance that should be used to perform the trials.
n (int) – The number of trials to generate.
ncpu (int | None) – The number of CPUs to use to generate the trials. If set to
None
the configured default value will be used.mean_n_sig (float) – The mean number of signal events that should be generated.
gamma_src (float) – The spectral index of the source.
gamma_min (float) – Lower bound of the gamma scan.
gamma_max (float) – Upper bound of the gamma scan.
n_gamma (int) – Number of steps of the gamma scan.
seed (int) – The seed for the random number generator.
gauss (dict | None) – Properties of the Gaussian time PDF. None or dictionary with {“mu”: float, “sigma”: float}.
box (dict | None) – Properties of the box time PDF. None or dictionary with {“start”: float, “stop”: float}.
tl (instance of TimeLord | None) – The optional instance of TimeLord to measure timing information.
ppbar (instance of ProgressBar | None) – The optional parent instance of ProgressBar.
- Returns:
trials (instance of numpy structured ndarray) – The numpy structured ndarray of length
n_trials
with the results for each trial. The array has the following fields:- seednumpy.int64
The seed value used to generate the trial.
- mean_n_signumpy.float64
The mean number of signal events of the trial.
- n_signumpy.int64
The actual number of signal events in the trial.
- gamma_srcnumpy.float64
The spectral index of the source.
- mu_signumpy.float64
The mean value of the Gaussian time PDF of the source.
- sigma_signumpy.float64
The sigma value of the Gaussian time PDF of the source.
- start_signumpy.float64
The start time of the box time PDF of the source.
- stop_signumpy.float64
The stop time of the box time PDF of the source.
- tsnumpy.float64
The test-statistic value of the trial.
- ns_fitnumpy.float64
The fitted number of signal events.
- ns_emnumpy.float64
The scaling factor of the flare.
- gamma_fitnumpy.float64
The fitted spectral index of the trial.
- gamma_emnumpy.float64
The spectral index of the best EM trial.
- mu_fitnumpy.float64
The fitted mean value of the Gaussian time PDF.
- sigma_fitnumpy.float64
The fitted sigma value of the Gaussian time PDF.
- skyllh.analyses.i3.publicdata_ps.time_dependent_ps.create_analysis(cfg, datasets, source, box=None, gauss=None, refplflux_Phi0=1, refplflux_E0=1000.0, refplflux_gamma=2.0, ns_seed=10.0, ns_min=0.0, ns_max=1000.0, gamma_seed=3.0, gamma_min=1.0, gamma_max=5.0, kde_smoothing=False, minimizer_impl='LBFGS', cut_sindec=None, spl_smooth=None, cap_ratio=False, compress_data=False, keep_data_fields=None, evt_sel_delta_angle_deg=10, construct_bkg_generator=True, construct_sig_generator=True, tl=None, ppbar=None, logger_name=None)
Creates the Analysis instance for this particular analysis.
- Parameters:
cfg (instance of Config) – The instance of Config holding the local configuration.
datasets (list of Dataset instances) – The list of Dataset instances, which should be used in the analysis.
source (PointLikeSource instance) – The PointLikeSource instance defining the point source position.
box (None or dictionary with start, end) – None if no box shaped time pdf, else dictionary of the format
{'start': float, 'stop': float}
.gauss (None or dictionary with mu, sigma) – None if no gaussian time pdf, else dictionary of the format
{'mu': float, 'sigma': float}
.refplflux_Phi0 (float) – The flux normalization to use for the reference power law flux model.
refplflux_E0 (float) – The reference energy to use for the reference power law flux model.
refplflux_gamma (float) – The spectral index to use for the reference power law flux model.
ns_seed (float) – Value to seed the minimizer with for the ns fit.
ns_min (float) – Lower bound for ns fit.
ns_max (float) – Upper bound for ns fit.
gamma_seed (float | None) – Value to seed the minimizer with for the gamma fit. If set to None, the refplflux_gamma value will be set as gamma_seed.
gamma_min (float) – Lower bound for gamma fit.
gamma_max (float) – Upper bound for gamma fit.
kde_smoothing (bool) – Apply a KDE-based smoothing to the data-driven background pdf. Default: False.
minimizer_impl (str | "LBFGS") – Minimizer implementation to be used. Supported options are “LBFGS” (L-BFG-S minimizer used from the
scipy.optimize
module), or “minuit” (Minuit minimizer used by theiminuit
module). Default: “LBFGS”.cut_sindec (list of float | None) – sin(dec) values at which the energy cut in the southern sky should start. If None, np.sin(np.radians([-2, 0, -3, 0, 0])) is used.
spl_smooth (list of float) – Smoothing parameters for the 1D spline for the energy cut. If None, [0., 0.005, 0.05, 0.2, 0.3] is used.
cap_ratio (bool) – If set to True, the energy PDF ratio will be capped to a finite value where no background energy PDF information is available. This will ensure that an energy PDF ratio is available for high energies where no background is available from the experimental data. If kde_smoothing is set to True, cap_ratio should be set to False! Default is False.
compress_data (bool) – Flag if the data should get converted from float64 into float32.
keep_data_fields (list of str | None) – List of additional data field names that should get kept when loading the data.
evt_sel_delta_angle_deg (float) – The delta angle in degrees for the event selection optimization methods.
construct_bkg_generator (bool) – Flag if the background generator should be constructed (
True
) or not (False
).construct_sig_generator (bool) – Flag if the signal generator should be constructed (
True
) or not (False
).tl (TimeLord instance | None) – The TimeLord instance to use to time the creation of the analysis.
ppbar (ProgressBar instance | None) – The instance of ProgressBar for the optional parent progress bar.
logger_name (str | None) – The name of the logger to be used. If set to
None
,__name__
will be used.
- Returns:
ana (instance of SingleSourceMultiDatasetLLHRatioAnalysis) – The Analysis instance for this analysis.
skyllh.analyses.i3.publicdata_ps.time_integrated_ps module
The time_integrated_ps analysis is a multi-dataset time-integrated single source analysis with a two-component likelihood function using a spacial and an energy event PDF.
- skyllh.analyses.i3.publicdata_ps.time_integrated_ps.create_analysis(cfg, datasets, source, refplflux_Phi0=1, refplflux_E0=1000.0, refplflux_gamma=2.0, ns_seed=10.0, ns_min=0.0, ns_max=1000.0, gamma_seed=3.0, gamma_min=1.0, gamma_max=5.0, kde_smoothing=False, minimizer_impl='LBFGS', cut_sindec=None, spl_smooth=None, cap_ratio=False, compress_data=False, keep_data_fields=None, evt_sel_delta_angle_deg=10, construct_sig_generator=True, tl=None, ppbar=None, logger_name=None)
Creates the Analysis instance for this particular analysis.
- Parameters:
cfg (instance of Config) – The instance of Config holding the local configuration.
datasets (list of Dataset instances) – The list of Dataset instances, which should be used in the analysis.
source (PointLikeSource instance) – The PointLikeSource instance defining the point source position.
refplflux_Phi0 (float) – The flux normalization to use for the reference power law flux model.
refplflux_E0 (float) – The reference energy to use for the reference power law flux model.
refplflux_gamma (float) – The spectral index to use for the reference power law flux model.
ns_seed (float) – Value to seed the minimizer with for the ns fit.
ns_min (float) – Lower bound for ns fit.
ns_max (float) – Upper bound for ns fit.
gamma_seed (float | None) – Value to seed the minimizer with for the gamma fit. If set to None, the refplflux_gamma value will be set as gamma_seed.
gamma_min (float) – Lower bound for gamma fit.
gamma_max (float) – Upper bound for gamma fit.
kde_smoothing (bool) – Apply a KDE-based smoothing to the data-driven background pdf. Default: False.
minimizer_impl (str) – Minimizer implementation to be used. Supported options are
"LBFGS"
(L-BFG-S minimizer used from thescipy.optimize
module), or"minuit"
(Minuit minimizer used by theiminuit
module). Default: “LBFGS”.cut_sindec (list of float | None) – sin(dec) values at which the energy cut in the southern sky should start. If None, np.sin(np.radians([-2, 0, -3, 0, 0])) is used.
spl_smooth (list of float) – Smoothing parameters for the 1D spline for the energy cut. If None, [0., 0.005, 0.05, 0.2, 0.3] is used.
cap_ratio (bool) – If set to True, the energy PDF ratio will be capped to a finite value where no background energy PDF information is available. This will ensure that an energy PDF ratio is available for high energies where no background is available from the experimental data. If kde_smoothing is set to True, cap_ratio should be set to False! Default is False.
compress_data (bool) – Flag if the data should get converted from float64 into float32.
keep_data_fields (list of str | None) – List of additional data field names that should get kept when loading the data.
evt_sel_delta_angle_deg (float) – The delta angle in degrees for the event selection optimization methods.
construct_sig_generator (bool) – Flag if the signal generator should be constructed (
True
) or not (False
).tl (TimeLord instance | None) – The TimeLord instance to use to time the creation of the analysis.
ppbar (ProgressBar instance | None) – The instance of ProgressBar for the optional parent progress bar.
logger_name (str | None) – The name of the logger to be used. If set to
None
,__name__
will be used.
- Returns:
ana (instance of SingleSourceMultiDatasetLLHRatioAnalysis) – The Analysis instance for this analysis.
skyllh.analyses.i3.publicdata_ps.utils module
- class skyllh.analyses.i3.publicdata_ps.utils.FctSpline1D(f, x_binedges, norm=False, **kwargs)
Bases:
object
Class to represent a 1D function spline using the PchipInterpolator class from scipy.
The evaluate the spline, use the
__call__
method.Creates a new 1D function spline using the PchipInterpolator class from scipy.
- Parameters:
f ((n_x,)-shaped 1D numpy ndarray) – The numpy ndarray holding the function values at the bin centers.
x_binedges ((n_x+1,)-shaped 1D numpy ndarray) – The numpy ndarray holding the bin edges of the x-axis.
norm (bool) – Whether to precalculate and save normalization internally.
- __call__(x, oor_value=0)
Evaluates the spline at the given x values. For x-values outside the spline’s range, the oor_value is returned.
- Parameters:
x ((n_x,)-shaped 1D numpy ndarray) – The numpy ndarray holding the x values at which the spline should get evaluated.
oor_value (float) – The value for out-of-range (oor) coordinates.
- Returns:
f ((n_x,)-shaped 1D numpy ndarray) – The numpy ndarray holding the evaluated values of the spline.
- evaluate(*args, **kwargs)
Alias for the __call__ method.
- class skyllh.analyses.i3.publicdata_ps.utils.FctSpline2D(f, x_binedges, y_binedges, **kwargs)
Bases:
object
Class to represent a 2D function spline using the RectBivariateSpline class from scipy.
The spline is constructed in the log10 space of the function value to ensure a smooth spline.
The evaluate the spline, use the
__call__
method.Creates a new 2D function spline using the RectBivariateSpline class from scipy.
- Parameters:
f ((n_x, n_y)-shaped 2D numpy ndarray) – he numpy ndarray holding the function values at the bin centers.
x_binedges ((n_x+1,)-shaped 1D numpy ndarray) – The numpy ndarray holding the bin edges of the x-axis.
y_binedges ((n_y+1,)-shaped 1D numpy ndarray) – The numpy ndarray holding the bin edges of the y-axis.
- __call__(x, y, oor_value=0)
Evaluates the spline at the given coordinates. For coordinates outside the spline’s range, the oor_value is returned.
- Parameters:
x ((n_x,)-shaped 1D numpy ndarray) – The numpy ndarray holding the x values at which the spline should get evaluated.
y ((n_y,)-shaped 1D numpy ndarray) – The numpy ndarray holding the y values at which the spline should get evaluated.
oor_value (float) – The value for out-of-range (oor) coordinates.
- Returns:
f ((n_x, n_y)-shaped 2D numpy ndarray) – The numpy ndarray holding the evaluated values of the spline.
- skyllh.analyses.i3.publicdata_ps.utils.clip_grl_start_times(grl_data)
Make sure that the start time of a run is not smaller than the stop time of the previous run.
- Parameters:
grl_data (instance of numpy structured ndarray) –
The numpy structured ndarray of length N_runs, with the following fields:
- startfloat
The start time of the run.
- stopfloat
The stop time of the run.
- skyllh.analyses.i3.publicdata_ps.utils.psi_to_dec_and_ra(rss, src_dec, src_ra, psi)
Generates random declinations and right-ascension coordinates for the given source location and opening angle psi.
- Parameters:
rss (instance of RandomStateService) – The instance of RandomStateService to use for drawing random numbers.
src_dec (float) – The declination of the source in radians.
src_ra (float) – The right-ascension of the source in radians.
psi (1d ndarray of float) – The opening-angle values in radians.
- Returns:
dec (1d ndarray of float) – The declination values.
ra (1d ndarray of float) – The right-ascension values.
- skyllh.analyses.i3.publicdata_ps.utils.create_energy_cut_spline(ds, exp_data, spl_smooth)
Create the spline for the declination-dependent energy cut that the signal generator needs for injection in the southern sky Some special conditions are needed for IC79 and IC86_I, because their experimental dataset shows events that should probably have been cut by the IceCube selection.