coords submodule

This submodule includes functions for computing celestial coordinates associated with time-ordered data, and projecting signals from maps to time domain and back again.

The main categories of functionality provided by this module are: the Projection Matrix, for quick projections between map and time-ordered domain; utilities for building useful map geometries.

The maps produced by this module are objects of class pixell.enmap; see pixell documentation.

Coordinate systems referred to here are defined in more detail in the coord_sys document, in the tod2maps_docs repository.

Conventions

Functions in this module that work with time-ordered data will work most easily if the data and pointing information are stored in an AxisManager, which we will denote tod, according to the following conventions:

  • Axes called “dets” and “samps”, with dets a LabelAxis.

  • tod['signal']: 2-d array with shape (dets, samps) that carries the detector measurements, often called “the data”.

  • tod['timestamps']: timestamps (as unix timestamps).

  • tod['boresight']: child AxisManager with boresight pointing in horizon coordinates.

  • tod['focal_plane']: child AxisManager, with shape (dets) and fields xi, eta, gamma representing detector position and orientation in the focal plane. See coord_sys documentation for exact details, but following fields:

    • 'xi': position of the detector, measured “to the right” of the boresight, i.e. parallel to increasing azimuth.

    • 'eta': position of the detector, measured “up” from the boresight, i.e. parallel to increasing elevation.

    • 'gamma': orientation of the detector, measured clockwise from North, i.e. 0 degrees is parall to the eta axis, and 90 degrees is parallel to the xi axis.

The following item is supported, as an alternative to tod['boresight'] and tod['timestamps']:

  • tod['boresight_equ']: child AxisManager giving the boresight pointing in equatorial coordinates ('ra', 'dec', 'psi', forming a lonlat triple, with 'psi' optional but recommended).

Many functions require an AxisManager tod as the main argument, and extract the above fields from that object by default. To override that, pass the field as a keyword argument.

Coordinate Conversions

To compute the RA and dec of all detectors at all times, use get_radec(). For example:

# Get coords for all detectors
radec = sotodlib.coords.get_radec(tod)

# Plot RA and dec for detector at index 10.
det_idx = 10
ra, dec = radec[det_idx,:,0], radec[det_idx,:,1]
plt.plot(ra, dec)

If want the horizon coordinates, use get_horizon().

Map Geometry Assistance

In the terminology of pixell, a “geometry” is the combination of a WCS describing a rectangular pixelization, and the number of pixels in each of the celestial axes. Supposing em is an enmap object, then the geometry is just:

geom = (em.shape[-2:], em.wcs)

So the geometry is a way to define a specific rectangular grid of pixels on the celestial sky, independent of a specific enmap.

Here are some common tasks related to geometries that sotodlib.coords helps with.

1. Get a map big enough to contain all samples from a TOD.

First call get_wcs_kernel() to define the pixelization you want (without specifying the size of the map). Then call get_footprint() to return a geometry. For example:

wcsk = coords.get_wcs_kernel('tan', 0., 0., 0.01*coords.DEG)
wcs = coords.get_footprint(tod, wcsk)

2. Given a list of geometries that are mutually pixel-compatible get a geometry that contains all of them.

Use get_supergeom().

For example, you might use get_wcs_kernel() to define a pixelization on the sky, then run get_footprint() on several TODs, with that same kernel, to get their individual geometries. Then call get_supergeom() on the set of footprints.

Projection Matrix

Usage

Assuming that you have an AxisManager “tod” with “signal”, “timestamps”, “boresight”, and “focal_plane” elements, then this should work:

from sotodlib import coords

wcsk = coords.get_wcs_kernel('car', 0., 0., 0.01*coords.DEG)
P = coords.P.for_tod(tod, wcs_kernel=wcsk)

map1 = P.remove_weights(tod)
map1.write('output_car.fits')

For more advice, see P in the module reference. See also the pwg-tutorials repository (this may be private to SO members), especially so_pipeline_pt2_20200623.ipynb.

Background

Mapmaking is based on the model

\[d = P m + n\]

In this representation, \(d\) is a vector of all time-ordered measurements, \(n\) is the vector of time-ordered noise, \(m\) is a vector of all map “pixels” (including all relevant spin-components, i.e. T, Q, U), and P is the (sparse) “Projection Matrix” (also called the “Pointing Matrix”) that transfers some combination of map elements into each sample.

In practice we want to solve the mapmaker’s equation:

\[(P^T N^{-1} P) m = P^T N^{-1} d\]

It’s assumed that we have measurements \(d\) and projection matrix \(P\), and we somehow know the noise covariance matrix \(N = \langle nn^T \rangle\); the solution of this equation is the maximum likelihood estimate of the map, \(m\).

Note that in the absence of fancy math formatting, that equation would be written as:

(P^T N^-1 P) m = P^T N^-1 d

This notation may be used in docstrings and comments.

Anyway, if the matrix acting on \(m\) in the LHS of the mapmaker’s equation, i.e. P^T N^-1 P, cannot be computed and inverted explicity, then the equation can be solved through iterative methods such as PCG.

One case in which the inversion of the matrix is tractable is if the noise matrix is perfectly diagonal (i.e. if the noise in each detector is white, and if the detectors are not correlated with each other). In that case, the pixel-pixel covariance matrix is nearly diagonal – the only off-diagonal components are ones that couple the T,Q,U elements of each sky pixel with each other. While P^T N^1 P in general has shape (n_pix*n_comp, n_pix*n_comp), it can be factored into n_pix separate matrices of shape (n_comp, n_comp). This is a much smaller structure, whose inverse is obtained by inverting each (n_comp, n_comp) sub-matrix.

The use of the present class for mapmaking is that it provides the operation P. In the case of uncorrelated white noise, it is also able to assist with the application of N^-1 and the inversion of P^T N^-1 P.

In “filter and bin” mapmaking, the signal d is filtered to produce d’ = F d, and the above equations are used to model d’ instead of d. The point of the filter is to suppress bright and/or correlated noise, so that the noise covariance matrix of the timestream d’ is more plausibly diagonal. The map m~ made from data d’ will be biased; but since the mapmaking equation can be solved directly, it is relatively cheap to make maps and the characterization of the bias involves similarly cheap operations on (noise-free) simulated signal.

Connecting to the methods of this class, for iterative map-making with a complicated noise covariance matrix:

  1. Compute the RHS, P^T N^-1 d:

    • Apply your inverse noise matrix to your data, to get d’ = N^-1 d.

    • Call .to_map() with signal set to d’. The result is P^T N^-1 d.

  2. During iteration, apply (P^T N^-1 P) to a map m:

    • Call from_map on map m, returning a timestream y.

    • Compute inverse-noise filtered timestream y’ = N^-1 y

    • Call to_map with signal y’. The result is the LHS.

