Preprocess

The preprocess module defines a standardized interface for TOD processing operations so that they can be easily implemented in automatic data analysis scripts. The core of the system is in two parts, the _Preprocess modules and the Pipeline object. The _Preprocess modules each define how a TOD operation is run on an AxisManager TOD and the Pipeline object is used to define the order of the operations and then run them. The site-pipeline.preprocess_tod script is used to run and save Pipelines on lists of observations, grouped by detset. The site-pipeline.preprocess_obs script is used for observation-level preprocessing. This module is similar to site-pipeline.preprocess_tod but removes grouping by detset so that the entire observation is loaded, without signal. For example, pipeline steps such as DetBiasFlags requires tod-level data including signal, whereas SSOFootprint does not and uses observation-level data.

Preprocessing Pipelines

A preprocessing pipeline is series of modules, each inheriting from _Preprocess, that are defined through a configuration file and intended to be run successively on an AxisManager containing time ordered data.

class sotodlib.preprocess.core._Preprocess(step_cfgs)[source]

The base class for Preprocessing modules which defines the required functions and keys required in the configurations.

Each preprocess module has four overwritable functions that are called by the processing scripts in site_pipeline. These four functions are each controlled by a specific key in a configuration dictionary passed to the module on creation.

The configuration dictionary has 6 special keys: name, process, calc, save, select, and plot. name is the name used to register the module with the PIPELINE registry. The other four keys are matched to functions in the module, if the key is not present then that function will be skipped when the preprocessing pipeline is run.

There are two special AxisManagers expected to be part of the preprocessing pipeline. aman is the “standard” time ordered data AxisManager that is loaded via our default styles. proc_aman is the preprocess AxisManager, this is carry the data products that will be saved to whatever Metadata Archive is connected to the preprocessing pipeline.

process(aman, proc_aman)[source]

This function makes changes to the time ordered data AxisManager. Ex: calibrating or detrending the timestreams. This function will use any configuration information under the process key of the configuration dictionary and is not expected to change or alter proc_aman.

Parameters:
  • aman (AxisManager) – The time ordered data

  • proc_aman (AxisManager) – Any information generated by previous elements in the preprocessing pipeline.

calc_and_save(aman, proc_aman)[source]

This function calculates data products of some sort off of the time ordered data AxisManager.

Ex: Calcuating the white noise of the timestream. This function will use any configuration information under the calc key of the configuration dictionary and can call the save function to make changes to proc_aman.

Parameters:
  • aman (AxisManager) – The time ordered data

  • proc_aman (AxisManager) – Any information generated by previous elements in the preprocessing pipeline.

save(proc_aman, *args)[source]

This function wraps new information into the proc_aman and will use any configuration information under the save key of the configuration dictionary.

Parameters:
  • proc_aman (AxisManager) – Any information generated by previous elements in the preprocessing pipeline.

  • args (any) – Any additional information calc_and_save needs to send to the save function.

select(meta, proc_aman=None)[source]

This function runs any desired data selection of the preprocessing pipeline results. Assumes the pipeline has already been run and that the resulting proc_aman is now saved under the preprocess key in the meta AxisManager loaded via context.

Ex: removing detectors with white noise above some limit. This function will use any configuration information under the select key.

Parameters:

meta (AxisManager) – Metadata related to the specific observation

Returns:

meta – Metadata where non-selected detectors have been removed

Return type:

AxisManager

plot(aman, proc_aman, filename)[source]

This function creates plots using results from calc_and_save.

Ex: Plotting det bias flags. This function will use any configuration information under the plot key of the configuration dictionary.

Parameters:
  • aman (AxisManager) – The time ordered data

  • proc_aman (AxisManager) – Any information generated by previous elements in the preprocessing pipeline.

  • filename (str) – Filename should be a concatenation of the global plot_dir config with a name with process step number and placeholder {name} as shown in Pipeline.run().

classmethod gen_metric(meta, proc_aman)[source]

Generate a QA metric from the output of this process.

Parameters:
  • meta (AxisManager) – Metadata related to the specific observation

  • proc_aman (AxisManager) – The output of the preprocessing pipeline.

Returns:

line – InfluxDB line entry elements to be fed to site_pipeline.monitor.Monitor.record

Return type:

dict

static register(process_class)[source]

Registers a new modules with the PIPELINE

