NuRadioMC.utilities.Veff module

NuRadioMC.utilities.Veff.remove_duplicate_triggers(triggered, gids)[source]

remove duplicate entried from triggered array

The hdf5 file contains a line per shower. One event can contain many showers, i.e. if we count all triggeres from all showers we overestimate the effective volume. This function modifies the triggered array such that it contains not more than one True value for each event group.

Parameters:
triggered: array of bools
gids: array of ints

the event group ids

Returns:
array of bools

the corrected triggered array

NuRadioMC.utilities.Veff.FC_limits(counts)[source]

Returns the 68% confidence belt for a number of counts, using the Feldman-Cousins method. Published in Phys. Rev. D 57, 3873, DOI 10.1103/PhysRevD.57.3873

Parameters:
counts: integer or float

Number of counts. Can be non-integer (weighted counts)

Returns:
(low_limit, upper_limit): float tuple

Lower and upper limits for the confidence belt.

NuRadioMC.utilities.Veff.get_Veff_water_equivalent(Veff, density_medium=5.723464435717068e+39, density_water=6.241509744511524e+39)[source]

convenience function to converte the effective volume of a medium with density density_medium to the water equivalent effective volume

Parameters:
Veff: float or array

the effective volume

dentity_medium: float (optional)

the density of the medium of the Veff simulation (default deep ice)

density water: float (optional)

the density of water

Returns:
water equivalen effective volume
NuRadioMC.utilities.Veff.get_veff_output(volume, counts, all_events)[source]
NuRadioMC.utilities.Veff.get_Veff_Aeff_single(filename, trigger_names, trigger_names_dict, trigger_combinations, deposited, station, veff_aeff='veff', bounds_theta=[0, 3.141592653589793])[source]

Calculates the effective volume or effective area from surface muons from a single NuRadioMC hdf5 file

The effective volume is NOT normalized to a water equivalent. It is also NOT multiplied with the solid angle (typically 4pi).

Parameters:
filename: string

filename of the hdf5 file

trigger_names: list of strings

list of the trigger names contained in the file

trigger_names_dict: dict

map from trigger name to index

trigger_combinations: dict, optional

keys are the names of triggers to calculate. Values are dicts again:

  • ‘triggers’: list of strings

    name of individual triggers that are combined with an OR

the following additional options are optional

  • ‘efficiency’: dict

    allows to apply an (analysis) efficiency cut for calculating effective volumes

    • ‘func’: function

      a function that paramaterized the efficiency as a function of SNR (=Vmax/Vrms)

    • ‘channel_ids’: array on ints

      the channels for which the maximum signal amplitude should be determined

    • ‘scale’: float

      rescaling of the efficiency curve by SNR’ = SNR * scale

  • ‘n_reflections’: int

    the number of bottom reflections of the ray tracing solution that likely triggered assuming that the solution with the shortest travel time caused the trigger, only considering channel 0

station: int

the station that should be considered

veff_aeff: string

specifiy if the effective volume or the effective area for surface muons is calculated can be

  • “veff” (default)

  • “aeff_surface_muons”

bounds_theta: list of floats

restrict theta to sub-range wrt. the simulated range in the file Note: assumes events were simulated uniformly in cos(theta) bounds_theta should be a (two-item) list, but will care only about the min/max values

Returns:
list of dictionaries

Each file is one entry. The dictionary keys store all relevant properties

NuRadioMC.utilities.Veff.get_Veff_Aeff_single_wrapper(args)[source]
NuRadioMC.utilities.Veff.get_Veff_Aeff(folder, trigger_combinations={}, station=101, veff_aeff='veff', n_cores=1, oversampling_theta=1)[source]

calculates the effective volume or effective area from surface muons from NuRadioMC hdf5 files

the effective volume is NOT normalized to a water equivalent. It is also NOT multiplied with the solid angle (typically 4pi).

Parameters:
folder: string

folder conaining the hdf5 files, one per energy OR filename

trigger_combinations: dict, optional

keys are the names of triggers to calculate. Values are dicts again:

  • ‘triggers’: list of strings

    name of individual triggers that are combined with an OR

the following additional options are optional

  • ‘efficiency’: dict

    allows to apply an (analysis) efficiency cut for calculating effective volumes

    • ‘func’: function

      a function that paramaterized the efficiency as a function of SNR (=Vmax/Vrms)

    • ‘channel_ids’: array on ints

      the channels for which the maximum signal amplitude should be determined

    • ‘scale’: float

      rescaling of the efficiency curve by SNR’ = SNR * scale

  • ‘n_reflections’: int

    the number of bottom reflections of the ray tracing solution that likely triggered assuming that the solution with the shortest travel time caused the trigger, only considering channel 0

