API Reference

This package contains tools for working with hardware information (as needed for analysis operations) as well as timestream processing, workflow management and other operations needed for time domain data.

Hardware Properties

For the purpose of this package, “hardware” refers to all properties of the telescopes, detectors, readout, etc that are needed to simulate and analyze the data. Initially this will be a fairly basic set of information, but that will expand as the instrument characterization progresses.

Data Format

In memory, the hardware configuration is stored as a set of nested dictionaries. This is wrapped by a simple “Hardware” class that has some methods for dumping / loading and selecting a subset of detectors. Some parameters may eventually reference external data files in a format and location scheme that is yet to be determined.

class sotodlib.sim_hardware.Hardware(path=None)[source]

Class representing a specific hardware configuration.

The data is stored in a dictionary, and can be loaded / dumped to disk as well as trimmed to include only a subset of detectors.

Parameters:

path (str, optional) – If specified, configuration is loaded from this file during construction.

dump(path, overwrite=False, compress=False)[source]

Write hardware config to a TOML file.

Dump data to a TOML format file, optionally compressing the contents with gzip and optionally overwriting the file.

Parameters:
  • path (str) – The file to write.

  • overwrite (bool) – If True, overwrite the file if it exists. If False, then existing files will cause an exception.

  • compress (bool) – If True, compress the data with gzip on write.

Returns:

None

load(path)[source]

Read data from a TOML file.

The file can either be regular text or a gzipped version of a TOML file.

Parameters:

path (str) – The file to read.

Returns:

None

wafer_map()[source]

Construct wafer mapping to other auxilliary data.

Given the current data state, build dictionaries to go from wafer_slots to all other non-detector info: telescopes, tube_slots, card_slots, crate_slots, and bands. This is a convenient mapping when pruning the hardware information or doing other kinds of lookups.

Returns:

Nested dictionaries from wafers to other properties.

Return type:

(dict)

select(telescopes=None, tube_slots=None, match={})[source]

Select a subset of detectors.

Select detectors whose properties match some criteria. A new Hardware object is created and returned. If a matching expression is not specified for a given property name, then this is equivalent to selecting all values of that property.

Before selecting on detector properties, any telescope / tube_slot filtering criteria are first applied.

Each key of the “match” dictionary should be the name of a detector property to be considered for selection (e.g. band, wafer, pol, pixel). The value is a matching expression which can be:

  • A list of explicit values to match.

  • A string containing a regex expression to apply.

Example

Imagine you wanted to select all 90GHz detectors on wafers 25 and 26 which have “A” polarization and are located in pixels 20-29 (recall the “.” matches a single character):

new = hw.select(match={"wafer_slot": ["w25", "w26"],
                "band": "SAT_f090",
                "pol": "A",
                "pixel": "02."})
Parameters:
  • telescopes (str) – A regex string to apply to telescope names or a list of explicit names.

  • tube_slots (str) – A regex string to apply to tube_slot names or a list of explicit names.

  • match (dict) – The dictionary of property names and their matching expressions.

Returns:

A new Hardware instance with the selected detectors.

Return type:

(Hardware)

To get an example hardware configuration as a starting point, you can use this function:

sotodlib.sim_hardware.sim_nominal()[source]

Return a simulated nominal hardware configuration.

This returns a simulated Hardware object with the nominal instrument properties / metadata, but with an empty set of detector locations.

This can then be passed to one of the detector simulation functions to build up the list of detectors.

Returns:

Hardware object with nominal metadata.

Return type:

(Hardware)

Simulated Detectors

Given a Hardware object loaded from disk or created in memory (for example using sotodlib.hardware.sim_nominal()), detector properties can be simulated for an entire telescope with:

sotodlib.sim_hardware.sim_detectors_toast(hw, tele, tube_slots=None)[source]

Update hardware model with simulated detector positions.