The preprocessing pipeline is defined in the Pipeline class. This class inherits from list so that you can easily find and interact with the various pipeline elements. Note that splicing a pipeline will return a list of process modules that can be used to make a new pipeline.

class sotodlib.preprocess.core.Pipeline(modules, plot_dir='./', logger=None, wrap_valid=True)[source]

This class is designed to create and run pipelines out of a series of different preprocessing modules (classes that inherent from _Preprocess). It inherits list object. It also contains the registration of all possible preprocess modules in Pipeline.PIPELINE

append(item)[source]

Append object to the end of the list.

insert(index, item)[source]

Insert object before index.

extend(index, other)[source]

Extend list by appending elements from the iterable.

run(aman, proc_aman=None, select=True)[source]

The main workhorse function for the pipeline class. This function takes an AxisManager TOD and successively runs the pipeline of preprocessing modules on the AxisManager. The order of operations called by run are:

for process in pipeline:
    process.process()
    process.calc_and_save()
        process.save() ## called by process.calc_and_save()
    process.select()
Parameters:
  • aman (AxisManager) – A TOD object. Generally expected to be raw, unprocessed data. This axismanager will be edited in place by the process and select functions of each preprocess module

  • proc_aman (AxisManager (Optional)) – A preprocess axismanager. If this is provided it is assumed that the pipeline has previously been run on this specific TOD and has returned this preprocess axismanager. In this case, calls to process.calc_and_save() are skipped as the information is expected to be present in this AxisManager.

  • select (boolean (Optional)) – if True, the aman detector axis is restricted as described in each preprocess module. Most pipelines are developed with select=True. Running select=False may produce unstable behavior

Returns:

proc_aman – A preprocess axismanager that contains all data products calculated throughout the running of the pipeline

Return type:

AxisManager

Processing Scripts

These scripts are designed to be the ones that interact with specific configuration files and specific manifest databases.

sotodlib.site_pipeline.preprocess_tod.preprocess_tod(obs_id, configs, group_list=None, overwrite=False, logger=None)[source]

Meant to be run as part of a batched script, this function calls the preprocessing pipeline a specific Observation ID and saves the results in the ManifestDb specified in the configs.

Parameters:
  • obs_id (string or ResultSet entry) – obs_id or obs entry that is passed to context.get_obs

  • configs (string or dictionary) – config file or loaded config directory

  • group_list (None or list) – list of groups to run if you only want to run a partial update

  • overwrite (bool) – if True, overwrite existing entries in ManifestDb

  • logger (logging instance) – the logger to print to

sotodlib.site_pipeline.preprocess_tod.load_preprocess_det_select(obs_id, configs, context=None, dets=None, meta=None)[source]

Loads the metadata information for the Observation and runs through any data selection specified by the Preprocessing Pipeline.

Parameters:
  • obs_id (multiple) – passed to context.get_obs to load AxisManager, see Notes for context.get_obs

  • configs (string or dictionary) – config file or loaded config directory

  • dets (dict) – dets to restrict on from info in det_info. See context.get_meta.

  • meta (AxisManager) – Contains supporting metadata to use for loading. Can be pre-restricted in any way. See context.get_meta.

sotodlib.site_pipeline.preprocess_tod.load_preprocess_tod(obs_id, configs='preprocess_configs.yaml', context=None, dets=None, meta=None)[source]

Loads the saved information from the preprocessing pipeline and runs the processing section of the pipeline.

Assumes preprocess_tod has already been run on the requested observation.

Parameters:
  • obs_id (multiple) – passed to context.get_obs to load AxisManager, see Notes for context.get_obs

  • configs (string or dictionary) – config file or loaded config directory

  • dets (dict) – dets to restrict on from info in det_info. See context.get_meta.

  • meta (AxisManager) – Contains supporting metadata to use for loading. Can be pre-restricted in any way. See context.get_meta.

sotodlib.site_pipeline.preprocess_obs.preprocess_obs(obs_id, configs, overwrite=False, logger=None)[source]

Meant to be run as part of a batched script, this function calls the preprocessing pipeline a specific Observation ID and saves the results in the ManifestDb specified in the configs.

Parameters:
  • obs_id (string or ResultSet entry) – obs_id or obs entry that is passed to context.get_obs

  • configs (string or dictionary) – config file or loaded config directory

  • overwrite (bool) – if True, overwrite existing entries in ManifestDb

  • logger (logging instance) – the logger to print to

