larndsim package¶
Quenching¶
Module to implement the quenching of the ionized electrons through the detector
- 
larndsim.quenching.quench(tracks, mode)¶ This CUDA kernel takes as input an array of track segments and calculates the number of electrons and photons that reach the anode plane after recombination. It is possible to pick among two models: Box (Baller, 2013 JINST 8 P08005) or Birks (Amoruso, et al NIM A 523 (2004) 275).
- Parameters
 tracks (
numpy.ndarray) – array containing the tracks segment informationmode (int) – recombination model (physics.BOX or physics.BIRKS).
Drifting¶
Module to implement the propagation of the electrons towards the anode.
- 
larndsim.drifting.drift(tracks)¶ This function takes as input an array of track segments and calculates the properties of the segments at the anode:
z coordinate at the anode
number of electrons taking into account electron lifetime
longitudinal diffusion
transverse diffusion
time of arrival at the anode
- Parameters
 tracks (
numpy.ndarray) – array containing the tracks segment information
Induced charge¶
Module that calculates the current induced by edep-sim track segments on the pixels
- 
larndsim.detsim.get_closest_waveform¶ This function, given a point on the pixel pad and a time, returns the closest tabulated waveformm from the response array.
- Parameters
 x (float) – x coordinate of the point
y (float) – y coordinate of the point
t (float) – time of the waveform
response (
numpy.ndarray) – array containing the tabulated waveforms
- Returns
 the value of the induced current at time t for a charge at (x,y)
- Return type
 float
- 
larndsim.detsim.get_pixel_coordinates¶ Returns the coordinates of the pixel center given the pixel ID
- 
larndsim.detsim.get_track_pixel_map(track_pixel_map, unique_pix, pixels)¶ This kernel fills a 2D array which contains, for each unique pixel, an array with the track indeces associated to that pixel.
- Parameters
 track_pixel_map (
numpy.ndarray) – 2D array that will contain the association between the unique pixels array and the track indecesunique_pix (
numpy.ndarray) – 1D array containing the unique pixelspixels (
numpy.ndarray) – 2D array containing the pixels for each track.
- 
larndsim.detsim.overlapping_segment¶ Calculates the segment of the track defined by start, end that overlaps with a circle centered at x,y
- 
larndsim.detsim.rho¶ Function that returns the amount of charge at a certain point in space
- Parameters
 point (tuple) – point coordinates
q (float) – total charge
start (tuple) – segment start coordinates
sigmas (tuple) – diffusion coefficients
segment (tuple) – segment sizes
- Returns
 the amount of charge at point.
- Return type
 float
- 
larndsim.detsim.sign¶ Sign function
- Parameters
 x (float) – input number
- Returns
 1 if x>=0 else -1
- Return type
 int
- 
larndsim.detsim.sum_pixel_signals(pixels_signals, signals, track_starts, pixel_index_map, track_pixel_map, pixels_tracks_signals)¶ This function sums the induced current signals on the same pixel.
- Parameters
 pixels_signals (
numpy.ndarray) – 2D array that will contain the summed signal for each pixel. First dimension is the pixel ID, second dimension is the time ticksignals (
numpy.ndarray) – 3D array with dimensions S x P x T, where S is the number of track segments, P is the number of pixels, and T is the number of time ticks.track_starts (
numpy.ndarray) – 1D array containing the starting time of each trackpixel_index_map (
numpy.ndarray) – 2D array containing the correspondence between the track index and the pixel ID index.track_pixel_map (
numpy.ndarray) – 2D array containing the association between the unique pixels array and the array containing the pixels for each track.pixels_tracks_signals (
numpy.ndarray) – 3D array that will contain the waveforms for each pixel and each track that induced current on the pixel.
- 
larndsim.detsim.time_intervals(track_starts, time_max, tracks)¶ Find the value of the longest signal time and stores the start time of each segment.
- Parameters
 track_starts (
numpy.ndarray) – array where we store the segments start timetime_max (
numpy.ndarray) – array where we store the longest signal timetracks (
numpy.ndarray) – array containing the segment information
- 
larndsim.detsim.track_point¶ This function returns the segment coordinates for a point along the z coordinate
- Parameters
 start (tuple) – start coordinates
