pueoAnalysisTools
Loading...
Searching...
No Matches
Phase Center Calibration

Applies SciPy's function minimize() to calibrate antenna phase centers. More...

Functions

pl.DataFrame generate_antenna_phase_center_pairs_with_placeholders (Path qrh_dot_dat)
 Creates antenna pairs with placeholder integers as antenna phase centers.
 
[int] __get_all_run_numbers (Path pueo_mc_data)
 A temporary helper function of load_dataset_and_measure_time_delays()
 
pl.DataFrame load_dataset_and_measure_time_delays (ROOT.pueo.Dataset dataset, pl.DataFrame phase_center_pairs)
 Uses SciPy's correlate() to "measure" the delays of waveforms between antenna pairs.
 
pl.DataFrame _scipy_cross_correlate_all_pairs_for_one_event (pl.DataFrame phase_center_pairs, pl.DataFrame waveforms)
 A helper function of load_dataset_and_measure_time_delays()
 
 retrieve_relevant_antennas_from_all_pairs (pl.DataFrame phase_center_pairs)
 Retrieves the (unique) antennas, given a bunch of antenna pairs.
 
int objective_function (pl.Series phase_center_guess, pl.Series phase_center_placeholders, pl.DataFrame pulser_directions_and_measured_time_delays, str axis)
 Computes the expected time delays and compares them with the measured time delays, designed to be run by SciPy's minimizer.
 
 _plot_calibration_result (str axis, pl.DataFrame before, np.ndarray after)
 

Detailed Description

Applies SciPy's function minimize() to calibrate antenna phase centers.

Function Documentation

◆ generate_antenna_phase_center_pairs_with_placeholders()

pl.DataFrame generate_antenna_phase_center_pairs_with_placeholders ( Path qrh_dot_dat)

Creates antenna pairs with placeholder integers as antenna phase centers.

Parameters
[in]qrh_dot_datPath to qrh.dat
Return values
nominal_phase_center_pairs
  • Output schema:
    $ A1_AntNum <enum>
    $ A1_Ring <u8>
    $ A1_Rho (placeholder) <u16>
    $ A1_Phi (placeholder) <u16>
    $ A1_Z (placeholder) <u16>
    $ A1_Rho <f64>
    $ A1_Phi <f64>
    $ A1_Z <f64>
    $ A2_AntNum <enum>
    $ A2_Ring <u8>
    $ A2_Rho (placeholder) <u16>
    $ A2_Phi (placeholder) <u16>
    $ A2_Z (placeholder) <u16>
    $ A2_Rho <f64>
    $ A2_Phi <f64>
    $ A2_Z <f64>
  • Details:
    1. Calls #antenna_attributes.read_MI_antenna_geometry() to read qrh.dat to obtain antenna face center positions
    2. Calls #antenna_attributes.get_MI_nominal_phase_center() to convert face centers to nominal phase centers
    3. Adds three columns to store placeholder values of the phase centers; these integer placeholders serve as dictionary keys.
      shape: (96, 9)
      ┌────────┬────────┬──────┬────────────────┬────────────────┬───────────────┬──────────┬───────────┬───────┐
      │ AntNum ┆ AntIdx ┆ Ring ┆ Rho ┆ Phi ┆ Z ┆ Rho ┆ Phi ┆ Z │
      │ --- ┆ --- ┆ --- ┆ (placeholder) ┆ (placeholder) ┆ (placeholder) ┆ --- ┆ --- ┆ --- │
      │ enum ┆ u32 ┆ u8 ┆ --- ┆ --- ┆ --- ┆ f64 ┆ f64 ┆ f64 │
      │ ┆ ┆ ┆ u16 ┆ u16 ┆ u16 ┆ ┆ ┆ │
      ╞════════╪════════╪══════╪════════════════╪════════════════╪═══════════════╪══════════╪═══════════╪═══════╡
      │ 101 ┆ 0 ┆ 1 ┆ 0 ┆ 96 ┆ 192 ┆ 1.188 ┆ 0.0 ┆ 5.114 │
      │ 201 ┆ 1 ┆ 2 ┆ 1 ┆ 97 ┆ 193 ┆ 2.296 ┆ 0.0 ┆ 1.457 │
      │ 301 ┆ 2 ┆ 3 ┆ 2 ┆ 98 ┆ 194 ┆ 2.296 ┆ 0.0 ┆ 0.728 │
      │ 401 ┆ 3 ┆ 4 ┆ 3 ┆ 99 ┆ 195 ┆ 2.296 ┆ 0.0 ┆ 0.0 │
      │ 102 ┆ 4 ┆ 1 ┆ 4 ┆ 100 ┆ 196 ┆ 1.076821 ┆ 0.261538 ┆ 4.476 │
      │ … ┆ … ┆ … ┆ … ┆ … ┆ … ┆ … ┆ … ┆ … │
      │ 423 ┆ 91 ┆ 4 ┆ 91 ┆ 187 ┆ 283 ┆ 2.295825 ┆ -0.523638 ┆ 0.0 │
      │ 124 ┆ 92 ┆ 1 ┆ 92 ┆ 188 ┆ 284 ┆ 1.07708 ┆ -0.26224 ┆ 4.476 │
      │ 224 ┆ 93 ┆ 2 ┆ 93 ┆ 189 ┆ 285 ┆ 2.295502 ┆ -0.261893 ┆ 1.457 │
      │ 324 ┆ 94 ┆ 3 ┆ 94 ┆ 190 ┆ 286 ┆ 2.295502 ┆ -0.261893 ┆ 0.728 │
      │ 424 ┆ 95 ┆ 4 ┆ 95 ┆ 191 ┆ 287 ┆ 2.295502 ┆ -0.261893 ┆ 0.0 │
      └────────┴────────┴──────┴────────────────┴────────────────┴───────────────┴──────────┴───────────┴───────┘
    4. Calls antenna_pairs.generate_MI_antenna_pairs with the above dataframe as the input.