sotodlib.site_pipeline.preprocess_obs.load_preprocess_obs(obs_id, configs='preprocess_obs_configs.yaml', context=None)[source]

Loads the saved information from the preprocessing pipeline and runs the processing section of the pipeline.

Assumes preprocess_tod has already been run on the requested observation.

Parameters:
  • obs_id (multiple) – passed to context.get_obs to load AxisManager, see Notes for context.get_obs

  • configs (string or dictionary) – config file or loaded config directory

Example TOD Pipeline Configuration File

Suppose we want to run a simple pipeline that runs the glitch calculator and estimates the white noise levels of the data. A configuration file for the processing pipeline would look like:

# Context for the data
context_file: 'context.yaml'

# Plot directory prefix
plot_dir: './plots'

# How to subdivide observations
subobs:
    use: detset
    label: detset

# Metadata index & archive filenaming
archive:
    index: 'preprocess_archive.sqlite'
    policy:
        type: 'simple'
        filename: 'preprocess_archive.h5'

process_pipe:
    - name : "fft_trim"
      process:
        axis: 'samps'
        prefer: 'right'

    - name: "trends"
      calc:
        max_trend: 30
        n_pieces: 5
      save: True
      select:
        kind: "any"

    - name: "glitches"
      calc:
        t_glitch: 0.002
        hp_fc: 0.5
        n_sig: 10
        buffer: 20
      save: True
      select:
        max_n_glitch: 20
        sig_glitch: 30

    - name: "detrend"
      process:
        method: "linear"
        count: 10

    - name: "calibrate"
      process:
        kind: "single_value"
        ## phase_to_pA: 9e6/(2*np.pi)
        val: 1432394.4878270582

    - name: "psd"
      process:
        detrend: False
        window: "hann"

    - name: "noise"
      calc:
        low_f: 5
        high_f: 10
      save: True
      select:
        max_noise: 2000

This pipeline can be run through the functions saved in site_pipeline. Each entry in “process_pipe” key will be used to generate a Preprocess module based on the name it is registered to. These entries will then be run in order through the processing pipe. The process function is always run before the calc_and_save function for each module. The plot function can be run after calc_and_save when plot: True for a module that supports it.

Example Obs Pipeline Configuration File

Suppose we want to run an observation-level pipeline that creates a SSO footprint. A configuration file for the processing pipeline would look like:

# Context for the data
context_file: 'context.yaml'

# Plot directory prefix
plot_dir: './plots'

# Metadata index & archive filenaming
archive:
    index: 'preprocess_archive.sqlite'
    policy:
        type: 'simple'
        filename: 'preprocess_archive.h5'

process_pipe:
    - name: "sso_footprint"
      calc:
        # If you want to search for nearby sources, exclude source_list
        source_list: ['jupiter']
        distance: 20
        nstep: 100
      save: True
      plot:
        wafer_offsets: {'ws0': [-2.5, -0.5],
                        'ws1': [-2.5, -13],
                        'ws2': [-13, -7],
                        'ws3': [-13, 5],
                        'ws4': [-2.5, 11.5],
                        'ws5': [8.5, 5],
                        'ws6': [8.5, -7]}
        focal_plane: 'focal_plane_positions.npz'

Processing Modules

TOD Operations

class sotodlib.preprocess.processes.FFTTrim(step_cfgs)[source]

Trim the AxisManager to optimize for faster FFTs later in the pipeline. All processing configs go to fft_trim

fft_trim(tod, axis='samps', prefer='right')

Restrict AxisManager sample range so that FFTs are efficient. This uses the find_inferior_integer function.

Parameters:
  • tod (AxisManager) – Target, which is modified in place.

  • axis (str) – Axis to target.

  • prefer (str) – One of [‘left’, ‘right’, ‘center’], indicating whether to trim away samples from the end, the beginning, or !equally at the beginning and end (respectively).

Returns:

The (start, stop) indices to use to slice an array and get these samples.

class sotodlib.preprocess.processes.Detrend(step_cfgs)[source]

Detrend the signal. All processing configs go to detrend_tod

detrend_tod(tod, method='linear', axis_name='samps', signal_name='signal', in_place=True, wrap_name=None, count=10)

Returns detrended data. Detrends data in place by default but pass in_place=False if you would like a copied array (such as if you’re just looking to use this in an FFT).