direction (tuple) – direction coordinates
z (float) – z coordinate corresponding to the x, y coordinates
- Returns
 the (x,y) pair of coordinates for the segment at z
- Return type
 tuple
- 
larndsim.detsim.tracks_current(signals, pixels, tracks, response)¶ This CUDA kernel calculates the charge induced on the pixels by the input tracks.
- Parameters
 signals (
numpy.ndarray) – empty 3D array with dimensions S x P x T, where S is the number of track segments, P is the number of pixels, and T is the number of time ticks. The output is stored here.pixels (
numpy.ndarray) – 2D array with dimensions S x P , where S is the number of track segments, P is the number of pixels and contains the pixel ID number.tracks (
numpy.ndarray) – 2D array containing the detector segments.response (
numpy.ndarray) – 3D array containing the tabulated response.
- 
larndsim.detsim.tracks_current_mc(signals, pixels, tracks, response, rng_states)¶ This CUDA kernel calculates the charge induced on the pixels by the input tracks using a MC method
- Parameters
 signals (
numpy.ndarray) – empty 3D array with dimensions S x P x T, where S is the number of track segments, P is the number of pixels, and T is the number of time ticks. The output is stored here.pixels (
numpy.ndarray) – 2D array with dimensions S x P , where S is the number of track segments, P is the number of pixels and contains the pixel ID number.tracks (
numpy.ndarray) – 2D array containing the detector segments.response (
numpy.ndarray) – 3D array containing the tabulated response.rng_states (
numpy.ndarray) – array of random states for noise generation
- 
larndsim.detsim.z_interval¶ Here we calculate the interval in the drift direction for the pixel pID using the impact factor
- Parameters
 start_point (tuple) – coordinates of the segment start
end_point (tuple) – coordinates of the segment end
x_p (float) – pixel center x coordinate
y_p (float) – pixel center y coordinate
tolerance (float) – maximum distance between the pixel center and the segment
- Returns
 z coordinate of the point of closest approach (POCA), z coordinate of the first slice, z coordinate of the last slice. (0,0,0) if POCA > tolerance.
- Return type
 tuple
Charge readout¶
Module that simulates the front-end electronics (triggering, ADC)
- 
larndsim.fee.ADC_BUSY_DELAY= 9¶ ADC busy delay in clock cycles
- 
larndsim.fee.ADC_COUNTS= 256¶ Number of ADC counts
- 
larndsim.fee.ADC_HOLD_DELAY= 15¶ ADC hold delay in clock cycles
- 
larndsim.fee.BUFFER_RISETIME= 0.1¶ Buffer risetime in \(\mu s\) (set >0 to include buffer response simulation)
- 
larndsim.fee.CLOCK_CYCLE= 0.1¶ Clock cycle time in \(\mu s\)
- 
larndsim.fee.DISCRIMINATION_THRESHOLD= 7000.0¶ Discrimination threshold in e-
- 
larndsim.fee.DISCRIMINATOR_NOISE= 650.0¶ Discriminator noise in e-
- 
larndsim.fee.EVENT_RATE= 100000¶ Average time between events in microseconds
- 
larndsim.fee.GAIN= 4e-12¶ Front-end gain in \(mV/e-\)
- 
larndsim.fee.MAX_ADC_VALUES= 20¶ Maximum number of ADC values stored per pixel
- 
larndsim.fee.RESET_CYCLES= 1¶ Reset time in clock cycles
- 
larndsim.fee.RESET_NOISE_CHARGE= 900.0¶ Reset noise in e-
- 
larndsim.fee.ROLLOVER_CYCLES= 2147483648¶ Clock rollover / reset time in larpix clock ticks
- 
larndsim.fee.UNCORRELATED_NOISE_CHARGE= 500.0¶ Uncorrelated noise in e-
- 
larndsim.fee.V_CM= 2.8800000000000004e-07¶ Common-mode voltage in \(mV\)
- 
larndsim.fee.V_PEDESTAL= 5.800000000000001e-07¶ Pedestal voltage in \(mV\)
- 
larndsim.fee.V_REF= 1.3e-06¶ Reference voltage in \(mV\)
- 
larndsim.fee.digitize(integral_list)¶ The function takes as input the integrated charge and returns the digitized ADC counts.
- Parameters
 integral_list (
numpy.ndarray) – list of charge collected by each pixel- Returns
 list of ADC values for each pixel