Definition at line 21 of file calibration.py.

◆ __get_all_run_numbers()

[int] __get_all_run_numbers ( Path pueo_mc_data)
private

A temporary helper function of load_dataset_and_measure_time_delays()

Parameters
[in]pueo_mc_dataSee the parameter of initialise.load_pueoEvent_Dataset
  • This function retrieves the run numbers (negative sign allowed) of the run/ subdirectories in pueo_mc_data
  • The runs are returned as a list of integers.
  • The function is here only for testing purposes, as it makes little sense to use every single run for antenna calibration.

Definition at line 103 of file calibration.py.

◆ load_dataset_and_measure_time_delays()

pl.DataFrame load_dataset_and_measure_time_delays ( ROOT.pueo.Dataset dataset,
pl.DataFrame phase_center_pairs )

Uses SciPy's correlate() to "measure" the delays of waveforms between antenna pairs.

Parameters
[in]datasetThe output of initialise.load_pueoEvent_Dataset
[in]phase_center_pairsThe output of generate_antenna_phase_center_pairs_with_placeholders
Return values
pulser_directions_and_measured_time_delays
  • The input dataset should have been initialized with a (valid) run number, but that is not the only run that gets loaded. For now, all runs are used.
  • Output schema:
    $ measured time delays (ns) <f64>
    $ max correlation <f64>
    $ A1_AntNum <enum>
    $ A1_Ring <u8>
    $ A1_Rho (placeholder) <u16>
    $ A1_Phi (placeholder) <u16>
    $ A1_Z (placeholder) <u16>
    $ A1_Rho <f64>
    $ A1_Phi <f64>
    $ A1_Z <f64>
    $ A2_AntNum <enum>
    $ A2_Ring <u8>
    $ A2_Rho (placeholder) <u16>
    $ A2_Phi (placeholder) <u16>
    $ A2_Z (placeholder) <u16>
    $ A2_Rho <f64>
    $ A2_Phi <f64>
    $ A2_Z <f64>
    $ incoming pulse direction (Phi Theta) <array[f64, 2]>
  • The dataset is also used to load the directions of the incoming pulses ( \(^\circ\phi\) and \(^\circ\theta\) in payload coordinates).
  • measured time delays (ns) is computed by the helper function _scipy_cross_correlate_all_pairs_for_one_event().