Using this with method =’mean’ and axis_name=’dets’ will remove a common mode from the detectors Using this with method =’median’ and axis_name=’dets’ will remove a common mode from the detectors with the median rather than the mean

Parameters:
  • tod (axis manager)

  • method (str) – method of detrending can be ‘linear’, ‘mean’, or median

  • axis_name (str) – the axis along which to detrend. default is ‘samps’

  • signal_name (str) – the name of the signal to detrend. defaults to ‘signal’. Can have any shape as long as axis_name can be resolved.

  • in_place (bool.) – If False it makes a copy of signal before detrending and returns the copy.

  • wrap_name (str or None.) – If not None, wrap the detrended data into tod with this name.

  • count (int) – Number of samples to use, on each end, when measuring mean level for ‘linear’ detrend. Values larger than 1 suppress the influence of white noise.

Returns:

signal – Detrended signal. Done in place or on a copy depend on in_place argument.

Return type:

array of type tod[signal_name]

class sotodlib.preprocess.processes.PSDCalc(step_cfgs)[source]

Calculate the PSD of the data and add it to the AxisManager under the “psd” field.

Example config block:

- "name : "psd"
  "signal: "signal" # optional
  "wrap": "psd" # optional
  "process":
    "psd_cfgs": # optional, kwargs to scipy.welch
      "nperseg": 1024
    "wrap_name": "psd" # optional
  "calc": True
  "save": True
calc_psd(aman, signal=None, timestamps=None, max_samples=262144, prefer='center', freq_spacing=None, merge=False, overwrite=True, **kwargs)

Calculates the power spectrum density of an input signal using signal.welch(). Data defaults to aman.signal and times defaults to aman.timestamps. By default the nperseg will be set to power of 2 closest to the 1/50th of the samples used, this can be overridden by providing nperseg or freq_spacing.

Parameters:
  • aman (AxisManager) – with (dets, samps) OR (channels, samps)axes.

  • signal (float ndarray) – data signal to pass to scipy.signal.welch().

  • timestamps (float ndarray) – timestamps associated with the data signal.

  • max_samples (int) – maximum samples along sample axis to send to welch.

  • prefer (str) – One of [‘left’, ‘right’, ‘center’], indicating what part of the array we would like to send to welch if cuts are required.

  • freq_spacing (float) – The approximate desired frequency spacing of the PSD. If None the default nperseg of ~1/50th the signal length is used. If an nperseg is explicitly passed then that will be used.

  • merge (bool) – if True merge results into axismanager.

  • overwrite (bool) – if true will overwrite f, pxx axes.

  • **kwargs – keyword args to be passed to signal.welch().

Returns:

array of frequencies corresponding to PSD calculated from welch. Pxx: array of PSD values.

Return type:

freqs

class sotodlib.preprocess.processes.Calibrate(step_cfgs)[source]

Calibrate the timestreams based on some provided information.

Type of calibration is decided by process[“kind”]

1. “single_value” : multiplies entire signal by the single value process[“val”]

2. “array” : takes the dot product of the array with the entire signal. The array is specified by process["cal_array"], which must exist in aman. The array can be nested within additional AxisManager objects, for instance det_cal.phase_to_pW.

Example config block(s):

- name: "calibrate"
  process:
    kind: "single_value"
    # phase_to_pA: 9e6/(2*np.pi)
    val: 1432394.4878270582
- name: "calibrate"
  process:
    kind: "array"
    cal_array: "cal.array"
class sotodlib.preprocess.processes.Apodize(step_cfgs)[source]

Apodize the edges of a signal. All process configs go to apodize_cosine

apodize_cosine(aman, signal='signal', apodize_samps=1600, in_place=True, apo_axis='apodized')

Function to smoothly filter the timestream to 0’s on the ends with a cosine function.

Parameters:
  • signal (str) – Axis to apodize

  • apodize_samps (int) – Number of samples on tod ends to apodize.

  • in_place (bool) – writes over signal with apodized version

  • apo_axis (str)

class sotodlib.preprocess.processes.SubPolyf(step_cfgs)[source]
Fit TOD in each subscan with polynominal of given order and subtract it.

All process configs go to sotodlib.tod_ops.sub_polyf.

subscan_polyfilter(aman, degree, signal=None, exclude_turnarounds=False, mask=None, in_place=True)