In the case of simple uncorrelated noise:

  1. Filter your timestream data somehow; call the result d.

  2. Compute (P^T N^-1 d), by calling to_map(). If your noise matrix includes per-detector weights, pass them through det_weights.

  3. Compute (P^T N-1 P)^-1, by calling .to_inverse_weights(). This doesn’t depend on the data, but only on the pointing information. If your noise matrix includes per-detector weights, pass them through det_weights.

  4. Apply the inverse pixel-pixel cov matrix to your binned map, using .remove_weights.

In fact, .remove_weights can perform steps (1), (2), and (3). But it is worth mentioning the full series, because although step (2) is the most expensive it does not depend on the input data and so its output map can be re-used for different filters or input signals.

Reference

Class and function references should be auto-generated here.

Projection Matrix

class sotodlib.coords.P(sight=None, fp=None, geom=None, comps='T', cuts=None, threads=None, det_weights=None, interpol=None)[source]

Projection Matrix.

This class provides functions to apply a Projection Matrix (or its transpose). The Projection Matrix, P, also sometimes called the Pointing Matrix, describes how a vector of time-ordered measurements d are determined by a vector of map pixels m:

d = P m

We are working in aspace, of course, where d is not a vector, but rather an array of vectors, and m is not a vector, it’s a multi-dimensional array with two of the dimensions corresponding to a rectangular pixelization of the sky.

If you are making filter+bin maps, you will want to use these functions:

  • to_map

  • to_inverse_weights

  • remove_weights

If you are solving for a map iteratively you will probably just need:

  • to_map

  • from_map

Important keyword arguments that are used in many functions:

  • tod: The AxisManager from which signal and pointing information should be taken.

  • signal: The array to use for time-ordered signal. If a string, it will be looked up in tod. Defaults to ‘signal’.

  • det_weights: A vector of floats representing the inverse variance of each detector, under the assumption of white noise.

  • cuts: A RangesMatrix that identifies samples that should be excluded from projection operations.

  • comps: the component code (e.g. ‘T’, ‘TQU’, ‘QU’, …) that specifies the spin-components being modeled in the map.

  • dest: the appropriately shaped array in which to place the computed result (as an alternative to allocating a new object to store the result).

Note that in the case of det_weights and cuts, the Projection Matrix may also have cached values for those. It is an error to pass either of these as a keyword argument to a projection routine if a value has been cached for it.

Objects of this class cache certain pre-computed information (such as the rotation taking boresight coordinates to celestial coordinates) and certain context-dependent settings (such as a map shape, WCS, and spin-component configuration). You may want to inspect or borrow these results, perhaps to reuse them when constructing new instances with slight modifications. The cached attributes are:

  • sight: A CelestialSightLine, representing the boresight pointing in celestial coordinates. [samps]

  • fp: G3VectorQuat representing the focal plane offsets of each detector. When constructing with a tod argument, fp can be automatically populated from tod.focal_plane; note that if hwp is passed as True then the gamma (polarization) angles will be reflected (gamma’ = -gamma). [dets]

  • geom: The target map geometry. This is a pixell.enmap.Geometry object, with attributes .shape and .wcs; or possibly (if tiled) a pixell.tilemap.TileGeometry.

  • comps: String indicating the spin-components to include in maps. E.g., ‘T’, ‘QU’, ‘TQU’.

  • rot: quat giving an additional fixed rotation to apply to get from boresight to celestial coordinates. Not for long…

  • cuts (optional): RangesMatrix indicating what samples to exclude from projection operations (the indicated samples have projection matrix element 0 in all components). [dets, samps]

  • threads (optional): list of RangesMatrix objects, to control the use of threads in TOD-to-map operations using OpenMP. See so3g documentation. [[threads, dets, samps], …]

  • det_weights (optional): weights (one per detector) to apply to time-ordered data when binning a map (and also when binning a weights matrix). [dets]

  • interpol (optional): How to interpolate the values for samples between pixel centers. Forwarded to Projectionist. Valid options are:

    • None, ‘nn’ or ‘nearest’: Standard nearest neighbor mapmaking.

    • ‘lin’ or ‘bilinear’: Linearly interpolate between the four closest pixels.

    Default: None

These things can be updated freely, with the following caveats:

  • If the number of “samples” or “detectors” is changed in one attribute, it will need to be changed in the others to match.

  • The threads attribute, if in use, needs to be recomputed if anything about the pointing changes (this includes map geometry but does not include map components).

Setting the “threads” argument to certain special values will activate different thread assignment algorithms:

  • False: do not use threading; to_map projections will be single-threaded.

  • True: use the default algorithm, ‘domdir’.

  • None: same as True.

  • ‘simple’: compute self.threads using simple map-stripe algorithm.

  • ‘domdir’: compute self.threads using dominant-direction algorithm (recommended).

  • ‘tiles’: for tiled geometries, design self.threads such that each tile is assigned to a single thread (each thread may be in charge of multiple tiles).

classmethod for_tod(tod, sight=None, geom=None, comps='T', rot=None, cuts=None, threads=None, det_weights=None, timestamps=None, focal_plane=None, boresight=None, boresight_equ=None, wcs_kernel=None, weather='typical', site='so', interpol=None, hwp=False)[source]

Set up a Projection Matrix for a TOD. This will ultimately call the main P constructor, but some missing arguments will be extracted from tod and computed along the way.

To determine the boresight pointing in celestial coordinates (ultimately passed to constructor as sight=), the first non-None item in the following list is used:

  • the sight= keyword argument.

  • the boresight_equ= keyword argument.

  • the boresight= keyword argument

  • tod.get(‘boresight_equ’)

  • tod.get(‘boresight’)

If the map geometry geom is not specified, but the wcs_kernel is provided, then get_footprint will be called to determine the geom.

The per-detector positions and polarization directions are extracted from focal_plane (xi, eta, gamma), or else tod.focal_plane. If hwp=True, then the rotation angles (gamma) are reflected (gamma’ = -gamma) before storing in self.fp.

classmethod for_geom(tod, geom, comps='TQU', timestamps=None, focal_plane=None, boresight=None, rot=None, cuts=None)[source]

Deprecated, use .for_tod.

zeros(super_shape=None, comps=None)[source]

Returns an enmap concordant with this object’s configured geometry and component count.

Parameters:
  • super_shape (tuple) – The leading dimensions of the array. If None, self._comp_count(comps) is used.

  • comps – The component list, to override self.comps.

Returns:

An enmap with shape super_shape + self.geom.shape.

to_map(tod=None, dest=None, comps=None, signal=None, det_weights=None, cuts=None, eigentol=None)[source]