Definition at line 124 of file calibration.py.

◆ _scipy_cross_correlate_all_pairs_for_one_event()

pl.DataFrame _scipy_cross_correlate_all_pairs_for_one_event ( pl.DataFrame phase_center_pairs,
pl.DataFrame waveforms )
protected

A helper function of load_dataset_and_measure_time_delays()

Parameters
[in]phase_center_pairsThe output of generate_antenna_phase_center_pairs_with_placeholders
[in]waveformsThe output of waveform_plots.load_waveforms
Return values
measured_time_delays
  • Output schema:
    $ measured time delays (ns) <f64>
    $ max correlation <f64>
    $ A1_AntNum <enum>
    $ A1_Ring <u8>
    $ A1_Rho (placeholder) <u16>
    $ A1_Phi (placeholder) <u16>
    $ A1_Z (placeholder) <u16>
    $ A1_Rho <f64>
    $ A1_Phi <f64>
    $ A1_Z <f64>
    $ A2_AntNum <enum>
    $ A2_Ring <u8>
    $ A2_Rho (placeholder) <u16>
    $ A2_Phi (placeholder) <u16>
    $ A2_Z (placeholder) <u16>
    $ A2_Rho <f64>
    $ A2_Phi <f64>
    $ A2_Z <f64>
  • Details:
    1. Load the waveforms of both antennas in the antenna pair
    2. Cross correlate the waveform pair
    3. Figure out by how much one needs to shift the second waveform against the first in order to achieve the maximum correlation.
    4. The above is called measured time delays (ns) in the output dataframe.
    5. For more details, see Explanation.

Definition at line 203 of file calibration.py.

◆ retrieve_relevant_antennas_from_all_pairs()

retrieve_relevant_antennas_from_all_pairs ( pl.DataFrame phase_center_pairs)

Retrieves the (unique) antennas, given a bunch of antenna pairs.

Parameters
[in]phase_center_pairsThe output of load_dataset_and_measure_time_delays
Return values
relevant_antennas
  • This step is not trivial, as explained below.
  • It is conceivable that not all antennas receive signals clear enough for calibration purposes (resulting in weak correlations).
  • As such, the weakly correlated antenna pairs should be filtered out. For example,
    pulser_directions_and_pairwise_time_delays = (
    load_dataset_and_measure_time_delays(dataset=ds, phase_center_pairs=nominal_pairs)
    .filter(pl.col("max correlation") > 0.8)
    )
  • It is possible that some antennas end up getting filtered out.
  • Sample output (note that only 76 out of 96 MI antennas survived the filter in this case):
    shape: (76, 8)
    ┌────────┬──────┬───────────────────┬───────────────────┬─────────────────┬──────────┬───────────┬───────┐
    │ AntNum ┆ Ring ┆ Rho (placeholder) ┆ Phi (placeholder) ┆ Z (placeholder) ┆ Rho ┆ Phi ┆ Z │
    │ --- ┆ --- ┆ --- ┆ --- ┆ --- ┆ --- ┆ --- ┆ --- │
    │ enum ┆ u8 ┆ u16 ┆ u16 ┆ u16 ┆ f64 ┆ f64 ┆ f64 │
    ╞════════╪══════╪═══════════════════╪═══════════════════╪═════════════════╪══════════╪═══════════╪═══════╡
    │ 101 ┆ 1 ┆ 0 ┆ 96 ┆ 192 ┆ 1.188 ┆ 0.0 ┆ 5.114 │
    │ 201 ┆ 2 ┆ 1 ┆ 97 ┆ 193 ┆ 2.296 ┆ 0.0 ┆ 1.457 │
    │ 301 ┆ 3 ┆ 2 ┆ 98 ┆ 194 ┆ 2.296 ┆ 0.0 ┆ 0.728 │
    │ 401 ┆ 4 ┆ 3 ┆ 99 ┆ 195 ┆ 2.296 ┆ 0.0 ┆ 0.0 │
    │ 102 ┆ 1 ┆ 4 ┆ 100 ┆ 196 ┆ 1.076821 ┆ 0.261538 ┆ 4.476 │
    │ … ┆ … ┆ … ┆ … ┆ … ┆ … ┆ … ┆ … │
    │ 423 ┆ 4 ┆ 91 ┆ 187 ┆ 283 ┆ 2.295825 ┆ -0.523638 ┆ 0.0 │
    │ 124 ┆ 1 ┆ 92 ┆ 188 ┆ 284 ┆ 1.07708 ┆ -0.26224 ┆ 4.476 │
    │ 224 ┆ 2 ┆ 93 ┆ 189 ┆ 285 ┆ 2.295502 ┆ -0.261893 ┆ 1.457 │
    │ 324 ┆ 3 ┆ 94 ┆ 190 ┆ 286 ┆ 2.295502 ┆ -0.261893 ┆ 0.728 │
    │ 424 ┆ 4 ┆ 95 ┆ 191 ┆ 287 ┆ 2.295502 ┆ -0.261893 ┆ 0.0 │
    └────────┴──────┴───────────────────┴───────────────────┴─────────────────┴──────────┴───────────┴───────┘
  • If we still pass these filtered antennas to the objective_function, SciPy's minimizer would act as if these antennas are relevant and try to still modify their phase centers during the minization.

