daq2lh5.buffer_processor package

This subpackage allows for preprocessing of data in daq and raw level files

Submodules

daq2lh5.buffer_processor.buffer_processor module

daq2lh5.buffer_processor.buffer_processor.buffer_processor(rb, db_dict=None)

Process raw data buffers.

Takes in a RawBuffer, performs any processes specified in the RawBuffer’s proc_spec attribute (a dictionary), and returns a Table with the processed buffer. This tmp_table shares columns with the RawBuffer’s LGDO (rb.lgdo), so no data is copied.

Currently implemented proc_spec processors:

"window": ["waveform", start_index, stop_index, "out_name"] (list)

Windows objects with a name specified by the first argument, the window start and stop indices are the next two arguments, and then updates the rb.lgdo with a name specified by the last argument. If the object is an WaveformTable, then the t0 attribute is updated accordingly. Although it is possible to use the DSP config to perform windowing, this hard-coded version avoids a conversion to float32.

"dsp_config": { <dsp_config> } (dict)

Performs DSP given by the <dsp_config> specification. See build_processing_chain() for more information on DSP configuration dictionaries. All fields in the output of the DSP are written to the rb.lgdo.

"drop": ["waveform" [, ...]] (list)

Drops any requested fields from the rb.lgdo.

"dtype_conv": {"lgdo": "dtype" [, ...]} (dict)

Casts lgdo to the requested data type.

"compression": {"lgdo": "codec_name" [, ...]} (dict)

Updates the compression attribute of lgdo to codec_name. The attribute sets the compression algorithm applied by read() before writing lgdo to disk. Can be used to apply custom waveform compression algorithms from lh5.compression.

"hdf5_settings": {"lgdo": { <HDF5 settings> }} (dict)

Updates the hdf5_settings attribute of lgdo. The attribute sets the HDF5 dataset options applied by read() before writing lgdo to disk.

Parameters:

rb (RawBuffer) – A RawBuffer to be processed, must contain a proc_spec attribute.

Return type:

Table

Notes

The original waveforms column in the table is not written to file if request! All updates are done on the tmp_table, which shares the fields with rb.lgdo and are done in place. The tmp_table is necessary so that the rb.lgdo keeps arrays needed by the table in the buffer. An example proc_spec in an build_raw() out_spec is below.

{
  "FCEventDecoder" : {
    "g{key:0>3d}" : {
      "key_list" : [[24, 64]],
      "out_stream" : "$DATADIR/{file_key}_geds.lh5:/geds",
      "proc_spec": {
        "window": ["waveform", 100, -100, "windowed_waveform"],
        "dsp_config": {
          "outputs": ["presummed_waveform", "t_sat_lo", "t_sat_hi"],
          "processors": {
            "presummed_waveform": {
              "function": "presum",
              "module": "dspeed.processors",
              "args": [
                "waveform",
                "presummed_waveform(shape=len(waveform)//16, period=waveform.period*16, offset=waveform.offset, 'f')"
              ],
              "unit": "ADC"
            },
            "t_sat_lo, t_sat_hi": {
              "function": "saturation",
              "module": "dspeed.processors",
              "args": ["waveform", 16, "t_sat_lo", "t_sat_hi"],
              "unit": "ADC"
            }
          }
        },
        "drop": ["waveform", "packet_ids"],
        "dtype_conv": {
          "presummed_waveform/values": "uint32",
          "t_sat_lo": "uint16",
          "t_sat_hi": "uint16",
        ,}
        "compression": {
          "windowed_waveform/values": RadwareSigcompress(codec_shift=-32768),
        }
        "hdf5_settings": {
          "presummed_waveform/values": {"shuffle": True, "compression": "lzf"},
        }
      }
    },
    "spms" : {
      "key_list" : [ [6,23] ],
      "out_stream" : "$DATADIR/{file_key}_spms.lh5:/spms"
    }
  }
}
daq2lh5.buffer_processor.buffer_processor.process_drop(rb, tmp_table)

Drops any requested fields from the rb.lgdo.

daq2lh5.buffer_processor.buffer_processor.process_dsp(rb, tmp_table, db_dict=None)

Run a DSP processing chain with dspeed.

Run a provided DSP config from rb.proc_spec using dsp.build_processing_chain(), and add specified outputs to the rb.lgdo.

Parameters:
  • rb (RawBuffer) – a RawBuffer that contains a proc_spec and an lgdo attribute.

  • tmp_table (Table) – a lgdo.Table that is temporarily created to be written to the raw file.

  • db_dict (dict) – a database dictionary storing parameters for each channel

Notes

rb.lgdo is assumed to be an :class`.Table` so that multiple DSP processor outputs can be written to it.

daq2lh5.buffer_processor.buffer_processor.process_dtype_conv(rb, tmp_table)

Change the types of fields in an rb.lgdo according to the values specified in the proc_spec’s dtype_conv list. It operates in place on tmp_table.

Notes

This assumes that name provided points to an object in the rb.lgdo that has an nda attribute.

daq2lh5.buffer_processor.buffer_processor.process_window(rb, tmp_table)

Window ArrayOfEqualSizedArrays.

Windows arrays of equal sized arrays according to specifications given in the rb.proc_spec window key.

First checks if the rb.lgdo is a Table or not. If it’s not a table, then we only process it if its rb.out_name is the same as the window in_name.

If rb.lgdo is a table, special processing is done if the window in_name field is an WaveformTable in order to update the t0s. Otherwise, windowing of the field is performed without updating any of the other attributes.

Parameters:
  • rb (RawBuffer) – a RawBuffer to be processed.

  • tmp_table (Table) – a Table that shares columns with the rb.lgdo.

Notes

This windowing hard-coded; it is done without calling dsp.build_dsp to avoid a conversion to float32.

daq2lh5.buffer_processor.buffer_processor.process_windowed_t0(t0s, dts, start_index)

In order for the processed data to work well with dsp.build_dsp, we need to keep t0 in its original units.

So we transform start_index to the units of t0 and add it to every t0 value.

Return type:

Array

daq2lh5.buffer_processor.buffer_processor.window_array_of_arrays(array_of_arrays, window_start_idx, window_end_idx)

Given an array of equal sized arrays, for each array it returns the view [window_start_idx:window_end_idx].

Return type:

ArrayOfEqualSizedArrays

daq2lh5.buffer_processor.lh5_buffer_processor module

daq2lh5.buffer_processor.lh5_buffer_processor.lh5_buffer_processor(lh5_raw_file_in, overwrite=False, out_spec=None, proc_file_name=None, db_dict=None)

Process raw buffers from an LH5 file.

This function performs data processing on an existing raw-level LH5 file according to out_spec. It iterates through any valid table in the lh5_raw_file_in file and checks that either the group_name or the out_name of this table are found in the out_spec. If there is a valid proc_spec entry in the out_spec, then the table is stored in a RawBuffer and the RawBuffer() is passed to the buffer processor. If no proc_spec is found, the table is written to the file without any modifications.

Parameters:
  • lh5_raw_file_in (str) – the path of a file created by build_raw().

  • overwrite (bool) – sets whether to overwrite the output file(s) if it (they) already exist.

  • out_spec (dict | None) – an out_spec for build_raw() also containing a dictionary named proc_spec used for data processing. See buffer_processor() for an example of such an out_spec.

  • proc_file_name (str | None) – The path to the processed output file created. If None, uses {lh5_raw_file_in}_proc.lh5 as the output filename.

Notes

This function assumes that all raw file data are stored in at most two h5py.Groups, as in e.g. ch000/raw or raw/geds