NuRadioReco.utilities package

Submodules

NuRadioReco.utilities.analytic_pulse module

NuRadioReco.utilities.analytic_pulse.amp_from_energy(energy)

energy is defined as the integral of squared voltage normalized to a time window of 128 ns

energy:

NuRadioReco.utilities.analytic_pulse.get_analytic_pulse(amp_p0, amp_p1, phase_p0, n_samples_time, sampling_rate, phase_p1=0, bandpass=None, quadratic_term=0, quadratic_term_offset=0)

Analytic pulse as described in PhD thesis Glaser and NuRadioReco paper in the time domain

amp_p0: float

amplitude parameter of analytic pulse

amp_p1:

slope parameter of analytic pulse

phase_p0:

phase parameter of analytic pulse

n_samples_time:

numer of samples in time-domain

sampling_rate:

sampling rate of trace

phase_p1:

default 0

bandpass:

default None

quadratic_term:

default 0

quadratic_term_offset:

default 0

NuRadioReco.utilities.analytic_pulse.get_analytic_pulse_freq(amp_p0, amp_p1, phase_p0, n_samples_time, sampling_rate, phase_p1=0, bandpass=None, quadratic_term=0, quadratic_term_offset=0)

Analytic pulse as described in PhD thesis Glaser and NuRadioReco paper in the frequency domain

amp_p0: float

amplitude parameter of analytic pulse

amp_p1:

slope parameter of analytic pulse

phase_p0:

phase parameter of analytic pulse

n_samples_time:

numer of samples in time-domain

sampling_rate:

sampling rate of trace

phase_p1:

default 0

bandpass:

default None

quadratic_term:

default 0

quadratic_term_offset:

default 0

NuRadioReco.utilities.fft module

NuRadioReco.utilities.fft.freq2time(spectrum, sampling_rate, n=None)

performs backward FFT with correct normalization that conserves the power

spectrum: complex np array

the frequency spectrum

sampling_rate: float

sampling rate of the spectrum

n: int

the number of sample in the time domain (relevant if time trace has an odd number of samples)

NuRadioReco.utilities.fft.time2freq(trace, sampling_rate)

performs forward FFT with correct normalization that conserves the power

trace: np.array

time trace to be transformed into frequency space

sampling_rate: float

sampling rate of the trace

NuRadioReco.utilities.geometryUtilities module

NuRadioReco.utilities.geometryUtilities.get_efield_in_spherical_coords(efield, theta, phi)

Get 3D electric field from cartesian coordinates in spherical coordinates, using the arrival directions theta and phi

NuRadioReco.utilities.geometryUtilities.get_fresnel_angle(zenith_incoming, n_2=1.3, n_1=1.0)

Apply Snell’s law for given zenith angle, when a signal travels from n1 to n2

NuRadioReco.utilities.geometryUtilities.get_fresnel_r_p(zenith_incoming, n_2=1.3, n_1=1.0)

returns the coefficient r which is the ratio of the reflected wave’s electric field amplitude to that of the incident wave for parallel polarization (p-wave) this polarization corresponds to the eTheta polarization

parallel and perpendicular refers to the signal’s polarization with respect to the ‘plane of incident’ which is defindes as: “the plane of incidence is the plane which contains the surface normal and the propagation vector of the incoming radiation.”

NuRadioReco.utilities.geometryUtilities.get_fresnel_r_s(zenith_incoming, n_2=1.3, n_1=1.0)

returns the coefficient r which is the ratio of the reflected wave’s electric field amplitude to that of the incident wave for perpendicular polarization (s-wave) this polarization corresponds to the ePhi polarization

parallel and perpendicular refers to the signal’s polarization with respect to the ‘plane of incident’ which is defindes as: “the plane of incidence is the plane which contains the surface normal and the propagation vector of the incoming radiation.”

NuRadioReco.utilities.geometryUtilities.get_fresnel_t_p(zenith_incoming, n_2=1.3, n_1=1.0)

returns the coefficient t which is the ratio of the transmitted wave’s electric field amplitude to that of the incident wave for parallel polarization (p-wave) this polarization corresponds to the eTheta polarization