Project time-ordered signal into a map. This performs the operation

m += P d

and returns m.

Parameters:
  • tod – AxisManager; possible source for ‘signal’, ‘det_weights’.

  • dest (enmap) – the map or array into which the data should be accumulated. (If None, a new enmap is created and initialized to zero.)

  • signal – The time-ordered data, d. If None, tod.signal is used.

  • det_weights – The per-detector weight vector. If None, self.det_weights will be used; if that is not set then uniform weights of 1 are applied.

  • cuts – Sample cuts to exclude from processing. If None, self.cuts is used.

  • eigentol – This is ignored.

to_weights(tod=None, dest=None, comps=None, signal=None, det_weights=None, cuts=None)[source]

Computes the weights matrix for the uncorrelated noise model and returns it. I.e.:

W += P N^-1 P^T

and returns W. Here the inverse noise covariance has shape (n_dets), and carries a single weight (1/var) value for each detector.

Parameters:
  • tod (AxisManager) – possible source for det_weights.

  • dest (enmap) – the map or array into which the weights should be accumulated. (If None, a new enmap is created and initialized to zero.)

  • det_weights – The per-detector weight vector. If None, tod.det_weights will be used; if that is not set then uniform weights of 1 are applied.

  • cuts – Sample cuts to exclude from processing. If None, self.cuts is used.

to_inverse_weights(weights_map=None, tod=None, dest=None, comps=None, signal=None, det_weights=None, cuts=None, eigentol=0.0001)[source]

Compute an inverse weights map, W^-1, from a weights map. If no weights_map is passed in, it will be computed by calling to_weights, passing through all other arguments.

remove_weights(signal_map=None, weights_map=None, inverse_weights_map=None, dest=None, **kwargs)[source]

Apply the inverse weights matrix to a signal map.

m’ = W^-1 m

If W or m are not fully specified, they will be computed by calling other routines inline, with relevant arguments passed through.

Parameters:
  • signal_map – The map m to filter.

  • inverse_weights_map – the matrix W^-1 to apply to the map. Shape should be (n_comp, n_comp, n_row, n_col), but only the upper diagonal in the first two dimensions needs to be populated. If this is None, then “weights_map” is taken as W, and it will be inverted and applied.

  • weights_map – the matrix W. Shape should be (n_comp, n_comp, n_row, n_col), but only the upper diagonal in the first two dimensions needs to be populated. If this is None, then W will be computed by a call to

from_map(signal_map, dest=None, comps=None, wrap=None, cuts=None, tod=None)[source]

Project from a map into the time-domain.

d += P m

Parameters:
  • signal_map – The map m. This can probably be just about anything supported by so3g.proj; it doesn’t have to match the internally configured geometry.

  • dest – Time-ordered data array, shape (dets, samps). If None, a new array will be created to hold the result. Otherwise, data are accumulated into d, so clear it manually if you are trying to do d = P m.

  • comps (str) – Projection components, if you want to override.

  • cuts – RangesMatrix, shape (dets, samps) flagging samples that should not be populated. Defaults to empty.

  • wrap (str) – If specified, wraps the result as tod[wrap] (after removing whatever was in there).

Returns:

The dest array.

Notes

Since this is a set of one-to-many operation, OpenMP can be used without carefully assigning samples to threads.

property tiled

Duck-typing to see if we’re tiled or not. Reload-safe, unlike isinstance

Coordinate Conversions

sotodlib.coords.get_radec(tod, wrap=False, dets=None, timestamps=None, focal_plane=None, boresight=None, sight=None)[source]

Get the celestial coordinates of all detectors at all times.

Parameters:
  • wrap (bool) – If True, the output is stored into tod[‘radec’]. This can also be a string, in which case the result is stored in tod[wrap].

  • dets (list of str) – If set, then coordinates are only computed for the requested detectors. Note you probably can’t wrap the result, in this case, as the dets axis is non-concordant.

Returns:

The returned array has shape (n_det, n_samp, 4). The four components in the last dimension correspond to (lon, lat, cos(psi), sin(psi)). The lon and lat are in radians and correspond to RA and dec of equatorial coordinates. Psi is is the parallactic rotation, measured from North towards West (opposite the direction of standard Position Angle).

sotodlib.coords.get_horiz(tod, wrap=False, dets=None, timestamps=None, focal_plane=None, boresight=None)[source]

Get the horizon coordinates of all detectors at all times.

Parameters:
  • wrap (bool) – If True, the output is stored into tod[‘horiz’]. This can also be a string, in which case the result is stored in tod[wrap].

  • dets (list of str) – If set, then coordinates are only computed for the requested detectors. Note you probably can’t wrap the result, in this case, as the dets axis is non-concordant.

Returns:

The returned array has shape (n_det, n_samp, 4). The four components in the last dimension correspond to (-lon, lat, cos(psi), sin(psi)). The -lon and lat are in radians and correspond to horizon azimuth and elevation. Psi is is the parallactic rotation, measured from North towards West (opposite the direction of standard Position Angle).

class sotodlib.coords.ScalarLastQuat(arr)[source]

Wrapper class for numpy arrays carrying quaternions with the ijk1 signature.

The practice in many codes, including TOAST, quaternionarray, and scipy is to store quaternions in numpy arrays with shape (…, 4), with the real part of the quaternion at index […, 3]. In contrast, spt3g_software inherits the 1ijk from boost.

This class serves two main purposes:

  • It can be used to annotate numpy arrays as carrying quaternions in the ijk1 signature (without disturbing their behavior as numpy arrays).

  • It provides conversion to and from G3 quaternion containers, including signature correction.

When instantiating an object of this class with a numpy array as an argument, the object will wrap a view of the array (not a copy).

When instantiating an object of this class with a G3 quat as an argument, the object will store a copy of the G3 quat, translated to ijk1 signature:

q_g3 = spt3g.core.quat(1.,2.,3.,4.)
q_tq = ScalarLastQuat(q_g3)
q_g3_again = q_tq.to_g3()

print(q_g3, q_tq, q_g3_again)
=>  (1,2,3,4) [2. 3. 4. 1.] (1,2,3,4)
to_g3()[source]

Return the ijk1-signature equivalent of the enclosed quaternion array, as a spt3g.core.quat (if 1-d) or a G3VectorQuat (if 2-d).

Map and Geometry helpers

sotodlib.coords.get_footprint(tod, wcs_kernel, dets=None, timestamps=None, boresight=None, focal_plane=None, sight=None, rot=None)[source]

Find a geometry (in the sense of enmap) based on wcs_kernel that is big enough to contain all data from tod. Returns (shape, wcs).

