SensitivityEstimator¶
-
class
gammapy.scripts.
SensitivityEstimator
(irf, livetime, slope=2.0, alpha=0.2, sigma=5.0, gamma_min=10.0, bkg_sys=0.05, random=0)[source]¶ Bases:
object
Estimate differential sensitivity for the ON/OFF method, ie 1D analysis
Parameters: irf :
CTAPerf
IRF object
livetime :
Quantity
Livetime (object with the units of time), e.g. 5*u.h
slope :
float
, optionalIndex of the spectral shape (Power-law), should be positive (>0)
alpha :
float
, optionalOn/OFF normalisation
sigma :
float
, optionalMinimum significance
gamma_min :
float
, optionalMinimum number of gamma-rays
bkg_sys :
float
, optionalFraction of Background systematics relative to the number of ON counts
random: :
int
, optionalNumber of random trial to derive the number of gamma-rays
Notes
For the moment, only the differential point-like sensitivity is computed at a fixed offset. This class allows to determine for each reconstructed energy bin the flux associated to the number of gamma-ray events for which the significance is ‘sigma’, and being larger than ‘gamma_min’ and ‘bkg_sys’% larger than the number of background events in the ON region
TODO:
- make the computation for any source size
- compute the integral sensitivity
- Add options to use different spectral shape?
Examples
Compute and plot a sensitivity curve for CTA:
from gammapy.scripts import CTAPerf, SensitivityEstimator filename = '$GAMMAPY_EXTRA/datasets/cta/perf_prod2/point_like_non_smoothed/South_5h.fits.gz' irf = CTAPerf.read(filename) sens = SensitivityEstimator(irf=irf, livetime='5h') sens.run() sens.print_results() sens.plot()
Further examples in cta_sensitivity.html .
Attributes Summary
diff_sensi_table
Differential sensitivity table ( Table
).Methods Summary
get_1TeV_differential_flux
(excess_counts, …)Compute the differential fluxes at 1 TeV. get_bkg
(bkg_rate)Return the Background rate for each energy bin get_excess
(bkg_counts)Compute the gamma-ray excess for each energy bin. plot
([ax])Plot the sensitivity curve. print_results
()Print results to the console. run
()Run the algorithm to compute the differential sensitivity as explained in the document of the class. Attributes Documentation
Methods Documentation
-
get_1TeV_differential_flux
(excess_counts, model, aeff, rmf)[source]¶ Compute the differential fluxes at 1 TeV.
Parameters: excess_counts :
ndarray
Array of gamma-ray excess (bins in reconstructed energy)
model :
SpectralModel
Spectral model
aeff :
EffectiveAreaTable
Effective area
rmf :
EnergyDispersion
RMF
Returns: flux :
Quantity
Array of 1TeV fluxes (bins in reconstructed energy)
-
get_bkg
(bkg_rate)[source]¶ Return the Background rate for each energy bin
Parameters: bkg_rate :
BgRateTable
Table of background rate
Returns: rate :
ndarray
Return an array of background counts (bins in reconstructed energy)
-
get_excess
(bkg_counts)[source]¶ Compute the gamma-ray excess for each energy bin.
Parameters: bkg_counts :
ndarray
Array of background counts (bins in reconstructed energy)
Returns: count :
ndarray
Array of gamma-ray excess (bins in reconstructed energy)
Notes
Find the number of needed gamma excess events using newtons method. Defines a function
significance_on_off(x, off, alpha) - self.sigma
and uses scipy.optimize.newton to find thex
for which this function is zero.