Given a Hardware model, generate all detector properties for the specified telescope and optionally a subset of optics tube slots (for the LAT). The detector dictionary of the hardware model is updated in place.

This function requires the toast subpackage (and hence toast) to be importable.

Parameters:
  • hw (Hardware) – The hardware object to update.

  • tele (str) – The telescope name.

  • tube_slots (list, optional) – The optional list of tube slots to include.

Returns:

None

OR with a more realistic method:

sotodlib.sim_hardware.sim_detectors_physical_optics(hw, tele, tube_slots=None)[source]

Update hardware model with simulated detector positions.

Given a Hardware model, generate all detector properties for the specified telescope and optionally a subset of optics tube slots (for the LAT). The detector dictionary of the hardware model is updated in place.

This function uses information from physical optics simulations to estimate the location of detectors.

Parameters:
  • hw (Hardware) – The hardware object to update.

  • tele (str) – The telescope name.

  • tube_slots (list, optional) – The optional list of tube slots to include.

Returns:

None

These functions update the Hardware object in place:

hw = sim_nominal()
# Adds detectors to input Hardware instance.
sim_detectors_toast(hw, "LAT")

Visualization

The detectors in a given Hardware model can be plotted with this function:

sotodlib.vis_hardware.plot_detectors(dets, outfile, width=None, height=None, labels=False, bandcolor=None, xieta=False, lat_corotate=True, lat_elevation=None, show_centers=False)[source]

Visualize a dictionary of detectors.

This makes a simple plot of the detector positions on the projected focalplane. The size of detector circles are controlled by the detector “fwhm” key, which is in arcminutes. If the bandcolor is specified it will override the defaults.

Parameters:
  • outfile (str) – Output PDF path.

  • dets (dict) – Dictionary of detector properties.

  • width (float) – Width of plot in degrees (None = autoscale).

  • height (float) – Height of plot in degrees (None = autoscale).

  • labels (bool) – If True, label each detector.

  • bandcolor (dict, optional) – Dictionary of color values for each band.

  • xieta (bool) – If True, plot in Xi / Eta / Gamma coordinates rather than focalplane X / Y / Z.

  • show_centers (bool) – If True, label pixel centers.

Returns:

None

To plot only a subset of detectors, first apply a selection to make a reduced hardware model and pass that to the plotting function. You can also dump out to the console a pretty formatted summary of the hardware configuration:

sotodlib.vis_hardware.summary_text(hw)[source]

Print a textual summary of a hardware configuration.

Parameters:

hw (Hardware) – A hardware dictionary.

Returns:

None

Data Processing

These modules are used to process detector timestream data within a G3Pipeline. The base class, DataG3Module, handles the translation between a G3TimestreamMap and any filtering/conditioning we want to do on the timestreams.

DataG3Module

class sotodlib.core.g3_core.DataG3Module(input='signal', output='signal_out')[source]

Base G3Module for processing or conditioning data.

Classes inheriting from DataG3Module only need to override the process() function to build a G3Module that will call process() on every detector in a given G3TimestreamMap.

input

key of G3TimestreamMap of source data

Type:

str

output

key of G3Timestream map of output data. if None, input will be overwritten with output

Type:

str or None

Todo

  • Add detector masking capabilities so we don’t waste CPU time on cut detectors

  • Make an option for input/output to be a list, so processing can be done on multiple timestreams

  • Make it possible for input/output to be a G3Timestream, instead of a G3TimestreamMap

  • Decide if “Modules that add members to frames” should be a default type

process(data, det_name)[source]
Parameters:
  • data (G3Timestream) – data for a single detector

  • det_name (str) – the detector name in the focal plane, in case it’s needed for accessing calibration info

Returns:

data (G3Timestream)

apply(f)[source]

Applies the G3Module on an individual G3Frame or G3TimestreamMap. Likely to be useful when debugging

Parameters:

f (G3Frame or G3TimestreamMap)

Returns:

if f is a G3Frame it returns the frame after the module has been applied G3TimestreamMap: if f is a G3TimestreamMap it returns self.output

Return type:

G3Frame

classmethod from_function(function, input='signal', output=None)[source]

Allows inline creation of DataG3Modules with just a function definition

Example

mean_sub = lambda data: data-np.mean(data) G3Pipeline.Add(DataG3Module.from_function(mean_sub) )

Parameters:

function (function) – function that takes a G3Timestream (or ndarray) and returns an ndarray

Filters

Data Filtering

This module contains code for filtering data in G3Frames

class sotodlib.g3_filter.Filter(input='signal', output='signal_filtered', filter_function=None)[source]

G3Module that takes the G3Timestream map and applies generic filter

input

the key to a G3Timestream map of the source

Type:

str

output

key of G3Timestream map of output data. if None, input will be overwritten with output

Type:

str

filter_function

function that takes frequency in Hz and returns a frequency filter

Type:

function

Todo

Get rid of numpy fft functions and get faster / better parallizable options.

process(data, det_name)[source]
Parameters:
  • data (G3Timestream) – data for a single detector

  • det_name (str) – the detector name in the focal plane in case it’s needed for accessing calibration info

Returns:

np.array filtered by filter_function

class sotodlib.g3_filter.LowPassButterworth(input='signal', output='signal_filtered', order=2, fc=1, gain=1)[source]

G3Module for a LowPassButterworth filter

order

order of butterworth

Type:

int

fc

cutoff frequency in Hertz

Type:

float

gain

filter gain

Type:

float

Conditioning

Data Conditioning

This module contains code for conditioning G3Timestream data in G3Frames

class sotodlib.g3_condition.MeanSubtract(input='signal', output='signal_out')[source]
process(data, det_name)[source]
Parameters:
  • data (G3Timestream) – data for a single detector

  • det_name (str) – the detector name in the focal plane, in case it’s needed for accessing calibration info

Returns:

data (G3Timestream)

class sotodlib.g3_condition.MedianSubtract(input='signal', output='signal_out')[source]
process(data, det_name)[source]
Parameters:
  • data (G3Timestream) – data for a single detector

  • det_name (str) – the detector name in the focal plane, in case it’s needed for accessing calibration info

Returns:

data (G3Timestream)

class sotodlib.g3_condition.Detrend(input='signal', output=None, info='detrend_values', type='linear')[source]

Module for Detrending data. Information is added to the frame so that the resulting data can be easily re-trended (if, for example the detrend is done just for filtering).

process(data, det_name)[source]
Parameters:
  • data (G3Timestream) – data for a single detector

  • det_name (str) – the detector name in the focal plane, in case it’s needed for accessing calibration info

Returns:

data (G3Timestream)

class sotodlib.g3_condition.Retrend(input='signal', output=None, detrend_info='detrend_values')[source]

Module for Retrending data that was Detrended with Detrend

process(data, det_name)[source]
Parameters:
  • data (G3Timestream) – data for a single detector

  • det_name (str) – the detector name in the focal plane, in case it’s needed for accessing calibration info

Returns:

data (G3Timestream)

class sotodlib.g3_condition.Decimate(input='signal', output=None, q=5, **kwargs)[source]

Module for decimating data. Uses scipy.signal.decimate()

process(data, det_name)[source]
Parameters:
  • data (G3Timestream) – data for a single detector

  • det_name (str) – the detector name in the focal plane, in case it’s needed for accessing calibration info

Returns:

data (G3Timestream)

class sotodlib.g3_condition.Resample(input='signal', output=None, num=3000, **kwargs)[source]

Module for resampling data. Uses scipy.signal.resample()

process(data, det_name)[source]
Parameters:
  • data (G3Timestream) – data for a single detector

  • det_name (str) – the detector name in the focal plane, in case it’s needed for accessing calibration info

Returns:

data (G3Timestream)