- Return type
 
- 
larndsim.fee.export_to_hdf5(event_id_list, adc_list, adc_ticks_list, unique_pix, current_fractions, track_ids, filename, event_start_times, is_first_event, light_trigger_times=None, light_trigger_event_id=None, light_trigger_modules=None, bad_channels=None)¶ Saves the ADC counts in the LArPix HDF5 format. :Parameters: * event_id_list (
numpy.ndarray) – list of event ids for each ADC value for each pixeladc_list (
numpy.ndarray) – list of ADC values for each pixeladc_ticks_list (
numpy.ndarray) – list of time ticks for each pixelunique_pix (
numpy.ndarray) – list of pixel IDscurrent_fractions (
numpy.ndarray) – array containing the fractional current induced by each track on each pixeltrack_ids (
numpy.ndarray) – 2D array containing the track IDs associated to each pixelfilename (str) – filename of HDF5 output file
event_times (
numpy.ndarray) – list of timestamps for start each unique event [in microseconds]is_first_event (bool) – True if this is the first event to save to the file
light_trigger_times (array) – 1D array of light trigger timestamps (relative to event t0) [in microseconds]
light_trigger_event_id (array) – 1D array of event id for each light trigger
light_trigger_modules (array) – 1D array of module id for each light trigger
bad_channels (dict) – dictionary containing as value a list of bad channels and as the chip key
- Returns
 a tuple containing the list of LArPix packets and the list of entries for the mc_packets_assn dataset
- Return type
 tuple
- 
larndsim.fee.gen_event_times(nevents, t0)¶ Generate sequential event times assuming events are uncorrelated
- Parameters
 nevents (int) – number of event times to generate
t0 (int) – offset to apply [microseconds]
- Returns
 shape (nevents,), sequential event times [microseconds]
- Return type
 array
- 
larndsim.fee.get_adc_values(pixels_signals, pixels_signals_tracks, time_ticks, adc_list, adc_ticks_list, time_padding, rng_states, current_fractions, pixel_thresholds)¶ Implementation of self-trigger logic
- Parameters
 pixels_signals (
numpy.ndarray) – list of induced currents for each pixelpixels_signals_tracks (
numpy.ndarray) – list of induced currents for each track that induces current on each pixeltime_ticks (
numpy.ndarray) – list of time ticks for each pixeladc_list (
numpy.ndarray) – list of integrated charges for each pixeladc_ticks_list (
numpy.ndarray) – list of the time ticks that correspond to each integrated chargetime_padding (float) – time interval to add to each time tick.
rng_states (
numpy.ndarray) – array of random states for noise generationcurrent_fractions (
numpy.ndarray) – 2D array that will contain the fraction of current induced on the pixel by each trackpixel_thresholds( – obj: numpy.ndarray): list of discriminator thresholds for each pixel
- 
larndsim.fee.rotate_tile(pixel_id, tile_id)¶ Returns the pixel ID of the rotated tile.
- Parameters
 pixel_id (int) – pixel ID
tile_id (int) – tile ID
- Returns
 pixel indeces
