km3io.tools

Module Contents

Classes

cached_property

A simple cache decorator for properties.

TimeConverter

Auxiliary class to convert Monte Carlo hit times to DAQ/triggered hit times.

Functions

unfold_indices(obj, indices)

Unfolds an index chain and returns the corresponding item

to_num(value)

Convert a value to a numerical one if possible

unique(array[, dtype])

Return the unique elements of an array with a given dtype.

uniquecount(array[, dtype])

Count the number of unique elements in a jagged Awkward1 array.

get_w2list_param(events, generator, param)

Get all the values of a specific parameter from the w2list

fitinf(fitparam, tracks)

Access fit parameters in tracks.fitinf.

count_nested(arr[, axis])

Count elements in a nested awkward Array.

get_multiplicity(tracks, rec_stages)

Tracks selection based on specific reconstruction stages.

best_track(tracks[, startend, minmax, stages])

Best track selection.

mask(arr[, sequence, startend, minmax, atleast])

Return a boolean mask which mask each nested sub-array for a condition.

has_jmuon(tracks)

Check if given tracks contain JMUON reconstruction.

has_jshower(tracks)

Check if given tracks contain JSHOWER reconstruction.

has_aashower(tracks)

Check if given tracks contain AASHOWER reconstruction.

has_dusjshower(tracks)

Check if given tracks contain AASHOWER reconstruction.

best_jmuon(tracks)

Select the best JMUON track.

best_jshower(tracks)

Select the best JSHOWER track.

best_aashower(tracks)

Select the best AASHOWER track.

best_dusjshower(tracks)

Select the best DISJSHOWER track.

is_cc(fobj)

Determin if events are a result of a Charged Curent interaction (CC)

usr(objects, field)

Return the usr-data for a given field.

is_bit_set(value, bit_position)

Returns true if a bit at the given position is 1.

is_3dshower(trigger_mask)

Returns True if the trigger mask contains the 3D shower flag.

is_mxshower(trigger_mask)

Returns True if the trigger mask contains the MX shower flag.

is_3dmuon(trigger_mask)

Returns True if the trigger mask contains the 3D muon flag.

get_w2list_idx(f)

Get the correct w2list_idx for the given file, or None if there is none.

is_nanobeacon(trigger_mask)

Returns True if the trigger mask contains the nano-beacon flag.

angle(v1, v2[, normalized])

Compute the unsigned angle between two vectors. For a stacked input, the

magnitude(v)

Calculates the magnitude of a vector or array of vectors.

theta(v)

Neutrino direction in polar coordinates.

theta_separg(dir_z)

phi(v)

Neutrino direction in polar coordinates.

phi_separg(dir_x, dir_y)

zenith(v)

Return the zenith angle in radians.

azimuth(v)

Return the azimuth angle in radians.

angle_between(v1, v2[, axis])

Returns the angle in radians between vectors 'v1' and 'v2'.

unit_vector(vector, **kwargs)

Returns the unit vector of the vector.

Attributes

BASKET_CACHE_SIZE

BASKET_CACHE

km3io.tools.BASKET_CACHE_SIZE[source]
km3io.tools.BASKET_CACHE[source]
class km3io.tools.cached_property(function)[source]

A simple cache decorator for properties.

km3io.tools.unfold_indices(obj, indices)[source]

Unfolds an index chain and returns the corresponding item

km3io.tools.to_num(value)[source]

Convert a value to a numerical one if possible

km3io.tools.unique(array, dtype=np.int64)[source]

Return the unique elements of an array with a given dtype.

The performance is better for pre-sorted input arrays.

km3io.tools.uniquecount(array, dtype=np.int64)[source]

Count the number of unique elements in a jagged Awkward1 array.

km3io.tools.get_w2list_param(events, generator, param)[source]

Get all the values of a specific parameter from the w2list in offline neutrino files.

Parameters:
eventskm3io.offline.OfflineBranch

events class in offline neutrino files.

generatorstr

the name of the software generating neutrinos, it is either ‘genhen’ or ‘gseagen’.

paramstr

the name of the parameters found in w2list as defined in the KM3NeT-Dataformat for both genhen and gseagen.

Returns:
awkward.Array

array of the values of interest.

km3io.tools.fitinf(fitparam, tracks)[source]

Access fit parameters in tracks.fitinf.

Parameters:
fitparamint

the fit parameter key according to fitparameters defined in KM3NeT-Dataformat (see km3io.definitions.fitparameters).

