pueoAnalysisTools
Loading...
Searching...
No Matches
Loading Waveforms

Load and plot waveforms from either the "IceFinal" files or the "dual output" files. More...

Functions

pl.DataFrame() load_waveforms (ROOT.pueo.Dataset dataset, str antennas=None)
 Loads time domain signals.
 
pl.DataFrame upsample_waveforms (pl.DataFrame waveforms, int upsample_factor)
 Calls SciPy's resample() to upsample the waveforms.
 
None plot_waveforms (pl.DataFrame waveforms, str plot_name)
 Uses Vega-Altair to make waveform plots.
 
tuple[pl.DataFrame(), np.array(1024)] load_waveforms_from_IceFinal (str filepath, [int] event_number=0, str antennas=None, int upsample_factor=None, bool want_plots=False, str plot_name=None)
 Loads time domain signals from an IceFinal file; users can optionally provide a list of antenna indices.
 

Detailed Description

Load and plot waveforms from either the "IceFinal" files or the "dual output" files.

Function Documentation

◆ load_waveforms()

pl.DataFrame() load_waveforms ( ROOT.pueo.Dataset dataset,
str antennas = None )

Loads time domain signals.

Parameters
[in]datasetThe output of initialise.load_pueoEvent_Dataset
[in]antennas(optional) whitespace separated list specifying AntNums.
Return values
waveformsThe desired channels and their waveforms [Volts]
Todo
Currently main instrument only, since oftentimes LF is not simulated
  • antenna:
    • Users can optionally provide a string containing a whitespace separated list of antenna indices (ie AntNum)
    • If no antennas is specified, all MI channels will be loaded into the DataFrame.
    • An exception is raised if the list contains invalid AntNums.
  • Output schema:
    $ GlobalChan <u32>
    $ AntNum <enum>
    $ AntIdx <u32>
    $ Pol <enum>
    $ AntType <enum>
    $ SurfNum <u8>
    $ SurfChanNum <u8>
    $ PhiSector <u8>
    $ Ring <u8>
    $ waveforms (volts) <array[f64, 1024]>
    $ step size (nanoseconds) <f64>

Definition at line 24 of file waveform_plots.py.

◆ upsample_waveforms()

pl.DataFrame upsample_waveforms ( pl.DataFrame waveforms,
int upsample_factor )

Calls SciPy's resample() to upsample the waveforms.

Parameters
[in]waveformsThe output of load_waveforms
[in]upsample_factor
Return values
upsampled_waveforms
  • The following columns are required in the input:
    1. waveforms (volts)
    2. step size (nanoseconds)
  • No additional columns will be attached to the output, waveforms (volts) and step size (nanoseconds) will be upsampled and replaced.

Definition at line 93 of file waveform_plots.py.

◆ plot_waveforms()

None plot_waveforms ( pl.DataFrame waveforms,
str plot_name )

Uses Vega-Altair to make waveform plots.

Parameters
[in]waveformsThe output of any other functions in load_waveforms
[in]plot_nameBe sure to specify the file format (eg. my_plot.png or my_plot.html)
  • The following columns are required in the input:
  • waveforms (volts)
  • step size (nanoseconds)
  • PhiSector
  • Ring
  • Pol
  • Example: Plot the vpol channels of three antennas using the last event in run0/.
    Place the figure in the /tmp directory.
    run_zero_data = load_pueoEvent_Dataset(pueo_mc_data=Path("/home/all_my_runs"), run_number=0)
    run_zero_data.last()
    wf = load_waveforms(dataset=run_zero_data, antennas="403 111 201")
    plot_waveforms(wf.filter(pl.col("Pol") == "V"), plot_name="/tmp/waveforms.html")
    Note
    The actual file name will be changed to /tmp/Vpol_waveforms.html. Had we chosen not to filter out the horizontally polarized channels, there will be another plot: /tmp/Hpol_waveforms.html.

Definition at line 122 of file waveform_plots.py.

◆ load_waveforms_from_IceFinal()

tuple[pl.DataFrame(), np.array(1024)] load_waveforms_from_IceFinal ( str filepath,
[int] event_number = 0,
str antennas = None,
int upsample_factor = None,
bool want_plots = False,
str plot_name = None )

Loads time domain signals from an IceFinal file; users can optionally provide a list of antenna indices.

Parameters
[in]filepathpath to IceFinal_*_allTree.root
[in]event_number(optional) an integer or a list of integers, defaults to event 0
[in]antennas(optional) whitespace separated string specifying the antenna names
[in]upsample_factor(optional)
[in]want_plots(optional) defaults to False
[in]plot_name(optional)
Return values
waveforms_framewaveforms recorded by the antennas [Volts]
time_meshtime mesh corresponding to the waveforms [seconds].
  • The function will error out if the passTree of the IceFinal file (see filepath) is empty. (ie. no triggered events)
  • If no antennas are specified, all waveforms will be loaded into the DataFrame, which presents them in ascending AntNum.
  • Example:
    load_waveforms_from_IceFinal("IceFinal_420_allTree.root", antennas="403 111 201", event_number=[0, 1])
    ┌──────────────┬────────┬─────────────┬─────────────────────┬───────────────────────┐
    │ event number ┆ AntNum ┆ pueoSim idx ┆ VPol waveforms (V) ┆ truth (signal origin) │
    │ --- ┆ --- ┆ --- ┆ --- ┆ --- │
    │ i32 ┆ enum ┆ u8 ┆ array[f64, 1024] ┆ array[f64, 3] │
    ╞══════════════╪════════╪═════════════╪═════════════════════╪═══════════════════════╡
    │ 0 ┆ 201 ┆ 1 ┆ […] ┆ […] │
    │ 0 ┆ 403 ┆ 11 ┆ […] ┆ […] │
    │ 0 ┆ 111 ┆ 40 ┆ […] ┆ […] │
    │ 1 ┆ 201 ┆ 1 ┆ […] ┆ […] │
    │ 1 ┆ 403 ┆ 11 ┆ […] ┆ […] │
    │ 1 ┆ 111 ┆ 40 ┆ […] ┆ […] │
    └──────────────┴────────┴─────────────┴─────────────────────┴───────────────────────┘
  • The waveforms can be upsampled by upsample_factor.
  • truth (signal origin) is a unit vector that represents the true direction of the incoming signal.
    • Note that in the IceFinal files, the variable RFdir_payload is the direction the signal is travelling to.
    • Therefore, to obtain the direction that the signal is coming from, we need the opposite vector:
      signal_origin = - RFdir_payload
    • This is useful for verifying the correctness of the correlation maps.
  • The function also optionally creates waveform plots and saves to the working directory (ie. wherever the user calls this script) via the internal function plot_waveforms(). Below is an example of a three-antenna plot (actually two plots would be made, one for each event):

    load_waveforms_from_IceFinal(
    "IceFinal_420_allTree.root", antennas="403 111 201",
    want_plots=True, event_number=[0, 1])

    Todo
    currently main instrument only
    Warning
    Will abort if non-existent antenna indices are provided.
    load_waveforms_from_IceFinal("IceFinal_420_allTree.root", antennas="201 511")
    In this case, Antenna 511 does not exist, since there are only 4 antennas in each PhiSector.
    Thus, you would see

    InvalidOperationError: conversion from str to enum failed in column '' for 1 out of 2 values: ["511"]

    Ensure that all values in the input column are present in the categories of the enum datatype.

    Perhaps you included a non-existent antenna index? aborting operation at function: load_waveforms_from_IceFinal()

Definition at line 193 of file waveform_plots.py.