parallel and perpendicular refers to the signal’s polarization with respect to the ‘plane of incident’ which is defindes as: “the plane of incidence is the plane which contains the surface normal and the propagation vector of the incoming radiation.”

NuRadioReco.utilities.geometryUtilities.get_fresnel_t_s(zenith_incoming, n_2=1.3, n_1=1.0)

returns the coefficient t which is the ratio of the transmitted wave’s electric field amplitude to that of the incident wave for perpendicular polarization (s-wave) this polarization corresponds to the ePhi polarization

parallel and perpendicular refers to the signal’s polarization with respect to the ‘plane of incident’ which is defindes as: “the plane of incidence is the plane which contains the surface normal and the propagation vector of the incoming radiation.”

NuRadioReco.utilities.geometryUtilities.get_time_delay_from_direction(zenith, azimuth, positions, n=1.000293)

Calculate the time delay between given positions for an arrival direction

zenith: float [rad]

Zenith angle in convention up = 0

azimuth: float [rad]

Azimuth angle in convention East = 0, counter-clock-wise

positions: array[N x 3]

Positions on ground

n: float (default: 1.000293)

Index of reflection of propagation medium. By default, air is assumed

NuRadioReco.utilities.geometryUtilities.rot_x(angle)

Angle helper function

NuRadioReco.utilities.geometryUtilities.rot_y(angle)

Angle helper function

NuRadioReco.utilities.geometryUtilities.rot_z(angle)

Angle helper function

NuRadioReco.utilities.ice module

implementation of ice models

NuRadioReco.utilities.ice.get_refractive_index(depth, site='southpole')

NuRadioReco.utilities.io_utilities module

NuRadioReco.utilities.io_utilities.read_pickle(filename, encoding='latin1')

Read in a pickle file and return the result This utility is supposed to provide compatibility for pickles created with different python versions. If a simple pickle.load fails, it will try to load the file with a specific encoding.

filename: string

Name of the pickle file to be opened

encoding: string

Encoding to be used if the first attempt to open the pickle fails

NuRadioReco.utilities.noise module

class NuRadioReco.utilities.noise.thermalNoiseGenerator(n_samples, sampling_rate, Vrms, threshold, time_coincidence, n_majority, time_coincidence_majority, n_channels, trigger_time, filt, noise_type='rayleigh')

Bases: object

generate_noise()

generates noise traces for all channels that will cause a high/low majority logic trigger

Returns np.array of shape (n_channels, n_samples)

NuRadioReco.utilities.templates module

class NuRadioReco.utilities.templates.Templates(path)

Bases: object

get_cr_ref_template(station_id)

returns one cosmic ray template for the reference direction, the same for all channels

get_cr_ref_templates(station_id)

returns one cosmic ray template per channel for the reference direction

get_nu_ref_template(station_id)

returns one neutrino template for the reference direction and on the cherenkov cone, the same for all channels

get_nu_ref_templates(station_id)

returns one neutrino template per channel for the reference direction and on the cherenkov cone

get_set_of_cr_templates(station_id, n=100)

gets set of n templates to allow for the calculation of average templates

loops first over different coreas pulses with different frequency content and then over azimuth angles of 0, 22.5 and 45 degree and then over zenith angles of 60, 50 and 70 degree

get_set_of_cr_templates_full(station_id, n=100)

gets set of n templates to allow for the calculation of average templates

get_set_of_nu_templates(station_id, n=100)

gets set of n templates to allow for the calculation of average templates

loops first over different viewing angles and then over azimuth angles of 0, 22.5 and 45 degree and then over zenith angles of 100, 120 and 140 degree

set_template_directory(path)

NuRadioReco.utilities.timing module

NuRadioReco.utilities.timing.analyze_timing(module_list, t_tot=None)

NuRadioReco.utilities.trace_utilities module

NuRadioReco.utilities.trace_utilities.apply_butterworth(spectrum, frequencies, passband, order=8)

Calculates the response from a Butterworth filter and applies it to the input spectrum

spectrum: array of complex

Fourier spectrum to be filtere

frequencies: array of floats

Frequencies of the input spectrum

passband: (float, float) tuple