sotodlib.coords.get_wcs_kernel(proj, ra=None, dec=None, res=None)[source]

Construct a WCS. This fixes the projection type (e.g. CAR, TAN), centered with respect to a reference point (at ra,dec = (0,0)), and resolution of a pixelization, without specifying a particular grid of pixels.

This interface is subject to change.

Parameters:
  • proj (str) –

    Either the name of a FITS projection to use (e.g. ‘car’, ‘cea’, ‘tan’), in which case “res” must also be specified, or a string containing both the projection and the resolution, in this format:

    proj-res

    where proj is a projection type (e.g. ‘car’, ‘tan, ‘gnom’) and res is the resolution, in appropriate units (deg, arcmin or arcsec). Examples of acceptable args:

    ’car-0.5deg’ ‘tan-3arcmin’ ‘gnom-25arcsec’

  • ra – Right Ascension (longitude) of the reference position, in radians.

  • dec – Declination (latitude) of the reference position, in radians.

  • res – Resolution, in radians.

Returns a WCS object that captures the requested pixelization.

sotodlib.coords.get_supergeom(*geoms, tol=0.001)[source]

Given a set of compatible geometries [(shape0, wcs0), (shape1, wcs1), …], return a geometry (shape, wcs) that includes all of them as a subset.

Planet Mapmaking support

class sotodlib.coords.planets.SlowSource(timestamp, ra, dec, v_ra=0.0, v_dec=0.0, precision=1.7453292519943296e-06)[source]

Class to track the time-dependent position of a slow-moving source, such as a Solar System planet in equatorial coordinates. “Slow-moving” is in relation to the time range of interest. A linear approximation is good enough for Solar System, arcsecond accuracy, on time scales of a few hours.

Parameters:
  • timestamp (float) – reference time, as a unix timestamp.

  • ra (float) – Right Ascension at reference time, in radians.

  • dec (float) – Declination at reference time, in radians.

  • v_ra (float) – rate of change of RA, in radians.

  • v_dec (float) – rate of change of dec, in radians.

  • precision – not implemented, don’t worry about it.

classmethod for_named_source(name, timestamp)[source]

Returns a SlowSource for planet name, with position and peculiar velocity measured at time timestamp (float, unix).

pos(timestamps)[source]

Get the (approximate) source position at the times given by the array of unix timestamps.

Returns two arrays (RA, dec), the same shape as timestamps, both in radians. The RA are not corrected into any particular branch (the mapping from timestamp to RA is continuous).

sotodlib.coords.planets.get_scan_q(tod, planet, refq=None)[source]

Identify the point (in time and azimuth space) at which the specified planet crosses the boresight elevation for the observation in tod. The rotation taking boresight coordinates to celestial coordinates at that moment in time (and pointing) is computed, and its conjugate is returned.

The returned rotation is useful because it represents a fixed way to rotate the celestial such that the target source ends up at (xi, eta) = 0, with the telescope scan direction parallel to xi and elevation parallel to eta; this is a useful system for measuring beam and pointing parameters.

sotodlib.coords.planets.get_scan_P(tod, planet, refq=None, res=None, size=None, **kw)[source]

Get a standard Projection Matrix targeting a planet (or some interesting fixed position), in source-scan coordinates.

Returns a Projection Matrix and the output from get_scan_q.

sotodlib.coords.planets.filter_for_sources(tod=None, signal=None, source_flags=None, n_modes=10, low_pass=None, wrap=None)[source]

Mask and gap-fill the signal at samples flagged by source_flags. Then PCA the resulting time ordered data. Restore the flagged signal, remove the strongest modes from PCA.

If signal is not passed in tod.signal will be modified directly. To leave tod.signal intact, pass in signal=tod.signal.copy().

Parameters:
  • tod – AxisManager from which defaults will be drawn.

  • signal – Time-ordered data to operate on. Defaults to tod.signal.

  • source_flags – RangesMatrix to use for source flagging. Defaults to tod.source_flags.

  • low_pass – Frequency, in Hz, at which to apply low pass filter to signal. If None, no filtering is done. You can pass in a filter from tod_ops.filters if you want.

  • n_modes (int) – Number of eigenmodes to remove… interface subject to change.

  • wrap (str) – If specified, the result will be stored at tod[wrap].

Returns:

The filtered signal.

sotodlib.coords.planets.get_source_pos(source_name, timestamp, site='_default')[source]

Get the equatorial coordinates of a planet (or fixed-position source, see note) at some time. Returns the apparent position, accounting for geographical position on earth, but assuming no atmospheric refraction.

Parameters:
  • source_name – Planet name; in capitalized format, e.g. “Jupiter”, or fixed source specification.

  • timestamp – unix timestamp.

  • site (str or so3g.proj.EarthlySite) – if this is a string, the site will be looked up in so3g.proj.SITES dict.

Returns:

in radians. dec (float): in radians. distance (float): in AU.

Return type:

ra (float)

Note

Before checking in the ephemeris, the source_name will be matched against a regular expression and if it has the format ‘Jxxx[+-]yyy’, where xxx and yyy are decimal numbers, then a fixed-position source at RA,Dec = xxx,yyy in degrees will be processed. In that case, the distance is returned as Inf.

sotodlib.coords.planets.get_source_azel(source_name, timestamp, site='_default')[source]

Get the apparent azimuth and elevation of a celestial source at a specific timestamp and observing site. Returns the apparent position, accounting for geographical position on earth, but assuming no atmospheric refraction.

Parameters:
  • source_name – Planet name; in capitalized format, e.g. “Jupiter”

  • timestamp (float) – The Unix timestamp representing the time for which to calculate azimuth and elevation.

  • site (str or so3g.proj.EarthlySite) – if this is a string, the

  • dict. (site will be looked up in so3g.proj.SITES)

Returns:

in radians. el (float): in radians. distance (float): in AU.

Return type:

az (float)

sotodlib.coords.planets.get_nearby_sources(tod=None, source_list=None, distance=1.0)[source]

Identify solar system objects (especially “planets”) that might be within a TOD’s scan footprint.

Parameters:
  • tod (AxisManager) – The data to check. Needs to have focal_plane, boresight, timestamps.

  • source_list (list or None) – A list of source names or None to use a default list. Use simple planet names in lower case (e.g. [‘uranus’]), or tuples with source name, RA, and dec in degrees (e.g. [(‘tau_a’, 83.63, 22.01)]).

  • distance (float) – Maximum distance from the source center, in degrees, to consider as “within the footprint”. (This should be at least the sum of the beam radius and the planet radius … though there’s usually no harm in going a bit larger than that.)

Returns:

List of tuples (source_name, SlowSource) that satisfy the “nearby” condition.