- Return type
 tuple
Scintillation light¶
Module that simulates smearing effects of the light incident on each photodetector
- 
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)¶ 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.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)¶ 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.calc_stat_fluctuations(light_sample_inc, light_sample_inc_disc, rng_states)¶ Simulates Poisson fluctuations in the number of PE per time tick.
- Parameters
 light_sample_inc (array) – shape (ndet, ntick), light incident on each detector
light_sample_inc_disc (array) – output array, shape (ndet, ntick), number PE in each time interval
rng_states (array) – shape (>ndet*ntick,), random number states
- 
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)¶ 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.export_to_hdf5(event_id, start_times, trigger_idx, op_channel_idx, waveforms, output_filename, event_times, waveforms_true_track_id, waveforms_true_photons)¶ 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, ntruth), segment ids contributing to each sample
waveforms_true_photons (array) – shape (ntrigs, ndet, nsamples, ntruth), true photocurrent at each sample
- 
larndsim.light_sim.gen_light_detector_noise(shape, light_det_noise)¶ 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_active_op_channel(light_incidence)¶ 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.get_nticks(light_incidence)¶ 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
 tuple
- 
larndsim.light_sim.get_triggers(signal, group_threshold, op_channel_idx)¶ 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
- Returns
 array of tick indices at each trigger (shape (ntrigs,)) and array of op channel index (shape (ntrigs, ndet_module))
- Return type
 tuple
- 
larndsim.light_sim.interp¶ Performs a simple linear interpolation of an array at a given floating point index
- Parameters
 idx (float) – index into array to interpolate
arr (array) – 1D array of values to interpolate
low (float) – value to return if index is less than 0
high (float) – value to return if index is above len(arr)-1
- Returns
 interpolated array value
- Return type
 float
- 
larndsim.light_sim.scintillation_model¶ Calculates the fraction of scintillation photons emitted during time interval time_tick to time_tick + 1
- Parameters
 time_tick (int) – time tick relative to t0
- Returns
 fraction of scintillation photons
- Return type
 float
- 
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)¶ 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.sipm_response_model¶ Calculates the SiPM response from a PE at time_tick relative to the PE time
- Parameters
 idet (int) – SiPM index
time_tick (int) – time tick relative to t0
- Returns
 response
- Return type
 float
- 
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)¶ 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
- 
larndsim.light_sim.xoroshiro128p_poisson_int32¶ 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]
- Parameters
 mean (float) – mean of poisson distribution
states (array) – array of RNG states
index (int) – offset in states to update
Constants¶
Set global variables with detector and physics properties
- 
larndsim.consts.load_properties(detprop_file, pixel_file, sim_file)¶ The function loads the detector properties, the pixel geometry, and the simulation YAML files and stores the constants as global variables
- Parameters
 detprop_file (str) – detector properties YAML filename
