skyllh.core.utils package
Submodules
skyllh.core.utils.analysis module
- skyllh.core.utils.analysis.pointlikesource_to_data_field_array(tdm, shg_mgr, pmm)
Function to transform a list of PointLikeSource sources into a numpy record ndarray. The resulting numpy record ndarray contains the following fields:
- ra: float
The right-ascention of the point-like source.
- dec: float
The declination of the point-like source.
- weight: float
The weight of the point-like source.
- Parameters:
tdm (instance of TrialDataManager) – The TrialDataManager instance.
shg_mgr (instance of SourceHypoGroupManager) – The instance of SourceHypoGroupManager that defines the sources.
pmm (instance of ParameterModelMapper) – The instance of ParameterModelMapper that defines the mapping of global parameters to local model parameters.
- Returns:
arr ((N_sources,)-shaped numpy record ndarray) – The numpy record ndarray holding the source parameters.
- skyllh.core.utils.analysis.calculate_pval_from_trials(ts_vals, ts_threshold, comp_operator='greater')
Calculates the percentage (p-value) of test-statistic trials that are above the given test-statistic critical value. In addition it calculates the standard deviation of the p-value assuming binomial statistics.
- Parameters:
ts_vals ((n_trials,)-shaped 1D ndarray of float) – The ndarray holding the test-statistic values of the trials.
ts_threshold (float) – The critical test-statistic value.
comp_operator (string, optional) – The comparison operator for p-value calculation. It can be set to one of the following options: ‘greater’ or ‘greater_equal’.
- Returns:
p, p_sigma (tuple(float, float))
- skyllh.core.utils.analysis.calculate_pval_from_gammafit_to_trials(ts_vals, ts_threshold, eta=3.0, n_max=500000)
Calculates the probability (p-value) of test-statistic exceeding the given test-statistic threshold. This calculation relies on fitting a gamma distribution to a list of ts values.
- Parameters:
ts_vals ((n_trials,)-shaped 1D ndarray of float) – The ndarray holding the test-statistic values of the trials.
ts_threshold (float) – The critical test-statistic value.
eta (float, optional) – Test-statistic value at which the gamma function is truncated from below. Default = 3.0.
n_max (int, optional) – The maximum number of trials that should be used during fitting. Default = 500,000
- Returns:
p, p_sigma (tuple(float, float))
- skyllh.core.utils.analysis.calculate_pval_from_trials_mixed(ts_vals, ts_threshold, switch_at_ts=3.0, eta=None, n_max=500000, comp_operator='greater_equal')
Calculates the probability (p-value) of test-statistic exceeding the given test-statistic threshold. This calculation relies on fitting a gamma distribution to a list of ts values if ts_threshold is larger than switch_at_ts. If ts_threshold is smaller then the pvalue will be taken from the trials directly.
- Parameters:
ts_vals ((n_trials,)-shaped 1D ndarray of float) – The ndarray holding the test-statistic values of the trials.
ts_threshold (float) – The critical test-statistic value.
switch_at_ts (float, optional) – Test-statistic value below which p-value is computed from trials directly. For thresholds greater than switch_at_ts the pvalue is calculated using a gamma fit.
eta (float, optional) – Test-statistic value at which the gamma function is truncated from below. Default is None.
n_max (int, optional) – The maximum number of trials that should be used during fitting. Default = 500,000
comp_operator (string, optional) – The comparison operator for p-value calculation. It can be set to one of the following options: ‘greater’ or ‘greater_equal’.
- Returns:
p, p_sigma (tuple(float, float))
- skyllh.core.utils.analysis.truncated_gamma_logpdf(a, scale, eta, ts_above_eta, N_above_eta)
Calculates the -log(likelihood) of a sample of random numbers generated from a gamma pdf truncated from below at x=eta.
- Parameters:
a (float) – Shape parameter.
scale (float) – Scale parameter.
eta (float) – Test-statistic value at which the gamma function is truncated from below.
ts_above_eta ((n_trials,)-shaped 1D ndarray) – The ndarray holding the test-statistic values falling in the truncated gamma pdf.
N_above_eta (int) – Number of test-statistic values falling in the truncated gamma pdf.
- Returns:
-logl (float)
- skyllh.core.utils.analysis.calculate_critical_ts_from_gamma(ts, h0_ts_quantile, eta=3.0)
Calculates the critical test-statistic value corresponding to h0_ts_quantile by fitting the ts distribution with a truncated gamma function.
- Parameters:
ts ((n_trials,)-shaped 1D ndarray) – The ndarray holding the test-statistic values of the trials.
h0_ts_quantile (float) – Null-hypothesis test statistic quantile.
eta (float, optional) – Test-statistic value at which the gamma function is truncated from below.
- Returns:
critical_ts (float)
- skyllh.core.utils.analysis.polynomial_fit(ns, p, p_weight, deg, p_thr)
Performs a polynomial fit on the p-values of test-statistic trials associated to each ns.. Using the fitted parameters it computes the number of signal events correponding to the given p-value critical value.
- Parameters:
ns (1D array_like object) – x-coordinates of the sample.
p (1D array_like object) – y-coordinates of the sample.
p_weight (1D array_like object) – Weights to apply to the y-coordinates of the sample points. For gaussian uncertainties, use 1/sigma.
deg (int) – Degree of the fitting polynomial function.
p_thr (float within [0,1]) – The critical p-value.
- Returns:
ns (float)
- skyllh.core.utils.analysis.estimate_mean_nsignal_for_ts_quantile(ana, rss, p, eps_p, mu_range, critical_ts=None, h0_trials=None, h0_ts_quantile=None, min_dmu=0.5, bkg_kwargs=None, sig_kwargs=None, ppbar=None, tl=None, pathfilename=None)
Calculates the mean number of signal events needed to be injected to reach a test statistic distribution with defined properties for the given analysis.
- Parameters:
ana (Analysis instance) – The Analysis instance to use for the calculation.
rss (instance of RandomStateService) – The RandomStateService instance to use for generating random numbers.
p (float) – Desired probability of signal test statistic for exceeding h0_ts_quantile part of null-hypothesis test statistic threshold.
eps_p (float) – Precision in p as stopping condition for the calculation.
mu_range (2-element sequence) – The range of mu (lower,upper) to search for mean number of signal events.
critical_ts (float | None) – The critical test-statistic value that should be overcome by the signal distribution. If set to None, the null-hypothesis test-statistic distribution will be used to compute the critical TS value.
h0_trials ((n_h0_trials,)-shaped ndarray | None) – The structured ndarray holding the trials for the null-hypothesis. If set to None, the number of trials is calculated from binomial statistics via h0_ts_quantile*(1-h0_ts_quantile)/eps**2, where eps is min(5e-3, h0_ts_quantile/10).
h0_ts_quantile (float | None) – Null-hypothesis test statistic quantile. If set to None, the critical test-statistic value that should be overcome by the signal distribution MUST be given.
min_dmu (float) – The minimum delta mu to use for calculating the derivative dmu/dp. The default is
0.5
.bkg_kwargs (dict | None) – Additional keyword arguments for the generate_events method of the background generation method class. An usual keyword argument is poisson.
sig_kwargs (dict | None) – Additional keyword arguments for the generate_signal_events method of the SignalGenerator class. An usual keyword argument is poisson. If poisson is set to True, the actual number of generated signal events will be drawn from a Poisson distribution with the mean number of signal events, mu.
ppbar (instance of ProgressBar | None) – The possible parent ProgressBar instance.
tl (instance of TimeLord | None) – The optional TimeLord instance that should be used to collect timing information about this function.
pathfilename (string | None) – Trial data file path including the filename. If set to None, generatedtrials won’t be saved.
- Returns:
mu (float) – Estimated mean number of signal events.
mu_err (None) – Error estimate needs to be implemented.
- skyllh.core.utils.analysis.estimate_sensitivity(ana, rss, h0_trials=None, h0_ts_quantile=0.5, p=0.9, eps_p=0.005, mu_range=None, min_dmu=0.5, bkg_kwargs=None, sig_kwargs=None, ppbar=None, tl=None, pathfilename=None)
Estimates the mean number of signal events that whould have to be injected into the data such that the test-statistic value of p*100% of all trials are larger than the critical test-statistic value c, which corresponds to the test-statistic value where h0_ts_quantile*100% of all null hypothesis test-statistic values are larger than c.
For sensitivity h0_ts_quantile, and p are usually set to 0.5, and 0.9, respectively.
- Parameters:
ana (Analysis) – The Analysis instance to use for sensitivity estimation.
rss (RandomStateService) – The RandomStateService instance to use for generating random numbers.
h0_trials ((n_h0_ts_vals,)-shaped ndarray | None) – The strutured ndarray holding the trials for the null-hypothesis. If set to None, the number of trials is calculated from binomial statistics via h0_ts_quantile*(1-h0_ts_quantile)/eps**2, where eps is min(5e-3, h0_ts_quantile/10).
h0_ts_quantile (float, optional) – Null-hypothesis test statistic quantile that defines the critical value.
p (float, optional) – Desired probability of the signal test statistic value to exceed the null-hypothesis test statistic value threshold, which is defined through the h0_ts_quantile value.
eps_p (float, optional) – Precision in p for execution to break.
mu_range (2-element sequence | None) – Range to search for the mean number of signal events. If set to None, the range (0, 10) will be used.
min_dmu (float) – The minimum delta mu to use for calculating the derivative dmu/dp. The default is
0.5
.bkg_kwargs (dict | None) – Additional keyword arguments for the generate_events method of the background generation method class. An usual keyword argument is poisson.
sig_kwargs (dict | None) – Additional keyword arguments for the generate_signal_events method of the SignalGenerator class. An usual keyword argument is poisson. If poisson is set to True, the actual number of generated signal events will be drawn from a Poisson distribution with the mean number of signal events, mu.
ppbar (instance of ProgressBar | None) – The possible parent ProgressBar instance.
tl (instance of TimeLord | None) – The optional TimeLord instance that should be used to collect timing information about this function.
pathfilename (string | None) – Trial data file path including the filename. If set to None, generated trials won’t be saved.
- Returns:
mu (float) – Estimated median number of signal events to reach desired sensitivity.
mu_err (float) – The uncertainty of the estimated mean number of signal events.
- skyllh.core.utils.analysis.estimate_discovery_potential(ana, rss, h0_trials=None, h0_ts_quantile=2.8665e-07, p=0.5, eps_p=0.005, mu_range=None, min_dmu=0.5, bkg_kwargs=None, sig_kwargs=None, ppbar=None, tl=None, pathfilename=None)
Estimates the mean number of signal events that whould have to be injected into the data such that the test-statistic value of p*100% of all trials are larger than the critical test-statistic value c, which corresponds to the test-statistic value where h0_ts_quantile*100% of all null hypothesis test-statistic values are larger than c.
For the 5 sigma discovery potential h0_ts_quantile, and p are usually set to 2.8665e-7, and 0.5, respectively.
- Parameters:
ana (Analysis) – The Analysis instance to use for discovery potential estimation.
rss (RandomStateService) – The RandomStateService instance to use for generating random numbers.
h0_trials ((n_h0_ts_vals,)-shaped ndarray | None) – The structured ndarray holding the trials for the null-hypothesis. If set to None, the number of trials is calculated from binomial statistics via h0_ts_quantile*(1-h0_ts_quantile)/eps**2, where eps is min(5e-3, h0_ts_quantile/10).
h0_ts_quantile (float, optional) – Null-hypothesis test statistic quantile that defines the critical value.
p (float, optional) – Desired probability of the signal test statistic value to exceed the critical value.
eps_p (float, optional) – Precision in p for execution to break.
mu_range (2-element sequence | None) – Range to search for the mean number of signal events. If set to None, the range (0, 10) will be used.
min_dmu (float) – The minimum delta mu to use for calculating the derivative dmu/dp. The default is
0.5
.bkg_kwargs (dict | None) – Additional keyword arguments for the generate_events method of the background generation method class. An usual keyword argument is poisson.
sig_kwargs (dict | None) – Additional keyword arguments for the generate_signal_events method of the SignalGenerator class. An usual keyword argument is poisson. If poisson is set to True, the actual number of generated signal events will be drawn from a Poisson distribution with the mean number of signal events, mu.
ppbar (instance of ProgressBar | None) – The possible parent ProgressBar instance.
tl (instance of TimeLord | None) – The optional TimeLord instance that should be used to collect timing information about this function.
pathfilename (string | None) – Trial data file path including the filename. If set to None, generated trials won’t be saved.
- Returns:
mu (float) – Estimated mean number of injected signal events to reach the desired discovery potential.
mu_err (float) – Estimated error of mu.
- skyllh.core.utils.analysis.generate_mu_of_p_spline_interpolation(ana, rss, h0_ts_vals, h0_ts_quantile, eps_p, mu_range, mu_step, kind='cubic', bkg_kwargs=None, sig_kwargs=None, ppbar=None, tl=None)
Generates a spline interpolation for mu(p) function for a pre-defined range of mu, where mu is the mean number of injected signal events and p the probability for the ts value larger than the ts value corresponding to the given quantile, h0_ts_quantile, of the null hypothesis test-statistic distribution.
- Parameters:
ana (instance of Analysis) – The Analysis instance to use for the calculation.
rss (instance of RandomStateService) – The RandomStateService instance to use for generating random numbers.
h0_ts_vals ((n_h0_ts_vals,)-shaped 1D ndarray | None) – The 1D ndarray holding the test-statistic values for the null-hypothesis. If set to None, 100/(1-h0_ts_quantile) null-hypothesis trials will be generated.
h0_ts_quantile (float) – Null-hypothesis test statistic quantile, which should be exceeded by the alternative hypothesis ts value.
eps_p (float) – The one sigma precision in p as stopping condition for the calculation for a single mu value.
mu_range (2-element sequence) – The range (lower,upper) of mean number of injected signal events to create the interpolation spline for.
mu_step (float) – The step size of the mean number of signal events.
kind (str) – The kind of spline to generate. Possble values are ‘linear’ and ‘cubic’ (default).
bkg_kwargs (dict | None) – Additional keyword arguments for the generate_events method of the background generation method class. An usual keyword argument is poisson.
sig_kwargs (dict | None) – Additional keyword arguments for the generate_signal_events method of the SignalGenerator class. An usual keyword argument is poisson. If poisson is set to True, the actual number of generated signal events will be drawn from a Poisson distribution with the mean number of signal events, mu.
ppbar (instance of ProgressBar | None) – The possible parent ProgressBar instance.
tl (instance of TimeLord | None) – The optional TimeLord instance that should be used to collect timing information about this function.
- Returns:
spline (callable) – The spline function mu(p).
- skyllh.core.utils.analysis.create_trial_data_file(ana, rss, n_trials, mean_n_sig=0, mean_n_sig_null=0, mean_n_bkg_list=None, minimizer_rss=None, bkg_kwargs=None, sig_kwargs=None, pathfilename=None, ncpu=None, ppbar=None, tl=None)
Creates and fills a trial data file with n_trials generated trials for each mean number of injected signal events from ns_min up to ns_max for a given analysis.
- Parameters:
ana (instance of Analysis) – The Analysis instance to use for the trial generation.
rss (instance of RandomStateService) – The RandomStateService instance to use for generating random numbers.
n_trials (int) – The number of trials to perform for each hypothesis test.
mean_n_sig (ndarray of float | float | 2- or 3-element sequence of float) – The array of mean number of injected signal events (MNOISEs) for which to generate trials. If this argument is not a ndarray, an array of MNOISEs is generated based on this argument. If a single float is given, only this given MNOISEs are injected. If a 2-element sequence of floats is given, it specifies the range of MNOISEs with a step size of one. If a 3-element sequence of floats is given, it specifies the range plus the step size of the MNOISEs.
mean_n_sig_null (ndarray of float | float | 2- or 3-element sequence of float) – The array of the fixed mean number of signal events (FMNOSEs) for the null-hypothesis for which to generate trials. If this argument is not a ndarray, an array of FMNOSEs is generated based on this argument. If a single float is given, only this given FMNOSEs are used. If a 2-element sequence of floats is given, it specifies the range of FMNOSEs with a step size of one. If a 3-element sequence of floats is given, it specifies the range plus the step size of the FMNOSEs.
mean_n_bkg_list (list of float | None) – The mean number of background events that should be generated for each dataset. This parameter is passed to the
do_trials
method of theAnalysis
class. If set to None (the default), the background generation method needs to obtain this number itself.minimizer_rss (instance of RandomStateService | None) – The instance of RandomStateService to use for generating random numbers for the minimizer, e.g. new initial fit parameter values. If set to
None
, a rss with the same seed asrss
will be initialized.bkg_kwargs (dict | None) – Additional keyword arguments for the generate_events method of the background generation method class. An usual keyword argument is poisson.
sig_kwargs (dict | None) – Additional keyword arguments for the generate_signal_events method of the SignalGenerator class. An usual keyword argument is poisson.
pathfilename (string | None) – Trial data file path including the filename. If set to None generated trials won’t be saved.
ncpu (int | None) – The number of CPUs to use.
ppbar (instance of ProgressBar | None) – The optional instance of the parent progress bar.
tl (instance of TimeLord | None) – The instance of TimeLord that should be used to measure individual tasks.
- Returns:
seed (int) – The seed used to generate the trials.
mean_n_sig (1d ndarray) – The array holding the mean number of signal events used to generate the trials.
mean_n_sig_null (1d ndarray) – The array holding the fixed mean number of signal events for the null-hypothesis used to generate the trials.
trial_data (structured numpy ndarray) – The generated trial data.
- skyllh.core.utils.analysis.extend_trial_data_file(ana, rss, n_trials, trial_data, mean_n_sig=0, mean_n_sig_null=0, mean_n_bkg_list=None, bkg_kwargs=None, sig_kwargs=None, pathfilename=None, **kwargs)
Appends to the trial data file n_trials generated trials for each mean number of injected signal events up to ns_max for a given analysis.
- Parameters:
ana (instance of Analysis) – The Analysis instance to use for sensitivity estimation.
rss (instance of RandomStateService) – The RandomStateService instance to use for generating random numbers.
n_trials (int) – The number of trials the trial data file needs to be extended by.
trial_data (structured numpy ndarray) – The structured numpy ndarray holding the trials.
mean_n_sig (ndarray of float | float | 2- or 3-element sequence of float) – The array of mean number of injected signal events (MNOISEs) for which to generate trials. If this argument is not a ndarray, an array of MNOISEs is generated based on this argument. If a single float is given, only this given MNOISEs are injected. If a 2-element sequence of floats is given, it specifies the range of MNOISEs with a step size of one. If a 3-element sequence of floats is given, it specifies the range plus the step size of the MNOISEs.
mean_n_sig_null (ndarray of float | float | 2- or 3-element sequence of float) – The array of the fixed mean number of signal events (FMNOSEs) for the null-hypothesis for which to generate trials. If this argument is not a ndarray, an array of FMNOSEs is generated based on this argument. If a single float is given, only this given FMNOSEs are used. If a 2-element sequence of floats is given, it specifies the range of FMNOSEs with a step size of one. If a 3-element sequence of floats is given, it specifies the range plus the step size of the FMNOSEs.
bkg_kwargs (dict | None) – Additional keyword arguments for the generate_events method of the background generation method class. An usual keyword argument is poisson.
sig_kwargs (dict | None) – Additional keyword arguments for the generate_signal_events method of the SignalGenerator class. An usual keyword argument is poisson.
pathfilename (string | None) – Trial data file path including the filename.
arguments (Additional keyword) –
---------------------------- –
create_trial_data_file (Additional keyword arguments are passed-on to the) –
function. –
- Returns:
trial_data – Trial data file extended by the required number of trials for each mean number of injected signal events..
- skyllh.core.utils.analysis.calculate_upper_limit_distribution(ana, rss, pathfilename, n_bkg=5000, n_bins=100)
Function to calculate upper limit distribution. It loads the trial data file containing test statistic distribution and calculates 10 percentile value for each mean number of injected signal event. Then it finds upper limit values which correspond to generated background trials test statistic values by linearly interpolated curve of 10 percentile values distribution.
- Parameters:
ana (instance of Analysis) – The Analysis instance to use for sensitivity estimation.
rss (instance of RandomStateService) – The RandomStateService instance to use for generating random numbers.
pathfilename (string) – Trial data file path including the filename.
n_bkg (int, optional) – Number of times to perform background analysis trial.
n_bins (int, optional) – Number of returned test statistic histograms bins.
- Returns:
result (dict) – Result dictionary which contains the following fields:
- ullist of float
List of upper limit values.
- meanfloat
Mean of upper limit values.
- medianfloat
Median of upper limit values.
- varfloat
Variance of upper limit values.
- ts_histnumpy ndarray
2D array of test statistic histograms calculated by axis 1.
- extentlist of float
Test statistic histogram boundaries.
- q_valueslist of float
q percentile values of test statistic for different injected events means.
skyllh.core.utils.coords module
- skyllh.core.utils.coords.rotate_spherical_vector(ra1, dec1, ra2, dec2, ra3, dec3)
Calculates the rotation matrix R to rotate the spherical vector (ra1,dec1) onto the direction (ra2,dec2), and performs this rotation on the spherical vector (ra3,dec3).
In practice (ra1,dec1) refers to the true location of a MC event, (ra2,dec2) the true location of the signal source, and (ra3,dec3) the reconstructed location of the MC event, which should get rotated according to the rotation of the two true directions.
- skyllh.core.utils.coords.rotate_signal_events_on_sphere(src_ra, src_dec, evt_true_ra, evt_true_dec, evt_reco_ra, evt_reco_dec)
Rotate signal events on a sphere to a given source position preserving position angle and separation (great circle distance) between the event’s true and reco directions.
- Parameters:
src_ra (instance of numpy.ndarray) – The (N_events,)-shaped 1D numpy.ndarray holding the true right-ascension of the source.
src_dec (instance of numpy.ndarray) – The (N_events,)-shaped 1D numpy.ndarray holding the true declination of the source.
evt_true_ra (instance of numpy.ndarray) – The (N_events,)-shaped 1D numpy.ndarray holding the true right-ascension of the MC event.
evt_true_dec (instance of numpy.ndarray) – The (N_events,)-shaped 1D numpy.ndarray holding the true declination of the MC event.
evt_reco_ra (instance of numpy.ndarray) – The (N_events,)-shaped 1D numpy.ndarray holding the reconstructed right-ascension of the MC event.
evt_reco_dec (instance of numpy.ndarray) – The (N_events,)-shaped 1D numpy.ndarray holding the reconstructed declination of the MC event.
- Returns:
rot_evt_reco_ra (instance of numpy.ndarray) – The (N_events,)-shaped 1D numpy.ndarray holding the rotated reconstructed event right-ascension.
rot_evt_reco_dec (instance of numpy.ndarray) – The (N_events,)-shaped 1D numpy.ndarray holding the rotated reconstructed event declination.
- skyllh.core.utils.coords.angular_separation(ra1, dec1, ra2, dec2, psi_floor=None)
Calculates the angular separation on the sphere between two vectors on the sphere.
- Parameters:
ra1 (instance of numpy.ndarray) – The (N_events,)-shaped numpy.ndarray holding the right-ascension or longitude coordinate of the first vector in radians.
dec1 (instance of numpy.ndarray) – The (N_events,)-shaped numpy.ndarray holding declination or latitude coordinate of the first vector in radians.
ra2 (instance of numpy.ndarray) – The (N_events,)-shaped numpy.ndarray holding the right-ascension or longitude coordinate of the second vector in radians.
dec2 (instance of numpy.ndarray) – The (N_events,)-shaped numpy.ndarray holding declination coordinate of the second vector in radians.
psi_floor (float | None) – If not
None
, specifies the floor value of psi.
- Returns:
psi (instance of numpy.ndarray) – The (N_events,)-shaped numpy.ndarray holding the calculated angular separation value of each event.
skyllh.core.utils.flux_model module
- skyllh.core.utils.flux_model.create_scipy_stats_rv_continuous_from_TimeFluxProfile(profile)
This function builds a scipy.stats.rv_continuous instance for a given
TimeFluxProfile
instance.It can be used to generate random numbers according to the given time flux profile function.
- Parameters:
profile (instance of TimeFluxProfile) – The instance of TimeFluxProfile providing the function of the time flux profile.
- Returns:
rv (instance of rv_continuous_frozen) – The instance of rv_continuous_frozen representing the time flux profile as a continuous random variate instance.
skyllh.core.utils.multidimgridpdf module
This module contains utility functions for creating and managing MultiDimGridPDF instances.
- skyllh.core.utils.multidimgridpdf.get_kde_pdf_sig_spatial_norm_factor_func(log10_psi_name='log10_psi')
Returns the standard normalization factor function for the spatial signal MultiDimGridPDF, which is created from KDE PDF values. It can be used for the
norm_factor_func
argument of thecreate_MultiDimGridPDF_from_photosplinetable
andcreate_MultiDimGridPDF_from_kde_pdf
function.- Parameters:
log10_psi_name (str) – The name of the event data field for the log10(psi) values.
- skyllh.core.utils.multidimgridpdf.get_kde_pdf_bkg_norm_factor_func()
Returns the standard normalization factor function for the background MultiDimGridPDF, which is created from KDE PDF values. It can be used for the
norm_factor_func
argument of thecreate_MultiDimGridPDF_from_photosplinetable
andcreate_MultiDimGridPDF_from_kde_pdf
function.
- skyllh.core.utils.multidimgridpdf.create_MultiDimGridPDF_from_photosplinetable(multidimgridpdf_cls, pmm, ds, data, info_key, splinetable_key, kde_pdf_axis_name_map_key='KDE_PDF_axis_name_map', norm_factor_func=None, cache_pd_values=False, tl=None, **kwargs)
Creates a MultiDimGridPDF instance with pdf values taken from a photospline pdf, i.e. a spline interpolation of KDE PDF values stored in a splinetable on disk.
- Parameters:
multidimgridpdf_cls (subclass of MultiDimGridPDF) – The MultiDimGridPDF class, which should be used.
pmm (instance of ParameterModelMapper) – The instance of ParameterModelMapper, which defines the mapping of global parameters to local model parameters.
ds (instance of Dataset) – The instance of Dataset the PDF applies to.
data (instance of DatasetData) – The instance of DatasetData that holds the experimental and monte-carlo data of the dataset.
info_key (str) – The auxiliary data name for the file containing PDF information.
splinetable_key (str) – The auxiliary data name for the name of the file containing the photospline spline table.
kde_pdf_axis_name_map_key (str) – The auxiliary data name for the KDE PDF axis name map.
norm_factor_func (callable | None) – The function that calculates a possible required normalization factor for the PDF value based on the event properties. For more information about this argument see the documentation of the
skyllh.core.pdf.MultiDimGridPDF.__init__()
method.cache_pd_values (bool) – Flag if the probability density values should get cached by the MultiDimGridPDF class.
tl (instance of TimeLord | None) – The optional instance of TimeLord to use for measuring timing information.
- Returns:
pdf (instance of
multidimgridpdf_cls
) – The created PDF instance of MultiDimGridPDF.
- skyllh.core.utils.multidimgridpdf.create_MultiDimGridPDF_from_kde_pdf(multidimgridpdf_cls, pmm, ds, data, numerator_key, denumerator_key=None, kde_pdf_axis_name_map_key='KDE_PDF_axis_name_map', norm_factor_func=None, cache_pd_values=False, tl=None, **kwargs)
Creates a MultiDimGridPDF instance with pdf values taken from KDE PDF values stored in the dataset’s auxiliary data.
- Parameters:
multidimgridpdf_cls (subclass of MultiDimGridPDF) – The MultiDimGridPDF class, which should be used.
pmm (instance of ParameterModelMapper) – The instance of ParameterModelMapper, which defines the mapping of global parameters to local model parameters.
ds (instance of Dataset) – The instance of Dataset the PDF applies to.
data (instance of DatasetData) – The instance of DatasetData that holds the auxiliary data of the dataset.
numerator_key (str) – The auxiliary data name for the PDF numerator array.
denumerator_key (str | None) – The auxiliary data name for the PDF denumerator array. This can be
None
, if no denumerator array is required.kde_pdf_axis_name_map_key (str) – The auxiliary data name for the KDE PDF axis name map.
norm_factor_func (callable | None) – The function that calculates a possible required normalization factor for the PDF value based on the event properties. For more information about this argument see the documentation of the
skyllh.core.pdf.MultiDimGridPDF.__init__()
method.cache_pd_values (bool) – Flag if the probability density values should get cached by the MultiDimGridPDF class.
tl (instance of TimeLord | None) – The optional instance of TimeLord to use for measuring timing information.
- Returns:
pdf (instance of
multidimgridpdf_cls
) – The created PDF instance of MultiDimGridPDF.
skyllh.core.utils.spline module
- skyllh.core.utils.spline.make_spline_1d(x, y, kind='linear', **kwargs)
Creates a 1D spline for the function y(x) using
scipy.interpolate.interp1d
.- Parameters:
x (array_like) – The x values.
y (array_like) – The y values.
kind (str) – The kind of the spline. See the
scipy.interpolate.interp1d
documentation for possible values. Default is'linear'
.**kwargs – Additional keyword arguments are passed to the
interp1d
function.
- Returns:
spline – The created 1D spline instance.
- class skyllh.core.utils.spline.CatmullRomRegular1DSpline(x, y, **kwargs)
Bases:
object
This class provides a one-dimensional Catmull-Rom spline which is a C^1 continous spline, where the control points coincide with the data points. The x data points need to be equal distant.
Note
The first and last data point are not part of the splined curve!
Creates a new CatmullRom1DSpline instance.
- Parameters:
x (instance of ndarray) – The x values of the data points.
y (instance of ndarray) – The y values of the data points.
- __call__(x, oor_value=nan)
Evaluates the spline given x-values in data coordinates.
- Parameters:
x (instance of ndarray) – The instance of ndarray holding the values for which the spline should get evaluated.
- Returns:
y (instance of ndarray) – The instance of ndarray with the spline values at the given x values.
skyllh.core.utils.tdm module
- skyllh.core.utils.tdm.get_tdm_field_func_psi(psi_floor=None)
Returns the TrialDataManager (TDM) field function for psi with an optional psi value floor.
- Parameters:
psi_floor (float | None) – The optional floor value for psi. This should be
None
for a standard point-source analysis that uses an analytic function for the detector’s point-spread-function (PSF).- Returns:
tdm_field_func_psi (function) – TrialDataManager (TDM) field function for psi.
skyllh.core.utils.trials module
This module contains utility functions related analysis trials.
- skyllh.core.utils.trials.create_pseudo_data_file(ana, rss, filename, mean_n_bkg_list=None, mean_n_sig=0, bkg_kwargs=None, sig_kwargs=None, tl=None)
Creates a pickle file that contains the pseudo data for a single trial by generating background and signal events.
- Parameters:
ana (Analysis) – The Analysis instance that should be used to generate the pseudo data.
rss (RandomStateService) – The RandomStateService instance to use for generating random numbers.
filename (str) – The data file name into which the generated pseudo data should get written to.
mean_n_bkg_list (list of float | None) – The mean number of background events that should be generated for each dataset. If set to None (the default), the background generation method needs to obtain this number itself.
mean_n_sig (float) – The mean number of signal events that should be generated for the trial. The actual number of generated events will be drawn from a Poisson distribution with this given signal mean as mean.
bkg_kwargs (dict | None) – Additional keyword arguments for the generate_events method of the background generation method class. An usual keyword argument is poisson.
sig_kwargs (dict | None) – Additional keyword arguments for the generate_signal_events method of the SignalGenerator class. An usual keyword argument is poisson.
tl (TimeLord | None) – The instance of TimeLord that should be used to time individual tasks.
- skyllh.core.utils.trials.load_pseudo_data(filename, tl=None)
Loads the pseudo data for a single trial from the given file name.
- Parameters:
filename (str) – The name of the file that contains the pseudo data.
tl (TimeLord | None) – The instance of TimeLord that should be used to time individual tasks.
- Returns:
mean_n_sig (float) – The mean number of signal events that was used to generate the pseudo data.
n_sig (int) – The actual total number of signal events in the pseudo data.
n_bkg_events_list (list of int) – The total number of background events for each data set of the pseudo data.
n_sig_events_list (list of int) – The total number of signal events for each data set of the pseudo data.
bkg_events_list (list of DataFieldRecordArray instances) – The list of DataFieldRecordArray instances containing the background pseudo data events for each data set.
sig_events_list (list of DataFieldRecordArray instances | None) – The list of DataFieldRecordArray instances containing the signal pseudo data events for each data set. If a particular dataset has no signal events, the entry for that dataset can be None.