+[docs]
+defeffective_area_3d_polar(
+ selected_events,
+ simulation_info,
+ true_energy_bins,
+ fov_offset_bins,
+ fov_position_angle_bins,
+):
+"""
+ Calculate effective area in bins of true energy, field of view offset, and field of view position angle.
+
+ Parameters
+ ----------
+ selected_events: astropy.table.QTable
+ DL2 events table, required columns for this function:
+ - `true_energy`
+ - `true_source_fov_offset`
+ - `true_source_fov_position_angle`
+ simulation_info: pyirf.simulations.SimulatedEventsInfo
+ The overall statistics of the simulated events
+ true_energy_bins: astropy.units.Quantity[energy]
+ The true energy bin edges in which to calculate effective area.
+ fov_offset_bins: astropy.units.Quantity[angle]
+ The field of view radial bin edges in which to calculate effective area.
+ fov_position_angle_bins: astropy.units.Quantity[radian]
+ The field of view azimuthal bin edges in which to calculate effective area.
+ """
+ area=np.pi*simulation_info.max_impact**2
+
+ hist_simulated=simulation_info.calculate_n_showers_3d_polar(
+ true_energy_bins,fov_offset_bins,fov_position_angle_bins
+ )
+
+ hist_selected,_=np.histogramdd(
+ np.column_stack(
+ [
+ selected_events["true_energy"].to_value(u.TeV),
+ selected_events["true_source_fov_offset"].to_value(u.deg),
+ selected_events["true_source_fov_position_angle"].to_value(u.rad),
+ ]
+ ),
+ bins=(
+ true_energy_bins.to_value(u.TeV),
+ fov_offset_bins.to_value(u.deg),
+ fov_position_angle_bins.to_value(u.rad),
+ ),
+ )
+
+ returneffective_area(hist_selected,hist_simulated,area)
+
+
+
+
+[docs]
+defeffective_area_3d_lonlat(
+ selected_events,
+ simulation_info,
+ true_energy_bins,
+ fov_longitude_bins,
+ fov_latitude_bins,
+ subpixels=20,
+):
+"""
+ Calculate effective area in bins of true energy, field of view longitude, and field of view latitude.
+
+ Parameters
+ ----------
+ selected_events: astropy.table.QTable
+ DL2 events table, required columns for this function:
+ - `true_energy`
+ - `true_source_fov_lon`
+ - `true_source_fov_lat`
+ simulation_info: pyirf.simulations.SimulatedEventsInfo
+ The overall statistics of the simulated events
+ true_energy_bins: astropy.units.Quantity[energy]
+ The true energy bin edges in which to calculate effective area.
+ fov_longitude_bins: astropy.units.Quantity[angle]
+ The field of view longitude bin edges in which to calculate effective area.
+ fov_latitude_bins: astropy.units.Quantity[angle]
+ The field of view latitude bin edges in which to calculate effective area.
+ """
+ area=np.pi*simulation_info.max_impact**2
+
+ hist_simulated=simulation_info.calculate_n_showers_3d_lonlat(
+ true_energy_bins,fov_longitude_bins,fov_latitude_bins,subpixels=subpixels
+ )
+
+ selected_columns=np.column_stack(
+ [
+ selected_events["true_energy"].to_value(u.TeV),
+ selected_events["true_source_fov_lon"].to_value(u.deg),
+ selected_events["true_source_fov_lat"].to_value(u.deg),
+ ]
+ )
+ bins=(
+ true_energy_bins.to_value(u.TeV),
+ fov_longitude_bins.to_value(u.deg),
+ fov_latitude_bins.to_value(u.deg),
+ )
+
+ hist_selected,_=np.histogramdd(selected_columns,bins=bins)
+
+ returneffective_area(hist_selected,hist_simulated,area)
+[docs]
+ @u.quantity_input(
+ energy_bins=u.TeV,fov_offset_bins=u.deg,fov_position_angle_bins=u.rad
+ )
+ defcalculate_n_showers_3d_polar(
+ self,energy_bins,fov_offset_bins,fov_position_angle_bins
+ ):
+"""
+ Calculate number of showers that were simulated in the given
+ energy and 2D fov bins in polar coordinates.
+
+ This assumes the events were generated uniformly distributed per solid angle,
+ and from a powerlaw in energy like CORSIKA simulates events.
+
+ Parameters
+ ----------
+ energy_bins: astropy.units.Quantity[energy]
+ The energy bin edges for which to calculate the number of simulated showers
+ fov_offset_bins: astropy.units.Quantity[angle]
+ The FOV radial bin edges for which to calculate the number of simulated showers
+ fov_position_angle_bins: astropy.units.Quantity[radian]
+ The FOV azimuthal bin edges for which to calculate the number of simulated showers
+
+ Returns
+ -------
+ n_showers: numpy.ndarray(ndim=3)
+ The expected number of events inside each of the
+ ``energy_bins``, ``fov_offset_bins`` and ``fov_position_angle_bins``.
+ Dimension (n_energy_bins, n_fov_offset_bins, n_fov_position_angle_bins)
+ This is a floating point number.
+ The actual numbers will follow a poissionian distribution around this
+ expected value.
+ """
+ e_fov_offset_integral=self.calculate_n_showers_per_energy_and_fov(
+ energy_bins,fov_offset_bins
+ )
+ viewcone_integral=self.calculate_n_showers_per_fov(
+ [self.viewcone_min,self.viewcone_max]*u.deg
+ )
+
+ n_bins_pa=len(fov_position_angle_bins)-1
+ position_angle_integral=np.full(n_bins_pa,viewcone_integral/n_bins_pa)
+
+ total_integral=e_fov_offset_integral[:,:,np.newaxis]*position_angle_integral
+
+ returntotal_integral/self.n_showers
+
+
+
+[docs]
+ @u.quantity_input(
+ energy_bins=u.TeV,fov_longitude_bins=u.deg,fov_latitude_bins=u.rad
+ )
+ defcalculate_n_showers_3d_lonlat(
+ self,energy_bins,fov_longitude_bins,fov_latitude_bins,subpixels=20
+ ):
+"""
+ Calculate number of showers that were simulated in the given
+ energy and 2D fov bins in nominal coordinates.
+
+ This assumes the events were generated uniformly distributed per solid angle,
+ and from a powerlaw in energy like CORSIKA simulates events.
+
+ Parameters
+ ----------
+ energy_bins: astropy.units.Quantity[energy]
+ The energy bin edges for which to calculate the number of simulated showers
+ fov_longitude_bins: astropy.units.Quantity[angle]
+ The FOV longitude bin edges for which to calculate the number of simulated showers
+ fov_latitude_bins: astropy.units.Quantity[angle]
+ The FOV latitude bin edges for which to calculate the number of simulated showers
+
+ Returns
+ -------
+ n_showers: numpy.ndarray(ndim=3)
+ The expected number of events inside each of the
+ ``energy_bins``, ``fov_longitude_bins`` and ``fov_latitude_bins``.
+ Dimension (n_energy_bins, n_fov_longitude_bins, n_fov_latitude_bins)
+ This is a floating point number.
+ The actual numbers will follow a poissionian distribution around this
+ expected value.
+ """
+
+ fov_overlap=_fov_lonlat_grid_overlap(
+ self,fov_longitude_bins,fov_latitude_bins,subpixels=subpixels
+ )
+
+ bin_grid_lon,bin_grid_lat=np.meshgrid(fov_longitude_bins,fov_latitude_bins)
+ bin_area=rectangle_solid_angle(
+ bin_grid_lon[:-1,:-1],
+ bin_grid_lon[1:,1:],
+ bin_grid_lat[:-1,:-1],
+ bin_grid_lat[1:,1:],
+ )
+ viewcone_area=cone_solid_angle(self.viewcone_max)-cone_solid_angle(self.viewcone_min)
+
+ shower_density=self.n_showers/viewcone_area
+
+ fov_integral=shower_density*bin_area
+ e_integral=self.calculate_n_showers_per_energy(energy_bins)
+
+ fov_integral=fov_overlap*fov_integral
+
+ return(e_integral[:,np.newaxis,np.newaxis]*fov_integral)/self.n_showers
Calculate number of showers that were simulated in the given
+energy and 2D fov bins in nominal coordinates.
+
This assumes the events were generated uniformly distributed per solid angle,
+and from a powerlaw in energy like CORSIKA simulates events.
+
+
Parameters:
+
+
energy_bins: astropy.units.Quantity[energy]
The energy bin edges for which to calculate the number of simulated showers
+
+
fov_longitude_bins: astropy.units.Quantity[angle]
The FOV longitude bin edges for which to calculate the number of simulated showers
+
+
fov_latitude_bins: astropy.units.Quantity[angle]
The FOV latitude bin edges for which to calculate the number of simulated showers
+
+
+
+
Returns:
+
+
n_showers: numpy.ndarray(ndim=3)
The expected number of events inside each of the
+energy_bins, fov_longitude_bins and fov_latitude_bins.
+Dimension (n_energy_bins, n_fov_longitude_bins, n_fov_latitude_bins)
+This is a floating point number.
+The actual numbers will follow a poissionian distribution around this
+expected value.
The FOV azimuthal bin edges for which to calculate the number of simulated showers
+
+
+
+
Returns:
+
+
n_showers: numpy.ndarray(ndim=3)
The expected number of events inside each of the
+energy_bins, fov_offset_bins and fov_position_angle_bins.
+Dimension (n_energy_bins, n_fov_offset_bins, n_fov_position_angle_bins)
+This is a floating point number.
+The actual numbers will follow a poissionian distribution around this
+expected value.