dspeed.processors package#
Contains a list of DSP processors, implemented using Numba’s
numba.guvectorize()
to implement NumPy’s numpy.ufunc
interface.
In other words, all of the functions are void functions whose outputs are given
as parameters. The ufunc
interface provides additional
information about the function signatures that enables broadcasting the arrays
and SIMD processing. Thanks to the ufunc
interface, they can
also be called to return a NumPy array, but if this is done, memory will be
allocated for this array, slowing things down.
The dspeed processors use the ufunc
framework, which is
designed to encourage highly performant python practices. These functions have
several advantages:
They work with
numpy.array
. NumPy arrays are arrays of same-typed objects that are stored adjacently in memory (equivalent to a dynamically allocated C-array or a C++ vector). Compared to standard Python lists, they perform computations much faster. They are ideal for representing waveforms and transformed waveforms.They perform vectorized operations. Vectorization causes commands to perform the same operation on all components of an array. For example, the
ufunc
np.add(a, b, out=c)
(equivalentlyc=a+b
) is equivalent to:for i in range(len(c)): c[i] = a[i] + b[i]
Loops are slow in python since it is an interpreted language; vectorized commands remove the loop and only call the Python interpreter once.
Furthermore,
ufunc
s are capable of broadcasting their dimensions. This involves a safety check to ensure the dimensions ofa
andb
are compatible sizes. It also will automatically replicate asize-1
dimension over another one, enabling the addition of a scalar to a vector quantity. This is useful, as it allows us to process multiple waveforms at once.One of the biggest advantages of vectorized
ufunc
s is that many of them will take advantage of SIMD (same input-multiple data) vectorization on a vector-CPU. Modern CPUs typically have 256- or 512-bit wide processing units, which can accommodate multiple 32- or 64-bit numbers. Programming with these, however, is quite difficult and requires specialized commands to be called. Luckily for us, many NumPyufunc
s will automatically use these for us, speeding up our code!ufunc
s are capable of calculating their output in place, meaning they can place calculated values in pre-allocated memory rather than allocating and returning new values. This is important because memory allocation is one of the slowest processes computers can perform, and should be avoided. Withufunc
s, this can be done using the out keyword in arguments (exnumpy.add(a, b, out=c)
, or more succinctly,numpy.add(a, b, c)
). While this may seem counterintuitive at first, the alternative (c = np.add(a,b)
orc = a+b
) causes an entirely new array to be allocated, with c pointing at that. These array allocations can add up very quickly:e = a*b + c*d
, for example, would allocate 3 different arrays: one fora*b
, one forc*d
, and one for the sum of those two. As we writeufunc
s, it is important that we try to use functions that operate in place as much as possible!
Submodules#
dspeed.processors.bl_subtract module#
- @numba.guvectorize dspeed.processors.bl_subtract.bl_subtract(w_in: ndarray, a_baseline: float, w_out: ndarray) None #
Options:
boundscheck=False
,cache=True
Precompiled signatures:
ff->f
,dd->d
Subtract the constant baseline from the entire waveform.
- Parameters:
JSON Configuration Example
"wf_bl": { "function": "bl_subtract", "module": "dspeed.processors", "args": ["waveform", "baseline", "wf_bl"], "unit": "ADC" }
dspeed.processors.convolutions module#
- @numba.guvectorize dspeed.processors.convolutions.convolve_wf(w_in: ndarray, kernel: array, mode_in: int8, w_out: ndarray) None #
Options:
boundscheck=False
,forceobj=True
,cache=True
Precompiled signatures:
ffbf->
,ddbd->
- Parameters:
w_in (ndarray) – the input waveform.
kernel (array) – the kernel to convolve with
mode – mode of convolution options are f : full, v : valid or s : same, explained here: https://numpy.org/doc/stable/reference/generated/numpy.convolve.html
w_out (ndarray) – the filtered waveform.
- @numba.guvectorize dspeed.processors.convolutions.fft_convolve_wf(w_in: ndarray, kernel: array, mode_in: int8, w_out: ndarray) None #
Options:
boundscheck=False
,forceobj=True
,cache=True
Precompiled signatures:
ffbf->
,ddbd->
- Parameters:
w_in (ndarray) – the input waveform.
kernel (array) – the kernel to convolve with
mode – mode of convolution options are f : full, v : valid or s : same, explained here: https://numpy.org/doc/stable/reference/generated/numpy.convolve.html
w_out (ndarray) – the filtered waveform.
dspeed.processors.dwt module#
- @numba.guvectorize dspeed.processors.dwt.discrete_wavelet_transform(w_in: ndarray, level: int, wave_type: int, coeff: int, w_out: ndarray) None #
Options:
boundscheck=False
,forceobj=True
,cache=True
Precompiled signatures:
fibbf->
,dlbbd->
Apply a discrete wavelet transform to the waveform and return only the detailed or approximate coefficients.
- Parameters:
w_in (ndarray) – The input waveform
level (int) – The level of decompositions to be performed
(1, 2, ...)
wave_type (int) – The wavelet type for discrete convolution
('h' = 'haar', 'd' = 'db1')
.coeff (int) – The coefficients to be saved
('a', 'd')
w_out (ndarray) – The approximate coefficients. The dimension of this array can be calculated by
out_dim = len(w_in)/(filter_length^level)
, wherefilter_length
can be obtained viapywt.Wavelet(wave_type).dec_len
JSON Configuration Example
"dwt_haar":{ "function": "discrete_wavelet_transform", "module": "dspeed.processors", "args": ["wf_blsub", 5, "'h'", "'a'", "dwt_haar(256, 'f')"], "unit": "ADC", "prereqs": ["wf_blsub"], }
dspeed.processors.energy_kernels module#
- @numba.guvectorize dspeed.processors.energy_kernels.cusp_filter(sigma: float, flat: int, decay: int, kernel: array) None #
Options:
boundscheck=False
,forceobj=True
,cache=True
Precompiled signatures:
ffff->
,dddd->
Calculates CUSP kernel.
- Parameters:
JSON Configuration Example
"kern_cusp": { "function": "cusp_filter", "module": "dspeed.processors", "args": ["10*us", "3*us", "400*us", "kern_cusp"], "unit": "ADC" }
- @numba.guvectorize dspeed.processors.energy_kernels.dplms(noise_mat: list, reference: list, a1: float, a2: float, a3: float, ff: int, kernel: array) None #
Options:
boundscheck=False
,forceobj=True
,cache=True
Precompiled signatures:
fffffff->
,ddddddd->
Calculate and apply an optimum DPLMS filter to the waveform.
The processor takes the noise matrix and the reference signal as input and calculates the optimum filter according to the provided length and penalized coefficients [DPLMS].
[DPLMS]V. D’Andrea et al. “Optimum Filter Synthesis with DPLMS Method for Energy Reconstruction” Eur. Phys. J. C 83, 149 (2023). https://doi.org/10.1140/epjc/s10052-023-11299-z
- Parameters:
noise_mat (list) – noise matrix
reference (list) – reference signal
a1 (float) – penalized coefficient for the noise matrix
a2 (float) – penalized coefficient for the reference matrix
a3 (float) – penalized coefficient for the zero area matrix
ff (int) – flat top length for the reference signal
kernel (array) – output kernel
JSON Configuration Example
"kern_dplms": { "function": "dplms", "module": "dspeed.processors", "args": ["db.dplms.noise_matrix", "db.dplms.reference", "50", "0.1", "1", "1" "kern_dplms"], "unit": "ADC", }
- @numba.guvectorize dspeed.processors.energy_kernels.zac_filter(sigma: float, flat: int, decay: int, kernel: array) None #
Options:
boundscheck=False
,forceobj=True
,cache=True
Precompiled signatures:
ffff->
,dddd->
Calculates ZAC (Zero Area CUSP) kernel.
- Parameters:
JSON Configuration Example
"kern_zac": { "function": "zac_filter", "module": "dspeed.processors", "args": ["10*us", "3*us", "400*us", "kern_zac"], "unit": "ADC" }
dspeed.processors.fftw module#
- dspeed.processors.fftw.dft(w_in: ndarray | ProcChainVarBase, w_out: ndarray | ProcChainVarBase) Callable #
Perform discrete Fourier transforms using the FFTW library.
- Parameters:
w_in (ndarray | ProcChainVarBase) – the input waveform.
w_out (ndarray | ProcChainVarBase) – the output fourier transform.
- Return type:
JSON Configuration Example
"wf_dft": { "function": "dft", "module": "dspeed.processors", "args": ["wf", "wf_dft"], "init_args": ["wf", "wf_dft"] }
Note
FFTW optimizes the FFT algorithm based on the size of the arrays, with SIMD parallelized commands. This optimization requires initialization, so this is a factory function that returns a Numba gufunc that performs the FFT. FFTW works on fixed memory buffers, so you must tell it what memory to use ahead of time. When using this with
ProcessingChain
, the output waveform’s size, dtype and coordinate grid units can be set automatically. The possible dtypes for the input/outputs are:Size
Size
float32
/float
\(n\)
complex64
\(n/2+1\)
float64
/double
\(n\)
complex128
\(n/2+1\)
float128
/longdouble
\(n\)
complex256
/clongdouble
\(n/2+1\)
complex64
\(n\)
complex64
\(n\)
complex128
\(n\)
complex128
\(n\)
complex256
/clongdouble
\(n\)
complex256
/clongdouble
\(n\)
- dspeed.processors.fftw.inv_dft(w_in: ndarray, w_out: ndarray) Callable #
Perform inverse discrete Fourier transforms using the FFTW library.
- Parameters:
- Return type:
JSON Configuration Example
"wf_invdft": { "function": "inv_dft", "module": "dspeed.processors", "args": ["wf_dft", "wf_invdft"], "init_args": ["wf_dft", "wf_invdft"] }
Note
FFTW optimizes the FFT algorithm based on the size of the arrays, with SIMD parallelized commands. This optimization requires initialization, so this is a factory function that returns a Numba gufunc that performs the FFT. FFTW works on fixed memory buffers, so you must tell it what memory to use ahead of time. When using this with
ProcessingChain
, the output waveform’s size, dtype and coordinate grid units can be set automatically. The automated behavior will produce a real output by default, unless you specify a complex output. Possible dtypes for the input/outputs are:Size
Size
complex64
\(n/2+1\)
float32
/float
\(n\)
complex128
\(n/2+1\)
float64
/double
\(n\)
complex256
/clongdouble
\(n/2+1\)
float128
/longdouble
\(n\)
complex64
\(n\)
complex64
\(n\)
complex128
\(n\)
complex128
\(n\)
complex256
/clongdouble
\(n\)
complex256
/clongdouble
\(n\)
- dspeed.processors.fftw.psd(w_in: ndarray, w_out: ndarray) Callable #
Perform discrete Fourier transforms using the FFTW library, and use it to get the power spectral density.
- Parameters:
- Return type:
JSON Configuration Example
"wf_psd": { "function": "psd", "module": "dspeed.processors", "args": ["wf", "wf_psd"], "init_args": ["wf", "wf_psd"] }
Note
FFTW optimizes the FFT algorithm based on the size of the arrays, with SIMD parallelized commands. This optimization requires initialization, so this is a factory function that returns a Numba gufunc that performs the FFT. FFTW works on fixed memory buffers, so you must tell it what memory to use ahead of time. When using this with
ProcessingChain
, the output waveform’s size, dtype and coordinate grid units can be set automatically. The possible dtypes for the input/outputs are:Size
Size
complex64
\(n\)
float32
/float
\(n\)
complex128
\(n\)
float64
/double
\(n\)
complex256
/clongdouble
\(n\)
float128
/longdouble
\(n\)
float32
/float
\(n\)
float32
/float
\(n/2+1\)
float64
/double
\(n\)
float64
/double
\(n/2+1\)
float128
/longdouble
\(n\)
float128
/longdouble
\(n/2+1\)
dspeed.processors.fixed_time_pickoff module#
- @numba.guvectorize dspeed.processors.fixed_time_pickoff.fixed_time_pickoff(w_in: ndarray, t_in: float, mode_in: int8, a_out: float)#
Options:
boundscheck=False
,cache=True
Precompiled signatures:
ffb->f
,ddb->d
Pick off the waveform value at the provided time.
For non-integral times, interpolate between samples using the method selected using mode_in. If the provided index t_in is out of range, return
numpy.nan
.- Parameters:
w_in (ndarray) – the input waveform.
t_in (float) – the waveform index to pick off.
mode_in (int8) –
character selecting which interpolation method to use. Note this must be passed as a
int8
, e.g.ord('i')
. Options:i
– integer t_in; equivalent tofixed_sample_pickoff()
n
– nearest-neighbor interpolation; defined at all values, but not continuousf
– floor, or value at previous neighbor; defined at all values but not continuousc
– ceiling, or value at next neighbor; defined at all values, but not continuousl
– linear interpolation; continuous at all values, but not differentiableh
– Hermite cubic spline interpolation; continuous and differentiable at all values but not twice-differentiables
– natural cubic spline interpolation; continuous and twice-differentiable at all values. This method is much slower than the others because it utilizes the entire input waveform!
a_out (float) – the output pick-off value.
Examples
"trapEftp": { "function": "fixed_time_pickoff", "module": "dspeed.processors", "args": ["wf_trap", "tp_0+10*us", "'h'", "trapEftp"], "unit": "ADC", }
dspeed.processors.gaussian_filter1d module#
dspeed.processors.get_multi_local_extrema module#
- @numba.guvectorize dspeed.processors.get_multi_local_extrema.get_multi_local_extrema(w_in: ndarray, a_delta_max_in: float, a_delta_min_in: float, search_direction: int, a_abs_max_in: float, a_abs_min_in: float, vt_max_out: ndarray, vt_min_out: ndarray, n_max_out: int, n_min_out: int) None #
Options:
boundscheck=False
,cache=True
Precompiled signatures:
ffffffffII->
,ddddddddII->
Get lists of indices of the local maxima and minima of data.
Converted from a MATLAB script by E. Billauer. See the parameters description for details about search modes and further settings.
- Parameters:
w_in (ndarray) – the array of data within which extrema will be found.
a_delta_min_in (float) – the relative level by which data must vary (in the search direction) about a candidate minimum in order for it to be tagged.
a_delta_max_in (float) – the relative level by which data must vary (in the search direction) about a candidate maximum in order for it to be tagged.
search_direction (int) –
the direction in which the input waveform is processed.
0
: one sweep, first to last sample1
: one sweep, last to first sample2
: two sweeps, in both directions. Largest common set of found extrema is returned (logical AND)3
: two sweeps, in both directions. Union of found extrema is returned (logical OR)
a_abs_min_in (float) – a candidate minimum is tagged only if its value is smaller than this.
a_abs_max_in (float) – a candidate maximum is tagged only if its value is larger than this.
vt_min_out (ndarray) – arrays of fixed length (padded with
numpy.nan
) that hold the indices of the identified minima.vt_max_out (ndarray) – arrays of fixed length (padded with
numpy.nan
) that hold the indices of the identified maxima.n_min_out (int) – the number of minima found in a waveform.
n_max_out (int) – the number of maxima found in a waveform.
JSON Configuration Example
"vt_max_out, vt_min_out, n_max_out, n_min_out": { "function": "get_multi_local_extrema", "module": "dspeed.processors", "args": [ "waveform", 5, 0.1, 1, 20, 0, "vt_max_out(20)", "vt_min_out(20)", "n_max_out", "n_min_out" ], "unit": ["ns", "ns", "none", "none"] }
dspeed.processors.get_wf_centroid module#
- @numba.guvectorize dspeed.processors.get_wf_centroid.get_wf_centroid(w_in: ndarray, shift: int, centroid: int) None #
Options:
boundscheck=False
,cache=True
Precompiled signatures:
ff->f
,dd->d
Calculate waveform centroid.
Note
This processor calculate the centroid position when provided the convolution product with a step function.
- Parameters:
JSON Configuration Example
"centroid": { "function": "get_wf_centroid", "module": "dspeed.processors", "args": ["waveform", "shift", "centroid"], "unit": "ADC" }
dspeed.processors.histogram module#
- @numba.guvectorize dspeed.processors.histogram.histogram(w_in: ndarray, weights_out: ndarray, borders_out: ndarray) None #
Options:
boundscheck=False
,cache=True
Precompiled signatures:
fff->
,ddd->
Produces and returns an histogram of the waveform.
- Parameters:
JSON Configuration Example
"hist_weights, hist_borders": { "function": "histogram", "module": "dspeed.processors.histogram", "args": ["waveform", "hist_weights(100)", "hist_borders(101)"], "unit": ["none", "ADC"] }
Note
This implementation is significantly faster than just wrapping
numpy.histogram()
.See also
histogram_stats
- @numba.guvectorize dspeed.processors.histogram.histogram_stats(weights_in: ndarray, edges_in: ndarray, mode_out: int, max_out: float, fwhm_out: float, max_in: float) None #
Options:
boundscheck=False
,cache=True
Precompiled signatures:
ffffff->
,dddddd->
Compute useful histogram-related quantities.
- Parameters:
weights_in (ndarray) – histogram weights.
edges_in (ndarray) – histogram bin edges.
max_in (float) – if not
numpy.nan
, this value is used as the histogram bin content at the mode. Otherwise the mode is computed automatically (from left to right).mode_out (int) – the computed mode of the histogram. If max_in is not
numpy.nan
then the closest waveform index to max_in is returned.max_out (float) – the histogram bin content at the mode.
fwhm_out (float) – the FWHM of the histogram, calculated by starting from the mode and descending left and right.
JSON Configuration Example
"fwhm, idx_out, max_out": { "function": "histogram_stats", "module": "dspeed.processors.histogram", "args": ["hist_weights","hist_borders","idx_out","max_out","fwhm","np.nan"], "unit": ["ADC", "none", "ADC"] }
See also
dspeed.processors.kernels module#
- @numba.guvectorize dspeed.processors.kernels.moving_slope(kernel)#
Options:
boundscheck=False
,forceobj=True
,cache=True
Precompiled signatures:
f->
,d->
Calculates the linear slope of kernel
- Parameters:
length – the length of the section to calculate slope
kernel – the output kernel
JSON Configuration Example
"kern_slopes": { "function": "moving_slope", "module": "dspeed.processors", "args": ["12", "kern_slopes"], "unit": "ADC" }
- @numba.guvectorize dspeed.processors.kernels.step(weight_pos: int, kernel: array) None #
Options:
boundscheck=False
,forceobj=True
,cache=True
Precompiled signatures:
ff->
,dd->
Process waveforms with a step function.
- Parameters:
weight_pos (int) – relative weight of positive step side.
kernel (array) – output kernel
JSON Configuration Example
"kern_step": { "function": "step", "module": "dspeed.processors", "args": ["16", "kern_step"], "unit": "ADC" }
- @numba.guvectorize dspeed.processors.kernels.t0_filter(rise: int, fall: int, kernel: array) None #
Options:
boundscheck=False
,forceobj=True
,cache=True
Precompiled signatures:
fff->
,ddd->
Apply a modified, asymmetric trapezoidal filter to the waveform.
- Parameters:
JSON Configuration Example
"t0_filter": { "function": "t0_filter", "module": "dspeed.processors", "args": ["128*ns", "2*us", "t0_filter"], "unit": "ADC", "init_args": ["128*ns", "2*us"] }
dspeed.processors.linear_slope_fit module#
- @numba.guvectorize dspeed.processors.linear_slope_fit.linear_slope_diff(w_in: ndarray, slope: float, intercept: float, mean: float, rms: float) None #
Options:
boundscheck=False
,cache=True
Precompiled signatures:
fff->ff
,ddd->dd
Calculate the mean, standard deviation, slope and y-intercept of waveform.
Uses Welford’s method and linear regression.
- Parameters:
JSON Configuration Example
"bl_mean, bl_std, bl_slope, bl_intercept": { "function": "linear_slope_fit", "module": "dspeed.processors", "args": ["wf_blsub[0:1650]", "bl_mean", "bl_std", "bl_slope", "bl_intercept"], "unit": ["ADC", "ADC", "ADC", "ADC"], }
- @numba.guvectorize dspeed.processors.linear_slope_fit.linear_slope_fit(w_in: ndarray, mean: float, stdev: float, slope: float, intercept: float) None #
Options:
boundscheck=False
,cache=True
Precompiled signatures:
f->ffff
,d->dddd
Calculate the mean and standard deviation of the waveform using Welford’s method as well as the slope an intercept of the waveform using linear regression.
- Parameters:
JSON Configuration Example
"bl_mean, bl_std, bl_slope, bl_intercept": { "function": "linear_slope_fit", "module": "dspeed.processors", "args": ["wf_blsub[0:1650]", "bl_mean", "bl_std", "bl_slope", "bl_intercept"], "unit": ["ADC", "ADC", "ADC", "ADC"], }
dspeed.processors.log_check module#
- @numba.guvectorize dspeed.processors.log_check.log_check(w_in: ndarray, w_log: ndarray) None #
Options:
boundscheck=False
,cache=True
Precompiled signatures:
f->f
,d->d
Calculate the logarithm of the waveform if all its values are positive; otherwise, return NaN.
- Parameters:
JSON Configuration Example
"wf_logged": { "function": "log_check", "module": "dspeed.processors", "args": ["wf_blsub[2100:]", "wf_logged"], "unit": "ADC" }
dspeed.processors.min_max module#
- @numba.guvectorize dspeed.processors.min_max.min_max(w_in: ndarray, t_min: int, t_max: int, a_min: float, a_max: float) None #
Options:
boundscheck=False
,cache=True
Precompiled signatures:
f->ffff
,d->dddd
Find the value and index of the minimum and maximum values in the waveform.
Note
The first found instance of each extremum in the waveform will be returned.
- Parameters:
JSON Configuration Example
"tp_min, tp_max, wf_min, wf_max": { "function": "min_max", "module": "dspeed.processors", "args": ["waveform", "tp_min", "tp_max", "wf_min", "wf_max"], "unit": ["ns", "ns", "ADC", "ADC"] }
- @numba.guvectorize dspeed.processors.min_max.min_max_norm(w_in: ndarray, a_min: float, a_max: float, w_out: ndarray) None #
Options:
boundscheck=False
,cache=True
Precompiled signatures:
fff->f
,ddd->d
Normalize waveform by minimum or maximum value, whichever is greater in absolute value.
- Parameters:
JSON Configuration Example
"wf_norm": { "function": "min_max_norm", "module": "dspeed.processors", "args": ["wf_blsub", "wf_min", "wf_max", "wf_norm"], "unit": ["ADC"] }
dspeed.processors.moving_windows module#
- @numba.guvectorize dspeed.processors.moving_windows.avg_current(w_in: ndarray, length: float, w_out: ndarray) None #
Options:
boundscheck=False
,cache=True
Precompiled signatures:
fff->
,ddd->
Calculate the derivative of a waveform, averaged across length samples.
- Parameters:
JSON Configuration Example
"curr": { "function": "avg_current", "module": "dspeed.processors", "args": ["wf_pz", 1, "curr(len(wf_pz)-1, f)"], "unit": "ADC/sample" }
- @numba.guvectorize dspeed.processors.moving_windows.moving_window_left(w_in: ndarray, length: float, w_out: ndarray) None #
Options:
boundscheck=False
,cache=True
Precompiled signatures:
ff->f
,dd->d
Applies a moving average window to the waveform.
Note
Starts from the left and assumes that the baseline is at zero.
- Parameters:
JSON Configuration Example
"wf_mw": { "function": "moving_window_left", "module": "dspeed.processors", "args": ["wf_pz", "96*ns", "wf_mw"], "unit": "ADC" }
- @numba.guvectorize dspeed.processors.moving_windows.moving_window_multi(w_in: ndarray, length: float, num_mw: int, mw_type: int, w_out: ndarray) None #
Options:
boundscheck=False
,cache=True
Precompiled signatures:
fffi->f
,dddi->d
Apply a series of moving-average windows to the waveform, alternating its application between the left and the right.
- Parameters:
JSON Configuration Example
"curr_av": { "function": "moving_window_multi", "module": "dspeed.processors", "args": ["curr", "96*ns", "3", "0", "curr_av"], "unit": "ADC/sample" }
- @numba.guvectorize dspeed.processors.moving_windows.moving_window_right(w_in: ndarray, length: float, w_out: ndarray) None #
Options:
boundscheck=False
,cache=True
Precompiled signatures:
ff->f
,dd->d
Applies a moving average window to the waveform from the right.
- Parameters:
JSON Configuration Example
"wf_mw": { "function": "moving_window_right", "module": "dspeed.processors", "args": ["wf_pz", "96*ns", "wf_mw"], "unit": "ADC" }
dspeed.processors.multi_a_filter module#
- @numba.guvectorize dspeed.processors.multi_a_filter.multi_a_filter(w_in, vt_maxs_in, va_max_out)#
Options:
boundscheck=False
,forceobj=True
,cache=True
Precompiled signatures:
fff->
,ddd->
Finds the maximums in a waveform and returns the amplitude of the wave at those points.
- Parameters:
w_in – the array of data within which amplitudes of extrema will be found.
vt_maxs_in – the array of max positions for each waveform.
va_max_out – an array (in-place filled) of the amplitudes of the maximums of the waveform.
dspeed.processors.multi_t_filter module#
- @numba.guvectorize dspeed.processors.multi_t_filter.multi_t_filter(w_in: ndarray, a_threshold_in: float, vt_max_in: ndarray, vt_min_in: ndarray, t_out: ndarray) None #
Options:
boundscheck=False
,forceobj=True
,cache=True
Precompiled signatures:
fffff->
,ddddd->
Gets list of indices of the start of leading edges of multiple peaks within a waveform.
Is built to handle afterpulses/delayed cross-talk and trains of pulses. Works by calling the vectorized functions
get_multi_local_extrema()
which returns a list of the maxima and minima in a waveform, and then the list of maxima is fed intotime_point_thresh()
which returns the final times that waveform is less than a specified threshold.- Parameters:
w_in (ndarray) – the array of data within which the list of leading edge times will be found.
a_threshold_in (float) – threshold to search for using
time_point_thresh()
.vt_max_in (ndarray) – the array of maximum positions for each waveform.
vt_min_in (ndarray) – the array of minimum positions for each waveform.
t_out (ndarray) – output array of fixed length (padded with
numpy.nan
) that hold the indices of the identified initial rise times of peaks in the signal.
See also
time_point_thresh
- @numba.guvectorize dspeed.processors.multi_t_filter.remove_duplicates(t_in: ndarray, vt_min_in: ndarray, t_out: ndarray) None #
Options:
boundscheck=False
,forceobj=True
,cache=True
Precompiled signatures:
ff->f
,dd->d
Helper function to remove duplicate peak positions.
time_point_thresh()
has issues with afterpulsing in waveforms that causes an aferpulse peak position to be sent to zero or the same index as the first pulse. This only happens when the relative minimum between the first pulse and the afterpulse is greater than the threshold. So, we sweep through the array again to ensure there are no duplicate indices. If there are duplicate indices caused by a misidentified position of an afterpulse, we replace its index by that of the corresponding minimum found usingget_multi_local_extrema()
. It also checks to make sure that the maximum of a waveform isn’t right at index 0.- Parameters:
See also
dspeed.processors.optimize module#
- class dspeed.processors.optimize.Model(func: Callable, w_in: ndarray, baseline: float, beg: int, end: int)#
Bases:
object
A class representing the function to minimize.
- errordef = 1.0#
- @numba.guvectorize dspeed.processors.optimize.optimize_1pz(w_in: ndarray, a_baseline_in: float, t_beg_in: int, t_end_in: int, p0_in: float, val0_out: float) None #
Options:
boundscheck=False
,forceobj=True
,cache=True
Precompiled signatures:
fffff->f
,ddddd->d
Find the optimal, single pole-zero cancellation’s parameter by minimizing the slope in the waveform’s specified time range.
- Parameters:
w_in (ndarray) – the input waveform.
a_baseline_in (float) – the resting baseline.
t_beg_in (int) – the lower bound’s index for the time range over which to optimize the pole-zero cancellation.
t_end_in (int) – the upper bound’s index for the time range over which to optimize the pole-zero cancellation.
p0_in (float) – the initial guess of the optimal time constant.
val0_out (float) – the output value of the best-fit time constant.
JSON Configuration Example
"tau0": { "function": "optimize_1pz", "module": "dspeed.processors", "args": ["waveform", "baseline", "0", "20*us", "500*us", "tau0"], "unit": "us" }
See also
pole_zero
,double_pole_zero
- @numba.guvectorize dspeed.processors.optimize.optimize_2pz(w_in: ndarray, a_baseline_in: float, t_beg_in: int, t_end_in: int, p0_in: float, p1_in: float, p2_in: float, val0_out: float, val1_out: float, val2_out: float) None #
Options:
boundscheck=False
,forceobj=True
,cache=True
Precompiled signatures:
fffffff->fff
,ddddddd->ddd
Find the optimal, double pole-zero cancellation’s parameters by minimizing the slope in the waveform’s specified time range.
- Parameters:
w_in (ndarray) – the input waveform.
a_baseline_in (float) – the resting baseline.
t_beg_in (int) – the lower bound’s index for the time range over which to optimize the pole-zero cancellation.
t_end_in (int) – the upper bound’s index for the time range over which to optimize the pole-zero cancellation.
p0_in (float) – the initial guess of the optimal, longer time constant.
p1_in (float) – the initial guess of the optimal, shorter time constant.
p2_in (float) – the initial guess of the optimal fraction.
val0_out (float) – the output value of the best-fit, longer time constant.
val1_out (float) – the output value of the best-fit, shorter time constant.
val2_out (float) – the output value of the best-fit fraction.
JSON Configuration Example
"tau1, tau2, frac": { "function": "optimize_2pz", "module": "dspeed.processors", "args": [ "waveform", "baseline", "0", "20*us", "500*us", "20*us", "0.02", "tau1", "tau2", "frac" ], "unit": "us" }
dspeed.processors.param_lookup module#
- dspeed.processors.param_lookup.param_lookup(param_dict: dict[int, float], default_val: float, dtype: str | dtype) ufunc #
Generate the
numpy.ufunc
lookup(channel, val)
, which returns a NumPy array of values corresponding to various channels that are looked up in the provided param_dict. If there is no key, use default_val instead.- Return type:
dspeed.processors.peak_snr_threshold module#
- @numba.guvectorize dspeed.processors.peak_snr_threshold.peak_snr_threshold(w_in: ndarray, idx_in: ndarray, ratio_in: float, width_in: int, idx_out: ndarray, n_idx_out: int) None #
Options:
boundscheck=False
,cache=True
Precompiled signatures:
ffff->fI
,dddd->dI
Search for local minima in a window around the provided indices.
If a minimum is found it is checked whether amplitude of the minimum divided by the amplitude of the waveform at the index is smaller then the given ratio. If this is the case the index is copied to the output.
- Parameters:
w_in (ndarray) – the input waveform.
idx_in (ndarray) – the array of indices of possible signal candidates.
ratio_in (float) – ratio threshold value.
width_in (int) – width of the local search window.
idx_out (ndarray) – indices of local minima.
n_idx_out (int) – number of non-
numpy.nan
indices in idx_out.
dspeed.processors.pole_zero module#
- @numba.guvectorize dspeed.processors.pole_zero.double_pole_zero(w_in: ndarray, t_tau1: float, t_tau2: float, frac: float, w_out: ndarray) ndarray #
Options:
boundscheck=False
,cache=True
Precompiled signatures:
ffff->f
,dddd->d
Apply a double pole-zero cancellation using the provided time constants to the waveform.
- Parameters:
w_in (ndarray) – the input waveform.
t_tau1 (float) – the time constant of the first exponential to be deconvolved.
t_tau2 (float) – the time constant of the second exponential to be deconvolved.
frac (float) – the fraction which the second exponential contributes.
w_out (ndarray) – the pole-zero cancelled waveform.
- Return type:
JSON Configuration Example
"wf_pz": { "function": "double_pole_zero", "module": "dspeed.processors", "args": ["wf_bl", "400*us", "20*us", "0.02", "wf_pz"], "unit": "ADC" }
Notes
This algorithm is an IIR filter to deconvolve the function
\[s(t) = A \left[ f \cdot \exp\left(-\frac{t}{\tau_2} \right) + (1-f) \cdot \exp\left(-\frac{t}{\tau_1}\right) \right]\](\(f\) = frac) into a single step function of amplitude \(A\). This filter is derived by \(z\)-transforming the input (\(s(t)\)) and output (step function) signals, and then determining the transfer function. For shorthand, define \(a=\exp(-1/\tau_1)\) and \(b=\exp(-1/\tau_2)\), the transfer function is then:
\[H(z) = \frac{1 - (a+b)z^{-1} + abz^{-2}} {1 + (fb - fa - b - 1)z^{-1}-(fb - fa - b)z^{-2}}\]By equating the transfer function to the ratio of output to input waveforms \(H(z) = w_\text{out}(z) / w_\text{in}(z)\) and then taking the inverse \(z\)-transform yields the filter as run in the function, where \(f\) is the frac parameter:
\[\begin{split}w_\text{out}[n] =& w_\text{in}[n] - (a+b)w_\text{in}[n-1] + abw_\text{in}[n-2] \\ & -(fb - fa - b - 1)w_\text{out}[n-1] + (fb - fa - b)w_\text{out}[n-2]\end{split}\]
- @numba.guvectorize dspeed.processors.pole_zero.pole_zero(w_in: ndarray, t_tau: float, w_out: ndarray) None #
Options:
boundscheck=False
,cache=True
Precompiled signatures:
ff->f
,dd->d
Apply a pole-zero cancellation using the provided time constant to the waveform.
- Parameters:
JSON Configuration Example
"wf_pz": { "function": "pole_zero", "module": "dspeed.processors", "args": ["wf_bl", "400*us", "wf_pz"], "unit": "ADC" }
dspeed.processors.presum module#
- @numba.guvectorize dspeed.processors.presum.presum(w_in: ndarray, do_norm: int, ps_fact: int, w_out: ndarray) None #
Options:
boundscheck=False
,cache=True
Precompiled signatures:
ffff->
,dddd->
Presum the waveform.
Combine bins in chunks of
len(w_in) / len(w_out)
, which is hopefully an integer. If it isn’t, then some samples at the end will be omitted.
dspeed.processors.pulse_injector module#
- @numba.guvectorize dspeed.processors.pulse_injector.inject_exp_pulse(wf_in: ndarray, t0: int, rt: float, a: float, decay: float, wf_out: ndarray) None #
Options:
boundscheck=False
,cache=True
Precompiled signatures:
fffff->f
,ddddd->d
Inject exponential pulse into existing waveform to simulate pileup.
- Parameters:
wf_in (ndarray) – the input waveform.
t0 (int) – the position of the injected waveform.
rt (float) – the rise time of the injected waveform.
a (float) – the amplitude of the injected waveform.
decay (float) – the exponential decay constant of the injected waveform.
wf_out (ndarray) – the output waveform.
- @numba.guvectorize dspeed.processors.pulse_injector.inject_sig_pulse(wf_in: ndarray, t0: int, rt: float, a: float, decay: float, wf_out: ndarray) None #
Options:
boundscheck=False
,cache=True
Precompiled signatures:
fffff->f
,ddddd->d
Inject sigmoid pulse into existing waveform to simulate pileup.
\[s(t) = \frac{A}{1 + \exp[-4 \log(99) (t - t_0 - t_r/2) / t_r]} e^{-(t-t_0)/\tau}\]- Parameters:
wf_in (ndarray) – the input waveform.
t0 (int) – the position \(t_0\) of the injected waveform.
rt (float) – the rise time \(t_r\) of the injected waveform.
a (float) – the amplitude \(A\) of the injected waveform.
decay (float) – the decay parameter \(\tau\) of the injected waveform.
wf_out (ndarray) – the output waveform.
dspeed.processors.round_to_nearest module#
- @numba.guvectorize dspeed.processors.round_to_nearest.round_to_nearest(val: ndarray, to_nearest: int | float, out: ndarray) None #
Options:
nopython=True
,boundscheck=False
,cache=True
Precompiled signatures:
BB->B
,HH->H
,II->I
,LL->L
,bb->b
,hh->h
,ii->i
,ll->l
,ff->f
,dd->d
Round value to nearest multiple of to_nearest.
- Parameters:
JSON Configuration Example
Note: this processor is aliased using round in ProcessingChain. The following two examples are equivalent.
"t_rounded": { "function": "round_to_nearest", "module": "dspeed.processors", "args": ["t_in", "1*us", "t_rounded"] "unit": ["ns"] }, "t_rounded": "round(t_in, 1*us)"
dspeed.processors.saturation module#
- @numba.guvectorize dspeed.processors.saturation.saturation(w_in: ndarray, bit_depth_in: int, n_lo_out: int, n_hi_out: int) None #
Options:
boundscheck=False
,cache=True
Precompiled signatures:
ff->ff
,dd->dd
Count the number of samples in the waveform that are saturated at the minimum and maximum possible values based on the bit depth.
- Parameters:
JSON Configuration Example
"sat_lo, sat_hi": { "function": "saturation", "module": "dspeed.processors", "args": ["waveform", "16", "sat_lo", "sat_hi"], "unit": "ADC" }
dspeed.processors.soft_pileup_corr module#
- @numba.guvectorize dspeed.processors.soft_pileup_corr.soft_pileup_corr(w_in: ndarray, n_in: int, tau_in: float, w_out: ndarray) None #
Options:
boundscheck=False
,cache=True
Precompiled signatures:
fff->f
,ddd->d
Fit the baseline to an exponential with the provided time constant and then subtract the best-fit function from the entire waveform.
- Parameters:
JSON Configuration Example
"wf_bl": { "function": "soft_pileup_corr", "module": "dspeed.processors", "args": ["waveform", "1000", "500*us", "wf_bl"], "unit": "ADC" }
- @numba.guvectorize dspeed.processors.soft_pileup_corr.soft_pileup_corr_bl(w_in: ndarray, n_in: int, tau_in: float, b_in: float, w_out: ndarray) None #
Options:
boundscheck=False
,cache=True
Precompiled signatures:
ffff->f
,dddd->d
Fit the baseline to an exponential with the provided time constant and then subtract the best-fit function from the entire waveform.
- Parameters:
w_in (ndarray) – the input waveform.
n_in (int) – the number of samples at the beginning of the waveform to fit to an exponential.
tau_in (float) – the fixed, exponential time constant.
b_in (float) – the fixed, exponential constant.
w_out (ndarray) – the output waveform with the exponential subtracted.
JSON Configuration Example
"wf_bl": { "function": "soft_pileup_corr_bl", "module": "dspeed.processors", "args": ["waveform", "1000", "500*us", "baseline", "wf_bl"], "unit": "ADC" }
dspeed.processors.svm module#
- dspeed.processors.svm.svm_predict(svm_file: str) Callable #
Apply a Support Vector Machine (SVM) to an input waveform to predict a data cleaning label.
Note
This processor is composed of a factory function that is called using the
init_args
argument. The input waveform and output label are passed usingargs
.- Parameters:
svm_file (str) – The name of the file with the trained SVM
svm_p*_r*_T***Z.sav
- Return type:
JSON Configuration Example
"svm_label":{ "function": "svm_predict", "module": "dspeed.processors", "args": ["dwt_norm", "svm_label"], "unit": "", "prereqs": ["dwt_norm"], "init_args": ["'svm_p*_r*_T***Z.sav'"] }
dspeed.processors.time_over_threshold module#
- @numba.guvectorize dspeed.processors.time_over_threshold.time_over_threshold(w_in: ndarray, a_threshold: float, n_samples: float) None #
Options:
boundscheck=False
,cache=True
Precompiled signatures:
ff->f
,dd->d
Calculates the number of samples in the input waveform over a_threshold
- Parameters:
JSON Configuration Example
"t_sat": { "function": "time_over_threshold", "module": "dspeed.processors", "args": ["wf_pz", "a_threshold", "t_sat"], "unit": "ns" }
dspeed.processors.time_point_thresh module#
- @numba.guvectorize dspeed.processors.time_point_thresh.interpolated_time_point_thresh(w_in: ndarray, a_threshold: float, t_start: int, walk_forward: int, mode_in: int8, t_out: float) None #
Options:
boundscheck=False
,cache=True
Precompiled signatures:
ffflb->f
,dddlb->d
Find the time where the waveform value crosses the threshold
Search performed walking either forward or backward from the starting index. Use interpolation to estimate a time between samples. Interpolation mode selected with mode_in.
- Parameters:
w_in (ndarray) – the input waveform.
a_threshold (float) – the threshold value.
t_start (int) – the starting index.
walk_forward (int) – the backward (
0
) or forward (1
) search direction.mode_in (int8) –
Character selecting which interpolation method to use. Note this must be passed as a
int8
, e.g.ord('i')
. Options:i
– integer t_in; equivalent tofixed_sample_pickoff()
f
– floor; interpolated values are at previous neighbor, so threshold crossing is at next neighborc
– ceiling, interpolated values are at next neighbor, so threshold crossing is at previous neighborn
– nearest-neighbor interpolation; threshold crossing is half-way between samplesl
– linear interpolationh
– Hermite cubic spline interpolation (not implemented)s
– natural cubic spline interpolation (not implemented)
t_out (float) – the index where the waveform value crosses the threshold.
JSON Configuration Example
"tp_0": { "function": "time_point_thresh", "module": "dspeed.processors", "args": ["wf_atrap", "bl_std", "tp_start", 0, "'l'", "tp_0"], "unit": "ns" }
- @numba.guvectorize dspeed.processors.time_point_thresh.time_point_thresh(w_in: ndarray, a_threshold: float, t_start: int, walk_forward: int, t_out: float) None #
Options:
boundscheck=False
,cache=True
Precompiled signatures:
ffff->f
,dddd->d
Find the index where the waveform value crosses the threshold, walking either forward or backward from the starting index.
- Parameters:
JSON Configuration Example
"tp_0": { "function": "time_point_thresh", "module": "dspeed.processors", "args": ["wf_atrap", "bl_std", "tp_start", 0, "tp_0"], "unit": "ns" }
dspeed.processors.transfer_function_convolver module#
- @numba.guvectorize dspeed.processors.transfer_function_convolver.transfer_function_convolver(w_in, a_in, b_in, w_out)#
Options:
boundscheck=False
,cache=True
Precompiled signatures:
fff->f
,ddd->d
Compute the difference equation of the form $a_0*w_out[i] = a_1*w_out[i-1] + ldots b_0*w_in[i] ldots$ which is equivalent of convolving a signal with a transfer function. a_in always needs at least one element The numba signature specifies these arrays as contiguous, so that way the dot product is as fast as possible
- Parameters:
w_in – An array of the values that this lti filter will be performed on
a_in – An array [a_0, …, a_n] of filter coefficients to apply recursively to the output
b_in – An array [b_0, …, b_m] of filter coefficients to apply to the input
w_out – The output array to store all the filtered values
JSON Configuration Example
"wf_pz": { "function": "transfer_function_convolver", "module": "dspeed.processors", "args": ["waveform", "np.array([-1,1])", "np.array([1,2])", "wf_pz"], "unit": "ADC" }
dspeed.processors.trap_filters module#
- @numba.guvectorize dspeed.processors.trap_filters.asym_trap_filter(w_in: ndarray, rise: int, flat: int, fall: int, w_out: ndarray) None #
Options:
boundscheck=False
,cache=True
Precompiled signatures:
fiii->f
,diii->d
Apply an asymmetric trapezoidal filter to the waveform, normalized by the number of samples averaged in the rise and fall sections.
- Parameters:
JSON Configuration Example
"wf_af": { "function": "asym_trap_filter", "module": "dspeed.processors", "args": ["wf_pz", "128*ns", "64*ns", "2*us", "wf_af"], "unit": "ADC" }
- @numba.guvectorize dspeed.processors.trap_filters.trap_filter(w_in: ndarray, rise: int, flat: int, w_out: ndarray) None #
Options:
boundscheck=False
,cache=True
Precompiled signatures:
fii->f
,dii->d
Apply a symmetric trapezoidal filter to the waveform.
- Parameters:
JSON Configuration Example
"wf_tf": { "function": "trap_filter", "module": "dspeed.processors", "args": ["wf_pz", "10*us", "3*us", "wf_tf"], "unit": "ADC" }
- @numba.guvectorize dspeed.processors.trap_filters.trap_norm(w_in: ndarray, rise: int, flat: int, w_out: ndarray) None #
Options:
boundscheck=False
,cache=True
Precompiled signatures:
fii->f
,dii->d
Apply a symmetric trapezoidal filter to the waveform, normalized by the number of samples averaged in the rise and fall sections.
- Parameters:
JSON Configuration Example
"wf_tf": { "function": "trap_norm", "module": "dspeed.processors", "args": ["wf_pz", "10*us", "3*us", "wf_tf"], "unit": "ADC" }
- @numba.guvectorize dspeed.processors.trap_filters.trap_pickoff(w_in: ndarray, rise: int, flat: int, t_pickoff: float, a_out: float) None #
Options:
boundscheck=False
,cache=True
Precompiled signatures:
fiif->f
,diid->d
Pick off the value at the provided index of a symmetric trapezoidal filter to the input waveform, normalized by the number of samples averaged in the rise and fall sections.
- Parameters:
JSON Configuration Example
"ct_corr": { "function": "trap_pickoff", "module": "dspeed.processors", "args": ["wf_pz", "1.5*us", 0, "tp_0", "ct_corr"], "unit": "ADC" }
dspeed.processors.upsampler module#
- @numba.guvectorize dspeed.processors.upsampler.interpolating_upsampler(w_in: ndarray, mode_in: int8, w_out: ndarray) None #
Options:
boundscheck=False
,cache=True
Precompiled signatures:
fbf->
,dbd->
Upsamples the waveform.
Resampling ratio is set by size of w_out and w_in. Using the interpolation technique specified by mode_in to fill newly added samples.
- Parameters:
w_in (ndarray) – waveform to upsample.
mode_in (int8) –
character selecting which interpolation method to use. Note this must be passed as a
int8
, e.g.ord('i')
. Options:i
– only set values at original samples; fill in between with zeros. Requires integer resampling ration
– nearest-neighbor interpolation; defined at all values, but not continuousf
– floor, or value at previous neighbor; defined at all values but not continuousc
– ceiling, or value at next neighbor; defined at all values, but not continuousl
– linear interpolation; continuous at all values, but not differentiableh
– Hermite cubic spline interpolation; continuous and differentiable at all values but not twice-differentiables
– natural cubic spline interpolation; continuous and twice- differentiable at all values. This method is much slower than the others because it utilizes the entire input waveform!
w_out (ndarray) – output array for upsampled waveform.
JSON Configuration Example
"wf_up": { "function": "interpolating_upsampler", "module": "dspeed.processors", "args": ["wf", "'s'", "wf_up(len(wf)*10, period=wf.period/10)"], "unit": "ADC" }
- @numba.guvectorize dspeed.processors.upsampler.upsampler(w_in: ndarray, upsample: float, w_out: ndarray) None #
Options:
boundscheck=False
,cache=True
Precompiled signatures:
fff->
,ddd->
Upsamples the waveform by the number specified.
Note
A series of moving windows should be applied afterwards for smoothing.
dspeed.processors.wf_alignment module#
- @numba.guvectorize dspeed.processors.wf_alignment.wf_alignment(w_in: ndarray, centroid: int, shift: int, size: int, w_out: ndarray) None #
Options:
boundscheck=False
,cache=True
Precompiled signatures:
fffff->
,ddddd->
Align waveform.
Note
This processor align the input waveform by setting the centroid position at the center of the output waveform.
- Parameters:
JSON Configuration Example
"wf_align": { "function": "wf_alignment", "module": "dspeed.processors", "args": ["waveform", "centroid", "shift", "size", "wf_align"], "unit": "ADC" }
dspeed.processors.wiener_filter module#
- dspeed.processors.wiener_filter.wiener_filter(file_name_array: list[str]) ndarray #
Apply a Wiener filter to the waveform.
Note
The convolution is performed in the frequency domain. This processor is composed of a factory function that is called using the init_args argument. The input and output waveforms are passed using args. The input must be the Fourier transform of the waveform. The output is the filtered waveform in the frequency domain.
- Parameters:
file_name_array (list[str]) – Array with path to an HDF5 file containing the time domain version of the superpulse in one column and noise waveform in another, the superpulse HDF5 group must be titled
spms/processed/superpulse
and the noise waveform must be calledspms/processed/noise_wf
.- Return type:
JSON Configuration Example
"wf_wiener": { "function": "wiener_filter", "module": "dspeed.processors", "args": ["wf_bl_fft", "wf_wiener(2000,f)"], "unit": "dB", "init_args": ["/path/to/file/wiener.lh5"] }
dspeed.processors.windower module#
- @numba.guvectorize dspeed.processors.windower.windower(w_in: ndarray, t0_in: int, w_out: ndarray) None #
Options:
boundscheck=False
,cache=True
Precompiled signatures:
fff->
,ddd->
Return a shorter sample of the waveform, starting at the specified index.
Note
The length of the output waveform is determined by the length of w_out rather than an input parameter. If the the length of w_out plus t0_in extends past the end of w_in or if t0_in is negative, remaining values are padded with
numpy.nan
.