Apply polynomial filtering to subscan segments in a data array. This function applies polynomial filtering to subscan segments within signal for each detector. Subscan segments are defined based on the presence of flags such as ‘left_scan’ and ‘right_scan’. Polynomial filtering is used to remove low-degree polynomial trends within each subscan segment.

Parameters:
  • aman (AxisManager)

  • degree (int) – The degree of the polynomial to be removed.

  • signal (array-like, optional) – The TOD signal to use. If not provided, aman.signal will be used.

  • exclude_turnarounds (bool) – Optional. If True, turnarounds are excluded from subscan identification. Default is False.

  • mask (str or RangesMatrix) – Optional. A mask used to select specific data points for filtering. Default is None. If None, no mask is applied. If the mask is given in str, aman.flags['mask'] is used as mask. Arbitrary mask can be specified in the style of RangesMatrix.

  • in_place (bool) – Optional. If True, aman.signal is overwritten with the processed signal.

Returns:

signal – The processed signal.

Return type:

array-like

class sotodlib.preprocess.processes.Jumps(step_cfgs)[source]

Run generic jump finding and fixing algorithm.

calc_cfgs should have ‘function’ defined as one of ‘find_jumps’, ‘twopi_jumps’ or ‘slow_jumps’. Any additional configs to the jump function goes in ‘jump_configs’.

Saves results in proc_aman under the “jumps” field.

Data section should define a maximum number of jumps “max_n_jumps”.

Example config block:

- name: "jumps"
  signal: "hwpss_remove"
  calc:
    function: "twopi_jumps"
  save:
    jumps_name: "jumps_2pi"
find_jumps(aman, signal=None, max_iters=1, min_sigma=None, min_size=None, win_size=20, nsigma=25, fix: Literal[False] = False, inplace=False, merge=True, overwrite=False, name='jumps', **filter_pars) Tuple[RangesMatrix, csr_array]
find_jumps(aman, signal=None, max_iters=1, min_sigma=None, min_size=None, win_size=20, nsigma=25, fix: Literal[True] = True, inplace=False, merge=True, overwrite=False, name='jumps', **filter_pars) Tuple[RangesMatrix, csr_array, ndarray[Any, dtype[floating]]]

Find jumps in aman.signal_name with a matched filter for edge detection. Expects aman.signal_name to be 1D of 2D.

Parameters:
  • aman – axis manager.

  • signal – Signal to jumpfind on. If None than aman.signal is used.

  • min_sigma – Number of standard deviations to count as a jump, note that the standard deviation here is computed by std_est and is the white noise standard deviation, so it doesn’t include contributions from jumps or 1/f. If min_size is provided it will be used instead of this.

  • min_size – The smallest jump size counted as a jump. By default this is set to None and min_sigma is used instead, if set this will override min_sigma. If both min_sigma and min_size are None then the IQR is used as min_size.

  • win_size – Size of window used by SG filter when peak finding. Also used for height estimation, should be of order jump width.

  • nsigma – Number of sigma above the mean for something to be a peak.

  • fix – Set to True to fix.

  • inplace – Whether of not signal should be fixed inplace.

  • merge – If True will wrap ranges matrix into aman.flags.<name>

  • overwrite – If True will overwrite existing content of aman.flags.<name>

  • name – String used to populate field in flagmanager if merge is True.

  • **filter_pars – Parameters to pass to _filter

Returns:

RangesMatrix containing jumps in signal,

if signal is 1D Ranges in returned instead. There is some uncertainty on order of a few samples. Jumps within a few samples of each other may not be distinguished.

heights: csr_array of jump heights.

fixed: signal with jump fixed. Only returned if fix is set.

Return type:

jumps

class sotodlib.preprocess.processes.FixJumps(step_cfgs)[source]

Repairs the jump heights given a set of jump flags and heights.

Example config block:

- name: "fix_jumps"
  signal: "signal" # optional
  process:
  jumps_aman: "jumps_2pi"
jumpfix_subtract_heights(x: ndarray[Any, dtype[floating]], jumps: RangesInt32 | RangesMatrix | ndarray[Any, dtype[bool_]], inplace: bool = False, heights: ndarray[Any, dtype[floating]] | csr_array | None = None, **kwargs) ndarray[Any, dtype[floating]]

