skyllh.i3 package
Subpackages
Submodules
skyllh.i3.background_generation module
- class skyllh.i3.background_generation.FixedScrambledExpDataI3BkgGenMethod(data_scrambler, **kwargs)
Bases:
BackgroundGenerationMethod
This class implements the background event generation method for the IceCube detector using scrambled experimental data as background hypothesis with a fixed number of background events equal to the number of events in the dataset. This background generation method is the one used in SkyLab.
Creates a new background generation method instance to generate background events from scrambled experimental data with a fixed number of events equal to the number of events in the dataset.
- Parameters:
data_scrambler (instance of DataScrambler) – The DataScrambler instance to use to generate scrambled experimental data.
- property data_scrambler
The DataScrambler instance that implements the data scrambling.
- generate_events(rss, dataset, data, **kwargs)
Generates background events from the given data, by scrambling the experimental data. The number of events is equal to the size of the given dataset.
- Parameters:
rss (instance of RandomStateService) – The instance of RandomStateService that should be used to generate random numbers from. It is used to scramble the experimental data.
dataset (instance of Dataset) – The Dataset instance describing the dataset for which background events should get generated.
data (instance of DatasetData) – The DatasetData instance holding the data of the dataset for which background events should get generated.
- Returns:
n_bkg (int) – The number of generated background events.
bkg_events (instance of DataFieldRecordArray) – The instance of DataFieldRecordArray holding the generated background events.
skyllh.i3.backgroundpdf module
- class skyllh.i3.backgroundpdf.BackgroundI3SpatialPDF(data_sin_dec, data_weights, sin_dec_binning, spline_order_sin_dec, **kwargs)
Bases:
SpatialPDF
,UsesBinning
,IsBackgroundPDF
This is the base class for all IceCube specific spatial background PDF models. IceCube spatial background PDFs depend solely on the zenith angle, and hence, on the declination of the event.
The IceCube spatial background PDF is modeled as a 1d spline function in sin(declination).
Creates a new IceCube spatial background PDF object.
- Parameters:
data_sin_dec (1d ndarray) – The array holding the sin(dec) values of the events.
data_weights (1d ndarray) – The array holding the weight of each event used for histogramming.
sin_dec_binning (BinningDefinition) – The binning definition for the sin(declination) axis.
spline_order_sin_dec (int) – The order of the spline function for the logarithmic values of the spatial background PDF along the sin(dec) axis.
- property spline_order_sin_dec
The order (int) of the logarithmic spline function, that splines the background PDF, along the sin(dec) axis.
- add_events(events)
Add events to spatial background PDF object and recalculate logarithmic spline function.
- Parameters:
events (numpy record ndarray) –
The array holding the event data. The following data fields must exist:
- sin_decfloat
The sin(declination) value of the event.
- reset()
Reset the logarithmic spline to the original function, which was calculated when the object was initialized.
- initialize_for_new_trial(tdm, tl=None, **kwargs)
Pre-cumputes the probability density values when new trial data is available.
- get_pd(tdm, params_recarray=None, tl=None)
Calculates the spatial background probability on the sphere of each event.
- Parameters:
tdm (instance of TrialDataManager) –
The TrialDataManager instance holding the trial event data for which to calculate the PDF values. The following data fields must exist:
- 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 numpy ndarray) – The (N_events,)-shaped numpy ndarray holding the background probability density value for each event.
grads (dict) – The dictionary holding the gradients of the probability density w.r.t. each global fit parameter. The background PDF does not depend on any global fit parameter, hence, this is an empty dictionary.
- class skyllh.i3.backgroundpdf.DataBackgroundI3SpatialPDF(data_exp, sin_dec_binning, spline_order_sin_dec=2, **kwargs)
Bases:
BackgroundI3SpatialPDF
This is the IceCube spatial background PDF, which gets constructed from experimental data.
Constructs a new IceCube spatial background PDF from experimental data.
- Parameters:
data_exp (instance of DataFieldRecordArray) –
The instance of DataFieldRecordArray holding the experimental data. The following data fields must exist:
- sin_decfloat
The sin(declination) of the data event.
sin_dec_binning (instance of BinningDefinition) – The binning definition for the sin(declination).
spline_order_sin_dec (int) – The order of the spline function for the logarithmic values of the spatial background PDF along the sin(dec) axis. The default is 2.
- class skyllh.i3.backgroundpdf.MCBackgroundI3SpatialPDF(data_mc, physics_weight_field_names, sin_dec_binning, spline_order_sin_dec=2, **kwargs)
Bases:
BackgroundI3SpatialPDF
This is the IceCube spatial background PDF, which gets constructed from monte-carlo data.
Constructs a new IceCube spatial background PDF from monte-carlo data.
- Parameters:
data_mc (instance of DataFieldRecordArray) –
The array holding the monte-carlo data. The following data fields must exist:
- sin_decfloat
The sin(reco_dec) of the trial data event.
physics_weight_field_names (str | list of str) – The name or the list of names of the monte-carlo data fields, which should be used as event weights. If a list is given, the weight values of all the fields will be summed to construct the final event weight.
sin_dec_binning (instance of BinningDefinition) – The binning definition for the sin(reco_dec).
spline_order_sin_dec (int) – The order of the spline function for the logarithmic values of the spatial background PDF along the sin(reco_dec) axis. The default is 2.
- class skyllh.i3.backgroundpdf.DataBackgroundI3EnergyPDF(data_exp, log10_energy_binning, sin_dec_binning, smoothing_filter=None, **kwargs)
Bases:
SingleConditionalEnergyPDF
,IsBackgroundPDF
This is the IceCube energy background PDF, which gets constructed from experimental data. This class is derived from
core.pdf.SingleConditionalEnergyPDF
.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:
log10_energy_binning.name
floatThe logarithm of the reconstructed energy value of the data event.
sin_dec_binning.name
floatThe sine of the reconstructed declination of the data event.
log10_energy_binning (instance of BinningDefinition) – The binning definition for the binning in log10(E).
sin_dec_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.i3.backgroundpdf.MCBackgroundI3EnergyPDF(data_mc, physics_weight_field_names, log10_energy_binning, sin_dec_binning, smoothing_filter=None, **kwargs)
Bases:
SingleConditionalEnergyPDF
,IsBackgroundPDF
This is the IceCube energy background PDF, which gets constructed from monte-carlo data. This class is derived from SingleConditionalEnergyPDF.
Constructs a new IceCube energy background PDF from monte-carlo data.
- Parameters:
data_mc (instance of DataFieldRecordArray) –
The array holding the monte-carlo data. The following data fields must exist:
log10_energy_binning.name
floatThe logarithm of the reconstructed energy value of the data event.
sin_dec_binning.name
floatThe sine of the reconstructed declination of the data event.
- mcweight: float
The monte-carlo weight of the event.
physics_weight_field_names (str | list of str) – The name or the list of names of the monte-carlo data fields, which should be used as physics event weights. If a list is given, the weight values of all the fields will be summed to construct the final event physics weight.
log10_energy_binning (instance of BinningDefinition) – The binning definition for the binning in log10(E_reco). The name of this binning definition defines the field name in the MC and trial data.
sin_dec_binning (instance of BinningDefinition) – The binning definition for the sin(declination). The name of this binning definition defines the field name in the MC and trial data.
smoothing_filter (instance of SmoothingFilter | None) – The smoothing filter to use for smoothing the energy histogram. If None, no smoothing will be applied.
skyllh.i3.config module
This file defines IceCube specific global configuration.
skyllh.i3.dataset module
- class skyllh.i3.dataset.I3Dataset(livetime=None, grl_pathfilenames=None, **kwargs)
Bases:
Dataset
The I3Dataset class is an IceCube specific Dataset class that adds IceCube specific properties to the Dataset class. These additional properties are:
good-run-list (GRL)
Creates a new IceCube specific dataset, that also can hold a list of GRL data files.
- Parameters:
livetime (float | None) – The live-time of the dataset in days. It can be
None
, if good-run-list data files are provided.grl_pathfilenames (str | sequence of str) – The sequence of pathfilenames pointing to the good-run-list (GRL) data files.
- static get_combined_grl_pathfilenames(datasets)
Creates the combined list of grl pathfilenames of all the given datasets.
- Parameters:
datasets (sequence of I3Dataset) – The sequence of I3Dataset instances.
- Returns:
grl_pathfilenames (list) – The combined list of grl pathfilenames.
- property grl_pathfilename_list
The list of file names of the good-run-list (GRL) data files for this dataset. If a file name is given with a relative path, it will be relative to the
root_dir
property of this Dataset instance.
- property grl_abs_pathfilename_list
(read-only) The list of absolute path file names of the good-run-list data files.
- property grl_field_name_renaming_dict
The dictionary specifying the field names of the good-run-list data which need to get renamed just after loading the data. The dictionary keys are the old names and their values are the new names.
- property exists
(read-only) Flag if all the data files of this data set exists. It is
True
if all data files exist andFalse
otherwise.
- create_file_list()
Creates the list of files of this dataset. The file paths are relative to the dataset’s root directory.
- Returns:
file_list (list of str) – The list of files of this dataset.
- load_grl(efficiency_mode=None, tl=None)
Loads the good-run-list and returns a DataFieldRecordArray instance which should contain the following data fields:
- runint
The run number.
- startfloat
The MJD start time of the run.
- stopfloat
The MJD stop time of the run.
- livetimefloat
The livetime in days of the run.
- eventsint
The number of experimental events in the run.
- Parameters:
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 data loading procedure.
- Returns:
grl_data (instance of DataFieldRecordArray) – The DataFieldRecordArray instance holding the good-run-list information of the dataset.
- load_data(livetime=None, keep_fields=None, dtc_dict=None, dtc_except_fields=None, efficiency_mode=None, tl=None)
Loads the data, which is described by the dataset. If a good-run-list (GRL) is provided for this dataset, only experimental data will be selected which matches the GRL.
- Parameters:
livetime (instance of Livetime | float | None) – If not None, uses this livetime (if float livetime in days) as livetime for the DatasetData instance, otherwise uses the live time from the Dataset instance or, if available, the livetime from the good-run-list (GRL).
keep_fields (list of str | None) – The list of user-defined data fields that should get loaded and kept in addition to the analysis required data fields.
dtc_dict (dict | None) – This dictionary defines how data fields of specific data types (key) should get converted into other data types (value). This can be used to use less memory. If set to None, no data conversion is performed.
dtc_except_fields (str | sequence of str | None) – The sequence of field names whose data type should not get converted.
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 that should be used to time the data load operation.
- Returns:
data (instance of DatasetData) – A DatasetData instance holding the experimental and monte-carlo data of this data set.
- prepare_data(data, tl=None)
Prepares the data for IceCube by pre-calculating the following experimental data fields:
- sin_dec: float
The sin value of the declination coordinate.
and monte-carlo data fields:
- sin_true_dec: float
The sin value of the true declination coordinate.
- Parameters:
data (DatasetData instance) – The DatasetData instance holding the data as numpy record ndarray.
tl (TimeLord instance | None) – The TimeLord instance that should be used to time the data preparation.
- class skyllh.i3.dataset.I3DatasetData(data, grl_data)
Bases:
DatasetData
The class provides the container for the loaded experimental and monte-carlo data of a data set. It’s the IceCube specific class that also holds the good-run-list (GRL) data.
Constructs a new I3DatasetData instance.
- Parameters:
data (instance of DatasetData) – The instance of DatasetData holding the experimental and monte-carlo data.
grl_data (instance of DataFieldRecordArray | None) – The instance of DataFieldRecordArray holding the good-run-list data of the dataset. This can be None, if no GRL data is available.
- property grl
The DataFieldRecordArray instance holding the good-run-list (GRL) data of the IceCube data set. It is None, if there is no GRL data available for this IceCube data set.
skyllh.i3.detector_model module
This module defines the IceCube detector model class.
- class skyllh.i3.detector_model.IceCubeDetectorModel(**kwargs)
Bases:
DetectorModel
Creates a new DetectorModel instance.
- Parameters:
name (str) – The name of the detector model.
location (instance of astropy.coordinates.EarthLocation) – The Earth location of the detector.
skyllh.i3.detsigyield module
This module contains classes for IceCube specific detector signal yields, for a variation of source model and flux model combinations.
- class skyllh.i3.detsigyield.I3DetSigYield(param_names, detector_model, dataset, fluxmodel, livetime, sin_dec_binning, **kwargs)
Bases:
DetSigYield
Abstract base class for all IceCube specific detector signal yield classes. It assumes that sin(dec) binning is required for calculating the detector effective area and hence the detector signal yield.
Constructor of the IceCube specific detector signal yield base class.
- Parameters:
param_names (sequence of str) – The sequence of parameter names this detector signal yield depends on. These are either fixed or floating parameters.
detector_model (instance of DetectorModel) – The instance of DetectorModel defining the detector for this detector signal yield.
dataset (Dataset instance) – The Dataset instance holding the monte-carlo event data.
fluxmodel (FluxModel) – The flux model instance. Must be an instance of FluxModel.
livetime (float | Livetime instance) – The live-time.
sin_dec_binning (BinningDefinition instance) – The BinningDefinition instance defining the sin(dec) binning.
- property sin_dec_binning
The BinningDefinition instance defining the sin(dec) binning.
- class skyllh.i3.detsigyield.I3DetSigYieldBuilder(sin_dec_binning=None, **kwargs)
Bases:
DetSigYieldBuilder
Abstract base class for an IceCube specific detector signal yield builder class.
Constructor of the IceCube specific detector signal yield builder class.
- Parameters:
sin_dec_binning (BinningDefinition instance) – The instance of BinningDefinition defining the sin(dec) binning.
- property sin_dec_binning
The BinningDefinition instance for the sin(dec) binning that should be used for computing the sin(dec) dependency of the detector signal yield. If None, the binning is supposed to be taken from the Dataset’s binning definitions.
- get_sin_dec_binning(dataset)
Gets the sin(dec) binning definition either as setting from this detector signal yield builder itself, or from the given dataset.
- class skyllh.i3.detsigyield.PointLikeSourceI3DetSigYield(param_names, detector_model, dataset, fluxmodel, livetime, sin_dec_binning, **kwargs)
Bases:
I3DetSigYield
Abstract base class for all IceCube specific detector signal yield classes for point-like sources.
Constructor of the IceCube specific detector signal yield base class for point-like sources.
- Parameters:
param_names (sequence of str) – The sequence of parameter names this detector signal yield depends on. These are either fixed or floating parameters.
detector_model (instance of DetectorModel) – The instance of DetectorModel defining the detector for this detector signal yield.
dataset (instance of Dataset) – The instance of Dataset holding the monte-carlo event data.
fluxmodel (instance of FluxModel) – The instance of FluxModel defining the source’s flux model.
sin_dec_binning (instance of BinningDefinition) – The instance of BinningDefinition defining the sin(dec) binning.
- sources_to_recarray(sources)
Converts the sequence of PointLikeSource sources into a numpy record array holding the information of the sources needed for the detector signal yield calculation.
- Parameters:
sources (SourceModel | sequence of SourceModel) – The source model(s) containing the information of the source(s).
- Returns:
recarr (numpy record ndarray) – The generated (N_sources,)-shaped 1D numpy record ndarray holding the information for each source.
- class skyllh.i3.detsigyield.PointLikeSourceI3DetSigYieldBuilder(sin_dec_binning=None, **kwargs)
Bases:
I3DetSigYieldBuilder
Abstract base class for all IceCube specific detector signal yield builders for point-like sources. All IceCube detector signal yield builders require a sin(dec) binning definition for the effective area. By default it is taken from the binning definitions stored in the dataset, but a user-defined sin(dec) binning can be specified if needed.
Initializes a new detector signal yield builder object.
- Parameters:
sin_dec_binning (BinningDefinition | None) – The BinningDefinition instance defining the sin(dec) binning that should be used to compute the sin(dec) dependency of the detector effective area. If set to None, the binning will be taken from the Dataset binning definitions.
- class skyllh.i3.detsigyield.FixedFluxPointLikeSourceI3DetSigYield(param_names, detector_model, dataset, fluxmodel, livetime, sin_dec_binning, log_spl_sinDec, **kwargs)
Bases:
PointLikeSourceI3DetSigYield
The detector signal yield class for a point-source with a fixed flux.
Constructs an IceCube detector signal yield instance for a point-like source with a fixed flux.
- Parameters:
param_names (sequence of str) – The sequence of parameter names this detector signal yield depends on. These are either fixed or floating parameters.
detector_model (instance of DetectorModel) – The instance of DetectorModel defining the detector for this detector signal yield.
dataset (Dataset instance) – The instance of Dataset holding the monte-carlo data this detector signal yield is made for.
fluxmodel (FluxModel instance) – The instance of FluxModel with fixed parameters this detector signal yield is made for.
livetime (float | Livetime instance) – The livetime in days or an instance of Livetime.
sin_dec_binning (BinningDefinition instance) – The binning definition for sin(dec).
log_spl_sinDec (scipy.interpolate.InterpolatedUnivariateSpline) – The spline instance representing the log value of the detector signal yield as a function of sin(dec).
- property log_spl_sinDec
The
scipy.interpolate.InterpolatedUnivariateSpline
instance representing the spline for the log value of the detector signal yield as a function of sin(dec).
- __call__(src_recarray, src_params_recarray=None)
Retrieves the detector signal yield for the list of given sources.
- Parameters:
src_recarray (numpy record ndarray) – The numpy record ndarray with the field
dec
holding the declination of the source.src_params_recarray (None) – Unused interface argument, because this detector signal yield does not depend on any source parameters.
- Returns:
values (numpy 1d ndarray) – The array with the detector signal yield for each source.
grads (dict) – This detector signal yield does not depend on any parameters. So there are no gradients and the dictionary is empty.
- class skyllh.i3.detsigyield.FixedFluxPointLikeSourceI3DetSigYieldBuilder(sin_dec_binning=None, spline_order_sinDec=2, **kwargs)
Bases:
PointLikeSourceI3DetSigYieldBuilder
,IsParallelizable
This detector signal yield builder constructs a detector signal yield for a fixed flux model, assuming a point-like source. This means that the detector signal yield does not depend on any source parameters, hence it is only dependent on the detector effective area. It constructs a one-dimensional spline function in sin(dec), using a
scipy.interpolate.InterpolatedUnivariateSpline
.This detector signal yield builder works with all flux models.
It is tailored to the IceCube detector at the South Pole, where the effective area depends soley on the zenith angle, and hence on the declination, of the source.
Creates a new IceCube detector signal yield builder object for a fixed flux model. It requires a sinDec binning definition to compute the sin(dec) dependency of the detector effective area. The construct_detsigyield class method of this builder will create a spline function of a given order in logarithmic space of the effective area.
- Parameters:
sin_dec_binning (BinningDefinition | None) – The BinningDefinition instance which defines the sin(dec) binning. If set to None, the binning will be taken from the Dataset binning definitions.
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.
- property spline_order_sinDec
The order (int) of the logarithmic spline function, that splines the detector signal yield, along the sin(dec) axis.
- construct_detsigyields(detector_model, dataset, data, shgs, ppbar=None)
Constructs a set of FixedFluxPointLikeSourceI3DetSigYield instances, one for each provided fluxmodel.
- Parameters:
detector_model (instance of DetectorModel) – The instance of DetectorModel defining the detector for the detector signal yield.
dataset (instance of Dataset) – The instance of Dataset holding meta information about the data.
data (instance of DatasetData) –
The instance of DatasetData holding the monte-carlo event data. The numpy record ndarray holding the monte-carlo event data must contain the following data fields:
- ’true_dec’float
The true declination of the data event.
- ’true_energy’float
The true energy value of the data event.
- ’mcweight’float
The monte-carlo weight of the data event in the unit GeV cm^2 sr.
shgs (sequence of instance of SourceHypoGroup) – The sequence of instance of SourceHypoGroup specifying the source hypothesis groups (i.e. flux model) for which the detector signal yields should get constructed.
ppbar (instance of ProgressBar | None) – The optional instance of ProgressBar of the parent progress bar.
- Returns:
detsigyields (list of instance of FixedFluxPointLikeSourceI3DetSigYield) – The list of instance of FixedFluxPointLikeSourceI3DetSigYield providing the detector signal yield function for a point-like source with each of the given fixed flux models.
- construct_detsigyield(detector_model, dataset, data, shg, ppbar=None)
Constructs a detector signal yield log spline function for the given fixed flux model.
This method calls the
construct_detsigyiels()
method of this class.- Parameters:
detector_model (instance of DetectorModel) – The instance of DetectorModel defining the detector for the detector signal yield.
dataset (instance of Dataset) – The instance of Dataset holding meta information about the data.
data (instance of DatasetData) –
The instance of DatasetData holding the monte-carlo event data. The numpy record ndarray holding the monte-carlo event data must contain the following data fields:
- ’true_dec’float
The true declination of the data event.
- ’true_energy’float
The true energy value of the data event.
- ’mcweight’float
The monte-carlo weight of the data event in the unit GeV cm^2 sr.
shg (instance of SourceHypoGroup) – The instance of SourceHypoGroup (i.e. sources and flux model) for which the detector signal yield should get constructed.
ppbar (instance of ProgressBar | None) – The optional instance of ProgressBar of the parent progress bar.
- Returns:
detsigyield (instance of FixedFluxPointLikeSourceI3DetSigYield) – The instance of FixedFluxPointLikeSourceI3DetSigYield providing the detector signal yield function for a point-like source with a fixed flux.
- get_detsigyield_construction_factory()
Returns the factory callable for constructing a set of instance of FixedFluxPointLikeSourceI3DetSigYield.
- Returns:
factory (callable) – The factory callable for constructing a set of instance of FixedFluxPointLikeSourceI3DetSigYield.
- class skyllh.i3.detsigyield.SingleParamFluxPointLikeSourceI3DetSigYield(param_name, detector_model, dataset, fluxmodel, livetime, sin_dec_binning, log_spl_sinDec_param, **kwargs)
Bases:
PointLikeSourceI3DetSigYield
The detector signal yield class for a point-like source with a flux model that depends on a single parameter, tailored for the IceCube detector.
Constructs the detector signal yield instance.
- Parameters:
param_name (str) – The parameter name this detector signal yield depends on. This is either a fixed or floating parameter.
detector_model (instance of DetectorModel) – The instance of DetectorModel defining the detector for this detector signal yield.
dataset (instance of Dataset) – The instance of Dataset holding the monte-carlo event data.
fluxmodel (instance of FluxModel) – The instance of FluxModel for which the detector signal yield is constructed.
livetime (float | Livetime instance) – The live-time of the dataset.
sin_dec_binning (instance of BinningDefinition) – The instance of BinningDefinition defining the sin(dec) binning.
log_spl_sinDec_param (instance of scipy.interpolate.RectBivariateSpline) – The 2D spline in sin(dec) and the flux model’s parameter this detector signal yield depends on.
- property log_spl_sinDec_param
The
scipy.interpolate.RectBivariateSpline
instance representing the spline for the log value of the detector signal yield as a function of sin(dec) and the flux model’s parameter.
- __call__(src_recarray, src_params_recarray)
Retrieves the detector signal yield for the given list of sources and their flux parameters.
- Parameters:
src_recarray (numpy record ndarray) – The numpy record ndarray with the field
dec
holding the declination of the source.src_params_recarray ((N_sources,)-shaped numpy record ndarray) – The numpy record ndarray containing the parameter values of the sources. The parameter values can be different for the different sources. The record array needs to contain two fields for each source parameter, one named <name> with the source’s local parameter name holding the source’s local parameter value, and one named <name:gpidx> holding the global parameter index plus one for each source value. For values mapping to non-fit parameters, the index should be negative.
- Returns:
values (numpy (N_sources,)-shaped 1D ndarray) – The array with the detector signal yield for each source.
grads (dict) – The dictionary holding the gradient values for each global floating parameter. The key is the global floating parameter index and the value is the (N_sources,)-shaped numpy ndarray holding the gradient value dY_k/dp_s.
- class skyllh.i3.detsigyield.SingleParamFluxPointLikeSourceI3DetSigYieldBuilder(param_grid, sin_dec_binning=None, spline_order_sinDec=2, spline_order_param=2, ncpu=None, **kwargs)
Bases:
PointLikeSourceI3DetSigYieldBuilder
,IsParallelizable
This detector signal yield builder constructs a detector signal yield for a variable flux model with 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 soley on the zenith angle, and hence on the declination, of the source.
Creates a new IceCube detector signal yield builder instance for a flux model with a single parameter. It requires a sinDec binning definition to compute the sin(dec) dependency of the detector effective area, 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 the parameter values. The name of the parameter is defined via the name property of the ParameterGrid instance.
sin_dec_binning (instance of BinningDefinition | None) – The instance of BinningDefinition which defines the sin(dec) binning. If set to None, the sin(dec) binning will be taken from the dataset’s binning definitions.
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 flux model’s parameter axis. The default is 2.
ncpu (int | None) – The number of CPUs to utilize. If set to
None
, configuration setting will take place.
- property param_grid
The ParameterGrid instance for the parameter grid that should be used for computing the parameter dependency of the detector signal yield.
- property spline_order_sinDec
The order (int) of the logarithmic spline function, that splines the detector signal yield, along the sin(dec) axis.
- property spline_order_param
The order (int) of the logarithmic spline function, that splines the detector signal yield, along the parameter axis.
- construct_detsigyield(detector_model, 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:
detector_model (instance of DetectorModel) – The instance of DetectorModel defining the detector for the detector signal yield.
dataset (instance of Dataset) – The instance of Dataset holding the sin(dec) binning definition.
data (instance of DatasetData) –
The instance of DatasetData holding the monte-carlo event data. The numpy record array for the monte-carlo data of the dataset must contain the following data fields:
'true_dec'
floatThe true declination of the data event.
'mcweight'
floatThe monte-carlo weight of the data event in the unit GeV cm^2 sr.
'true_energy'
floatThe true energy value of the data event.
shg (instance of SourceHypoGroup) – The instance of SourceHypoGroup for which the detector signal yield should get constructed.
ppbar (instance of ProgressBar | None) – The instance of ProgressBar of the optional parent progress bar.
- Returns:
detsigyield (instance of SingleParamFluxPointLikeSourceI3DetSigYield) – The I3DetSigYield instance for a point-like source with a flux model of a single parameter.
skyllh.i3.livetime module
- class skyllh.i3.livetime.I3Livetime(*args, **kwargs)
Bases:
Livetime
The I3Livetime class provides the functionality to load a Livetime object from a good-run-list data file.
Creates a new instance of I3Livetime.
- classmethod from_grl_data(grl_data)
Creates an I3LiveTime instance from the given good-run-list (GRL) data.
- Parameters:
grl_data (instance of numpy structured ndarray.) –
The numpy structured ndarray of length N_runs holding the start end end times of the good runs. The following fields need to exist:
- startfloat
The MJD of the run start.
- endfloat
The MJD of the run stop.
- Returns:
livetime (instance of I3Livetime) – The created instance of I3Livetime for the provided GRL data.
- static from_grl_files(pathfilenames)
Loads an I3Livetime instance from the given good-run-list (GRL) data file. The data file needs to contain the following data fields:
- startfloat
The MJD of the run start.
- stopfloat
The MJD of the run stop.
- Parameters:
pathfilenames (str | list of str) – The list of fully qualified file names of the GRL data files.
- Returns:
livetime (instance of I3Livetime) – The created instance of I3Livetime for the provided GRL data.
skyllh.i3.scrambling module
- class skyllh.i3.scrambling.I3TimeScramblingMethod(time_generator, **kwargs)
Bases:
TimeScramblingMethod
The I3TimeScramblingMethod class provides a data scrambling method to perform time scrambling of the data, by drawing a MJD time from a given time generator.
Initializes a new I3 time scrambling instance.
- Parameters:
time_generator (instance of TimeGenerator) – The time generator that should be used to generate random MJD times.
- scramble(rss, dataset, data)
Draws a time from the time generator and calculates the right ascension coordinate from the azimuth angle according to the time. Sets the values of the
time
andra
keys of data.- Parameters:
rss (RandomStateService) – The random state service providing the random number generator (RNG).
dataset (instance of Dataset) – The instance of Dataset for which the data should get scrambled.
data (DataFieldRecordArray instance) – The DataFieldRecordArray instance containing the to be scrambled data.
- Returns:
data (numpy record ndarray) – The given numpy record ndarray holding the scrambled data.
- class skyllh.i3.scrambling.I3SeasonalVariationTimeScramblingMethod(data, **kwargs)
Bases:
DataScramblingMethod
The I3SeasonalVariationTimeScramblingMethod class provides a data scrambling method to perform data coordinate scrambling based on a generated time, which follows seasonal variations within the experimental data.
Initializes a new seasonal time scrambling instance.
- Parameters:
data (instance of I3DatasetData) – The instance of I3DatasetData holding the experimental data and good-run-list information.
- scramble(rss, dataset, data)
Scrambles the given data based on random MJD times, which are generated uniformly within the data runs, where the data runs are weighted based on their amount of events compared to the total events.
- Parameters:
rss (instance of RandomStateService) – The random state service providing the random number generator (RNG).
dataset (instance of Dataset) – The instance of Dataset for which the data should get scrambled.
data (instance of DataFieldRecordArray) – The DataFieldRecordArray instance containing the to be scrambled data.
- Returns:
data (instance of DataFieldRecordArray) – The given DataFieldRecordArray holding the scrambled data.
skyllh.i3.signal_generation module
- skyllh.i3.signal_generation.source_sin_dec_shift_linear(x, w, L, U)
Calculates the shift of the sine of the source declination, in order to allow the construction of the source sine declination band with sin(dec_src) +/- w. This shift function, S(x), is implemented as a line with the following points:
S(L) = w S((L+U)/2) = 0 S(U) = -w
- Parameters:
x (1D numpy ndarray) – The sine of the source declination for each source.
w (float) – The half size of the sin(dec)-window.
L (float) – The lower value of the allowed sin(dec) range.
U (float) – The upper value of the allowed sin(dec) range.
- Returns:
S (1D numpy ndarray) – The sin(dec) shift of the sin(dec) values of the given sources, such that
sin(dec_src) + S
is the new sin(dec) of the source, andsin(dec_src) + S +/- w
is always within the sin(dec) range [L, U].
- skyllh.i3.signal_generation.source_sin_dec_shift_cubic(x, w, L, U)
Calculates the shift of the sine of the source declination, in order to allow the construction of the source sine declination band with sin(dec_src) +/- w. This shift function, S(x), is implemented as a cubic function with the following points:
S(L) = w S((L+U)/2) = 0 S(U) = -w
- Parameters:
x (1D numpy ndarray) – The sine of the source declination for each source.
w (float) – The half size of the sin(dec)-window.
L (float) – The lower value of the allowed sin(dec) range.
U (float) – The upper value of the allowed sin(dec) range.
- Returns:
S (1D numpy ndarray) – The sin(dec) shift of the sin(dec) values of the given sources, such that
sin(dec_src) + S
is the new sin(dec) of the source, andsin(dec_src) + S +/- w
is always within the sin(dec) range [L, U].
- class skyllh.i3.signal_generation.PointLikeSourceI3SignalGenerationMethod(src_sin_dec_half_bandwidth=0.01745240643728351, src_sin_dec_shift_func=None, energy_range=None, src_batch_size=128, **kwargs)
Bases:
SignalGenerationMethod
This class provides a signal generation method for point-like sources seen in the IceCube detector.
Constructs a new signal generation method instance for a point-like source detected with IceCube.
- Parameters:
src_sin_dec_half_bandwidth (float) – The half-width of the sin(dec) band to take MC events from around a source. The default is sin(1deg), i.e. a 1deg half-bandwidth.
src_sin_dec_shift_func (callable | None) – The function that provides the source sin(dec) shift needed for constructing the source declination bands from where to draw monte-carlo events from. If set to None, the default function
source_sin_dec_shift_linear
will be used.energy_range (2-element tuple of float | None) – The energy range from which to take MC events into account for signal event generation. If set to None, the entire energy range [0, +inf] is used.
src_batch_size (int) – The source processing batch size used for the signal event flux calculation.
- property src_sin_dec_half_bandwidth
The half-width of the sin(dec) band to take MC events from around a source.
- property src_sin_dec_shift_func
The function that provides the source sin(dec) shift needed for constructing the source declination bands from where to draw monte-carlo events from.
- property src_batch_size
The source processing batch size used for the signal event flux calculation.
- calc_source_signal_mc_event_flux(data_mc, shg)
Calculates the signal flux of each given MC event for each source hypothesis of the given source hypothesis group.
- Parameters:
data_mc (numpy record ndarray) – The numpy record array holding the MC events of a dataset.
shg (SourceHypoGroup instance) – The source hypothesis group, which defines the list of sources, and their flux model.
- Returns:
ev_idx_arr (ndarray) – The (N_selected_signal_events,)-shaped 1D ndarray holding the index of the MC event.
shg_src_idx_arr (ndarray) – The (N_selected_signal_events,)-shaped 1D ndarray holding the index of the source within the given source hypothesis group for each signal candidate event.
flux_arr (ndarray) – The (N_selected_signal_events,)-shaped 1D ndarray holding the flux value of each signal candidate event.
- signal_event_post_sampling_processing(shg, shg_sig_events_meta, shg_sig_events)
Rotates the generated signal events to their source location for a given source hypothesis group.
- Parameters:
shg (SourceHypoGroup instance) – The source hypothesis group instance holding the sources and their locations.
shg_sig_events_meta (numpy record ndarray) –
The numpy record ndarray holding meta information about the generated signal events for the given source hypothesis group. The length of this array must be the same as shg_sig_events. It needs to contain the following data fields:
- ’shg_src_idx’: int
The source index within the source hypothesis group.
shg_sig_events (numpy record ndarray) – The numpy record ndarray holding the generated signal events for the given source hypothesis group and in the format of the original MC events.
- Returns:
shg_sig_events (numpy record ndarray) – The numpy record ndarray with the processed MC signal events.
skyllh.i3.signalpdf module
- class skyllh.i3.signalpdf.SignalI3EnergyPDFSet(cfg, data_mc, log10_energy_binning, sin_dec_binning, fluxmodel, param_grid_set, smoothing_filter=None, ncpu=None, ppbar=None, **kwargs)
Bases:
SignalSingleConditionalEnergyPDFSet
This is the signal energy PDF for IceCube. It creates a set of SingleConditionalEnergyPDF objects for a discrete set of energy signal parameters. Energy signal parameters influence the source’s flux model.
Creates a new IceCube energy signal PDF for a given flux model and a set of parameter grids for the flux model. It creates a set of SingleConditionalEnergyPDF objects for each signal parameter value permutation and stores it in an internal dictionary, where the hash of the parameters dictionary is the key.
- Parameters:
cfg (instance of Config) – The instance of Config holding the local configuration.
data_mc (instance of DataFieldRecordArray) –
The instance of DataFieldRecordArray holding the monte-carlo data. The following data fields must exist:
- true_energyfloat
The true energy value of the data event.
log10_energy_binning.name
floatThe base10-logarithm of the reconstructed energy value of the data event.
sin_dec_binning.name
floatThe declination of the data event.
- mcweightfloat
The monte-carlo weight value of the data events in unit GeV cm^2 sr.
log10_energy_binning (instance of BinningDefinition) – The binning definition for the reconstructed energy binning in log10(E_reco). The name of this binning definition defines the field name in the MC and trial data.
sin_dec_binning (instance of BinningDefinition) – The binning definition for the binning in sin(declination). The name of this binning definition defines the field name in the MC and trial data.
fluxmodel (instance of FluxModel) – The flux model to use to create the signal energy PDF.
param_grid_set (instance of ParameterGridSet |) – instance of ParameterGrid The set of parameter grids. A ParameterGrid instance for each energy parameter, for which an SingleConditionalEnergyPDF object needs to be created.
smoothing_filter (instance of SmoothingFilter | None) – The smoothing filter to use for smoothing the energy histogram. If
None
, no smoothing will be applied.ncpu (int | None) – The number of CPUs to use to create the different SingleConditionalEnergyPDF instances for the different parameter grid values. If set to
None
, the configured default number of CPUs will be used.ppbar (instance of ProgressBar | None) – The instance of ProgressBar of the optional parent progress bar.