DetMatch
The sotodlib.coords.det_match
module allows us to map resonators from one
source to another, using information such as resonator frequency, bias-line
assignments, and pointing information. This is particularly useful to create
a map from resonators in a SMuRF tune-file, to real detectors from either a
design-file, or a handmade solutions file based on a separate tune-file.
This works by translating the detector matching problem into an instance
of the well-studied assignment problem, in which, given a bipartite graph and
edge weights, one can efficiently find the minimum-cost matching.
Here, the two sets of the bipartite graph are the two resonator sets, with
additional nodes added to represent the possibility of a resonator being
unmatched. The edge weights computed using a cost-function that uses
resonator properties to determine resonator-to-resonator costs, and
resonator-to-unmatched costs. The function
scipy.optimize.linear_sum_assignment
is then used to find the minimum-cost
match.
Usage
Below is an example of how to match a resonator-set based on a smurf tune-file and sodetlib bgmap to a resonator-set based on a handmade solution file, and save the output.
from sotodlib.coords import det_match as dm
tunefile = <path_to_tune_file>
bgmap_file = <path_to_bgmap_file>
sol_file = <path_to_solution_file>
src = dm.ResSet.from_tunefile(
tune_file, north_is_highband=False, bgmap_file=bgmap_file,
name=f'SAT Tuning'
)
dst = dm.ResSet.from_solutions(
sol_file, north_is_highband=True, name='solution'
)
match = dm.Match(src, dst)
match.save('match.h5')
To validate, you can check match.stats
to see information such as how many
resonators in each set were matched or left unassigned, and to see how many
matches had bias-line mismaps, etc. If the majority of resonators have bias-line
mismaps, then it is likely the north_is_highband
flag was set incorrectly,
and you are attempting to match opposite sides of the UFM.
To validate, you can also use the dm.plot_match_freqs
to view how well the
resonators match in frequency space. With a proper solution file, this should
be rather well, as is seen below
Tuning Matching Params
Depending on what data you have available, you may want to tune the matching algorithm and the cost function. For instance, once we have accurate pointing data we’ll want to rely less solely on the frequency pairings, and more on pointing information. We can do this by passing in a MatchParams object to set various parameters. For instance, below is an example of how to tell the match function to be stricter with the pointing cost penalties, and more lenient with the frequency penalty:
from sotodlib.io.coords import det_match as dm
tunefile = <path_to_tune_file>
bgmap_file = <path_to_bgmap_file>
sol_file = <path_to_solution_file>
src = dm.ResSet.from_tunefile(
tune_file, north_is_highband=False, bgmap_file=bgmap_file,
name=f'SAT Tuning'
)
dst = dm.ResSet.from_solutions(
sol_file, north_is_highband=True, name='solution'
)
mpars = dm.MatchParams(
freq_width=5 #MHz
dist_width=np.deg2rad(0.3), # radians
)
match = dm.Match(src, dst, match_params=mpars)
API
- sotodlib.coords.det_match.map_band_chans(b1, c1, b2, c2, chans_per_band=512)[source]
Returns an index mapping of length nchans1 from one set of bands and channels to another. Note that unmapped indices are returned as -1, so this must be handled before (or after) indexing or else you’ll get weird results. :param b1: Array of length nchans1 containing the smurf band of each channel :type b1: np.ndarray :param c1: Array of length nchans1 containing the smurf channel of each channel :type c1: np.ndarray :param b2: Array of length nchans2 containing the smurf band of each channel :type b2: np.ndarray :param c2: Array of length nchans2 containing the smurf channel of each channel :type c2: np.ndarray :param chans_per_band: Lets just hope this never changes. :type chans_per_band: int
- sotodlib.coords.det_match.get_north_is_highband(bands, bgs)[source]
Checks if north is highband based on bgmapping. This will tell you if the majority of dets on the north side of the ufm (bgs 0-5) belong to highband (bands 4-7).
- class sotodlib.coords.det_match.PointingConfig(fp_file: str, wafer_slot: str, tel_type: str, zemax_path: str | None = None)[source]
Bases:
object
Helper class for getting pointing info from an optics model.
- Parameters:
fp_file (str) – Path to focal-plane file that is used by the optics module.
wafer_slot (str) – Wafer slot of the UFM. For example: “ws0”
tel_type (str) – Tel type for the optics model. Either “SAT” or “LAT”
zemax_path (str) – If running for a “LAT” tel_type, the path to the zemax file must be specified.
- class sotodlib.coords.det_match.Resonator(idx: int, is_north: int, res_freq: float, res_qi: float = nan, smurf_res_idx: int = -1, smurf_band: int = -1, smurf_channel: int = -1, smurf_subband: int = -1, readout_id: str = '', xi: float = nan, eta: float = nan, gamma: float = nan, bg: int = -1, det_x: float = nan, det_y: float = nan, det_row: int = 0, det_col: int = 0, pixel_num: int = 0, det_rhomb: str = '', det_pol: str = '', det_freq: int = 0, det_bandpass: str = '', det_angle_raw_deg: float = nan, det_angle_actual_deg: float = nan, det_type: str = '', det_id: str = 'NO_MATCH', is_optical: int = 1, mux_bondpad: int = 0, mux_subband: str = '', mux_band: int = -1, mux_channel: int = -1, mux_layout_pos: int = -1, matched: int = 0, match_idx: int = -1)[source]
Bases:
object
Data structure to hold any resonator information.
- sotodlib.coords.det_match.apply_design_properties(smurf_res, design_res, in_place=False, apply_pointing=True)[source]
Combines two resonators into one, taking smurf-properties such as res-idx, smurf-band, and smurf-channel one, and design properties such as det position and polarization from the other.
- class sotodlib.coords.det_match.ResSet(resonances: List[Resonator], name=None)[source]
Bases:
object
Class to hold a group of resonances. This provides easy interfaces for accessing Resonance fields as np arrays for fast computations and provides initialization functions from different data sources.
- classmethod from_array(arr, name=None, ignore_extra_fields=True)[source]
Creates a ResSet from a numpy structured array (resulting from
as_array
method).- Parameters:
- classmethod from_aman(aman, stream_id, det_cal=None, name=None)[source]
Load a resonator set from a Context object based on an obs_id
- Parameters:
aman (AxisManager) – Axis manager containing metadata
stream_id (str) – Stream id for ResSet to load
det_cal (AxisManager) – Detector calibration metadata. If not specified, will default to
aman.det_cal
- classmethod from_tunefile(tunefile, name=None, north_is_highband=True, resfit_file=None, bgmap_file=None)[source]
Creates an instance based on a smurf-tune file. If a resfit or bgmap file is included, that data will be added to the Resonance objects as well.
- classmethod from_wafer_info_file(wafer_info_file, array_name, name=None, pt_cfg: PointingConfig | None = None)[source]
Initialize a ResSet from a wafer info file. This is a file that contains detector design information.
- Parameters:
wafer_info_file (str) – Path to wafer info file
array_name (str) – Array name, which is the key in the wafer-info-file. For example: “mv7”.
name (str) – Name to assign to the ResSet
pt_cfg (PointingConfig) – If set, this will be used to get pointing info based on the optics model. If not set, pointing info will not be included in the ResSet.
- classmethod from_solutions(sol_file, north_is_highband=True, name=None, fp_pars=None, platform='SAT', zemax_path=None)[source]
Creates an instance from an input-solution file. This will include both design data, along with smurf-band and smurf-channel info. Resonance frequencies used here are the VNA freqs measured by Kaiwen.
- Parameters:
sol_file (str) – Path to solutions file
name (str) – Name to label this ResSet
north_is_highband (bool) – True if the north-side of the array corresponds to bands 4-8
fp_pars (dict) – Result of the function
sotododlib.coords.optics.get_ufm_to_fp_pars
. If this is None, detector positions will not be mapped to pointing angles.platform (str) – ‘SAT’ or ‘LAT’. Used to determine which focal plane function to use for pointing
zemax_path (str) – zemax path, required to get pointing for LAT optics
- class sotodlib.coords.det_match.MatchParams(unassigned_slots: int = 1000, freq_offset_mhz: float = 0.0, freq_width: float = 2.0, dist_width: float = 0.01, unmatched_good_res_pen: float = 10.0, good_res_qi_thresh: float = 100000.0, force_src_pointing: bool = False, assigned_bg_unmatched_pen: float = 100000, unassigned_bg_unmatched_pen: float = 10000, assigned_bg_mismatch_pen: float = 100000, unassigned_bg_mismatch_pen: float = 1)[source]
Bases:
object
Any constants / hardcoded values go here to be fiddled with
- Parameters:
unassigned_slots (int) – Number of additional “unassigned” node to use per-side
freq_offset_mhz (float) – constant offset between resonator-set frequencies to use in matching.
freq_width (float) – width of exponential to use in the frequency cost function (MHz).
dist_width (float) – width of exponential to use in the pointing cost function (rad)
unmatched_good_res_pen (float) – penalty to apply to leaving a resonator with a good qi unassigned
good_res_qi_thresh (float) – qi threshold that is considered “good”
force_src_pointing (bool) – If true, will assign a np.inf penalty to leaving a src resonator with a provided pointing unmatched.
assigned_bg_unmatched_pen (float) – Penalty to apply to leaving a resonator with an assigned bg unmatched
unassigned_bg_unmatched_pen (float) – Penalty to apply to leaving a resonator with an unassigned bg unmatched
assigned_bg_mismatch_pen (float) – Penalty to apply for matching a resonator with an assigned bias line to another one with a mis-matched bias line.
unassigned_bg_mismatch_pen (float) – Penalty to apply for matching a resonator with no bias line to another one with a mis-matched bias line.
- class sotodlib.coords.det_match.MatchingStats(unmatched_src: int = 0, unmatched_dst: int = 0, unmatched_src_with_pointing: int = 0, matched_chans: int = 0, mismatched_bg: int = 0, freq_diff_avg: float = 0.0, freq_err_avg: float = 0.0, pointing_err_avg: float = 0.0)[source]
Bases:
object
- class sotodlib.coords.det_match.Match(src: ResSet, dst: ResSet, match_pars: MatchParams | None = None, apply_dst_pointing=True)[source]
Bases:
object
Class for performing a Resonance Matching between two sets of resonators, labeled src and dst. In the matching algorithm there is basically no difference between src and dst res-sets, except:
When merged, smurf-data such as band, channel, and res-idx will be taken from the
src
res-setThe
force_src_pointing
param can be used to assign a very high penalty to leaving any src resonator that has pointing info unassigned.
- Parameters:
src (ResSet) – The source resonator set
dst (ResSet) – The dest resonator set
match_pars (MatchParams) – MatchParams object used in the matching algorithm. This can be used to tune the cost-function and matching.
apply_dst_pointing (bool) – If True, the
merged
res-set will take its pointing information fromdst
instead ofsrc
.
- match_pars
MatchParams object used in the matching algorithm. This can be used to tune the cost-function and matching.
- Type:
- matching
A 2xN array of indices, where the first row corresponds to the indices of the src resonators, and the second row corresponds to the indices of the dst resonators. If the index is larger than the size of the corresponding res-set, that means the paired resonator was not matched.
- Type:
np.ndarray
- merged
A ResSet containing the merged resonators. This is created by applying the “design” properties of the dst resonators to the source resonators. If a source resonator is not matched to any dest resonator, it will be copied as-is.
- Type:
- stats
A MatchingStats object containing some statistics about the matching.
- Type:
- get_match_iter(include_unmatched=True) Iterator[Tuple[Resonator | None, Resonator | None]] [source]
Returns an iterator over matched resonators (r1, r2).
- Parameters:
include_unmatched (bool) – If True, will include unmatched resonators, with the pair set to None.
- get_stats() MatchingStats [source]
Gets stats associated with current matching.