sotodlib.coords.planets.compute_source_flags(tod=None, P=None, mask=None, wrap=None, center_on=None, res=None, max_pix=4000000.0)[source]

Process masking instructions and create RangesMatrix that flags samples in the TOD that are within the masked region. This masking makes use of a map with the footprint encoded in P, so flagging boundaries will correspond to pixel edges.

The interface for “mask” is subject to change – use of a simple text file is interim.

Parameters:
  • tod (AxisManager) – the observation.

  • P (Projection Matrix) – if passed in, must include a map geom (shape and WCS). If None, will be created from get_scan_P using center_on and res parameters.

  • mask – source masking instructions (see note).

  • wrap – key in tod at which to store the result.

  • res – If P is None, sets the target mask map resolution (radians).

  • max_pix – If P is None, this sets the maximum acceptable number of pixels for the mask map. This is to catch cases where an incorrect source has been passed in, for example, leading to a weird map footprint

Returns:

RangesMatrix marking the samples inside the masked region.

Notes

The mask can be a dict or a list of dicts. Each dict must be of the form:

{'shape': 'circle',
 'xyr': (XI, ETA, R)}

where R is the radius of the circular mask, and XI and ETA are the center of the circle, all in degrees.

sotodlib.coords.planets.load_detector_splits(tod=None, filename=None, dataset=None, source=None, wrap=None)[source]

Convert a partition of detectors into a dict of disjoint RangesMatrix objects; such an object can be passed to make_map() to efficiently make detector-split maps.

The “detector split” data can be read from an HDF5 dataset, or passed in directly as an AxisManager.

Parameters:
  • tod (AxisManager) – This is required, to get the list of dets and the samps count.

  • filename (str) – The HDF filename, or filename:dataset.

  • dataset (str) – The HDF dataset (if not passed in with filename).

  • source (array, ResultSet or AxisManager) – If not None, then filename and dataset are ignored and this object is processed (as though it had just been loaded from HDF).

  • wrap (str) – If not None, the address in tod where to store the loaded split data.

The format of detector splits, in an HDF5 dataset, is as one would write from a ResultSet with columns [‘dets:name’, ‘group’] (both str). All detectors sharing a value in the group column will grouped together and the group label will be that value.

If passing in “source” directly as a ResultSet, it should have columns ‘dets:name’ and ‘group’; if as an AxisManager then it should have a ‘dets’ axis and a vector ‘group’ with shape (‘dets’,) providing the group name for each detector. If it is a numpy array, it is assumed to correspond one-to-one with the .dets axis of TOD and the array gives the group name for each detector.

Returns:

Each entry of the dict is a

RangesMatrix that can be interpreted as cuts to apply during mapmaking. In this case the RangesMatrix will simply mark each detector as either fully cut (flagged) or fully uncut.

Return type:

data_splits (dict of RangesMatrix)

sotodlib.coords.planets.get_det_weights(tod, signal=None, wrap=None, outlier_clip=None)[source]

Compute detector weights, based on variance of signal. If outlier_clip is set, it is used to trim outliers. See code for details, but think of it as the number of stdev away from the mean weight that will be kept, and try 2.5.

Returns:

A det_weights array compatible with projection matrix functions.

sotodlib.coords.planets.write_det_weights(tod, filename, dataset, det_weights=None)[source]

Save detector weights to an HDF file dataset.

Parameters:
  • tod (AxisManager) – provides the dets axis, and det_weights if not passed explicitly.

  • filename (str or h5py.File) – destination filename (or open File)

  • dataset (str) – address in the HDF5 file to put the result (it will be overwritten if exists already).

  • det_weights (array) – the weights to write; must be in correspondence with tod.dets. Defaults to tod.det_weights.

Returns:

ResultSet with the info that was written.

sotodlib.coords.planets.make_map(tod, center_on=None, scan_coords=True, thread_algo=False, res=0.00017453292519943296, size=None, wcs_kernel=None, comps='TQU', signal=None, det_weights=None, filename=None, source_flags=None, cuts=None, data_splits=None, low_pass=None, n_modes=10, eigentol=0.001, info={})[source]

Make a compact source map from the TOD. Specify filename to write things to disk; this should be a format string, for example ‘{obs_id}_{map}.fits’, where ‘map’ will be given values of [‘binned’, ‘solved’, ‘weights’] and any other keys (obs_id in this case) must be passed through info.

Mapmaking for HWP data

sotodlib.coords.demod.make_map(tod, P=None, wcs_kernel=None, res=0.0017453292519943296, dsT=None, demodQ=None, demodU=None, cuts=None, det_weights=None, det_weights_demod=None, center_on=None, flip_gamma=True)[source]

Generates maps of temperature and polarization from demodulated HWP data. It is assumed that each detector contributes a T signal (which has been low-pass filtered to avoid the modulated signal) stored in dsT, as well as separate Q and U timestreams, corresponding to the cosine-phase (demodQ = Q cos 4 chi) and sine-phase (demodU = U sin 4 chi) response of the HWP.

The demodQ and demodU signals are assumed to have been computed without regard for the polarization angle of the detector, nor the on-sky parallactic angle. The impact of these is handled by the projection routines in this function.

Parameters:
  • tod (sotodlib.core.AxisManager) – An AxisManager object

  • P (sotodlib.coords.pmat) – Projection Matrix to be used for mapmaking. If None, it is generated from wcs_kernel or center_on.

  • wcs_kernel (enlib.wcs.WCS or None, optional) – The WCS object used to generate the output map. If None, a new WCS object with a Cartesian projection and a resolution of res will be created.

  • center_on (None or str) – Name of point source. If specified, Source-centered map with a tangent projection will be made. If None, wcs_kernel will be used to generate a projection matrix.

  • size (float or None) – If specified, trim source-centored map to a local square with size, in radian. If None, not trimming will be applied. Only valid when centor_on is specified.

  • res (float, optional) – The resolution of the output map, in radian.

  • dsT (array-like or None, optional) – The input dsT timestream data. If None, the ‘dsT’ field of tod will be used.

  • demodQ (array-like or None, optional) – The input demodulated Q timestream data. If None, the ‘demodQ’ field of tod will be used.

  • demodU (array-like or None, optional) – The input demodulated U timestream data. If None, the ‘demodU’ field of tod will be used.

  • cuts (RangesMatrix or None, optional) – A RangesMatrix that identifies samples that should be excluded from projection operations. If None, no cuts will be applied.

  • det_weights (array-like or None, optional) – The detector weights to use in the map-making for the dsT timestream.

  • det_weights_demod (array-like or None, optional) – The detector weights to use in the map-making for the demodulated Q and U timestreams. If both of det_weights and det_weights_demod are None, uniform detector weights will be used. If only one of two are provided, the other weight is provided by det_weights = 2 * det_weights_demod.

  • flip_gamma (bool or None) – Defaults to True. When constructing the pointing matrix, reflect the nominal focal plane polarization angles about zero (gamma -> -gamma). If you pass in your own P, make sure it was constructed with hwp=True.