pixel_file (str) – pixel layout YAML filename
sim_file (str) – simulation properties YAML filename
Set detector constants
- 
larndsim.consts.detector.DEFAULT_PLANE_INDEX= 48879¶ Default value for pixel_plane, to indicate out-of-bounds edep
- 
larndsim.consts.detector.DRIFT_LENGTH= 0¶ TPC drift length in \(cm\)
- 
larndsim.consts.detector.ELECTRON_LIFETIME= 2200.0¶ Electron lifetime in \(\mu s\)
- 
larndsim.consts.detector.E_FIELD= 0.5¶ Electric field magnitude in \(kV/cm\)
- 
larndsim.consts.detector.LAR_DENSITY= 1.38¶ Liquid argon density in \(g/cm^3\)
- 
larndsim.consts.detector.LONG_DIFF= 4e-06¶ Longitudinal diffusion coefficient in \(cm^2/\mu s\)
- 
larndsim.consts.detector.MODULE_TO_IO_GROUPS= {}¶ Association between modules and io groups
- 
larndsim.consts.detector.MODULE_TO_TPCS= {}¶ Association between modules and tpcs
- 
larndsim.consts.detector.N_PIXELS= (0, 0)¶ Total number of pixels
- 
larndsim.consts.detector.N_PIXELS_PER_TILE= (0, 0)¶ Number of pixels in each tile
- 
larndsim.consts.detector.PIXEL_CONNECTION_DICT= {}¶ Dictionary between pixel ID and its position in the pixel array
- 
larndsim.consts.detector.PIXEL_PITCH= 0.4434¶ Pixel pitch in \(cm\)
- 
larndsim.consts.detector.RESPONSE_BIN_SIZE= 0.04434¶ Spatial sampling in the pixel reponse file in \(cm\)
- 
larndsim.consts.detector.RESPONSE_SAMPLING= 0.1¶ Time sampling in the pixel response file in \(\mu s\)
- 
larndsim.consts.detector.SAMPLED_POINTS= 40¶ Number of sampled points for each segment slice
- 
larndsim.consts.detector.TEMPERATURE= 87.17¶ Detector temperature in K
- 
larndsim.consts.detector.TILE_BORDERS= array([[0., 0.], [0., 0.]])¶ Pixel tile borders in \(cm\)
- 
larndsim.consts.detector.TILE_CHIP_TO_IO= {}¶ Association between chips and io channels
- 
larndsim.consts.detector.TILE_MAP= ()¶ Map of tiles in each anode
- 
larndsim.consts.detector.TILE_ORIENTATIONS= {}¶ Tile orientations in each anode
- 
larndsim.consts.detector.TILE_POSITIONS= {}¶ Tile position wrt the center of the anode in \(cm\)
- 
larndsim.consts.detector.TIME_INTERVAL= (0, 200.0)¶ Drift time window in \(\mu s\)
- 
larndsim.consts.detector.TIME_PADDING= 10¶ Signal time window padding in \(\mu s\)
- 
larndsim.consts.detector.TIME_SAMPLING= 0.1¶ Time sampling in \(\mu s\)
- 
larndsim.consts.detector.TIME_TICKS= array([0.000e+00, 1.000e-01, 2.000e-01, ..., 1.998e+02, 1.999e+02, 2.000e+02])¶ Numpy array containing all the time ticks in the drift time window
- 
larndsim.consts.detector.TIME_WINDOW= 8.9¶ Time window of current response in \(\mu s\)
- 
larndsim.consts.detector.TPC_BORDERS= array([], shape=(0, 3, 2), dtype=float64)¶ Borders of each TPC volume in \(cm\)
- 
larndsim.consts.detector.TPC_OFFSETS= array([], shape=(0, 3, 2), dtype=float64)¶ TPC offsets wrt the origin in \(cm\)
- 
larndsim.consts.detector.TRAN_DIFF= 8.8e-06¶ Transverse diffusion coefficient in \(cm^2/\mu s\)
- 
larndsim.consts.detector.V_DRIFT= 0.1648¶ Drift velocity in \(cm/\mu s\)
- 
larndsim.consts.detector.electron_mobility(efield, temperature)¶ Calculation of the electron mobility w.r.t temperature and electric field. .. rubric:: References
https://doi.org/10.1016/j.nima.2016.01.073 (parameterization)
- Parameters
 efield (float) – electric field in kV/cm
temperature (float) – temperature
- Returns
 electron mobility in cm^2/kV/us
- Return type
 float
- 
larndsim.consts.detector.set_detector_properties(detprop_file, pixel_file)¶ The function loads the detector properties and the pixel geometry YAML files and stores the constants as global variables
- Parameters
 detprop_file (str) – detector properties YAML filename