station: int

the station that should be considered

veff_aeff: string

specifiy if the effective volume or the effective area for surface muons is calculated can be

  • “veff” (default)

  • “aeff_surface_muons”

n_cores: int

the number of cores to use

oversampling_theta: int

calculate the effective volume for finer binning (<oversampling_theta> data points per file):

  • 1: no oversampling

  • >1: oversampling with <oversampling_theta> equal-size cos(theta) bins within thetamin/max of the input file

Note

oversampling assumes that events were simulated uniformly in cos(theta)

Returns:
list of dictionaries.

Each file is one entry. The dictionary keys store all relevant properties

NuRadioMC.utilities.Veff.get_Veff_Aeff_array(data)[source]

calculates a multi dimensional array of effective volume or effective area from surface muons calculations for fast slicing

the array dimensions are (energy, zenith bin, triggername, 5) where the last tuple is the effective volume, its uncertainty, the weighted sum of triggered events, lower 68% uncertainty, upper 68% uncertainty

Parameters:
data: dict

the result of the get_Veff function

Returns:
  • (n_energy, n_zenith_bins, n_triggernames, 5) dimensional array of floats
  • array of unique mean energies (the mean is calculated in the logarithm of the energy)
  • array of unique lower bin edges of energies
  • array of unique upper bin edges of energies
  • array of unique zenith bins
  • array of unique trigger names

Examples

To plot the full sky effective volume for ‘all_triggers’ do

output, uenergies, uzenith_bins, utrigger_names, zenith_weights = get_Veff_array(data)

fig, ax = plt.subplots(1, 1)
tname = "all_triggers"
Veff = np.average(output[:,:,get_index(tname, utrigger_names),0], axis=1, weights=zenith_weights)
Vefferror = Veff / np.sum(output[:,:,get_index(tname, utrigger_names),2], axis=1)**0.5
ax.errorbar(uenergies/units.eV, Veff/units.km**3 * 4 * np.pi, yerr=Vefferror/units.km**3 * 4 * np.pi, fmt='-o', label=tname)

ax.legend()
ax.semilogy(True)
ax.semilogx(True)
fig.tight_layout()
plt.show()

To plot the effective volume for different declination bands do

fig, ax = plt.subplots(1, 1)
tname = "LPDA_2of4_100Hz"
iZ = 9
Veff = output[:,iZ,get_index(tname, utrigger_names)]
ax.errorbar(uenergies/units.eV, Veff[:,0]/units.km**3, yerr=Veff[:,1]/units.km**3,
            label=f"zenith bin {uzenith_bins[iZ][0]/units.deg:.0f} - {uzenith_bins[iZ][1]/units.deg:.0f}")

iZ = 8
Veff = output[:,iZ,get_index(tname, utrigger_names)]
ax.errorbar(uenergies/units.eV, Veff[:,0]/units.km**3, yerr=Veff[:,1]/units.km**3,
            label=f"zenith bin {uzenith_bins[iZ][0]/units.deg:.0f} - {uzenith_bins[iZ][1]/units.deg:.0f}")
iZ = 7
Veff = output[:,iZ,get_index(tname, utrigger_names)]
ax.errorbar(uenergies/units.eV, Veff[:,0]/units.km**3, yerr=Veff[:,1]/units.km**3,
            label=f"zenith bin {uzenith_bins[iZ][0]/units.deg:.0f} - {uzenith_bins[iZ][1]/units.deg:.0f}")
iZ = 10
Veff = output[:,iZ,get_index(tname, utrigger_names)]
ax.errorbar(uenergies/units.eV, Veff[:,0]/units.km**3, yerr=Veff[:,1]/units.km**3,
            label=f"zenith bin {uzenith_bins[iZ][0]/units.deg:.0f} - {uzenith_bins[iZ][1]/units.deg:.0f}")

ax.legend()
ax.semilogy(True)
ax.semilogx(True)
fig.tight_layout()
plt.show()
NuRadioMC.utilities.Veff.get_index(value, array)[source]
NuRadioMC.utilities.Veff.export(filename, data, trigger_names=None, export_format='yaml')[source]

export effective volumes (or effective areas) into a human readable JSON or YAML file

Parameters:
filename: string

the output filename of the JSON file

data: array

the output of the getVeff function

trigger_names: list of strings (optional, default None)

save only specific trigger names, if None all triggers are exported

export_format: string (default “yaml”)

specify output format, choose

  • “yaml”

  • “json”