tracksak.Array or km3io.rootio.Branch

reconstructed tracks with .fitinf attribute

Returns:
awkward1.Array

awkward array of the values of the fit parameter requested. Missing values are set to NaN.

km3io.tools.count_nested(arr, axis=0)[source]

Count elements in a nested awkward Array.

Parameters:
arrawkward1.Array

Array of data. Example tracks.fitinf or tracks.rec_stages.

axisint, optional

axis = 0: to count elements in the outmost level of nesting. axis = 1: to count elements in the first level of nesting. axis = 2: to count elements in the second level of nesting.

Returns:
awkward1.Array or int

counts of elements found in a nested awkward1 Array.

km3io.tools.get_multiplicity(tracks, rec_stages)[source]

Tracks selection based on specific reconstruction stages.

Counts how many tracks with the specific reconstructions stages are found per event.

Parameters:
trackskm3io.offline.OfflineBranch

tracks or a subste of tracks.

rec_stageslist

Reconstruction stages (the ordering is respected) e.g. [1, 2, 3, 4, 5].

Returns:
awkward1.Array

tracks multiplicty.

km3io.tools.best_track(tracks, startend=None, minmax=None, stages=None)[source]

Best track selection.

Parameters:
tracksawkward.Array

A list of tracks or doubly nested tracks, usually from OfflineReader.events.tracks or subarrays of that, containing recunstructed tracks.

startend: tuple(int, int), optional

The required first and last stage in tracks.rec_stages.

minmax: tuple(int, int), optional

The range (minimum and maximum) value of rec_stages to take into account.

stageslist or set, optional
  • list: the order of the rec_stages is respected.

  • set: a subset of required stages; the order is irrelevant.

Returns:
awkward.Array or namedtuple

Be aware that the dimensions are kept, which means that the final track attributes are nested when multiple events are passed in. If a single event (just a list of tracks) is provided, a named tuple with a single track and flat attributes is created.

Raises:
ValueError

When invalid inputs are specified.

km3io.tools.mask(arr, sequence=None, startend=None, minmax=None, atleast=None)[source]

Return a boolean mask which mask each nested sub-array for a condition.

Parameters:
arrawkward.Array with ndim>=2

The array to mask.

startend: tuple(int, int), optional

True for entries where the first and last element are matching the tuple.

minmax: tuple(int, int), optional

True for entries where each element is within the min-max-range.

sequencelist(int), optional

True for entries which contain the exact same elements (in that specific order)

atleastlist(int), optional

True for entries where at least the provided elements are present.

An extensive discussion about this implementation can be found at:
https://github.com/scikit-hep/awkward-1.0/issues/580
Many thanks for Jim for the fruitful discussion and the final implementation.
km3io.tools.has_jmuon(tracks)[source]

Check if given tracks contain JMUON reconstruction.

km3io.tools.has_jshower(tracks)[source]

Check if given tracks contain JSHOWER reconstruction.

km3io.tools.has_aashower(tracks)[source]

Check if given tracks contain AASHOWER reconstruction.

km3io.tools.has_dusjshower(tracks)[source]

Check if given tracks contain AASHOWER reconstruction.

km3io.tools.best_jmuon(tracks)[source]

Select the best JMUON track.

km3io.tools.best_jshower(tracks)[source]

Select the best JSHOWER track.

km3io.tools.best_aashower(tracks)[source]

Select the best AASHOWER track.

km3io.tools.best_dusjshower(tracks)[source]

Select the best DISJSHOWER track.

km3io.tools.is_cc(fobj)[source]

Determin if events are a result of a Charged Curent interaction (CC) or a Neutral Curent interaction (NC).

Parameters:
fobjkm3io.offline.OfflineReader

offline neutrino file object.

Returns:
awkward1.highlevel.Array

a mask where True corresponds to an event resulting from a CC interaction.

Raises:
ValueError

if the simulations program used to generate the neutrino file is neither gseagen nor genhen.

km3io.tools.usr(objects, field)[source]

Return the usr-data for a given field.

Parameters:
objectsawkward.Array

Events, tracks, hits or whatever objects which have usr and usr_names fields (e.g. OfflineReader().events).

km3io.tools.is_bit_set(value, bit_position)[source]

Returns true if a bit at the given position is 1.

value: int or array([u]int64)

The value to check, can be a single value or an array of values.

bit_position: int

0 for the first position, 1 for the second etc.