pixel_file (str) – pixel layout YAML filename
Set physics constants
- 
larndsim.consts.physics.BIRKS_Ab= 0.8¶ Recombination \(A_b\) value for the Birks Model
- 
larndsim.consts.physics.BIRKS_kb= 0.0486¶ Recombination \(k_b\) value for the Birks Model in \((kV/cm)(g/cm^2)/MeV\)
- 
larndsim.consts.physics.BOX_ALPHA= 0.93¶ Recombination \(\alpha\) constant for the Box model
- 
larndsim.consts.physics.BOX_BETA= 0.207¶ Recombination \(\beta\) value for the Box model in \((kV/cm)(g/cm^2)/MeV\)
- 
larndsim.consts.physics.E_CHARGE= 1.602e-19¶ Electron charge in Coulomb
- 
larndsim.consts.physics.W_ION= 2.36e-05¶ Average energy expended per ion pair in LAr in \(MeV\) from Phys. Rev. A 10, 1452
Sets ligth-related constants
- 
larndsim.consts.light.DEFAULT_LIGHT_GAIN= -2.3¶ Conversion from PE/microsecond to ADC
- 
larndsim.consts.light.IMPULSE_MODEL= array([1, 0])¶ Arbitrary input model (normalized to sum of 1)
- 
larndsim.consts.light.IMPULSE_TICK_SIZE= 0.001¶ Arbitrary input model tick size [microseconds]
- 
larndsim.consts.light.LIGHT_DET_NOISE_SAMPLE_SPACING= 0.01¶ Sample rate for input noise spectrum [microseconds]
- 
larndsim.consts.light.LIGHT_DIGIT_SAMPLE_SPACING= 0.01¶ Light waveform sample rate [microseconds]
- 
larndsim.consts.light.LIGHT_NBIT= 10¶ Light digitizer bits
- 
larndsim.consts.light.LIGHT_OSCILLATION_PERIOD= 0.095¶ Reponse oscillation period [microseconds]
- 
larndsim.consts.light.LIGHT_RESPONSE_TIME= 0.055¶ Response RC time [microseconds]
- 
larndsim.consts.light.LIGHT_TICK_SIZE= 0.001¶ Step size for light simulation [microseconds]
- 
larndsim.consts.light.LIGHT_TRIG_THRESHOLD= array([], dtype=float64)¶ Total detector light threshold [ADC] (one value for every OP_CHANNEL_PER_TRIG detector sum)
- 
larndsim.consts.light.LIGHT_TRIG_WINDOW= (0.9, 1.66)¶ Light digitization window [microseconds]
- 
larndsim.consts.light.LIGHT_WINDOW= (1, 10)¶ Pre- and post-window for light simulation [microseconds]
- 
larndsim.consts.light.MAX_MC_TRUTH_IDS= 0¶ Number of true segments to track for each time tick (MAX_MC_TRUTH_IDS=0 to disable complete truth tracking)
- 
larndsim.consts.light.MC_TRUTH_THRESHOLD= 0.1¶ Threshold for propogating truth information on a given SiPM
- 
larndsim.consts.light.OP_CHANNEL_PER_TRIG= 6¶ Number of SiPMs per detector (used by trigger)
- 
larndsim.consts.light.SCINT_PRESCALE= 1¶ Prescale factor analogous to ScintPreScale in LArSoft FIXME
- 
larndsim.consts.light.SINGLET_FRACTION= 0.3¶ Fraction of total light emitted from singlet state
- 
larndsim.consts.light.SIPM_RESPONSE_MODEL= 0¶ Set response model type (0=RLC response, 1=arbitrary input)
- 
larndsim.consts.light.TAU_S= 0.001¶ Singlet decay time [microseconds]
- 
larndsim.consts.light.TAU_T= 1.53¶ Triplet decay time [microseconds]
- 
larndsim.consts.light.W_PH= 1.95e-05¶ Ion + excitation work function in MeV
- 
larndsim.consts.light.set_light_properties(detprop_file)¶ The function loads the detector properties YAML file and stores the light-related constants as global variables
- Parameters
 detprop_file (str) – detector properties YAML filename