Returns:

  • A dictionary which contains

  • ’map’ (enmap.ndmap) – map of temperature and polarization

  • ’weighted_map’ (enmap.ndmap) – The inverse variance weighted map of temperature and polarization

  • ’weight’ (enmap.ndmap) – The map of inverse variance weights used in the map-making process.

sotodlib.coords.demod.from_map(tod, signal_map, cuts=None, flip_gamma=True, wrap=False, modulated=False)[source]

Generate simulated TOD with HWP from a given signal map.

Parameters:
  • tod – an axisManager object

  • signal_map – pixell.enmap.ndmap containing (Tmap, Qmap, Umap) representing the signal.

  • cuts (RangesMatrix, optional) – Cuts to apply to the data. Default is None.

  • flip_gamma (bool, optional) – Whether to flip detector coordinate. If you use the HWP, keep it True. Default is True.

  • wrap (bool, optional) – Whether to wrap the simulated data. Default is False.

  • modulated (bool, optional) – If True, return modulated signal. If False, return the demodulated signal

  • (dsT

  • demodQ

  • False. (and demodU). Default is)

Returns:

A tuple containing the TOD (np.array) of dsT, demodQ and demodU. modulate==True : The modulated TOD (np.array)

Return type:

modulate==False

Local coordinate systems

sotodlib.coords.local.ellipsoid(model='WGS84')[source]

Return the major and minor semiaxis given an ellipsoid model.

Returns:

  • a (astropy.Quantity, semiaxis major in meter)

  • b (astropy.Quantity, semiaxis minor in meter)

sotodlib.coords.local.ecef2lonlat(x, y, z, ell='WGS84')[source]

Convert ECEF to geodetic coordinates

Based on: You, Rey-Jer. (2000). ‘Transformation of Cartesian to Geodetic Coordinates without Iterations’. Journal of Surveying Engineering. doi: 10.1061/(ASCE)0733-9453

Parameters:
  • x (float, array (numpy or Quantity)) – target x ECEF coordinate. If numpy float or array, it is considered in meters

  • y (float, array (numpy or Quantity)) – target y ECEF coordinate. If numpy float or array, it is considered in meters

  • z (float, array (numpy or Quantity)) – target z ECEF coordinate. If numpy float or array, it is considered in meters

  • ell (string, optional) – reference ellipsoid

Returns:

  • lat (float, array (Quantity)) – target geodetic latitude

  • lon (float, array (Quantity)) – target geodetic longitude

  • alt (float, array (Quantity)) – target altitude above geodetic ellipsoid

sotodlib.coords.local.hor2enu(az, el, srange, deg=True)[source]

Azimuth, Elevation, Slant range to target to East, North, Up

Parameters:
  • azimuth (float, array (numpy or Quantity)) – azimuth, if numpy float or array the unit is determined by the deg parameter

  • elevation (float) – elevation, if numpy float or array the unit is determined by the deg parameter

  • srange (float, array (numpy or Quantity)) – slant range. If numpy float or array, it is considered in meters

  • deg (bool, optional) – if azimuth and elevation are not astropy quantities, if True set them to degrees

Returns:

  • E (Quantity) – East ENU coordinate (meters)

  • N (float) – North ENU coordinate (meters)

  • U (Quantity) – Up ENU coordinate (meters)

sotodlib.coords.local.enu2ecef(E, N, U, lon, lat, alt, ell='WGS84', deg=True)[source]

East, North, Up to ECEF.

Parameters:
  • E (float, array (numpy or Quantity)) – target east ENU coordinate. If numpy float or array, it is considered in meters

  • N (float, array (numpy or Quantity)) – target north ENU coordinate. If numpy float or array, it is considered in meters

  • U (float, array (numpy or Quantity)) – target up ENU coordinate. If numpy float or array, it is considered in meters

  • lon (float, array (numpy or Quantity)) – geodetic longitude of the observer, if numpy float or array the unit is determined by the deg parameter

  • lat (float, array (numpy or Quantity)) – geodetic latitude of the observer, if numpy float or array the unit is determined by the deg parameter

  • alt (float, array (numpy or Quantity)) – altitude above geodetic ellipsoid of the observer, If numpy float or array, it is considered in meters

  • ell (string, optional) – reference ellipsoid

  • deg (bool, optional) – if lon and lat are not quantities, if True set them to degrees

Returns:

  • x (Quantity) – target x ECEF coordinate

  • y (Quantity) – target y ECEF coordinate

  • z (Quantity) – target z ECEF coordinate

sotodlib.coords.local.lonlat2ecef(lon, lat, alt, ell='WGS84', deg=True, uncertainty=False, delta_lon=0, delta_lat=0, delta_alt=0)[source]

Convert geodetic coordinates to ECEF

Parameters:
  • lon (float, array (numpy or Quantity)) – geodetic longitude, if numpy float or array the unit is determined by the deg parameter

  • lat (float, array (numpy or Quantity)) – geodetic latitude, if numpy float or array the unit is determined by the deg parameter

  • alt (float, array (numpy or Quantity)) – altitude above geodetic ellipsoid, If numpy float or array, it is considered in meters

  • ell (string, optional) – reference ellipsoid

  • deg (bool, optional) – if azimuth and elevation are not astropy quantities, if True set them to degrees

  • ell – reference ellipsoid

  • uncertainty (bool, optional) – if True computes the uncertainties associated in ECEF coordinates

  • delta_lon (float, array (numpy or Quantity)) – delta geodetic longitude, if numpy float or array the unit is determined by the deg parameter

  • delta_lat (float, array (numpy or Quantity)) – delta geodetic latitude, if numpy float or array the unit is determined by the deg parameter

  • delta_alt (float, array (numpy or Quantity)) – delta altitude above geodetic ellipsoid, if numpy float or array the unit is determined by the deg parameter

Returns:

  • x (Quantity) – target x ECEF coordinate

  • y (Quantity) – target y ECEF coordinate

  • z (Quantity) – target z ECEF coordinate

sotodlib.coords.local.ecef2enu(ecef_obs_x, ecef_obs_y, ecef_obs_z, ecef_target_x, ecef_target_y, ecef_target_z, delta_obs_x, delta_obs_y, delta_obs_z, delta_target_x, delta_target_y, delta_target_z, lon, lat, deg=True)[source]

Return the position relative to a refence point in a ENU system.

If one of delta_obs or delta_target are different from zeros then uncertainties is calculated automatically. The array returns as EAST, NORTH, UP

Parameters:
  • ecef_obs_x (float or array, numpy or Quantity) – x ECEF coordinates array of the observer, if numpy float or array the unit is meter

  • ecef_obs_y (float or array, numpy or Quantity) – y ECEF coordinates array of the observer, if numpy float or array the unit is meter

  • ecef_obs_z (float or array, numpy or Quantity) – z ECEF coordinates array of the observer, if numpy float or array the unit is meter

  • ecef_target_x (float or array, numpy or Quantity) – x ECEF coordinates array of the target, if numpy float or array the unit is meter

  • ecef_target_y (float or array, numpy or Quantity) – y ECEF coordinates array of the target, if numpy float or array the unit is meter

  • ecef_target_z (float or array, numpy or Quantity) – z ECEF coordinates array of the target, if numpy float or array the unit is meter

  • delta_obs_x (numpy or Quantity, numpy or Quantity) – x ECEF Coordinates error of the observer, if numpy array the unit is meter

  • delta_obs_y (numpy or Quantity, numpy or Quantity) – x ECEF Coordinates error of the observer, if numpy array the unit is meter

  • delta_obs_z (numpy or Quantity, numpy or Quantity) – x ECEF Coordinates error of the observer, if numpy array the unit is meter

  • delta_target_x (numpy or Quantity, numpy or Quantity) – x ECEF Coordinates error of the target, if numpy array the unit is meter

  • delta_target_y (numpy or Quantity, numpy or Quantity) – y ECEF Coordinates error of the target, if numpy array the unit is meter

  • delta_target_z (numpy or Quantity, numpy or Quantity) – z ECEF Coordinates error of the target, if numpy array the unit is meter

  • lon (float, array (numpy or Quantity)) – geodetic longitude, if numpy float or array the unit is determined by the deg parameter

  • lat (float, array (numpy or Quantity)) – geodetic latitude, if numpy float or array the unit is determined by the deg parameter

  • deg (bool, optional) – if azimuth and elevation are not astropy quantities, if True set them to degrees

Returns:

  • E (Quantity) – East ENU coordinate error

  • N (float) – North ENU coordinate error

  • U (Quantity) – Up ENU coordinate error

  • delta_E (Quantity) – East ENU coordinate error

  • delta_N (float) – North ENU coordinate error

  • delta_U (Quantity) – Up ENU coordinate error

sotodlib.coords.local.enu2hor(E, N, U, delta_E, delta_N, delta_U, degrees=True)[source]

Compute Horizontal coordinates given the coordinates in the East-North-Up system.

Parameters:
  • E (float, array (numpy or Quantity)) – target east ENU coordinate. If numpy float or array, it is considered in meters

  • N (float, array (numpy or Quantity)) – target north ENU coordinate. If numpy float or array, it is considered in meters

  • U (float, array (numpy or Quantity)) – target up ENU coordinate. If numpy float or array, it is considered in meters

Focal Plane from Physical Optics

Map positions on physical focal plane to sky using physical optics. Also includes tools for computing transforms between point clouds.

LAT code adapted from code provided by Simon Dicker.

sotodlib.coords.optics.gen_pol_endpoints(x, y, pol)[source]

Get end points of unit vectors that will be centered on the provided xy positions and have the angles specified by pol.

Parameters:
  • x – X positions nominally in mm.

  • y – Y positions nominally in mm.

  • pol – Angles in degrees.

Returns:

X postions of endpoints where the even indices are starts and the odds are ends.

y_pol: Y postions of endpoints where the even indices are starts and the odds are ends.

Return type:

x_pol

sotodlib.coords.optics.get_gamma(pol_xi, pol_eta)[source]

Convert xi, eta endpoints to angles that correspond to gamma.

Parameters:
  • pol_xi – 2n xi values where the ones with even indices are starting points and odd are ending points where each pair forms a vector whose angle is gamma.

  • pol_eta – 2n eta values where the ones with even indices are starting points and odd are ending points where each pair forms a vector whose angle is gamma.

Returns:

Array of n gamma values in radians.

Return type:

gamma

sotodlib.coords.optics.load_ufm_to_fp_config(config_path)[source]

Load and cache config file with the parameters to transform from UFM to focal_plane coordinates.

Parameters:

config_path – Path to the yaml config file.

Returns:

Dictionairy containing the config information.

Return type:

config

sotodlib.coords.optics.get_ufm_to_fp_pars(telescope_flavor, wafer_slot, config_path)[source]

Get (and cache) the parameters to transform from UFM to focal plane coordinates for a specific slot of a given telescope’s focal plane.

Parameters:
  • telescope_flavor – The telescope, should be LAT or SAT.

  • wafer_slot – The wafer slot to get parameters for.

  • config_path – Path to the yaml with the parameters.

Returns:

Dict of transformation parameters that can be passed to ufm_to_fp.

Return type:

transform_pars

sotodlib.coords.optics.ufm_to_fp(aman, x=None, y=None, pol=None, theta=0, dx=0, dy=0)[source]

Transform from coords internal to wafer to focal plane coordinates.

Parameters:
  • aman – AxisManager assumed to contain aman.det_info.wafer. If provided outputs will be wrapped into aman.focal_plane.

  • x – X position in wafer’s internal coordinate system in mm. If provided overrides the value from aman.

  • y – Y position in wafer’s internal coordinate system in mm. If provided overrides the value from aman.

  • pol – Polarization angle in wafer’s internal coordinate system in deg. If provided overrides the value from aman.

  • theta – Internal rotation of the UFM in degrees.

  • dx – X offset in mm.

  • dy – Y offset in mm.

Returns:

X position on focal plane.

y_fp: Y position on focal plane.

pol_fp: Pol angle on focal plane.

Return type:

x_fp

sotodlib.coords.optics.LAT_pix2sky(x, y, sec2xi, sec2eta, array2secx, array2secy, rot=0, opt2cryo=30.0)[source]

Routine to map pixels from arrays to sky.

Parameters:
  • x – X position on focal plane (currently zemax coord).

  • y – Y position on focal plane (currently zemax coord).

  • sec2xi – Function that maps positions on secondary to on sky xi.

  • sex2eta – Function that maps positions on secondary to on sky eta.

  • array2secx – Function that maps positions on tube’s focal plane to x position on secondary.

  • array2secy – Function that maps positions on tube’s focal plane to y position on secondary.

  • rot – Rotation about the line of site = xi - 60 - corotator.

  • opt2cryo – The rotation to get from cryostat coordinates to zemax coordinates (TBD, prob 30 deg).

Returns:

The on sky elevation in radians.

eta: The on sky eta in radians.

Return type:

xi

sotodlib.coords.optics.load_zemax(path)[source]

Load zemax_dat from path

Parameters:

path – Path to zemax data

Returns:

Dictionairy with data from zemax

Return type:

zemax_dat

sotodlib.coords.optics.LAT_optics(zemax_path)[source]

Compute mapping from LAT secondary to sky.

Parameters:

zemax_path – Path to LAT optics data from zemax.

Returns:

Function that maps positions on secondary to on sky elevation

sex2eta: Function that maps positions on secondary to on sky eta.

Return type:

sec2xi

sotodlib.coords.optics.LATR_optics(zemax_path, tube_slot)[source]

Compute mapping from LAT secondary to sky.

Parameters:
  • zemax_path – Path to LATR optics data from zemax.

  • tube_slot – Either the tube name as a string or the tube number as an int.

Returns:

Function that maps positions on tube’s focal plane to x position on secondary.

array2secy: Function that maps positions on tube’s focal plane to y position on secondary.

Return type:

array2secx

sotodlib.coords.optics.LAT_focal_plane(aman, zemax_path, x=None, y=None, pol=None, roll=0, tube_slot='c1')[source]

Compute focal plane for a wafer in the LAT.

Parameters:
  • aman – AxisManager nominally containing aman.focal_plane.x_fp and aman.focal_plane.y_fp. If provided focal plane will be stored in aman.focal_plane.

  • zemax_path – Path to LATR optics data from zemax.

  • x – Detector x positions, if provided will override positions loaded from aman.

  • y – Detector y positions, if provided will override positions loaded from aman.

  • pol – Detector polarization angle, if provided will override positions loaded from aman.

  • roll – Rotation about the line of site = elev - 60 - corotator.

  • tube_slot – Either the tube name as a string or the tube number as an int.

Returns:

Detector xi on sky from physical optics in radians.

If aman is provided then will be wrapped as aman.focal_plane.xi.

eta: Detector eta on sky from physical optics in radians.

If aman is provided then will be wrapped as aman.focal_plane.eta.

gamma: Detector gamma on sky from physical optics in radians.

If aman is provided then will be wrapped as aman.focal_plane.eta.

Return type:

xi

sotodlib.coords.optics.sat_to_sky(x, theta)[source]

Interpolate x and theta values to create mapping from SAT focal plane to sky. This function is a wrapper whose main purpose is the cache this mapping.

Parameters:
  • x – X values in mm, should be all positive.

  • theta – Theta values in radians, should be all positive. Theta is defined by ISO coordinates.

Returns:

Interp object with the mapping from the focal plane to sky.

Return type:

sat_to_sky

sotodlib.coords.optics.SAT_focal_plane(aman, x=None, y=None, pol=None, roll=0, mapping_data=None)[source]

Compute focal plane for a wafer in the SAT.

Parameters:
  • aman – AxisManager nominally containing aman.focal_plane.x_fp and aman.focal_plane.y_fp. If provided focal plane will be stored in aman.focal_plane.

  • x – Detector x positions, if provided will override positions loaded from aman.

  • y – Detector y positions, if provided will override positions loaded from aman.

  • pol – Detector polarization angle, if provided will override positions loaded from aman.

  • roll – Rotation about the line of site = -1*boresight.

  • mapping_data – Tuple of (x, theta) that can be interpolated to map the focal plane to the sky. Leave as None to use the default mapping.

Returns:

Detector xi on sky from physical optics in radians.

If aman is provided then will be wrapped as aman.focal_plane.xi.

eta: Detector eta on sky from physical optics in radians.

If aman is provided then will be wrapped as aman.focal_plane.eta.

gamma: Detector gamma on sky from physical optics in radians.

If aman is provided then will be wrapped as aman.focal_plane.eta.

Return type:

xi

sotodlib.coords.optics.get_focal_plane(aman, x=None, y=None, pol=None, roll=0, telescope_flavor=None, tube_slot=None, wafer_slot='ws0', config_path=None, zemax_path=None, mapping_data=None, return_fp=False)[source]

Unified interface for LAT_focal_plane and SAT_focal_plane including the step of applying ufm_to_fp.

Parameters:
  • aman – AxisManager nominally containing aman.det_info.wafer.det_x and aman.det_info.wafer.det_y. If provided focal plane will be stored in aman.focal_plane.

  • x – Detector x positions, if provided will override positions loaded from aman.

  • y – Detector y positions, if provided will override positions loaded from aman.

  • pol – Detector polarization angle, if provided will override positions loaded from aman.

  • roll – Rotation about the line of sight. For the LAT this is elev - 60 - corotator. For the SAT this is -1*boresight.

  • telescope_flavor – What the telescope flavor is (‘LAT’ or ‘SAT’). Note that this is case insenstive. If None and aman contains obs_info, will be extracted from there.

  • tube_slot – Which optics tube of the telescope to use. Only used by the LAT. If None and aman contains obs_info, will be extracted from there.

  • wafer_slot – Which wafer slot to use.

  • config_path – Path to the ufm_to_fp config file.

  • zemax_path – Path to the data file from Zemax. Only used by the LAT.

  • mapping_data – Tuple of (x, theta) that can be interpolated to map the focal plane to the sky. Leave as None to use the default mapping. Only used by the SAT.

  • return_fp – If True also return the transformed focal plane in mm.

Returns:

Detector xi on sky from physical optics in radians.

If aman is provided then will be wrapped as aman.focal_plane.xi.

eta: Detector eta on sky from physical optics in radians.

If aman is provided then will be wrapped as aman.focal_plane.eta.

gamma: Detector gamma on sky from physical optics in radians.

If aman is provided then will be wrapped as aman.focal_plane.eta.

x_fp: X position on focal plane.

Only returned if return_fp is True.

y_fp: Y position on focal plane.

Only returned if return_fp is True.

pol_fp: Pol angle on focal plane.

Only returned if return_fp is True.

Return type:

xi

sotodlib.coords.optics.gen_template(wafer_info, stream_id, **pointing_cfg)[source]

Generate a pointing template from a wafer info ResultSet.

Parameters:
  • wafer_info – Either the path to a wafer_info ResultSet or an open h5py File object.

  • stream_id – The UFM to generate the template for (ie: ufm_mv5).

  • **pointing_cfg – Arguments to pass to get_focal_plane, should be arguments from rot onwards.

Returns:

The detector ids of the detectors in the template.

template: (ndet, 3) array of (xi, eta, gamma) for all dets in template.

is_optical: (ndet) mask that is True for all optical dets.

Return type:

template_det_ids