Naive jump fixing routine where we subtract known heights between jumps. Note that you should exepect a glitch at the jump locations. Works best if you buffer the jumps mask by a bit.

Parameters:
  • x – Data to jumpfix on, expects 1D or 2D.

  • jumps – Boolean mask or Ranges(Matrix) of jump locations. Should be the same shape at x.

  • inplace – Whether of not x should be fixed inplace.

  • heights – Array of jump heights, can be sparse. If None will be computed.

  • **kwargs – Additional arguments to pass to estimate_heights if heights is None.

Returns:

x with jumps removed.

If inplace is True this is just a reference to x.

Return type:

x_fixed

Flagging and Products

class sotodlib.preprocess.processes.Trends(step_cfgs)[source]

Calculate the trends in the data to look for unlocked detectors. All calculation configs go to get_trending_flags.

Saves results in proc_aman under the “trend” field.

Data selection can have key “kind” equal to “any” or “all.”

Example config block:

- name : "trends"
  signal: "signal" # optional
  calc:
    max_trend: 2.5
    n_pieces: 10
  save: True
  select:
    kind: "any"

Flag Detectors with trends larger than max_trend. This function can be used to find unlocked detectors. Note that this is a rough cut and unflagged detectors can still have poor tracking.

Parameters:
  • aman (AxisManager) – The tod

  • max_trend (float) – Slope at which detectors are unlocked. The default is for use with phase units.

  • n_pieces (int) – Number of pieces to cut the timestream in to to look for trends.

  • max_samples (int) – Maximum samples to compute the slope with.

  • signal (array) – (Optional). Signal to use to generate flags, if None default is aman.signal.

  • timestamps (array) – (Optional). Timestamps to use to generate flags, default is aman.timestamps.

  • merge (bool) – If true, merges the generated flag into aman.

  • overwrite (bool) – If true, write over flag. If false, don’t.

  • name (str) – Name of flag to add to aman.flags if merge is True.

  • full_output (bool) – If true, returns calculated slope sizes

Returns:

  • cut (RangesMatrix) – RangesMatrix of trending regions

  • trends (AxisManager) – If full_output is true, calculated slopes and the sample edges where they were calculated.

class sotodlib.preprocess.processes.GlitchDetection(step_cfgs)[source]

Run glitch detection algorithm to find glitches. All calculation configs go to get_glitch_flags

Saves retsults in proc_aman under the “glitches” field.

Data section should define a glitch significant “sig_glitch” and a maximum number of glitches “max_n_glitch.”

Example configuration block:

- name: "glitches"
  calc:
    signal_name: "hwpss_remove"
    t_glitch: 0.00001
    buffer: 10
    hp_fc: 1
    n_sig: 10
  save: True
  select:
    max_n_glitch: 10
    sig_glitch: 10
get_glitch_flags(aman, t_glitch=0.002, hp_fc=0.5, n_sig=10, buffer=200, detrend=None, signal_name=None, merge=True, overwrite=False, name='glitches', full_output=False, edge_guard=2000)

Find glitches with fourier filtering. Translation from moby2 as starting point

Parameters:
  • aman (AxisManager) – The tod.

  • t_glitch (float) – Gaussian filter width.

  • hp_fc (float) – High pass filter cutoff.

  • n_sig (int or float) – Significance of detection.

  • buffer (int) – Amount to buffer flags around found location

  • detrend (str) – Detrend method to pass to fourier_filter

  • signal_name (str) – Field name in aman to detect glitches on if None, defaults to signal

  • merge (bool)) – If true, add to aman.flags

  • name (string) – Name of flag to add to aman.flags

  • overwrite (bool) – If true, write over flag. If false, raise ValueError if name already exists in AxisManager

  • full_output (bool) – If true, return sparse matrix with the significance of the detected glitches

  • edge_guard (int) – Number of samples at the beginning and end of the tod to exclude from the returned glitch RangesMatrix. Defaults to 2000 samples (10 sec).

Returns:

flag – RangesMatrix object containing glitch mask.

Return type:

RangesMatrix

class sotodlib.preprocess.processes.GlitchFill(step_cfgs)[source]

Fill glitches. All process configs go to fill_glitches.

Example configuration block:

- name: "glitchfill"
  signal: "hwpss_remove"
  flag_aman: "jumps_2pi"
  flag: "jump_flag"
  process:
    nbuf: 10
    use_pca: False
    modes: 1
