reboost.hpge package¶
Submodules¶
reboost.hpge.psd module¶
- reboost.hpge.psd._drift_time_heuristic_impl(dt, edep)¶
Low-level implementation of the HPGe drift-time-based pulse-shape heuristic.
Accepts Awkward arrays and uses Numba to speed up the computation.
For each hit (collection of steps), the drift times and corresponding energies are sorted in ascending order. The function finds the optimal split point \(m\) that maximizes the identification metric:
\[I = \frac{|T_1 - T_2|}{E_\text{s}(E_1, E_2)}\]where:
\[T_1 = \frac{\sum_{i < m} t_i E_i}{\sum_{i < m} E_i} \quad \text{and} \quad T_2 = \frac{\sum_{i \geq m} t_i E_i}{\sum_{i \geq m} E_i}\]are the energy-weighted mean drift times of the two groups.
\[E_\text{scale}(E_1, E_2) = \frac{1}{\sqrt{E_1 E_2}}\]is the scaling factor.
The function iterates over all possible values of \(m\) and selects the maximum I as the drift time heuristic value.
- Parameters:
dt (ak.Array)
edep (ak.Array)
- Return type:
NDArray
- reboost.hpge.psd.drift_time(xloc, yloc, zloc, dt_map, coord_offset=<Quantity([0 0 0], 'meter')>)¶
Calculates drift times for each step (cluster) in an HPGe detector.
- Parameters:
xloc (_SupportsArray[dtype[Any]] | _NestedSequence[_SupportsArray[dtype[Any]]] | bool | int | float | complex | str | bytes | _NestedSequence[bool | int | float | complex | str | bytes]) – awkward array of x coordinate position.
yloc (_SupportsArray[dtype[Any]] | _NestedSequence[_SupportsArray[dtype[Any]]] | bool | int | float | complex | str | bytes | _NestedSequence[bool | int | float | complex | str | bytes]) – awkward array of y coordinate position.
zloc (_SupportsArray[dtype[Any]] | _NestedSequence[_SupportsArray[dtype[Any]]] | bool | int | float | complex | str | bytes | _NestedSequence[bool | int | float | complex | str | bytes]) – awkward array of z coordinate position.
dt_map (HPGeScalarRZField) – the drift time map.
coord_offset (Quantity | Position) – this (x, y, z) coordinates will be subtracted to (xloc, yloc, zloc)` before drift time computation. The length units must be the same as xloc, yloc and zloc.
- Return type:
- reboost.hpge.psd.drift_time_heuristic(drift_time, edep)¶
HPGe drift-time-based pulse-shape heuristic.
See
_drift_time_heuristic_impl()for a description of the algorithm.- Parameters:
drift_time (_SupportsArray[dtype[Any]] | _NestedSequence[_SupportsArray[dtype[Any]]] | bool | int | float | complex | str | bytes | _NestedSequence[bool | int | float | complex | str | bytes]) – drift time of charges originating from steps/clusters. Can be calculated with
drift_time().edep (_SupportsArray[dtype[Any]] | _NestedSequence[_SupportsArray[dtype[Any]]] | bool | int | float | complex | str | bytes | _NestedSequence[bool | int | float | complex | str | bytes]) – energy deposited in step/cluster (same shape as drift_time).
- Return type:
reboost.hpge.surface module¶
reboost.hpge.utils module¶
- class reboost.hpge.utils.HPGeScalarRZField(φ, r_units, z_units, φ_units)¶
Bases:
NamedTupleA scalar field defined in the cylindrical-like (r, z) HPGe plane.
Create new instance of HPGeScalarRZField(φ, r_units, z_units, φ_units)
- _asdict()¶
Return a new dict which maps field names to their values.
- _field_defaults = {}¶
- _fields = ('φ', 'r_units', 'z_units', 'φ_units')¶
- classmethod _make(iterable)¶
Make a new HPGeScalarRZField object from a sequence or iterable
- _replace(**kwds)¶
Return a new HPGeScalarRZField object replacing specified fields with new values
- reboost.hpge.utils.get_hpge_scalar_rz_field(filename, obj, field, out_of_bounds_val=nan, **kwargs)¶
Create an interpolator for a gridded scalar HPGe field defined on (r, z).
Reads from disk the following data structure:
FILENAME/ └── OBJ · struct{r,z,FIELD} ├── r · array<1>{real} ── {'units': 'UNITS'} ├── z · array<1>{real} ── {'units': 'UNITS'} └── FIELD · array<2>{real} ── {'units': 'UNITS'}where
FILENAME,OBJandFIELDare provided as arguments to this function. obj is aStruct, r and z are one dimensional arrays specifying the radial and z coordinates of the rectangular grid — not the coordinates of each single grid point. In this coordinate system, the center of the p+ contact surface is at (0, 0), with the p+ contact facing downwards. field is instead a two-dimensional array specifying the field value at each grid point. The first and second dimensions are r and z, respectively. NaN values are interpreted as points outside the detector profile in the (r, z) plane.Before returning a
HPGeScalarRZField, the gridded field is fed toscipy.interpolate.RegularGridInterpolator.- Parameters:
- Return type: