Scintillation Light (light_sim)¶
Simulates scintillation-photon propagation (via a visibility LUT), SiPM response, Poisson fluctuations, trigger logic, and waveform digitisation.
Module that simulates smearing effects of the light incident on each photodetector
- larndsim.light_sim.get_nticks(light_incidence)[source][source]¶
Calculates the number of time ticks needed to simulate light signals of the event (plus the desired pre- and post-intervals)
- Parameters:
light_incidence (array) – shape (ntracks, ndet), containing first hit time and number of photons on each detector
- Returns:
number of time ticks (int), time of first tick (float) [in microseconds]
- Return type:
- larndsim.light_sim.get_active_op_channel(light_incidence)[source][source]¶
Returns an array of optical channels that need to be simulated
- Parameters:
light_incidence (array) – shape (ntracks, ndet), containing first hit time and number of photons on each detector
- Returns:
shape (ndet_active,) op detector index of each active channel (int)
- Return type:
array
- larndsim.light_sim.sum_light_signals(segments, segment_voxel, segment_track_id, light_inc, op_channel, lut, start_time, light_sample_inc, light_sample_inc_true_track_id, light_sample_inc_true_photons, sorted_indices, t0_profile_length)[source]¶
Sums the number of photons observed by each light detector at each time tick
- Parameters:
segments (array) – shape (ntracks,), edep-sim tracks to simulate
segment_voxel (array) – shape (ntracks, 3), LUT voxel for eack edep-sim track
segment_track_id (array) – shape (ntracks,), unique id for each track segment (for MC truth backtracking)
light_inc (array) – shape (ntracks, ndet), number of photons incident on each detector and voxel id
op_channel (array) – shape (ntracks, ndet_active), optical channel index, will use lut[:,:,:,op_channel%lut.shape[3]] to look up timing information
lut (array) – shape (nx,ny,nz,ndet_tpc), light look up table
start_time (float) – start time of light simulation in microseconds
light_sample_inc (array) – output array, shape (ndet, nticks), number of photons incident on each detector at each time tick (propogation delay only)
light_sample_inc_true_track_id (array) – output array, shape (ndet, nticks, maxtracks), true track ids on each detector at each time tick (propogation delay only)
light_sample_inc_true_photons (array) – output array, shape (ndet, nticks, maxtracks), number of photons incident on each detector at each time tick from each track
sorted_indices (array) – shape (maxtracks,), indices of segments sorted by how much light they contribute
- larndsim.light_sim.scintillation_model(time_tick)[source]¶
Calculates the fraction of scintillation photons emitted during time interval time_tick to time_tick + 1
- larndsim.light_sim.scintillation_array(scint_model)[source]¶
Calculates the fraction of scintillation photons emitted during time interval time_tick to time_tick + 1 for the entire input array.
- Parameters:
scint_model – array to store result
- larndsim.light_sim.calc_scintillation_effect(light_sample_inc, light_sample_inc_true_track_id, light_sample_inc_true_photons, light_sample_inc_scint, light_sample_inc_scint_true_track_id, light_sample_inc_scint_true_photons, scint_model)[source]¶
Applies a smearing effect due to the liquid argon scintillation time profile using a two decay component scintillation model.
- Parameters:
light_sample_inc (array) – shape (ndet, ntick), light incident on each detector
light_sample_inc_scint (array) – output array, shape (ndet, ntick), light incident on each detector after accounting for scintillation time
- larndsim.light_sim.xoroshiro128p_poisson_int32(mean, states, index)[source]¶
Return poisson distributed int32 and advance states[index]. For efficiency, if mean > 30, returns a gaussian distributed int32 with mean == mean and std = sqrt(mean) truncated at 0 (approximately equivalent to a poisson- distributed number)
[DOI:10.1007/978-1-4613-8643-8_10]
- larndsim.light_sim.calc_stat_fluctuations(light_sample_inc, light_sample_inc_disc, rng_states)[source]¶
Simulates Poisson fluctuations in the number of PE per time tick.
- Parameters:
light_sample_inc (array) – shape (ndet, ntick), effective photocurrent on each detector
light_sample_inc_disc (array) – output array, shape (ndet, ntick), effective photocurrent on each detector (with stocastic fluctuations)
rng_states (array) – shape (>ndet*ntick,), random number states
- larndsim.light_sim.interp(idx, arr, low, high)[source]¶
Performs a simple linear interpolation of an array at a given floating point index
- larndsim.light_sim.sipm_response_model(time_tick)[source]¶
Calculates the SiPM response from a PE at time_tick relative to the PE time
- larndsim.light_sim.sipm_response_array(sipm_response)[source]¶
Calculates the SiPM response from a PE at every time_tick relative to the PE time for the given array
- Parameters:
sipm_response – array to store response
- larndsim.light_sim.calc_light_detector_response(light_sample_inc, light_sample_inc_true_track_id, light_sample_inc_true_photons, light_response, light_response_true_track_id, light_response_true_photons, light_gain, sipm_response)[source]¶
Simulates the SiPM reponse and digit
- Parameters:
light_sample_inc (array) – shape (ndet, ntick), PE produced on each SiPM at each time tick
light_response (array) – shape (ndet, ntick), ADC value at each time tick
- larndsim.light_sim.gen_light_detector_noise(shape, light_det_noise)[source][source]¶
Generates uncorrelated noise with a defined frequency spectrum
- Parameters:
shape (tuple) – desired shape of output noise, shape[0] must equal light_det_noise.shape[0]
light_det_noise (array) – FFT of noise, light_det_noise.ndim == 2
- Returns:
shape (shape[0], shape[1]), randomly generated sample noise
- Return type:
array
- larndsim.light_sim.get_triggers(signal, group_threshold, op_channel_idx, i_subbatch)[source][source]¶
Identifies each simulated ticks that would initiate a trigger taking into account the ADC digitization window
- Parameters:
signal (array) – shape (ndet, nticks), simulated signal on each channel
group_threshold (array) – shape (ngrp,), threshold on group sum (requires ndet/ngrp == OP_CHANNEL_PER_TRIG)
op_channel_idx (array) – shape (ndet,), optical channel index for each signal
i_subbatch (int) – index of the sub_batch numbering (“itrk in the batch for loop”)
- Returns:
array of tick indices at each trigger (shape (ntrigs,)) and array of op channel index (shape (ntrigs, ndet_module))
- Return type:
- larndsim.light_sim.digitize_signal(signal, signal_op_channel_idx, trigger_idx, trigger_op_channel_idx, signal_true_track_id, signal_true_photons, digit_signal, digit_signal_true_track_id, digit_signal_true_photons)[source]¶
Interpolate signal to the appropriate sampling frequency
- Parameters:
signal (array) – shape (ndet, nticks), simulated signal on each channel
signal_op_channel_idx (array) – shape (ndet,), optical channel index for each simulated signal
trigger_idx (array) – shape (ntrigs,), tick index for each trigger
trigger_op_channel_idx (array) – shape (ntrigs, ndet_module), optical channel index for each trigger
digit_signal (array) – output array, shape (ntrigs, ndet_module, nsamples), digitized signal
- larndsim.light_sim.sim_triggers(bpg, tpb, signal, signal_op_channel_idx, signal_true_track_id, signal_true_photons, trigger_idx, op_channel_idx, digit_samples, light_det_noise)[source][source]¶
Generates digitized waveforms at specified simulation tick indices
- Parameters:
bpg (tuple) – blocks per grid used to generate digitized waveforms, len(bpg) == 3, prod(bpg) * prod(tpb) >= digit_samples.size
tpb (tuple) – threads per grid used to generate digitized waveforms, len(bpg) == 3, bpg[i] * tpb[i] >= digit_samples.shape[i]
signal (array) – shape (ndet, nticks), simulated signal on each channel
signal_op_channel_idx (array) – shape (ndet,), optical channel index for each simulated signal
signal_true_track_id (array) – shape (ndet, nticks, ntruth), true segments associated with each tick
signal_true_photons (array) – shape (ndet, nticks, ntruth), true photons associated with each tick from each track
trigger_idx (array) – shape (ntrigs,), tick index for each trigger to digitize
op_channel_idx (array) – shape (ntrigs, ndet_module), optical channel indices for each trigger
digit_samples (int) – number of digitizations per waveform
light_det_noise (array) – shape (ndet, nnoise_bins), noise spectrum for each channel (only used if waveforms extend past simulated signal)
- Returns:
shape (ntrigs, ndet_module, digit_samples), digitized waveform on each channel for each trigger
- Return type:
array
- larndsim.light_sim.zero_suppress_waveform_truth(waveforms_true_track_id, waveforms_true_photons, i_evt, i_trig, i_mod=-1)[source][source]¶
Filter empty light waveform backtracking which is filled with the default value ‘-1’, and flatten the backtracking info to 1D array.
- Parameters:
i_evt (int) – event id
waveforms_true_track_id (array) – shape (ntrigs, ndet, nsamples), segment ids contributing to each sample
waveforms_true_photons (array) – shape (ntrigs, ndet, nsamples), true photocurrent at each sample
i_evt – true event id
i_trig (int) – light trigger or light event id
i_mod (int) – module id. The default value is -1 which indicates that there is no modular variation activated.
- Returns:
1D array which logs [‘trigger_id’, ‘op_channel_id’, ‘tick’, ‘event_id’, ‘segment_id’, ‘pe_current’]. The first three locate a unique data point on a recorded light waveform, and the last three provide the truth information associated to it. event_id can be infered from segment_id but the information helps searching the corresponding reco information from a true event.
- larndsim.light_sim.export_light_wvfm_to_hdf5(event_id, waveforms, output_filename, waveforms_true_track_id, waveforms_true_photons, i_trig, i_mod=-1, compression=None)[source][source]¶
Saves waveforms to output file
- Parameters:
event_id (array) – shape (ntrigs,), event id for each trigger
waveforms (array) – shape (ntrigs, ndet_module, nsamples), simulated waveforms to save
output_filename (str) – output hdf5 file path
waveforms_true_track_id (array) – shape (ntrigs, ndet, nsamples), segment ids contributing to each sample
waveforms_true_photons (array) – shape (ntrigs, ndet, nsamples), true photocurrent at each sample
i_mod (int) – module id. The default value is -1 which indicates that there is no modular variation activated.
- larndsim.light_sim.export_light_trig_to_hdf5(event_id, start_times, trigger_idx, op_channel_idx, output_filename, event_times, compression=None)[source][source]¶
Saves light trigger to output file
- Parameters:
event_id (array) – shape (ntrigs,), event id for each trigger
start_times (array) – shape (ntrigs,), simulation time offset for each trigger [microseconds]
trigger_idx (array) – shape (ntrigs,), simulation time tick of each trigger
op_channel_idx (array) – shape (ntrigs, ndet_module), optical channel index for each trigger
output_filename (str) – output hdf5 file path
event_times (array) – shape (nevents,), global event t0 for each unique event [microseconds]
- larndsim.light_sim.export_to_hdf5(event_id, start_times, trigger_idx, op_channel_idx, waveforms, output_filename, event_times, waveforms_true_track_id, waveforms_true_photons, i_trig, i_mod, compression=None)[source][source]¶
Saves waveforms to output file
- Parameters:
event_id (array) – shape (ntrigs,), event id for each trigger
start_times (array) – shape (ntrigs,), simulation time offset for each trigger [microseconds]
trigger_idx (array) – shape (ntrigs,), simulation time tick of each trigger
op_channel_idx (array) – shape (ntrigs, ndet_module), optical channel index for each trigger
waveforms (array) – shape (ntrigs, ndet_module, nsamples), simulated waveforms to save
output_filename (str) – output hdf5 file path
event_times (array) – shape (nevents,), global event t0 for each unique event [microseconds]
waveforms_true_track_id (array) – shape (ntrigs, ndet, nsamples), segment ids contributing to each sample
waveforms_true_photons (array) – shape (ntrigs, ndet, nsamples), true photocurrent at each sample