toolbox_scs.detectors.pes

Beam Arrival Monitor related sub-routines

Copyright (2021) SCS Team.

(contributions preferrably comply with pep8 code structure guidelines.)

Module Contents

Functions

get_pes_tof(proposal, runNB, mnemonic[, start, ...])

Extracts time-of-flight spectra from raw digitizer traces. The spectra

get_pes_params(run[, channel])

Extract PES parameters for a given extra_data DataCollection.

save_pes_avg_traces(proposal, runNB[, channels, subdir])

Save average traces of PES into an h5 file.

load_pes_avg_traces(proposal, runNB[, channels, subdir])

Load existing PES average traces.

toolbox_scs.detectors.pes.get_pes_tof(proposal, runNB, mnemonic, start=0, origin=None, width=None, subtract_baseline=False, baseStart=None, baseWidth=40, merge_with=None)[source]

Extracts time-of-flight spectra from raw digitizer traces. The spectra are aligned by pulse Id using the SASE 3 bunch pattern. If origin is not None, a time coordinate in nanoseconds ‘time_ns’ is computed and added to the DataArray.

Parameters:
  • proposal (int) – The proposal number.

  • runNB (int) – The run number.

  • mnemonic (str) – mnemonic for PES, e.g. “PES_2Araw”.

  • start (int) – starting sample of the first spectrum in the raw trace.

  • origin (int) – sample of the trace that corresponds to time-of-flight origin, also called prompt. Used to compute the ‘time_ns’ coordinates. If None, computation of ‘time_ns’ is skipped.

  • width (int) – number of samples per spectra. If None, the number of samples for 4.5 MHz repetition rate is used.

  • subtract_baseline (bool) – If True, subtract baseline defined by baseStart and baseWidth to each spectrum.

  • baseStart (int) – starting sample of the baseline.

  • baseWidth (int) – number of samples to average (starting from baseStart) for baseline calculation.

  • merge_with (xarray Dataset) – If provided, the resulting Dataset will be merged with this one.

Returns:

pes – DataArray containing the PES time-of-flight spectra.

Return type:

xarray DataArray

Example

>>> import toolbox_scs as tb
>>> import toolbox_scs.detectors as tbdet
>>> proposal, runNB = 900447, 12
>>> pes = tbdet.get_pes_tof(proposal, runNB, 'PES_2Araw',
>>>                         start=2557, origin=76)
toolbox_scs.detectors.pes.get_pes_params(run, channel=None)[source]

Extract PES parameters for a given extra_data DataCollection. Parameters are gas, binding energy, retardation voltages or all voltages of the MPOD.

Parameters:
  • run (extra_data.DataCollection) – DataCollection containing the digitizer data

  • channel (str) – Channel name or PES mnemonic, e.g. ‘2A’ or ‘PES_1Craw’. If None, or if the channel is not found in the data, the retardation voltage for all channels is retrieved.

Returns:

params – dictionnary of PES parameters

Return type:

dict

toolbox_scs.detectors.pes.save_pes_avg_traces(proposal, runNB, channels=None, subdir='usr/processed_runs')[source]

Save average traces of PES into an h5 file.

Parameters:
  • proposal (int) – The proposal number.

  • runNB (int) – The run number.

  • channels (str or list) – The PES channels or mnemonics, e.g. ‘2A’, [‘2A’, ‘3C’], [‘PES_1Araw’, ‘PES_4Draw’, ‘3B’]

  • subdir (str) – subdirectory. The data is stored in <proposal path>/<subdir>/r{runNB:04d}/f’r{runNB:04d}-pes-data.h5’

  • Output

  • ------

  • traces. (xarray Dataset saved in a h5 file containing the PES average) –

toolbox_scs.detectors.pes.load_pes_avg_traces(proposal, runNB, channels=None, subdir='usr/processed_runs')[source]

Load existing PES average traces.

Parameters:
  • proposal (int) – The proposal number.

  • runNB (int) – The run number.

  • channels (str or list) – The PES channels or mnemonics, e.g. ‘2A’, [‘2A’, ‘3C’], [‘PES_1Araw’, ‘PES_4Draw’, ‘3B’]

  • subdir (str) – subdirectory. The data is stored in <proposal path>/<subdir>/r{runNB:04d}/f’r{runNB:04d}-pes-data.h5’

  • Output

  • ------

  • ds (xarray Dataset) – dataset containing the PES average traces.