Definition at line 285 of file calibration.py.

◆ objective_function()

int objective_function ( pl.Series phase_center_guess,
pl.Series phase_center_placeholders,
pl.DataFrame pulser_directions_and_measured_time_delays,
str axis )

Computes the expected time delays and compares them with the measured time delays, designed to be run by SciPy's minimizer.

Parameters
[in]phase_center_guessThe output of retrieve_relevant_antennas_from_all_pairs
[in]phase_center_placeholdersThe output of retrieve_relevant_antennas_from_all_pairs
[in]pulser_directions_and_measured_time_delaysThe output of load_dataset_and_measure_time_delays
[in]axisAxis to fit: \(\rho\), \(\phi\), or \(z\)
Return values
total_difference_bw_expected_time_delay_and_measured_time_delayA scalar
  • In load_dataset_and_measure_time_delays(), the measured_time_delay has already been computed.
  • The job of this function is to
    1. Given phase_center_guess, replace the placeholder antenna phase centers in pulser_directions_and_pairwise_time_delays
    2. Compute the pairwise antenna displacement vectors (call these \(\bf{d}\)), using the updated antenna positions.
    3. Compute the expected time delay based on the dot product between \(\bf{d}\) and incoming pulse direction (Phi Theta), see this figure for an illustration.
      ┌───────────────────────────┬───────────────────────────┐
      │ measured time delays (ns) ┆ expected time delays (ns) │
      │ --- ┆ --- │
      │ f64 ┆ f64 │
      ╞═══════════════════════════╪═══════════════════════════╡
      │ -0.866667 ┆ -0.853384 │
      │ -0.866667 ┆ -0.853384 │
      │ -7.466667 ┆ -7.487111 │
      │ 0.6 ┆ 0.611216 │
      │ -9.466667 ┆ -9.48138 │
      │ … ┆ … │
      │ 8.333333 ┆ 8.33698 │
      │ 0.2 ┆ 0.15054 │
      │ -0.866667 ┆ -0.855732 │
      │ 0.8 ┆ 0.760589 │
      │ 1.466667 ┆ 1.463433 │
      └───────────────────────────┴───────────────────────────┘
    4. Quantify the difference between the expectation and the measurement; this is done by summing the difference over all rows above.
  • SciPy's minimizer tries to minimize the total difference by updating phase_center_guess.

Definition at line 345 of file calibration.py.

◆ _plot_calibration_result()

_plot_calibration_result ( str axis,
pl.DataFrame before,
np.ndarray after )
protected
Parameters
[in]axisAxis to fit: \(\rho\), \(\phi\), or \(z\)
[in]beforeThe output of of retrieve_relevant_antennas_from_all_pairs
[in]afterThe restul.x of running scipy.optimize.minimize on objective_function

Nothing crazy here!

Definition at line 427 of file calibration.py.