pueoAnalysisTools
Loading...
Searching...
No Matches
waveform_plots.py
Go to the documentation of this file.
1r"""!
2@file waveform_plots.py
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.
7 ```
8 pip3 install altair[save]
9 ```
10"""
11
12# Examples are in the main function at the bottom of the script
13
14import ROOT
15import numpy as np
16import polars as pl
17from pathlib import Path
18from initialise import load_pueoEvent_Dataset
19
20
21
23
24def load_waveforms(dataset: ROOT.pueo.Dataset, antennas: str = None) -> pl.DataFrame():
25 r"""! Loads time domain signals
26 @ingroup lwf
27
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]
33
34@todo Currently main instrument only, since oftentimes LF is not simulated
35
36* `antenna`:
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.
41
42* Output schema:
43 ```
44 $ GlobalChan <u32>
45 $ AntNum <enum>
46 $ AntIdx <u32>
47 $ Pol <enum>
48 $ AntType <enum>
49 $ SurfNum <u8>
50 $ SurfChanNum <u8>
51 $ PhiSector <u8>
52 $ Ring <u8>
53 $ waveforms (volts) <array[f64, 1024]>
54 $ step size (nanoseconds) <f64>
55 ```
56 """
57 from antenna_attributes import read_global_channel_mapping, ALL_PUEO_ANTENNA_NAMES
58
59 if antennas is None:
60 channel_mapping: pl.DataFrame = read_global_channel_mapping()
61 else: # if user provides a list of antennas, filter out unused antennas
62 try: # see if the user had provided a bad antenna index
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()")
69 exit(1)
70
71 channel_mapping = (
72 read_global_channel_mapping().filter(pl.col("AntNum").is_in(antennas))
73 )
74
75 # Filter out LF channels
76 channel_mapping = channel_mapping.filter(pl.col("AntType") != "LF")
77
78 waveforms = channel_mapping.select(
79 pl.all().exclude("Comment"),
80
81 pl.col("GlobalChan").map_batches(
82 lambda s: np.array([dataset.useful().volts[i] for i in s])
83 ).alias("waveforms (volts)"),
84
85 pl.col("GlobalChan").map_batches(
86 lambda s: np.array([dataset.useful().dt[i] for i in s])
87 ).alias("step size (nanoseconds)")
88 )
89
90 return waveforms
91
92
93def upsample_waveforms(waveforms: pl.DataFrame, upsample_factor: int) -> pl.DataFrame:
94 r"""!Calls SciPy's
95 [resample()](https://docs.scipy.org/doc/scipy/reference/generated/scipy.signal.resample.html)
96 to upsample the waveforms.
97 @ingroup lwf
98 @param[in] waveforms The output of #load_waveforms
99 @param[in] upsample_factor
100 @retval upsampled_waveforms
101
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.
107
108 """
109 from scipy.signal import resample
110 signal_length = len(waveforms["waveforms (volts)"][0])
111
112 upsampled_waveforms = waveforms.with_columns(
113 pl.col("waveforms (volts)").map_batches(
114 lambda wf: resample(wf, signal_length * upsample_factor, axis=1)
115 ),
116 pl.col("step size (nanoseconds)") / upsample_factor
117 )
118
119 return upsampled_waveforms
120
121
122def plot_waveforms(waveforms: pl.DataFrame, plot_name: str) -> None:
123 r"""!Uses Vega-Altair to make waveform plots.
124 @ingroup lwf
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`)
127
128* The following columns are required in the input:
129 * `waveforms (volts)`
130 * `step size (nanoseconds)`
131 * `PhiSector`
132 * `Ring`
133 * `Pol`
134
135* Example: Plot the vpol channels of three antennas using the last event in `run0/`.\n
136 Place the figure in the `/tmp` directory.
137```python
138run_zero_data = load_pueoEvent_Dataset(pueo_mc_data=Path("/home/all_my_runs"), run_number=0)
139run_zero_data.last()
140wf = load_waveforms(dataset=run_zero_data, antennas="403 111 201")
141plot_waveforms(wf.filter(pl.col("Pol") == "V"), plot_name="/tmp/waveforms.html")
142```
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`.
146
147 \htmlonly
148 <iframe src="example_three_antenna_waveforms_from_pueoEvent.html" width="950" height="1000"></iframe>
149 \endhtmlonly
150 """
151 import altair as alt
152
153 # this forces function to error out instead of failing silently
154 waveforms.select("PhiSector", "Ring")
155
156 time = np.zeros_like(waveforms["waveforms (volts)"])
157 signal_length = time.shape[1]
158
159 for i, dt in enumerate(waveforms["step size (nanoseconds)"]):
160 time[i] = dt * np.arange(signal_length)
161
162 waveforms_list = (
163 waveforms.with_columns(pl.Series(time).alias("time (ns)"))
164 .explode("waveforms (volts)", "time (ns)")
165 .partition_by("Pol")
166 )
167
168 plot_path = Path(plot_name)
169 destination = plot_path.parent
170
171 # make the plots, let Altair take care of subplot positions automatically
172 for frame in waveforms_list:
173 pol = frame["Pol"][0]
174
175 chart = (
176 alt.Chart(frame).mark_line()
177 .encode(
178 x='time (ns):Q',
179 y=alt.Y('waveforms (volts)').title("[Volts]"),
180 )
181 .facet(
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))
184 )
185 .properties(title=f"{pol}Pol Waveforms")
186 .configure_title(fontSize=16)
187 )
188
189 pn = destination / f"{pol}pol_{plot_path.name}"
190 chart.save(f"{str(pn)}")
191
192
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.
199 @ingroup lwf
200
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].
209
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`.
214* Example:
215 ```
216 load_waveforms_from_IceFinal("IceFinal_420_allTree.root", antennas="403 111 201", event_number=[0, 1])
217 ```
218 ```
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 └──────────────┴────────┴─────────────┴─────────────────────┴───────────────────────┘
231 ```
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
234 signal.
235 * Note that in the IceFinal files, the variable `RFdir_payload` is the direction the signal is
236 travelling **to**.
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).
241
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):
245 ```python
246 load_waveforms_from_IceFinal(
247 "IceFinal_420_allTree.root", antennas="403 111 201",
248 want_plots=True, event_number=[0, 1])
249 ```
250 \anchor exampleThreeAntennaPlot
251 \htmlonly
252 <iframe src="example_three_antenna_waveforms.html" width="950" height="1000"></iframe>
253 \endhtmlonly
254 @todo **currently main instrument only**
255
256 @warning Will abort if non-existent antenna indices are provided.
257 ```
258 load_waveforms_from_IceFinal("IceFinal_420_allTree.root", antennas="201 511")
259 ```
260 In this case, Antenna `511` does not exist, since there are only 4 antennas in each `PhiSector`.\n
261 Thus, you would see
262 > InvalidOperationError: conversion from `str` to `enum` failed in column '' for 1 out of 2 values: ["511"]
263 >
264 > Ensure that all values in the input column are present in the categories of the enum datatype.
265 >
266 > Perhaps you included a non-existent antenna index?
267 > aborting operation at function: load_waveforms_from_IceFinal()
268 """
269 from antenna_attributes import generate_dataframe_from_pueodata_qrh, ALL_PUEO_ANTENNA_NAMES
270
271 # turn int into a list[int] so that it's iterable
272 if isinstance(event_number, int):
273 event_number = [event_number]
274
275 # if user provides a list of antennas, filter out unused antennas
276 if antennas is not None:
277 # see if the user had provided a bad antenna index
278 try:
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()")
285 exit(1)
286
287 _index_conversion = (
288 generate_dataframe_from_pueodata_qrh().filter(pl.col("AntNum").is_in(antennas))
289 )
290 else:
291 _index_conversion: pl.DataFrame = generate_dataframe_from_pueodata_qrh()
292
293 ROOT.gSystem.Load("libNiceMC.so")
294 ROOT.gSystem.Load("libPueoSim.so")
295
296 _rootFile: ROOT.TFile = ROOT.TFile.Open(filepath)
297 _passTree: ROOT.TTree = _rootFile.Get('passTree')
298
299 waveforms: [pl.DataFrame] = []
300 for ev in event_number:
301
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.")
305 exit(1)
306
307 _passTree.GetEvent(ev)
308 _df: pl.DataFrame = (
309 _index_conversion
310 .select(
311 pl.lit(ev).alias("event number"),
312 pl.col("AntNum", "pueoSim idx"),
313 pl.col("pueoSim idx").map_batches(
314 lambda s: np.array([
315 _passTree.detectorEvents[0].waveformsV[i].GetY()
316 for i in s
317 ])
318 ).alias("VPol waveforms (V)"),
319 pl.col("pueoSim idx").map_batches(
320 lambda s: np.array([
321 -(_passTree.detectorEvents[0].RFdir_payload)
322 for _ in s
323 ])
324 ).alias("truth (signal origin)"),
325 )
326 )
327 waveforms.append(_df)
328 # time mesh in seconds
329 time_mesh = np.array(_passTree.detectorEvents[0].waveformsV[0].GetX())
330
331 # turn the list of waveform dataframes into one dataframe and apply SciPy upsampling if needed
332 waveforms = pl.concat(waveforms)
333
334 # if want_plots:
335 # _plot_waveforms_frame( # plot in nanoseconds
336 # _explode_waveforms_frame(waveforms, time_mesh * 1e9), plot_name=plot_name
337 # )
338
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)
344 )
345 )
346 _, upsampled_time_mesh = resample(waveforms["VPol waveforms (V)"][0],
347 num=len(time_mesh) * upsample_factor, t=time_mesh)
348
349 return upsampled_waveforms, upsampled_time_mesh
350
351 return waveforms, time_mesh
352
353
354if __name__ == "__main__":
355 # Below we illustrate the basic features of the Dataset class,
356 # and how to pass it to the analysis functions.
357 # For more pueo::Dataset class methods, see its header file:
358 # https://github.com/PUEOCollaboration/pueoEvent/blob/main/src/pueo/Dataset.h
359
360 # Use run zero (for example) to initialize the Dataset
361 # Warning: make sure that the run is not empty (has events),
362 # otherwise pueo::Dataset's constructor will crash
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()}")
365
366 # suppose we now want to use another run, maybe run one. Warning: `/tmp/run1/` has to exist!
367 if not myDataset.loadRun(1, ROOT.pueo.Dataset.PUEO_MC_DATA):
368 exit(1)
369 else:
370 print(f"The current run is {myDataset.getCurrRun()}")
371
372 # Suppose we want to analyze the 15th entry
373 if not myDataset.getEntry(14):
374 exit(1)
375 else:
376 print(f"Loaded event {myDataset.current()}")
377
378 # We could also analyze the last event
379 # myDataset.last()
380
381 # ---------------------------- Examples ---------------------------- #
382
383 # Example 1: plot the VPol waveforms of just three antennas:
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")
386
387 # Example 2: plot upsampled VPol and HPol waveforms of all antennas
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")