3@brief `main()` shows how to use
4[load_waveforms()](#waveform_plots.load_waveforms) and
5[load_waveforms_from_IceFinal()](#waveform_plots.load_waveforms_from_IceFinal).
6@note Altair is needed to make the plot.
8 pip3 install altair[save]
17from pathlib
import Path
18from initialise
import load_pueoEvent_Dataset
24def load_waveforms(dataset: ROOT.pueo.Dataset, antennas: str =
None) -> pl.DataFrame():
25 r"""! Loads time domain signals
28 @param[in] dataset The output of #initialise.load_pueoEvent_Dataset
29 @param[in] antennas (optional) **whitespace separated list** specifying
30 [AntNum](#antenna_attributes.ALL_PUEO_ANTENNA_NAMES)s.
31 @retval waveforms The desired [channels](#antenna_attributes.read_global_channel_mapping)
32 and their waveforms [Volts]
34@todo Currently main instrument only, since oftentimes LF is not simulated
37 * Users can optionally provide a string containing a whitespace separated list of antenna
38 indices (ie [AntNum](#antenna_attributes.read_global_channel_mapping))
39 * If no antennas is specified, all MI channels will be loaded into the DataFrame.
40 * An exception is raised if the list contains invalid `AntNum`s.
53 $ waveforms (volts) <array[f64, 1024]>
54 $ step size (nanoseconds) <f64>
57 from antenna_attributes
import read_global_channel_mapping, ALL_PUEO_ANTENNA_NAMES
60 channel_mapping: pl.DataFrame = read_global_channel_mapping()
63 antennas = antennas.split()
64 pl.Series(antennas, dtype=ALL_PUEO_ANTENNA_NAMES)
65 except pl.exceptions.InvalidOperationError
as exc:
66 print(
"InvalidOperationError:", exc)
67 print(
"\n\nPerhaps you included a non-existent antenna index?")
68 print(
"aborting operation at function: load_waveforms()")
72 read_global_channel_mapping().filter(pl.col(
"AntNum").is_in(antennas))
76 channel_mapping = channel_mapping.filter(pl.col(
"AntType") !=
"LF")
78 waveforms = channel_mapping.select(
79 pl.all().exclude(
"Comment"),
81 pl.col(
"GlobalChan").map_batches(
82 lambda s: np.array([dataset.useful().volts[i]
for i
in s])
83 ).alias(
"waveforms (volts)"),
85 pl.col(
"GlobalChan").map_batches(
86 lambda s: np.array([dataset.useful().dt[i]
for i
in s])
87 ).alias(
"step size (nanoseconds)")
93def upsample_waveforms(waveforms: pl.DataFrame, upsample_factor: int) -> pl.DataFrame:
95 [resample()](https://docs.scipy.org/doc/scipy/reference/generated/scipy.signal.resample.html)
96 to upsample the waveforms.
98 @param[in] waveforms The output of #load_waveforms
99 @param[in] upsample_factor
100 @retval upsampled_waveforms
102 * The following columns are required in the input:
103 -# `waveforms (volts)`
104 -# `step size (nanoseconds)`
105 * No additional columns will be attached to the output,
106 `waveforms (volts)` and `step size (nanoseconds)` will be upsampled and replaced.
109 from scipy.signal
import resample
110 signal_length = len(waveforms[
"waveforms (volts)"][0])
112 upsampled_waveforms = waveforms.with_columns(
113 pl.col(
"waveforms (volts)").map_batches(
114 lambda wf: resample(wf, signal_length * upsample_factor, axis=1)
116 pl.col(
"step size (nanoseconds)") / upsample_factor
119 return upsampled_waveforms
122def plot_waveforms(waveforms: pl.DataFrame, plot_name: str) ->
None:
123 r"""!Uses Vega-Altair to make waveform plots.
125 @param[in] waveforms The output of any other functions in @ref load_waveforms
126 @param[in] plot_name Be sure to specify the file format (eg. `my_plot.png` or `my_plot.html`)
128* The following columns are required in the input:
129 * `waveforms (volts)`
130 * `step size (nanoseconds)`
135* Example: Plot the vpol channels of three antennas using the last event in `run0/`.\n
136 Place the figure in the `/tmp` directory.
138run_zero_data = load_pueoEvent_Dataset(pueo_mc_data=Path("/home/all_my_runs"), run_number=0)
140wf = load_waveforms(dataset=run_zero_data, antennas="403 111 201")
141plot_waveforms(wf.filter(pl.col("Pol") == "V"), plot_name="/tmp/waveforms.html")
143 @note The actual file name will be changed to `/tmp/Vpol_waveforms.html`.
144 Had we chosen not to filter out the horizontally polarized channels,
145 there will be another plot: `/tmp/Hpol_waveforms.html`.
148 <iframe src="example_three_antenna_waveforms_from_pueoEvent.html" width="950" height="1000"></iframe>
154 waveforms.select(
"PhiSector",
"Ring")
156 time = np.zeros_like(waveforms[
"waveforms (volts)"])
157 signal_length = time.shape[1]
159 for i, dt
in enumerate(waveforms[
"step size (nanoseconds)"]):
160 time[i] = dt * np.arange(signal_length)
163 waveforms.with_columns(pl.Series(time).alias(
"time (ns)"))
164 .explode(
"waveforms (volts)",
"time (ns)")
168 plot_path = Path(plot_name)
169 destination = plot_path.parent
172 for frame
in waveforms_list:
173 pol = frame[
"Pol"][0]
176 alt.Chart(frame).mark_line()
179 y=alt.Y(
'waveforms (volts)').title(
"[Volts]"),
182 column=alt.Column(
"PhiSector:O", header=alt.Header(labelFontSize=16, titleFontSize=16)),
183 row=alt.Row(
"Ring:O", header=alt.Header(labelFontSize=16, titleFontSize=16))
185 .properties(title=f
"{pol}Pol Waveforms")
186 .configure_title(fontSize=16)
189 pn = destination / f
"{pol}pol_{plot_path.name}"
190 chart.save(f
"{str(pn)}")
193def load_waveforms_from_IceFinal(filepath: str, event_number: [int] = 0,
194 antennas: str =
None, upsample_factor: int =
None,
195 want_plots: bool =
False,
196 plot_name: str =
None) -> tuple[pl.DataFrame(), np.array(1024)]:
197 r"""! Loads time domain signals from an IceFinal file;
198 users can optionally provide a list of antenna indices.
201 @param[in] filepath path to `IceFinal_*_allTree.root`
202 @param[in] event_number (optional) an integer or a **list** of **integers**, defaults to event 0
203 @param[in] antennas (optional) **whitespace separated string** specifying the antenna names
204 @param[in] upsample_factor (optional)
205 @param[in] want_plots (optional) defaults to False
206 @param[in] plot_name (optional)
207 @retval waveforms_frame waveforms recorded by the antennas [Volts]
208 @retval time_mesh time mesh corresponding to the waveforms [seconds].
210* The function will error out if the `passTree` of the IceFinal file (see `filepath`) is empty.
211 (ie. no triggered events)
212* If no antennas are specified, all waveforms will be loaded into the DataFrame,
213 which presents them in [ascending](#antenna_attributes.ALL_PUEO_ANTENNA_NAMES) `AntNum`.
216 load_waveforms_from_IceFinal("IceFinal_420_allTree.root", antennas="403 111 201", event_number=[0, 1])
219 ┌──────────────┬────────┬─────────────┬─────────────────────┬───────────────────────┐
220 │ event number ┆ AntNum ┆ pueoSim idx ┆ VPol waveforms (V) ┆ truth (signal origin) │
221 │ --- ┆ --- ┆ --- ┆ --- ┆ --- │
222 │ i32 ┆ enum ┆ u8 ┆ array[f64, 1024] ┆ array[f64, 3] │
223 ╞══════════════╪════════╪═════════════╪═════════════════════╪═══════════════════════╡
224 │ 0 ┆ 201 ┆ 1 ┆ […] ┆ […] │
225 │ 0 ┆ 403 ┆ 11 ┆ […] ┆ […] │
226 │ 0 ┆ 111 ┆ 40 ┆ […] ┆ […] │
227 │ 1 ┆ 201 ┆ 1 ┆ […] ┆ […] │
228 │ 1 ┆ 403 ┆ 11 ┆ […] ┆ […] │
229 │ 1 ┆ 111 ┆ 40 ┆ […] ┆ […] │
230 └──────────────┴────────┴─────────────┴─────────────────────┴───────────────────────┘
232* The waveforms can be upsampled by @p upsample_factor.
233* `truth (signal origin)` is a unit vector that represents the true direction of the incoming
235 * Note that in the IceFinal files, the variable `RFdir_payload` is the direction the signal is
237 * Therefore, to obtain the direction that the signal is coming **from**,
238 we need the opposite vector: \n
239 `signal_origin` = - `RFdir_payload`
240 * This is useful for verifying the correctness of the [correlation maps](@ref example_corr_map).
242* The function also optionally creates waveform plots and saves to the working directory
243 (ie. wherever the user calls this script) via the internal function #plot_waveforms().
244 Below is an example of a three-antenna plot (actually two plots would be made, one for each event):
246 load_waveforms_from_IceFinal(
247 "IceFinal_420_allTree.root", antennas="403 111 201",
248 want_plots=True, event_number=[0, 1])
250 \anchor exampleThreeAntennaPlot
252 <iframe src="example_three_antenna_waveforms.html" width="950" height="1000"></iframe>
254 @todo **currently main instrument only**
256 @warning Will abort if non-existent antenna indices are provided.
258 load_waveforms_from_IceFinal("IceFinal_420_allTree.root", antennas="201 511")
260 In this case, Antenna `511` does not exist, since there are only 4 antennas in each `PhiSector`.\n
262 > InvalidOperationError: conversion from `str` to `enum` failed in column '' for 1 out of 2 values: ["511"]
264 > Ensure that all values in the input column are present in the categories of the enum datatype.
266 > Perhaps you included a non-existent antenna index?
267 > aborting operation at function: load_waveforms_from_IceFinal()
269 from antenna_attributes
import generate_dataframe_from_pueodata_qrh, ALL_PUEO_ANTENNA_NAMES
272 if isinstance(event_number, int):
273 event_number = [event_number]
276 if antennas
is not None:
279 antennas = antennas.split()
280 pl.Series(antennas, dtype=ALL_PUEO_ANTENNA_NAMES)
281 except pl.exceptions.InvalidOperationError
as exc:
282 print(
"InvalidOperationError:", exc)
283 print(
"\n\nPerhaps you included a non-existent antenna index?")
284 print(
"aborting operation at function: load_waveforms_from_IceFinal()")
287 _index_conversion = (
288 generate_dataframe_from_pueodata_qrh().filter(pl.col(
"AntNum").is_in(antennas))
291 _index_conversion: pl.DataFrame = generate_dataframe_from_pueodata_qrh()
293 ROOT.gSystem.Load(
"libNiceMC.so")
294 ROOT.gSystem.Load(
"libPueoSim.so")
296 _rootFile: ROOT.TFile = ROOT.TFile.Open(filepath)
297 _passTree: ROOT.TTree = _rootFile.Get(
'passTree')
299 waveforms: [pl.DataFrame] = []
300 for ev
in event_number:
302 if ev < 0
or ev >= _passTree.GetEntries():
303 print(f
"Invalid event index \"{ev}\" specified for {filepath}")
304 print(f
"There are only {_passTree.GetEntries()} events.")
307 _passTree.GetEvent(ev)
308 _df: pl.DataFrame = (
311 pl.lit(ev).alias(
"event number"),
312 pl.col(
"AntNum",
"pueoSim idx"),
313 pl.col(
"pueoSim idx").map_batches(
315 _passTree.detectorEvents[0].waveformsV[i].GetY()
318 ).alias(
"VPol waveforms (V)"),
319 pl.col(
"pueoSim idx").map_batches(
321 -(_passTree.detectorEvents[0].RFdir_payload)
324 ).alias(
"truth (signal origin)"),
327 waveforms.append(_df)
329 time_mesh = np.array(_passTree.detectorEvents[0].waveformsV[0].GetX())
332 waveforms = pl.concat(waveforms)
339 if upsample_factor
is not None:
340 from scipy.signal
import resample
341 upsampled_waveforms = waveforms.with_columns(
342 pl.col(
"VPol waveforms (V)").map_batches(
343 lambda s: resample(s, len(time_mesh) * upsample_factor, axis=1)
346 _, upsampled_time_mesh = resample(waveforms[
"VPol waveforms (V)"][0],
347 num=len(time_mesh) * upsample_factor, t=time_mesh)
349 return upsampled_waveforms, upsampled_time_mesh
351 return waveforms, time_mesh
354if __name__ ==
"__main__":
363 myDataset: ROOT.pueo.Dataset = load_pueoEvent_Dataset(pueo_mc_data=Path(
"/tmp"), run_number=0)
364 print(f
"The current run is {myDataset.getCurrRun()}")
367 if not myDataset.loadRun(1, ROOT.pueo.Dataset.PUEO_MC_DATA):
370 print(f
"The current run is {myDataset.getCurrRun()}")
373 if not myDataset.getEntry(14):
376 print(f
"Loaded event {myDataset.current()}")
384 wf: pl.DataFrame = load_waveforms(dataset=myDataset, antennas=
"403 111 201")
385 plot_waveforms(wf.filter(pl.col(
"Pol") ==
"V"), plot_name=
"./three_waveforms.html")
388 wf: pl.DataFrame = load_waveforms(dataset=myDataset)
389 up: pl.DataFrame = upsample_waveforms(waveforms=wf, upsample_factor=3)
390 plot_waveforms(up, plot_name=
"./all_waveforms.html")