km3io.tools.is_3dshower(trigger_mask)[source]

Returns True if the trigger mask contains the 3D shower flag.

Parameters:
trigger_maskint or array(int)

A value or an array of the trigger_mask, either of an event, or a hit.

km3io.tools.is_mxshower(trigger_mask)[source]

Returns True if the trigger mask contains the MX shower flag.

Parameters:
trigger_maskint or array(int)

A value or an array of the trigger_mask, either of an event, or a hit.

km3io.tools.is_3dmuon(trigger_mask)[source]

Returns True if the trigger mask contains the 3D muon flag.

Parameters:
trigger_maskint or array(int)

A value or an array of the trigger_mask, either of an event, or a hit.

km3io.tools.get_w2list_idx(f)[source]

Get the correct w2list_idx for the given file, or None if there is none.

Parameters:
fkm3io.OfflineReader

The file.

km3io.tools.is_nanobeacon(trigger_mask)[source]

Returns True if the trigger mask contains the nano-beacon flag.

Parameters:
trigger_maskint or array(int)

A value or an array of the trigger_mask, either of an event, or a hit.

class km3io.tools.TimeConverter(event)[source]

Auxiliary class to convert Monte Carlo hit times to DAQ/triggered hit times.

FRAME_TIME_NS = 100000000.0[source]
get_time_of_frame(frame_index)[source]

Get start time of frame in ns since start of run for a given frame index

Parameters:
frame_index: int

The index of the DAQ frame

get_DAQ_time(t0)[source]

Get DAQ/triggered hit time

Parameters:
t0: float or array(float)

Simulated time [ns]

get_MC_time(t0)[source]

Get Monte Carlo hit time

Parameters:
t0: float or array(float)

DAQ/trigger hit time [ns]

km3io.tools.angle(v1, v2, normalized=False)[source]

Compute the unsigned angle between two vectors. For a stacked input, the angle is computed pairwise. Inspired by the “vg” Python package.

Parameters:
v1np.array

With shape (3,) or a kx3 stack of vectors.

v2np.array

A vector or stack of vectors with the same shape as v1.

normalizedbool

By default, the vectors will be normalised unless normalized is True.

Returns:
A float or a vector of floats with the angle in radians.
km3io.tools.magnitude(v)[source]

Calculates the magnitude of a vector or array of vectors.

Parameters:
vnp.array

With shape (3,) or kx3 stack of vectors.

Returns:
A float or a vector of floats with the magnitudes.
km3io.tools.theta(v)[source]

Neutrino direction in polar coordinates.

Parameters:
varray (x, y, z)

Notes

Downgoing event: theta = 180deg Horizont: 90deg Upgoing: theta = 0

Angles in radians.

km3io.tools.theta_separg(dir_z)[source]
km3io.tools.phi(v)[source]

Neutrino direction in polar coordinates.

Parameters:
varray (x, y, z)

Notes

phi, theta is the opposite of zenith, azimuth.

Angles in radians.

km3io.tools.phi_separg(dir_x, dir_y)[source]
km3io.tools.zenith(v)[source]

Return the zenith angle in radians.

Parameters:
varray (x, y, z)

Notes

Defined as ‘Angle respective to downgoing’. Downgoing event: zenith = 0 Horizont: 90deg Upgoing: zenith = 180deg

km3io.tools.azimuth(v)[source]

Return the azimuth angle in radians.

Parameters:
varray (x, y, z)

Notes

phi, theta is the opposite of zenith, azimuth.

This is the ‘normal’ azimuth definition – beware of how you define your coordinates. KM3NeT defines azimuth differently than e.g. SLALIB, astropy, the AAS.org

km3io.tools.angle_between(v1, v2, axis=0)[source]

Returns the angle in radians between vectors ‘v1’ and ‘v2’.

If axis=1 it evaluates the angle between arrays of vectors.

Examples

>>> angle_between((1, 0, 0), (0, 1, 0))
1.5707963267948966
>>> angle_between((1, 0, 0), (1, 0, 0))
0.0
>>> angle_between((1, 0, 0), (-1, 0, 0))
3.141592653589793
>>> v1 = np.array([[1, 0, 0], [1, 0, 0]])
>>> v2 = np.array([[0, 1, 0], [-1, 0, 0]])
>>> angle_between(v1, v2, axis=1)
array([1.57079633, 3.14159265])
km3io.tools.unit_vector(vector, **kwargs)[source]

Returns the unit vector of the vector.