Tuple indicating the cutoff frequencies

order: integer

Filter order

filtered_spectrum: array of complex

The filtered spectrum

NuRadioReco.utilities.trace_utilities.butterworth_filter_trace(trace, sampling_frequency, passband, order=8)

Filters a trace using a Butterworth filter.

trace: array of floats

Trace to be filtered

sampling_frequency: float

Sampling frequency

passband: (float, float) tuple

Tuple indicating the cutoff frequencies

order: integer

Filter order

filtered_trace: array of floats

The filtered trace

NuRadioReco.utilities.trace_utilities.delay_trace(trace, sampling_frequency, time_delay, delayed_samples)

Delays a trace by transforming it to frequency and multiplying by phases. Since this method is cyclic, the trace has to be cropped. It only accepts positive delays, so some samples from the beginning are thrown away and then some samples from the end so that the total number of samples is equal to the argument delayed samples.

trace: array of floats

Array containing the trace

sampling_frequency: float

Sampling rate for the trace

time_delay: float

Time delay used for transforming the trace. Must be positive or 0

delayed_samples: integer

Number of samples that the delayed trace must contain

delayed_trace: array of floats

The delayed, cropped trace

NuRadioReco.utilities.trace_utilities.get_channel_voltage_from_efield(station, electric_field, channels, detector, zenith, azimuth, antenna_pattern_provider, return_spectrum=True)

Returns the voltage traces that would result in the channels from the station’s E-field.

station: Station electric_field: ElectricField channels: array of int

IDs of the channels for which the expected voltages should be calculated

detector: Detector zenith, azimuth: float

incoming direction of the signal. Note that reflection and refraction at the air/ice boundary are already being taken into account.

antenna_pattern_provider: AntennaPatternProvider return_spectrum: boolean

if True, returns the spectrum, if False return the time trace

NuRadioReco.utilities.trace_utilities.get_efield_antenna_factor(station, frequencies, channels, detector, zenith, azimuth, antenna_pattern_provider)

Returns the antenna response to a radio signal coming from a specific direction

station: Station frequencies: array of complex

frequencies of the radio signal for which the antenna response is needed

channels: array of int

IDs of the channels

detector: Detector zenith, azimuth: float, float

incoming direction of the signal. Note that refraction and reflection at the ice/air boundary are taken into account

antenna_pattern_provider: AntennaPatternProvider

NuRadioReco.utilities.trace_utilities.get_electric_field_energy_fluence(electric_field_trace, times, signal_window_mask=None, noise_window_mask=None)
NuRadioReco.utilities.trace_utilities.upsampling_fir(trace, original_sampling_frequency, int_factor=2, ntaps=128)

This function performs an upsampling by inserting a number of zeroes between samples and then applying a finite impulse response (FIR) filter.

trace: array of floats

Trace to be upsampled

original_sampling_frequency: float

Sampling frequency of the input trace

int_factor: integer

Upsampling factor. The resulting trace will have a sampling frequency int_factor times higher than the original one

ntaps: integer

Number of taps (order) of the FIR filter

upsampled_trace: array of floats

The upsampled trace

NuRadioReco.utilities.units module

standard system of units

You should use the units defined in this file whenever you have a dimensional quantity in your code. For example, write:

s = 1.5 * units.km

instead of:

s = 1.5   # don't forget this is in km!

The conversion factors defined in this file convert your data into Auger base units, so that all dimensional quantities in the code are in a single system of units! You can also use the conversions defined here to, for example, display data with the unit of your choice. For example:

print "s = " , s / units.mm, " mm"

The base units are:

  • meter (meter)

  • nanosecond (nanosecond)

  • electron Volt (eV)

  • positron charge (eplus)

  • degree Kelvin (kelvin)

  • the amount of substance (mole)

  • luminous intensity (candela)

  • radian (radian)

  • steradian (steradian)

The SI numerical value of the positron charge is defined here, as it is needed for conversion factor : positron charge = eSI (coulomb)

Adapted from Offline, the reconstruction Framework of the Pierre Auger Collaboration. This is a slightly modified version of the units definitions written by the Geant4 collaboration

Module contents