fill_glitches(aman, nbuf=10, use_pca=False, modes=3, signal=None, glitch_flags=None, wrap=True)

This function fills pre-computed glitches provided by the caller in time-ordered data using either a polynomial (default) or PCA-based approach. Wraps the other functions in the tod_ops.gapfill module.

Parameters:
  • aman (AxisManager) – AxisManager to fill glitches in

  • nbuf (int) – Number of buffer samples to use in polynomial gap filling.

  • use_pca (bool) – Whether or not to fill glitches using pca model. Default is False

  • modes (int) – Number of modes in the pca to use if pca=True. Default is 3.

  • signal (ndarray or None) – Array of data to fill glitches in. If None then uses aman.signal. Default is None.

  • glitch_flags (RangesMatrix or None) – RangesMatrix containing flags to use for gap filling. If None then uses aman.flags.glitches.

  • wrap (bool or str) – If True wraps new field called gap_filled, if False returns the gap filled array, if a string wraps new field with provided name.

Returns:

signal – Returns ndarray with gaps filled from input signal.

Return type:

ndarray

class sotodlib.preprocess.processes.Noise(step_cfgs)[source]

Estimate the white noise levels in the data. Assumes the PSD has been wrapped into the AxisManager. All calculation configs goes to calc_wn.

Saves the results into the “noise” field of proc_aman.

Can run data selection of a “max_noise” value.

Example config block:

- name: "noise"
  calc:
    low_f: 5
    high_f: 10
  save: True
  select:
    max_noise: 2000

If fit: True this operation will run sotodlib.tod_ops.fft_ops.fit_noise_model(), else it will run sotodlib.tod_ops.fft_ops.calc_wn().

class sotodlib.preprocess.processes.FlagTurnarounds(step_cfgs)[source]
From the Azimuth encoder data, flag turnarounds, left-going, and right-going.

All process configs go to get_turnaround_flags. If the method key is not included in the preprocess config file calc configs then it will default to ‘scanspeed’.

get_turnaround_flags(aman, az=None, method='scanspeed', name='turnarounds', merge=True, merge_lr=True, overwrite=True, t_buffer=2.0, kernel_size=400, peak_threshold=0.1, rel_distance_peaks=0.3, truncate=False, qlim=1)

Compute turnaround flags for a dataset.

Parameters:
  • aman (AxisManager) – Input axis manager.

  • az (Array) – (Optional). Azimuth data for turnaround flag computation. If not provided, it uses aman.boresight.az.

  • method (str) – (Optional). The method for computing turnaround flags. Options are az or scanspeed.

  • name (str) – (Optional). The name of the turnaround flag in aman.flags. Default is turnarounds

  • merge (bool) – (Optional). Merge the computed turnaround flags into aman.flags if True.

  • merge_lr (bool) – (Optional). Merge left and right scan flags as aman.flags.left_scan and aman.flags.right_scan if True.

  • overwrite (bool) – (Optional). Overwrite an existing flag in aman.flags with the same name.

  • t_buffer (float) – (Optional). Buffer time (in seconds) for flagging turnarounds in scanspeed method.

  • kernel_size (int) – (Optional). Size of the step-wise matched filter kernel used in scanspeed method.

  • peak_threshold (float) – (Optional). Peak threshold for identifying peaks in the matched filter response. It is a value used to determine the minimum peak height in the signal.

  • rel_distance_peaks (float) – (Optional). Relative distance between peaks. It specifies the minimum distance between peaks as a fraction of the approximate number of samples in one scan period.

  • truncate (bool) – (Optional). Truncate unstable scan segments if True in scanspeed method.

  • qlim (float) – (Optional). Azimuth threshold percentile for az method turnaround detection.

Returns:

Ranges – The turnaround flags as a Ranges object.

Return type:

RangesMatrix

Obs Operations

class sotodlib.preprocess.processes.SSOFootprint(step_cfgs)[source]

Find nearby sources within a given distance and get SSO footprint and plot each source on the focal plane.

get_sso(aman, sso, nstep=100)

Function for getting xi, eta position of given sso.

Parameters:
  • aman (AxisManager) – Input axis manager.

  • sso (str) – Name of input sso.

  • nstep (int) – Number of steps to downsample the TOD.

Returns:

  • xi (array) – Array of xi positions.

  • eta (array) – Array of eta positions.