maps

The maps project defines map projections, along with G3SkyMap subclasses that provide sky maps in those map projections and tools for format/projection conversions of these data types.

In addition, this library contains three pipeline modules (MapBinner, SingleDetectorMapBinner, and SingleDetectorBoresightBinner) to make binned maps from time-ordered data, as well as a module (MapMockObserver) to mock-observe a provided sky map, generating fake time-ordered-data from it corresponding to some stored instrument pointing itinerary. A few other utility pipeline modules (described below) are provided for some map manipulation tasks.

The key data types defined here are:

HealpixSkyMap

Implements Healpix over all or a fraction of the sky, either in nested or ring mode. The underlying sky map data are represented in one of three ways: as a dense 1-D array (full sky), as a locally-dense region surrounded by zeroes (ring mode only), or as a list of non-zero pixels and their values. The second two modes efficiently represent partial-sky maps.

FlatSkyMap

Implements a flat-sky map (similar to a 2D numpy array) in any of the supported projections. The stored map is either a dense 2D array, or a locally dense region (a set of neighboring columns containing non-zero values, each of which contains a single contiguous block of non-zero values).

Map Attributes

The following attributes are common to all G3SkyMap subclasses:

coord_ref

The coordinate system on the sky to which each map pixel is referenced, stored as an instance of the MapCoordReference enum. Currently supported coordinate systems are Equatorial (FK5 J2000), Galactic and Local (telescope azimuth and elevation).

pol_type

The Stokes polarization of the map object, which is an instance of the MapPolType enum, and can have the value T, Q, U or None.

pol_conv

The polarization convention used to encode the Q and U Stokes orientations relative to the coordinate axes. This attribute is an instance of the MapPolConv enum, which can have the value IAU, COSMO or None. Both IAU and COSMO polarization conventions are supported in polarization-aware functions (e.g. FlattenPol), but most default to using the IAU convention. Warnings will be raised when a polarized map is used without a polarization convention set. Use the SetPolConv pipeline module to change the polarization convention between IAU and COSMO for an entire map frame. This will result in flipping the sign of all pixels in the U map as well as the TU and QU weights.

units

The units system in which the map is computed, stored as an instance of the G3TimestreamUnits enum, typically Tcmb.

weighted

A boolean attribute indicating whether the data in the map have been normalized by the inverse of the appropriate Mueller matrix (weighted=False) or not (weighted=True). See more information on map weights below.

File Format Conversions

We support writing maps into FITS files that can be read with other tools (such as DS9), using the fitsio.save_skymap_fits function. FITS files with compatible headers can be read in using the fitsio.load_skymap_fits function.

T, Q, U and corresponding G3SkyMapWeights objects are written to a single file as a sequence of HDUs. FlatSkyMap objects are stored in dense format (see below) to CompImageHDU objects if compression is enabled, and otherwise stored in dense format to standard ImageHDU objects. The latter can be loaded using “old style” fits readers, such as the idlastro fits utilities.

HealpixSkyMap objects are stored in a sequence of BinTableHDU objects, a format that is compatible with the healpy.read_map function. Dense maps (see below) are stored using implicit indexing, and sparse maps are stored using explicit indexing with an additional pixel index column.

Indexing

Values in maps can be set and retrieved using the standard python (or C++) [] operator. Both flat and Healpix maps support a 1-D indexing convention. For flat-sky maps, this 1-D index follows C ordering; for Healpix maps, this is the normal 1-D Healpix pixel number. Flat-sky maps also accept 2-D indices, which have ordering following normal language conventions for 2-D indices ((y, x) in Python, (x, y) in C++).

Note that sky maps do not support numpy-style slicing operations, except for 2-D indexing of flat-sky maps (see below), which makes a copy of the underlying map data. To perform operations with other numpy arrays, use numpy.asarray, which will convert the map to its dense representation (see below) and provides read-write access the map’s internal buffer, which requires no meaningful CPU time or memory.

Sparsity

By default, both Healpix and flat-sky maps are initialized in sparse mode. This imposes a slight performance penalty but will result in the map storing only non-zero portions (with caveats, see details above), substantially reducing RAM usage. Some map operations, in particular casting to numpy arrays, will result in the implicit conversion of the map to dense storage, which can result in sudden increases in RAM usage. The current sparsity mode can be examined or changed with the sparse property (flat sky maps) or the dense, ringsparse, or indexedsparse properties (Healpix maps). Serialization to .g3 files will maintain the current sparsity scheme, as do arithmetic operators where possible. Serialization to .fits files implicitly converts flat sky maps to dense mode, but preserves the sparsity of Healpix maps. The current number of stored pixels can be obtained using the npix_allocated property, and the number of non-zero pixels can be obtained using the npix_nonzero property. Dense maps can be efficiently compactified in memory using the G3SkyMap.compact method, or the CompactMaps pipeline module.

Beyond paying attention to implicit conversions to dense storage and the performance impact of sparse storage (which is small, at least for FlatSkyMap objects), users of this code do not need to worry about the storage mode–all interfaces are identical in all modes.

Masking

Maps containing only boolean data for each pixel are stored as G3SkyMapMask objects. Such mask objects have a .parent attribute which is a shallow clone of the map object with which they are associated (to check for shape compatibility).

Masks are returned when using comparison operators with map objects, e.g. map1 > 5 or map1 == map2. The supported comparison operators are: >, >=, ==, !=, <=, <. Masks can also be combined together using logical operators, e.g. mask3 = mask1 & mask2 or mask1 ^= mask2. The supported comparison operators are: &, &=, |, |=, ^, ^=. Masks can also be checked for equality to other masks using == and != operators.

Mask objects can be clone’ed the same way as maps. A map can be converted to a boolean mask using G3SkyMap.to_mask(), which returns a mask which is True wherever the map is non-zero (optionally excluding nan or inf pixels). A mask can be converted back to a map object using G3SkyMapMask.to_map(), which returns a sparse, unit-less, unweighted, unpolarized map object of the same type as G3SkyMapMask.parent, containing double 1.0 wherever the mask is True.

Masks can also be applied to maps or masks using the appropriate .apply_mask method, with optional inversion; alternatively maps can also be directly multiplied by a compatible mask object. A list of non-zero pixels can be returned using .nonzero() (note that this returns a single vector of pixel positions), and mask contents can be checked using .all(), .any() and .sum(). Mask contents can be inverted in-place using .invert().

Mask objects cannot be accessed using numpy slicing, or converted directly to arrays, because numpy does not represent boolean values as single bits. To be able to use numpy tools with masks, you need to first convert the mask to a dense map using .to_map(). All associated methods of the parent map are accessible as attributes of the mask object in python, e.g. mask.angles_to_pixels() works as one would expect.

Mask Memory Usage

The current implementation of masks is to use a dense std::vector<bool> as the data storage backend, which uses 64x less memory than a dense map (std::vector<double>) of the same dimensions. This implementation is sufficient for FlatSkyMap objects, since these are typically O(50%) full populated in their sparse state; however, the memory savings for HealpixSkyMap objects is not as significant when observing sufficiently small patches of sky. Future work would enable a similar sparse storage backend for masks.

In general, when working with high-resolution maps of any sort, it is important to think carefully about doing the sorts of operations that can balloon memory usage, e.g. taking care to preserve the sparsity of maps by avoiding numpy operations if possible, or using in-place operations to avoid unintentionally creating extra maps or masks in memory.

Statistics

Most numpy.ufunc-like methods are defined for map objects, namely .all(), .any(), .sum(), .mean(), .median(), .var(), .std(), .min(), .max(), .argmin(), .argmax(). All methods take an optional where argument, which can be a compatible G3SkyMapMask object, or size-compatible 1-D numpy array that can be converted into one. In addition, these methods are called under the hood when using the numpy equivalent functions (numpy.all(), etc), in order to preserve the sparsity of the input map. Methods that ignore NaN values are also defined (.nansum(), etc), which behave much like the standard methods, except that calling numpy.nansum() and friends on a map object does not preserve sparsity.

Map values can be tested using .isnan(), .isinf(), .isfinite() methods as well; these return G3SkyMapMask objects.

Map Interpolation

Several interpolation and rebinning utilities are provided. The method G3SkyMap.get_interp_values can be used for extracting map values at arbitrary sky positions using bilinear interpolation. The method G3SkyMap.rebin can be used to downgrade the map resolution in a way that preserves the total power within each map pixel.

The functions healpix_to_flatsky and flatsky_to_healpix functions are provided to reproject maps between flat sky and curved sky systems, with options to use interpolation or rebinning to improve the accuracy of the reprojection.

The more general reproj_map function can also be used to convert between flat sky projections.

Note: The interpolation routine for healpix maps produces results that differ from those of the equivalent healpy.get_interp_val routine. The interpolation routine implemented here is area-preserving in the computation of bilinear weights, whereas the healpy routine is not.

Map Weights

The G3SkyMapWeights class combines the six unique components of the Mueller weight matrix into one object. The individual matrix terms can be accessed using the attributes G3SkyMapWeights.TT, etc, or as keyed elements (e.g. weights['TT']). The full matrix for an individual map pixel can be accessed using the standard [] operator. In python, this returns a symmetric 3x3 numpy array that is a copy of the values in the underlying maps, and in C++ this returns a MuellerMatrix object, with scalar attributes MuellerMatrix.tt, etc that are writable references to elements of the underlying map objects. The G3SkyMapWeights.polarized attribute determines whether the weight structure contains polarization information. For unpolarized weights, only the TT element is set, and the [] operator returns a scalar value in python, and a MuellerMatrix with just the TT element set in C++.

In C++ there is also a StokesVector object that is analogous to the MuellerMatrix object. It has scalar attributes StokesVector.t etc, that are writable references to elements of map objects. Matrix operations on the StokesVector and MuellerMatrix objects are well defined.

Weights are removed from or applied to a set of Stokes T/Q/U maps simultaneously, using the remove_weights or apply_weights functions, or their corresponding pipeline modules.

Map Frames and Pipelines

Maps and associated weights are generally stored in memory and on disk in G3Frames of type G3FrameType.Map, with keys 'T', 'Q', 'U', 'Wpol' defined for polarized maps, and 'T', 'Wunpol' defined for unpolarized maps. Map frames can be checked for validity using the ValidateFrames pipeline module, which raises errors or warnings for missing keys or inconsistent attributes.

Map frames can be manipulated in a pipeline using some memory-efficient pipeline modules. Weights can be applied or removed from their corresponding Stokes maps using the ApplyWeights or RemoveWeights pipeline modules. Maps can be converted to polarized or unpolarized versions using the MakeMapPolarized and MakeMapUnpolarized modules. They can also be compactified to their most sparse representation using the CompactMaps module.

Existing maps can be injected into a pipeline using the InjectMaps module, and map stubs can be injected using InjectMapStub or ReplicateMaps. Maps can also be extracted from a pipeline using the ExtractMaps module.

Flat Sky Map Projections

For flat-sky maps, we support the following map projections:

ProjSansonFlamsteed

Sanson-Flamsteed (also called the sinusoidal projection). It has equal-area pixels, defined by multiplying azimuth distances by cos(latitude). Mercator-esque in that lines of constant latitude are transformed to lines of constant y. Distances are not preserved. Also known as “proj 0”.

ProjPlateCarree

The Plate-Carree projection just plots latitude and longitude on a grid: latitude lines are at constant y and equally spaced, while longitude lines are at constant x and equally spaced. Pixels are not equal-area. Also known as “proj 1”. A variant of this projection, called ProjBICEP (or “proj 9”), adjusts the resolution along x to scale with the cosine of the latitude of the center of the map.

ProjOrthographic

The projection of the sphere onto a plane – the sky looks like a circle. Can only show one hemisphere. Lines drawn on the map do not correspond to latitude or longitude. Pixels are not equal-area. Also known as “proj 2”.

ProjStereographic

Another projection of the sphere onto a plane that makes it look like a circle. Differs from an orthographic projection in that it lets you see both hemispheres. Popularized in the form of the UN logo. Lines drawn on the map do not correspond to latitude or longitude. Pixels are not equal-area. Also known as “proj 4”.

ProjLambertAzimuthalEqualArea

Yet another mapping of the sphere to a circle, but this one has equal-area pixels. Largely distance-preserving, which makes it particularly useful for power-spectrum analyses. Also known as “proj 5”.

ProjGnomonic

Another projection of the sphere onto a circle. This one has the property that straight lines correspond to geodesics. Does not have equal-area pixels. Can show less than half a sphere. Also known as a “tangent projection” or “proj 6”.

ProjCylindricalEqualArea

The Lambert cylindrical equal-area projection (CEA) maps the sphere to a rectangle. Has equal-area pixels. Lines of constant x correspond to constant longitude; lines of constant y are constant latitude. Latitudes get closer together (by sin(latitude)) at the poles. Also known as “proj 7”.

Flat Sky Map Manipulation

Flat sky maps have additional functions defined for efficient manipulation in memory.

The FlattenPol pipeline module flattens the Q and U stokes parameters to align with the pixel coordinate grid, which is necessary for computing power spectra in the flat sky approximation.

Small patches can be extracted from and inserted into larger flat sky maps using the FlatSkyMap.extract_patch and FlatSkyMap.insert_patch methods, respectively. Also, maps can be padded and cropped using the FlatSkyMap.reshape method, which keeps the patch centered in the output map. All of these preserve the map pixelization and correspondence to angle on the sky.

As an equivalent and more Pythonic alternative, you can also extract portions of the map using numpy-style slicing operations (e.g. map[45:130,114:182]), which will produce a map with the same contents as the numpy operation but without converting it to a dense map and with all the coordinate information set appropriately (and is equivalent to extract_patch()). This also works with setting, but the coordinates have to match the sub-subcoordinates (as you would have gotten them from getting a slice or extract_patch()). Note that this slicing creates a copy of the underlying data, so in-place operations (e.g. map[45:130,114:182] += 5) will work, but are not necessarily memory efficient.

Map Pointing

This package also provides functions and pipeline modules for creating and manipulating the quaternions necessary for mapmaking. In general, there are two forms of quaternions that are used throughout the code: pointing quaternions and rotation quaternions.

Pointing Quaternions

Pointing quaternions encode the two-dimensional sky coordinate angles in their vector component. These quaternions can be created using the ang_to_quat function, and their sky coordinates extracted using the quat_to_ang function. The various methods of the G3SkyMap classes return or accept pointing quaternions. Note that local (horizon) coordinates have a different parity than sky coordinates (equatorial, galactic); the z vector coordinate encodes -sin(elevation) in local coordinates, but +sin(dec) in sky coordinates.

Rotation Quaternions

Conversion between coordinate systems is done by constructing rotation quaternions. A pointing quaternion q_p can be rotated to a new coordinate system by the rotation quaternion q_r by using quaternion multiplication: q_p_rot = q_r * q_p / q_r. For example, the module FillCoordTransRotations can be used to construct rotation quaternions for rotating detector offset coordinates into local or on-sky coordinate systems. Rotation quaternions can be rotated into Galactic coordinates using the EquatorialToGalacticTransRotations module.

Detector Pointing

Detector pointing timestreams are constructed by first using the offsets_to_quat function to construct the detector offset quaternion in boresight coordinates, then rotating that pointing quaternion onto the sky by applying a rotation quaternion constructed from the boresight pointing timestreams. This is done internally for each detector in each of the mapmaking pipeline modules (MapBinner, MapMockObserver, etc), which all require an input BolometerPropertiesMap object with offsets for each detector, and pre-computed timestreams of boresight rotation quaternions associated with each input Scan frame.

Frame Objects in spt3g.maps

spt3g.maps.FlatSkyMap

FlatSkyMap is a G3SkyMap with the extra meta information about the particular flat sky projection included. In practice it behaves (mostly) like a 2d numpy array. The pixels are normally indexed with an 1d pixel index. If you find that you need numpy functionality from a FlatSkyMap, e.g. for slicing across the two dimensions, you can access a numpy representation of the map using np.asarray(m). This does not copy the data, so any changes to the resulting array will affect the data stored in the map. Alternatively, you can use 2d slice indexing directly on the map to access a copy of the data with the coordinate representation intact. The latter method is most efficient for extracting small patches from sparse maps.

Constructors:

FlatSkyMap( (FlatSkyMap)arg2)

FlatSkyMap( (int)x_len, (int)y_len, (float)res [, (bool)weighted=True [, (MapProjection)proj=spt3g.maps.MapProjection.ProjNone [, (float)alpha_center=0 [, (float)delta_center=0 [, (MapCoordReference)coord_ref=spt3g.maps.MapCoordReference.Equatorial [, (G3TimestreamUnits)units=spt3g.core.G3TimestreamUnits.Tcmb [, (MapPolType)pol_type=spt3g.maps.MapPolType.none [, (float)x_res=nan [, (float)x_center=nan [, (float)y_center=nan [, (bool)flat_pol=False [, (MapPolConv)pol_conv=spt3g.maps.MapPolConv.none]]]]]]]]]]]])

FlatSkyMap( (object)obj, (float)res [, (bool)weighted=True [, (MapProjection)proj=spt3g.maps.MapProjection.ProjNone [, (float)alpha_center=0 [, (float)delta_center=0 [, (MapCoordReference)coord_ref=spt3g.maps.MapCoordReference.Equatorial [, (G3TimestreamUnits)units=spt3g.core.G3TimestreamUnits.Tcmb [, (MapPolType)pol_type=spt3g.maps.MapPolType.none [, (float)x_res=nan [, (float)x_center=nan [, (float)y_center=nan [, (bool)flat_pol=False [, (MapPolConv)pol_conv=spt3g.maps.MapPolConv.none]]]]]]]]]]]])

FlatSkyMap( [ (FlatSkyMap)flat_map])

Members:

  • proj: Map projection (one of maps.MapProjection)

  • alpha_center: Horizontal axis center position

  • delta_center: Vertical axis center position

  • x_center: Horizontal axis center pixel position

  • y_center: Vertical axis center pixel position

  • res: Map resolution in angular units for maps with square pixels

  • x_res: Resolution in X direction for maps with rectangular pixels

  • y_res: Resolution in Y direction for maps with rectangular pixels

  • sparse: True if the map is stored with column and row zero-suppression, False if every pixel is stored. Map sparsity can be changed by setting this to True (or False).

  • flat_pol: True if this map has been flattened using flatten_pol.

  • wcs: astropy.wcs.WCS instance containing projection information

Methods:

spt3g.maps.FlatSkyMap.array_clone

array_clone( (FlatSkyMap)arg1, (object)array) -> G3SkyMap :

Return a map of the same type, populated with a copy of the input numpy array

spt3g.maps.FlatSkyMap.pixel_to_angle

pixel_to_angle( (G3SkyMap)arg1, (int)pixel) -> tuple :

Compute the sky coordinates of the given 1D pixel

pixel_to_angle( (FlatSkyMap)arg1, (int)x, (int)y) -> DoubleVector :

Compute the sky coordinates of the given 2D pixel (also see xy_to_angle()).

spt3g.maps.FlatSkyMap.xy_to_angle

xy_to_angle( (FlatSkyMap)arg1, (float)x, (float)y) -> DoubleVector :

Compute the sky coordinates of the input flat 2D coordinates

xy_to_angle( (FlatSkyMap)arg1, (DoubleVector)x, (DoubleVector)y) -> tuple :

Compute the sky coordinates of the input flat 2D coordinates (vectorized)

spt3g.maps.FlatSkyMap.angle_to_xy

angle_to_xy( (FlatSkyMap)arg1, (float)alpha, (float)delta) -> DoubleVector :

Compute the flat 2D coordinates of the input sky coordinates

angle_to_xy( (FlatSkyMap)arg1, (DoubleVector)alpha, (DoubleVector)delta) -> tuple :

Compute the flat 2D coordinates of the input sky coordinates (vectorized)

spt3g.maps.FlatSkyMap.xy_to_pixel

xy_to_pixel( (FlatSkyMap)arg1, (float)x, (float)y) -> int :

Compute the pixel index of the input flat 2D coordinates

xy_to_pixel( (FlatSkyMap)arg1, (DoubleVector)x, (DoubleVector)y) -> UInt64Vector :

Compute the pixel indices of the input flat 2D coordinates (vectorized)

spt3g.maps.FlatSkyMap.pixel_to_xy

pixel_to_xy( (FlatSkyMap)arg1, (int)pixel) -> DoubleVector :

Compute the flat 2D coordinates of the input pixel index

pixel_to_xy( (FlatSkyMap)arg1, (UInt64Vector)pixel) -> tuple :

Compute the flat 2D coordinates of the input pixel indices (vectorized)

spt3g.maps.FlatSkyMap.nonzero_pixels

nonzero_pixels( (FlatSkyMap)arg1) -> tuple :

Returns a list of the indices of the non-zero pixels in the map and a list of the values of those non-zero pixels.

spt3g.maps.FlatSkyMap.extract_patch

extract_patch( (FlatSkyMap)arg1, (int)x0, (int)y0, (int)width, (int)height [, (float)fill=0]) -> G3SkyMap :

Returns a map of shape (width, height) containing a rectangular patch of the parent map. Pixel (width // 2, height // 2) of the output map corresponds to pixel (x0, y0) in the parent map, and the angular location of each pixel on the sky is maintained. Surrounding empty pixels can be optionally filled.

spt3g.maps.FlatSkyMap.insert_patch

insert_patch( (FlatSkyMap)arg1, (FlatSkyMap)patch [, (bool)ignore_zeros=False]) -> None :

Inserts a patch (e.g. as extracted using extract_patch) into the parent map. The coordinate system and angular center of the patch must match that of the parent map. If ignore_zeros is True, zero-valued pixels in the patch are not inserted. This is useful for preserving the sparsity of the parent map.

spt3g.maps.FlatSkyMap.reshape

reshape( (FlatSkyMap)arg1, (int)width, (int)height [, (float)fill=0]) -> G3SkyMap :

Returns a map of shape (width, height) containing the parent map centered within it. The angular location of each pixel on the sky is maintained. Surrounding empty pixels can be optionally filled.

spt3g.maps.G3SkyMap

Base class for 1- and 2-D skymaps of various projections. Usually you want a subclass of this (e.g. FlatSkyMap) rather than using it directly.

Constructors:

Raises an exception This class cannot be instantiated from Python

Members:

  • coord_ref: Coordinate system (maps.MapCoordReference) of the map (e.g. Galactic, Equatorial, etc.)

  • pol_type: Polarization type (maps.MapPolType) of the map (e.g. maps.MapPolType.Q).

  • pol_conv: Polarization convention (maps.MapPolConv) of the map (e.g. maps.MapPolConv.IAU or maps.MapPolConv.COSMO).

  • polarized: True if the pol_conv property is set to IAU or COSMO, False otherwise.

  • units: Unit class (core.G3TimestreamUnits) of the map (e.g. core.G3TimestreamUnits.Tcmb). Within each unit class, further conversions, for example from K to uK, should use core.G3Units.

  • weighted: True if map is multiplied by weights

  • size: Number of pixels in map

  • shape: Shape of map

  • npix_allocated: Number of pixels in map currently stored in memory

  • npix_nonzero: Number of nonzero pixels in map currently stored in memory

  • overflow: Combined value of data processed by the map maker but outside of the map area

Methods:

spt3g.maps.G3SkyMap.copy

copy( (G3SkyMap)arg1) -> G3SkyMap :

Return a copy of the map object

spt3g.maps.G3SkyMap.clone

clone( (G3SkyMap)arg1 [, (bool)copy_data=True]) -> G3SkyMap :

Return a map of the same type, populated with a copy of the data if the argument is true (default), empty otherwise.

spt3g.maps.G3SkyMap.array_clone

array_clone( (G3SkyMap)arg1, (object)array) -> G3SkyMap :

Return a map of the same type, populated with a copy of the input numpy array

spt3g.maps.G3SkyMap.compatible

compatible( (G3SkyMap)arg1, (G3SkyMap)arg2) -> bool :

Returns true if the input argument is a map with matching dimensions and boundaries on the sky.

spt3g.maps.G3SkyMap.nonzero

nonzero( (G3SkyMap)arg1) -> UInt64Vector :

Return indices of non-zero pixels in the map

spt3g.maps.G3SkyMap.angles_to_pixels

angles_to_pixels( (G3SkyMap)arg1, (DoubleVector)alphas, (DoubleVector)deltas) -> UInt64Vector :

Compute the 1D pixel location for each of the sky coordinates (vectorized)

spt3g.maps.G3SkyMap.pixels_to_angles

pixels_to_angles( (G3SkyMap)arg1, (UInt64Vector)pixels) -> tuple :

Compute the sky coordinates of each of the given 1D pixels (vectorized)

spt3g.maps.G3SkyMap.angle_to_pixel

angle_to_pixel( (G3SkyMap)arg1, (DoubleVector)alphas, (DoubleVector)deltas) -> UInt64Vector :

Compute the 1D pixel location for each of the sky coordinates (vectorized)

angle_to_pixel( (G3SkyMap)arg1, (float)alpha, (float)delta) -> int :

Compute the 1D pixel location of the given sky position.

spt3g.maps.G3SkyMap.pixel_to_angle

pixel_to_angle( (G3SkyMap)arg1, (UInt64Vector)pixels) -> tuple :

Compute the sky coordinates of each of the given 1D pixels (vectorized)

pixel_to_angle( (G3SkyMap)arg1, (int)pixel) -> tuple :

Compute the sky coordinates (alpha, delta) of the given 1D pixel

spt3g.maps.G3SkyMap.quats_to_pixels

quats_to_pixels( (G3SkyMap)arg1, (G3VectorQuat)arg2) -> UInt64Vector :

Compute the 1D pixel location for each of the sky coordinates (vectorized), expressed as quaternion rotations from the pole.

spt3g.maps.G3SkyMap.pixels_to_quats

pixels_to_quats( (G3SkyMap)arg1, (UInt64Vector)arg2) -> G3VectorQuat :

Compute the sky coordinates, expressed as quaternion rotations from the pole, for each of the given 1-D pixel coordinates.

spt3g.maps.G3SkyMap.quat_to_pixel

quat_to_pixel( (G3SkyMap)arg1, (quat)arg2) -> int :

Compute the 1D pixel location of the given sky position, expressed as a quaternion rotation from the pole.

quat_to_pixel( (G3SkyMap)arg1, (G3VectorQuat)arg2) -> UInt64Vector :

Compute the 1D pixel location for each of the sky coordinates (vectorized), expressed as quaternion rotations from the pole.

spt3g.maps.G3SkyMap.pixel_to_quat

pixel_to_quat( (G3SkyMap)arg1, (int)arg2) -> quat :

Compute the quaternion rotation from the pole corresponding to the given 1D pixel.

pixel_to_quat( (G3SkyMap)arg1, (UInt64Vector)arg2) -> G3VectorQuat :

Compute the sky coordinates, expressed as quaternion rotations from the pole, for each of the given 1-D pixel coordinates.

spt3g.maps.G3SkyMap.query_disc

query_disc( (G3SkyMap)arg1, (float)alpha, (float)delta, (float)radius) -> UInt64Vector :

Return a list of pixel indices whose centers are located within a disc of the given radius at the given sky coordinates.

query_disc( (G3SkyMap)arg1, (quat)quat, (float)radius) -> UInt64Vector :

Return a list of pixel indices whose centers are located within a disc of the given radius at the given sky coordinates.

spt3g.maps.G3SkyMap.query_alpha_ellipse

query_alpha_ellipse( (G3SkyMap)arg1, (float)alpha, (float)delta, (float)a, (float)b) -> UInt64Vector :

Return a list of pixel indices whose centers are located within an ellipse extended in the alpha direction, at the given alpha and delta sky coordinates, with semimajor and semiminor axes a and b.

query_alpha_ellipse( (G3SkyMap)arg1, (quat)quat, (float)a, (float)b) -> UInt64Vector :

Return a list of pixel indices whose centers are located within an ellipse extended in the alpha direction, at the given alpha and delta sky coordinates, with semimajor and semiminor axes a and b.

spt3g.maps.G3SkyMap.get_interp_values

get_interp_values( (G3SkyMap)arg1, (DoubleVector)alphas, (DoubleVector)deltas) -> DoubleVector :

Return the values at each of the input coordinate locations. Computes each value using bilinear interpolation over the map pixels.

get_interp_values( (G3SkyMap)arg1, (G3VectorQuat)quats) -> DoubleVector :

Return the values at each of the input coordinate locations. Computes each value using bilinear interpolation over the map pixels.

spt3g.maps.G3SkyMap.rebin

rebin( (G3SkyMap)arg1, (int)scale [, (bool)norm=True]) -> G3SkyMap :

Rebin the map into larger pixels by summing (if norm is false) or averaging (if norm is true) scale-x-scale blocks of pixels together. Returns a new map object. Map dimensions must be a multiple of the rebinning scale.

spt3g.maps.G3SkyMap.compact

compact( (G3SkyMap)arg1 [, (bool)zero_nans=False]) -> None :

Convert the map to its default sparse representation, removing empty pixels, and optionally also removing NaN values. A map that is already sparse will be compactified in place in its current representation without additional memory overhead.

spt3g.maps.G3SkyMap.to_mask

to_mask( (G3SkyMap)arg1 [, (bool)zero_nans=False [, (bool)zero_infs=False]]) -> G3SkyMapMask :

Create a G3SkyMapMask object from the parent map, with pixels set to true where the map is non-zero (and optionally non-nan and/or finite).

spt3g.maps.G3SkyMap.apply_mask

apply_mask( (G3SkyMap)arg1, (G3SkyMapMask)mask [, (bool)inverse=False]) -> None :

Apply a mask in-place to the map, optionally inverting which pixels are zeroed. If inverse = False, this is equivalent to in-place multiplication by the mask.

spt3g.maps.G3SkyMap.median

median( (G3SkyMap)arg1 [, (G3SkyMapMask)where=None]) -> float

spt3g.maps.G3SkyMap.nansum

nansum( (G3SkyMap)arg1 [, (G3SkyMapMask)where=None]) -> float

spt3g.maps.G3SkyMap.nanmean

nanmean( (G3SkyMap)arg1 [, (G3SkyMapMask)where=None]) -> float

spt3g.maps.G3SkyMap.nanmedian

nanmedian( (G3SkyMap)arg1 [, (G3SkyMapMask)where=None]) -> float

spt3g.maps.G3SkyMap.nanvar

nanvar( (G3SkyMap)arg1 [, (int)ddof=0 [, (G3SkyMapMask)where=None]]) -> float

spt3g.maps.G3SkyMap.nanstd

nanstd( (G3SkyMap)arg1 [, (int)ddof=0 [, (G3SkyMapMask)where=None]]) -> float

spt3g.maps.G3SkyMap.nanmin

nanmin( (G3SkyMap)arg1 [, (G3SkyMapMask)where=None]) -> float

spt3g.maps.G3SkyMap.nanmax

nanmax( (G3SkyMap)arg1 [, (G3SkyMapMask)where=None]) -> float

spt3g.maps.G3SkyMap.nanargmin

nanargmin( (G3SkyMap)arg1 [, (G3SkyMapMask)where=None]) -> int

spt3g.maps.G3SkyMap.nanargmax

nanargmax( (G3SkyMap)arg1 [, (G3SkyMapMask)where=None]) -> int

spt3g.maps.G3SkyMap.isinf

isinf( (G3SkyMap)arg1 [, (G3SkyMapMask)where=None]) -> G3SkyMapMask

spt3g.maps.G3SkyMap.isnan

isnan( (G3SkyMap)arg1 [, (G3SkyMapMask)where=None]) -> G3SkyMapMask

spt3g.maps.G3SkyMap.isfinite

isfinite( (G3SkyMap)arg1 [, (G3SkyMapMask)where=None]) -> G3SkyMapMask

spt3g.maps.G3SkyMapMask

Boolean mask of a sky map. Set pixels to use to true, pixels to ignore to false. If use_data set in contrast, mask initialized to true where input map is non-zero; otherwise, all elements are initialized to zero. Use zero_nans or zero_infs to exclude nan or inf elements from the mask.

Constructors:

G3SkyMapMask( (G3SkyMapMask)arg2)

G3SkyMapMask( (G3SkyMap)parent [, (bool)use_data=False [, (bool)zero_nans=False [, (bool)zero_infs=False]]]) :

Instantiate a G3SkyMapMask from a parent G3SkyMap

G3SkyMapMask( (G3SkyMap)parent, (object)data [, (bool)zero_nans=False [, (bool)zero_infs=False]]) :

Instantiate a G3SkyMapMask from a 1D numpy array

Members:

  • parent: “Parent” map which contains no data, but can be used to retrieve the parameters of the map to which this mask corresponds.

  • size: Number of pixels in mask

  • __array_interface__: No Doc (Shame!)

Methods:

spt3g.maps.G3SkyMapMask.clone

clone( (G3SkyMapMask)arg1 [, (bool)copy_data=True]) -> G3SkyMapMask :

Return a mask of the same type, populated with a copy of the data if the argument is true (default), empty otherwise.

spt3g.maps.G3SkyMapMask.array_clone

array_clone( (G3SkyMapMask)arg1, (object)data [, (bool)zero_nans=False [, (bool)zero_infs=False]]) -> G3SkyMapMask :

Return a mask of the same type, populated from the input numpy array

spt3g.maps.G3SkyMapMask.compatible

compatible( (G3SkyMapMask)arg1, (G3SkyMapMask)arg2) -> bool :

Returns true if the two masks can be applied to the same map.

compatible( (G3SkyMapMask)arg1, (G3SkyMap)arg2) -> bool :

Returns true if this mask can be applied to the given map.

spt3g.maps.G3SkyMapMask.invert

invert( (G3SkyMapMask)arg1) -> G3SkyMapMask :

Invert all elements in mask

spt3g.maps.G3SkyMapMask.nonzero

nonzero( (G3SkyMapMask)arg1) -> UInt64Vector :

Return a list of indices of non-zero pixels in the mask

spt3g.maps.G3SkyMapMask.apply_mask

apply_mask( (G3SkyMapMask)arg1, (G3SkyMapMask)mask [, (bool)inverse=False]) -> None :

Apply a mask in-place to the mask, optionally inverting which pixels are zeroed. If inverse = False, this is equivalent to in-place element-wise logical-and with the mask.

spt3g.maps.G3SkyMapMask.to_map

to_map( (G3SkyMapMask)arg1) -> G3SkyMap :

Create a skymap with data set to the contents of this mask (1.0 where True, 0.0 where False), which can be useful for plotting.

spt3g.maps.G3SkyMapWeights

Polarized (Mueller matrix) or unpolarized (scalar) map pixel weights.Weights are polarized if the pol_conv attribute of the reference map is set.

Constructors:

G3SkyMapWeights()

G3SkyMapWeights( (G3SkyMapWeights)arg2)

G3SkyMapWeights( (object)skymap)

Members:

  • TT: Mueller matrix component map

  • TQ: Mueller matrix component map

  • TU: Mueller matrix component map

  • QQ: Mueller matrix component map

  • QU: Mueller matrix component map

  • UU: Mueller matrix component map

  • size: Number of pixels in weights

  • shape: Shape of weights

  • polarized: True if all components are set, False if only the TT component is set

  • congruent: True if all components are internally compatible with each other

Methods:

spt3g.maps.G3SkyMapWeights.copy

copy( (G3SkyMapWeights)arg1) -> G3SkyMapWeights :

Return a copy of the weights object

spt3g.maps.G3SkyMapWeights.compatible

compatible( (G3SkyMapWeights)arg1, (G3SkyMap)arg2) -> bool :

Returns true if the input argument is a map with matching dimensions and boundaries on the sky.

spt3g.maps.G3SkyMapWeights.rebin

rebin( (G3SkyMapWeights)arg1, (int)scale) -> G3SkyMapWeights :

Rebin the weights into larger pixels by summing scale-x-scale blocks of pixels together. Returns a new weights object. Map dimensions must be a multiple of the rebinning scale.

spt3g.maps.G3SkyMapWeights.compact

compact( (G3SkyMapWeights)arg1 [, (bool)zero_nans=False]) -> None :

Convert the map to its default sparse representation, excluding empty pixels, and optionally converting NaN values to zeroes.

spt3g.maps.G3SkyMapWeights.det

det( (G3SkyMapWeights)arg1) -> G3SkyMap :

Return the determinant of the Mueller matrix for each pixel

spt3g.maps.G3SkyMapWeights.cond

cond( (G3SkyMapWeights)arg1) -> G3SkyMap :

Return the condition number of the Mueller matrix for each pixel

spt3g.maps.G3SkyMapWeights.inv

inv( (G3SkyMapWeights)arg1) -> G3SkyMapWeights :

Return the inverse of the Mueller matrix for each pixel

spt3g.maps.G3SkyMapWeights.apply_mask

apply_mask( (G3SkyMapWeights)arg1, (G3SkyMapMask)mask [, (bool)inverse=False]) -> None :

Apply a mask in-place to the weights, optionally inverting which pixels are zeroed. If inverse = False, this is equivalent to in-place multiplication by the mask.

spt3g.maps.G3SkyMapWeights.clone

clone( (G3SkyMapWeights)arg1 [, (bool)copy_data=True]) -> G3SkyMapWeights :

Return weights of the same type, populated with a copy of the data if the argument is true (default), empty otherwise.

spt3g.maps.HealpixSkyMap

HealpixSkyMap is a G3SkyMap with the extra meta information about the particular Healpix pixelization used. In practice it behaves (mostly) like a 1d numpy array. If you find that you need numpy functionality from a HealpixSkyMap, e.g. for slicing across the array, you can access a numpy representation of the map using np.asarray(m). This does not copy the data, so any changes to the resulting array will affect the data stored in the map.

Constructors:

HealpixSkyMap( (HealpixSkyMap)arg2)

HealpixSkyMap( (int)nside [, (bool)weighted=True [, (bool)nested=False [, (MapCoordReference)coord_ref=spt3g.maps.MapCoordReference.Equatorial [, (G3TimestreamUnits)units=spt3g.core.G3TimestreamUnits.Tcmb [, (MapPolType)pol_type=spt3g.maps.MapPolType.none [, (bool)shift_ra=False [, (MapPolConv)pol_conv=spt3g.maps.MapPolConv.none]]]]]]]) :

Instantiate a HealpixSkyMap with given nside

HealpixSkyMap( (object)data [, (bool)weighted=True [, (bool)nested=False [, (MapCoordReference)coord_ref=spt3g.maps.MapCoordReference.Equatorial [, (G3TimestreamUnits)units=spt3g.core.G3TimestreamUnits.Tcmb [, (MapPolType)pol_type=spt3g.maps.MapPolType.none [, (bool)shift_ra=False [, (MapPolConv)pol_conv=spt3g.maps.MapPolConv.none]]]]]]]) :

Instantiate a Healpix map from existing data. If the data are a single numpy array, assumes this is a dense map. Otherwise, pass an (indices, data, nside) tuple.

HealpixSkyMap( [ (HealpixSkyMap)healpix_map])

Members:

  • nside: Healpix resolution parameter

  • res: Map resolution in angular units

  • nested: True if pixel ordering is nested, False if ring-ordered

  • shift_ra: True if the ringsparse representation of the map is stored with the rings centered at ra = 0 deg, rather than ra = 180 deg.

  • dense: True if the map is stored with all elements, False otherwise. If set to True, converts the map to a dense representation.

  • ringsparse: True if the map is stored as a dense 2D region using ring ordering (analogous to FlatSkyMap’s sparse mode). Ring-sparsity is efficient for dense blocks on a ring-ordered map (e.g. a continuous sky region), but is inefficient otherwise (e.g. nested pixel ordering or discontinous coverage). If set to True, converts the map to this representation.

  • indexedsparse: True if the map is stored as a list of non-zero pixels and values. More efficient than ring-sparse for maps with holes or very small filling factors. If set to True, converts the map to this representation.

Methods:

spt3g.maps.HealpixSkyMap.nonzero_pixels

nonzero_pixels( (HealpixSkyMap)arg1) -> tuple :

Returns a list of the indices of the non-zero pixels in the map and a list of the values of those non-zero pixels.

Modules in spt3g.maps

spt3g.maps.HitsBinner

HitsBinner(map_id, stub_map, pointing, timestreams, bolo_properties_name=”BolometerProperties”, map_per_scan=False)

Bins up pointing for a set of channels into a hits map with properties (projection, etc.) specified by the stub_map argument.

Parameters:
map_idstring

Will be set as the “Id” element of the output frame from this module.

stub_mapG3SkyMap

Template of the map in which to accumulate timestream data. All parameters of the output map (projection, boundaries, pixel size, etc.) are copied from this map, which is not modified.

pointingstring (frame object G3TimestreamQuat)

Name of a frame object containing the boresight pointing quaternion timestream. Must be present in every Scan frame.

timestreamsstring (frame object G3TimestreamMap)

Name of a frame object containing the timestreams to be binned into the output map. Must exist in every Scan frame that is to be binned. This key is used only to determine the subset of channels to bin into the hits map.

bolo_properties_namestring, optional (frame object BolometerPropertiesMap)

Name of a bolometer properties map containing detector pointing offsets from boresight. These are usually named “BolometerProperties”, the default, and this parameter only need be set if you wish to use alternative values.

map_per_scan: boolean or callable, optional

Defaults to False. If set to True, will make an output map for every input Scan frame, which can be useful for fine-grained time-domain maps. Can also be set to a Python callable for complex situations. In this case, the callable will be called on each Scan frame and passed the frame as an argument. If it returns True, the map binner will emit a map frame with the data since the last emitted map frame but before the current Scan; if False, the binner will continue accumulating data.

Emits:

Emits a frame of type Map at the end of processing containing the output map:

Frame (Map) [ “Id” (spt3g.core.G3String) => “ra0hdec-57.5-150GHz” “H” (spt3g.maps.FlatSkyMap) => 2700 x 1500 (45 x 25 deg) ZEA centered at (1350, 749.5) = (0, -57.5 deg) in equatorial coordinates (Tcmb) ]

See Also:
MapBinner, SingleDetectorMapBinner, SingleDetectorBoresightBinner :

Alternative map binners.

FlatSkyMap, HealpixSkyMap :

Possible output types.

Examples:
pipe.Add(
    maps.HitsBinner,
    map_id="150GHz",
    stub_map=map_params,
    timestreams="DeflaggedTimestreams150GHz",
    pointing="OfflineRaDecRotation",
)

spt3g.maps.MapBinner

MapBinner(map_id, stub_map, pointing, timestreams, detector_weights, bolo_properties_name=”BolometerProperties”, store_weight_map=True, map_per_scan=False)

Bins up timestreams into a map with properties (projection, etc.) specified by the stub_map argument.

Parameters:
map_idstring

Will be set as the “Id” element of the output frame from this module.

stub_mapG3SkyMap

Template of the map in which to accumulate timestream data. All parameters of the output map (projection, boundaries, pixel size, etc.) are copied from this map, which is not modified. T-only maps are made if the pol_conv property of stub_map is None; polarized maps are made if it is set to IAU or Cosmo.

pointingstring (frame object G3VectorQuat or G3TimestreamQuat)

Name of a frame object containing the boresight pointing quaternion timestream. Must have the same number of elements as the data in timestreams and be present in every Scan frame.

timestreamsstring (frame object G3TimestreamMap)

Name of a frame object containing the timestreams to be binned into the output map. Must exist in every Scan frame, though may be empty if the frame should be ignored. Units in the output map are taken from the units of the detector timestreams.

detector_weightsstring, optional (frame object G3MapDouble)

Name of a frame object containing weights to assign data from each detector. This frame object maps detector names (as in timestreams) to a scalar weight. If unset or set to an empty string (“”), detectors will be weighted equally. If set, every detector in timestreams must occur in detector_weights, but detector_weights may contain additional elements not present in timestreams.

bolo_properties_namestring, optional (frame object BolometerPropertiesMap)

Name of a bolometer properties map containing detector pointing offsets from boresight. These are usually named “BolometerProperties”, the default, and this parameter only need be set if you wish to use alternative values.

store_weight_mapboolean, optional

Defaults to True. If set to False, the output weight map will not be stored, though the output maps will still be weighted. It should be set to False if and only if you know what the weights are a priori (e.g. mock observing, where they are the same as the data maps).

map_per_scan: boolean or callable, optional

Defaults to False. If set to True, will make an output map for every input Scan frame, which can be useful for fine-grained time-domain maps. Can also be set to a Python callable for complex situations. In this case, the callable will be called on each Scan frame and passed the frame as an argument. If it returns True, the map binner will emit a map frame with the data since the last emitted map frame but before the current Scan; if False, the binner will continue accumulating data.

Emits:

Emits a frame of type Map at the end of processing containing the output map:

Frame (Map) [ “Id” (spt3g.core.G3String) => “ra0hdec-57.5-150GHz” “T” (spt3g.maps.FlatSkyMap) => 2700 x 1500 (45 x 25 deg) ZEA centered at (1350, 749.5) = (0, -57.5 deg) in equatorial coordinates (Tcmb) “Q” (spt3g.maps.FlatSkyMap) => 2700 x 1500 (45 x 25 deg) ZEA centered at (1350, 749.5) = (0, -57.5 deg) in equatorial coordinates (Tcmb) “U” (spt3g.maps.FlatSkyMap) => 2700 x 1500 (45 x 25 deg) ZEA centered at (1350, 749.5) = (0, -57.5 deg) in equatorial coordinates (Tcmb) “Wpol” (spt3g.maps.G3SkyMapWeights) => G3SkyMapWeights ]

For unpolarized maps, the “Q” and “U” outputs will be absent and the “Wpol” frame object will be named “Wunpol”.

See Also:
SingleDetectorMapBinner, SingleDetectorBoresightBinner :

Alternative map binners for making single-detector maps.

MapMockObserver :

The inverse of this module. Produces timestreams from aninput map.

FlatSkyMap, HealpixSkyMap :

Possible output types.

Examples:
pipe.Add(
    maps.MapBinner,
    map_id="150GHz",
    stub_map=map_params,
    timestreams="DeflaggedTimestreams150GHz",
    pointing="OfflineRaDecRotation",
    detector_weights="TodWeights",
)

spt3g.maps.MapMockObserver

MapMockObserver(pointing, timestreams, band, T, Q=None, U=None, bolo_properties_name=”BolometerProperties”, interp=False)

Creates a new set of timestreams by sampling from an input map.

Parameters:
pointingstring (frame object G3TimestreamQuat)

Name of a frame object containing the boresight pointing quaternion timestream. Used to construct detector pointing with which to sample the input map(s). Must be present in every Scan frame.

timestreamsstring (frame object G3TimestreamMap)

Name of a frame object that will contain the output timestreams sampled from the input map(s). Units of the timestreams are taken from the input map. Will be inserted into every Scan frame.

bandfloat

Frequency band for which the input map was constructed. Only detectors whose band matches this value will be included in the output timestreams.

Tframe object G3SkyMap

The input temperature map to be sampled.

Q, Uframe object G3SkyMap, optional

Optional Q and U polarization maps to include in the output timestreams, using the polarization angle and efficiency to determine the appropriate coupling to each Stokes component.

bolo_properties_namestring, optional (frame object BolometerPropertiesMap)

Name of a bolometer properties map containing detector pointing offsets from boresight. These are usually named “BolometerProperties”, the default, and this parameter only need be set if you wish to use alternative values.

interpbool, optional

If True, use bilinear interpolation to sample the input map(s). If False (default), use the nearest-neighbor pixel value.

error_on_zerobool, optional

If True (default), complain loudly if the simulation scans across zeroes in the input map or out of the input map boundaries. If False, allow this to happen without comment.

See Also:
MapBinner :

The inverse of this module. Bins timestreams into an output map.

FlatSkyMap, HealpixSkyMap :

Possible input map types

Examples:
pipe.Add(
    maps.MapMockObserver,
    pointing="OfflineRaDecRotation",
    timestreams="CalTimestreams150GHz",
    band=150 * core.G3Units.GHz,
    T=map_frame["T"],
)

spt3g.maps.MapTODMasker

Creates a filter mask <tod_mask> to indicate to later TOD filtering that a detector has passed over a “bad” part of the sky. This is used, for example, to avoid fitting a polynomial to a timestream while the detector is looking at a bright point source. All timestream values for the detectors in <timestreams> are set to be ignored (true in the output <tod_mask>) if they point at pixels of the input <mask> map containing a non-zero value. Boresight pointing is obtained from the quaternion vector specified by <pointing>. Detector pointing offsets and polarization angles and efficiencies are obtained from the specified BolometerPropertiesMap, which can generally be left at its default value.

spt3g.maps.MapTODPointing

MapTODPointing(pointing, timestreams, stub_map, tod_pointing, bolo_properties_name=”BolometerProperties”)

Compute pixel pointing for timestreams for a map with properties (projection, etc.) specified by the stub_map argument.

The outputs of this module are identical to internally-calculated quantities in the map maker and can be used for debugging map-making or studying individual-detector pixel coverage. The outputs are not consumed by any stage of the map-making pipeline.

Parameters:
pointingstring (frame object G3VectorQuat or G3TimestreamQuat)

Name of a frame object containing the boresight pointing quaternion timestream. Must have the same number of elements as the data in timestreams and be present in every Scan frame.

timestreamsstring (frame object G3TimestreamMap)

Name of a frame object that contains the detector timestreams for the channels for which to compute pointing timestreams.

stub_mapG3SkyMap

Template of the map for which to compute pixel pointing.

tod_pointingstring

The key to which the pixel pointing G3MapVectorInt object is stored in the output Scan frames from this module.

bolo_properties_namestring, optional (frame object BolometerPropertiesMap)

Name of a bolometer properties map containing detector pointing offsets from boresight. These are usually named “BolometerProperties”, the default, and this parameter only need be set if you wish to use alternative values.

Emits:

Adds a new key to every Scan frame containing the output pointing data:

Frame (Scan) [ … “OnlineRaDecRotation” (spt3g.core.G3TimestreamQuat) => 15716 quaternions at 152.6 Hz “RawTimestreams_I” (spt3g.core.G3TimestreamMap) => Timestreams from 14112 detectors “TodPointing” (spt3g.core.G3MapVectorInt) => 14112 elements … ]

See Also:
MapTODMasker :

Projection of sky masks to detector timestreams.

MapBinner, SingleDetectorMapBinner, SingleDetectorBoresightBinner :

Map binners for making coadded and single-detector maps.

MapMockObserver :

The inverse of MapBinner. Produces timestreams from an input map.

FlatSkyMap, HealpixSkyMap :

Possible map types.

Examples:
pipe.Add(
    maps.MapTODPointing,
    pointing="OnlineRaDecRotation",
    timestreams="RawTimestreams_I",
    stub_map=map_params,
    tod_pointing="TodPointing",
)

spt3g.maps.SingleDetectorBoresightBinner

SingleDetectorBoresightBinner(stub_map, pointing, timestreams)

Makes simple binned maps of the sky, in boresight coordinates, for every detector present in the given timestreams. This module is intended for use when making maps for detector-pointing calibration.

Parameters:
stub_mapG3SkyMap
Template of the map in which to accumulate timestream data. All

parameters of the output map (projection, boundaries, pixel size, etc.) are copied from this map, which is not modified.

pointingstring (frame object G3VectorQuat or G3TimestreamQuat)

Name of a frame object containing the boresight pointing quaternion timestream. Must have the same number of elements as the data in timestreams and be present in every Scan frame.

timestreamsstring (frame object G3TimestreamMap)

Name of a frame object containing the timestreams to be binned into the output map. Must exist in every Scan frame, though may be empty if the frame should be ignored. Units in the output map are taken from the units of the detector timestreams. Because of the single stored weight/hit map, the set of detectors in every scan must be identical.

Emits:

At the end of processing, emits one frame of type Map for each detector containing the corresponding map. The “Id” property of the output frames gives the detector ID to which the map corresponds. Because the maps are unweighted (weights don’t make sense for single-detector maps) and the effective pointing for all detectors is identical, the weight/hit map is identical for all detectors and is stored only in the first output frame.

Frame (Map) [ “Id” (spt3g.core.G3String) => “2019.000” “T” (spt3g.maps.FlatSkyMap) => 360 x 360 (3 x 3 deg) ZEA centered at (179.5, 179.5) = (0, 0 deg) in local coordinates (Power) “Wunpol” (spt3g.maps.G3SkyMapWeights) => G3SkyMapWeights ] Frame (Map) [ “Id” (spt3g.core.G3String) => “2019.005” “T” (spt3g.maps.FlatSkyMap) => 360 x 360 (3 x 3 deg) ZEA centered at (179.5, 179.5) = (0, 0 deg) in local coordinates (Power) ]

See Also:
SingleDetectorMapBinner :

Produces single-detector maps, but in sky coordinates, taking into account individual-detector pointing offsets.

MapBinner :

Produces co-added maps from many detectors.

MapMockObserver :

The inverse of MapBinner. Produces timestreams from an input map.

FlatSkyMap, HealpixSkyMap :

Possible output types.

Examples:
pipe.Add(
    maps.SingleDetectorBoresightBinner,
    stub_map=smstub,
    timestreams='PolyFilteredTimestreams',
    pointing='OffsetRotation',
)

spt3g.maps.SingleDetectorMapBinner

Makes a simple binned map of the sky, in sky coordinates, for every detector present in the given <timestreams>. Boresight pointing is specified by the <pointing> argument, with per-detector offsets from boresight as stored in the given bolometer properties map. The map parameters are copied from <stub_map>. When processing ends, this module will emit one map frame per detector, including per-detector (unpolarized) weights.

spt3g.maps.azel.EquatorialToGalacticPointing

Converts a set of timestreams in Scan frames representing RA and Declination pointing of the telescope into Galactic longitude and latitude timestreams, stored in the frame under their respective names.

Definition:

EquatorialToGalacticPointing(frame, ra_timestream='BoresightRa', dec_timestream='BoresightDec', glon_timestream='BoresightGalLon', glat_timestream='BoresightGalLat')

spt3g.maps.azel.LocalToAstronomicalPointing

Converts a set of timestreams in Scan frames representing Az and El pointing of the telescope into RA and Declination timestreams, stored in the frame under their respective names.

Definition:

LocalToAstronomicalPointing(frame, az_timestream='BoresightAz', el_timestream='BoresightEl', ra_timestream='BoresightRa', dec_timestream='BoresightDec', Telescope=<EarthLocation (710.21505704, -701.59071905, -6359587.23641261) m>)

spt3g.maps.coordsysmodules.AddLocalTransRotations

Creates the transform for boresight pointing for az el based maps. This is equivalent to FillCoordTransRotations with end_coord_sys in Local coordinates.

Right now it’s using the coordinate system where delta = -el because of implementation details. All of the maps generated with this will be upside down in el.

Definition:

AddLocalTransRotations(frame, az_key='RawBoresightAz', el_key='RawBoresightEl', out_key='RawAzElRotation')

spt3g.maps.coordsysmodules.EquatorialToGalacticTransRotations

Takes a quaternion vector specifying the rotation to FK5 (Equatorial) boresight and converts it into a rotation to Galactic J2000 boresight coordinates.

Use this to convert the output of FillCoordTransRotations to Galactic coordinates.

Definition:

EquatorialToGalacticTransRotations(frame, eq_trans_key='OnlineRaDecRotation', out_key='OnlineGalacticRotation')

spt3g.maps.coordsysmodules.FillCoordTransRotations

Calculates the rotation quaternions that take the point (1,0,0) (so az=el=0) in local coordinates to the coordinates specified by end_coord_sys and stores them in transform_store_key. This encodes the boresight pointing and any rotations about this boresight pointing due to coordinate system changes, az/el bearing tilt, etc.

Arguments:
transform_store_keystring

The key where the output transformation quaternion will be stored. If already present in the frame, this calculation will be skipped.

end_coord_sysMapCoordReference

If Local, the transformation is computed using the negative of the detector delta angle. Otherwise the detector angle is not inverted.

do_bad_transformbool

If end_coord_sys is not Local and this argument is True, the offset keys are ignored and the coordinate transformation does not take into account rotation about the boresight.

bs_az_key, bs_el_keystring

Boresight coordinates in the local coordinate system. If end_coord_sys is Local, only these two keys are required.

bs_ra_key, bs_dec_keystring

Boresight coordinates in the output coordinate system. If do_bad_transform is True, only these two keys and the previous two keys are required.

offset_az_key, offset_el_key, offset_ra_key, offset_dec_keystring

Local and output coordinates computed at a small offset from boresight. These keys are required if do_bad_transform is False, and will be used to account for any rotation about boresight in the coordinate transformation.

Definition:

FillCoordTransRotations(frame, transform_store_key='OnlineRaDecRotation', end_coord_sys=spt3g.maps.MapCoordReference.Equatorial, do_bad_transform=False, bs_az_key='RawBoresightAz', bs_el_key='RawBoresightEl', bs_ra_key='OnlineBoresightRa', bs_dec_key='OnlineBoresightDec', offset_az_key='OffsetBoresightAz', offset_el_key='OffsetBoresightEl', offset_ra_key='OnlineOffsetRa', offset_dec_key='OnlineOffsetDec')

spt3g.maps.fitsio.SaveMapFrame

Save a map frame to a FITS file. See save_skymap_fits for details. The map frame should contain T maps and (optionally) unpolarized weights, or T/Q/U maps and (optionally) polarized weights to store in the output file.

Arguments:
output_filestr or callable

Fits filename to which the map will be written. Maybe a callable function that takes a frame object as its sole input argument and returns a string filename.

hdrdict

If defined, extra keywords to be appened to the FITS header. The dict can contain entries such as hdr['NEWKEY'] = 'New value' or hdr['NEWKEY'] = ('New value', "Comment for New value").

compressstr or bool

If defined, and if input maps are FlatSkyMap objects, store these in a series of compressed image HDUs, one per map. Otherwise, store input maps in a series of standard ImageHDUs, which are readable with older FITS readers (e.g. idlastro). If defined, the compression algorithm to be used by the Astropy class astropy.io.fits.CompImageHDU. Can be: ‘RICE_1’, ‘RICE_ONE’, ‘PLIO_1’, ‘GZIP_1’, ‘GZIP_2’ or ‘HCOMPRESS_1’. Only GZIP_1 and GZIP_2 are lossless, although only for integer data.

quantize_levelfloat

Floating point quantization level for compression. Higher values result in more accurate floating point representation, but worse compression ratio. See the astropy FITS image documention for details: https://docs.astropy.org/en/stable/io/fits/api/images.html

overwritebool

If True, any existing file with the same name will be ovewritten.

Definition:

SaveMapFrame(frame, output_file=None, hdr=None, compress=False, quantize_level=16.0, overwrite=False)

spt3g.maps.map_modules.ApplyWeights

Apply weights to the input maps. The operation is performed in place to minimize memory use.

Definition:

ApplyWeights(frame)

spt3g.maps.map_modules.CoaddMaps

Coadd maps and weights, optionally collating by map Id. This class can be used as an argument to pipe.Add() as a standard pipeline module, or instantiated as a standalone instance. In the latter case, the object is treated as a callable for input map frames, and the coadd_frame or coadd_frames attribute contains the running coadd(s).

The output coadd frames contain two additional keys: 'InputMapIds' and 'InputFiles', which are both lists of unique map Id’s and filenames that are associated with the frames that contribute to each coadd. When one coadd is added to another, these keys are updated recursively, so that the resulting coadd includes the Id’s and filenames that contributed to both itself and the other coadd. The list of filenames can be populated by combining this module with a G3Reader whose track_filename option is set to True; however, this feature is fragile and may not work as expected with complex pipelines.

Attributes:
coadd_frameG3Frame

Output coadd map frame, also injected into the pipeline on EndProcessing. This attribute is only populated if the collate option is set to False.

coadd_framesdict of G3Frames

Output coadd map frames, keyed by input map Id. Each frame is also injected into the pipeline on EndProcessing. This attribute is only populated if the collate option is set to True.

Methods:
get_map_id :

Takes a map frame as an argument and returns an identifier string for the coadd to which it should be added, or None if the map should be ignored. This method can be modified by subclassing the CoaddMaps module.

Arguments:
map_idslist of str

List of map Id’s to include in the coadd(s). If None, any maps in the pipeline are included. Otherwise, the output of the get_map_ids method is compared with this list, and the input frame is discarded if no match is found.

output_map_idstr

Id to assign to the output frame. If collate is True, this argument is required and treated as a prefix to which each input map Id is appended.

collatebool

If True, coadd unique map Id’s into separate output map frames.

weightedbool

If True (default), ensure that maps have had weights applied before coadding. Otherwise, coadd maps without checking the weights.

ignore_missing_weightsbool

If False (default), a warning is issued when the frame contains weighted Stokes maps without a weights map. Set this option to True when feeding single bolometer map frames with common weights through a pipeline.

drop_input_framesbool

If True, drop input map frames from the pipeline that are included in any coadds.

record_obs_idbool

If True, include source name and observation ID info in the output coadd frame InputMapIds key, along with the map ID for each input frame. If False, only the map frame ID is included.

Constructor:

CoaddMaps(self, map_ids=None, output_map_id='Coadd', collate=False, weighted=True, ignore_missing_weights=False, drop_input_frames=False, record_obs_id=False)

spt3g.maps.map_modules.CompactMaps

Compact all maps in a frame to their default sparse representation. Optionally remove NaN values as well. Removing NaN values will reduce memory use, but will remove the distinction in unweighted (or weight-removed) maps between unobserved regions and regions with zero temperature.

Definition:

CompactMaps(frame, zero_nans=False)

spt3g.maps.map_modules.ExtractMaps

Cache maps that come through the pipeline. Initialize an instance of this module before adding to a pipeline.. Any maps that pass through the pipe are stored in the .maps attribute of the object after the pipeline is run.

Arguments:
map_idstring

If supplied, select only map frames that match this ID.

copybool

If True, make a copy of the map on extraction.

ignore_missing_weightsbool

If False (default), a warning is issued when the frame contains weighted Stokes maps without a weights map. Set this option to True when feeding single bolometer map frames with common weights through a pipeline.

Constructor:

ExtractMaps(self, map_id=None, copy=False, ignore_missing_weights=False)

spt3g.maps.map_modules.FlattenPol

For maps defined on the sphere the direction of the polarization angle is is defined relative to the direction of North. When making maps we follow this definition.

For any flat sky estimators, the polarization angle is defined relative to the vertical axis. For some map projections the direction of north is not the same as the vertical axis. This function applies a rotation to the Q and U values to switch the curved sky Q/U definition to the flat sky Q/U definition.

If for whatever reason you want to reverse the process set the invert argument to True.

Definition:

FlattenPol(frame, invert=False)

spt3g.maps.map_modules.InjectMapStub

Inject a new map frame from a map stub.

Arguments:
map_idstring

Id to assign to the new map frame

map_stubG3SkyMap instance

Map stub from which to clone the Stokes maps and weights. If the weighted attribute of the stub is True, the output frame will include weights. If the pol_conv attribute of the stub is not None, the output frame will include Q and U maps (and polarized weights).

Constructor:

InjectMapStub(self, map_id, map_stub)

spt3g.maps.map_modules.InjectMaps

Inject a set of maps into a new map frame.

Arguments:
map_idstring

Id to assign to the new map frame

maps_inlist or dict

Maps to add to the frame. If a list, contains Stokes maps with valid pol_type and weights. If a dict, contains Stokes and weights maps keyed by the standard map frame names.

ignore_missing_weights [False]bool

Skip warning about missing weights. Useful for masks.

Constructor:

InjectMaps(self, map_id, maps_in, ignore_missing_weights=False)

spt3g.maps.map_modules.MakeMapsPolarized

Converts individual unpolarized maps to polarized versions of the same map, with the given polarization convention

This module is only a shim that creates null Q and U maps and populates a properly invertible Wpol array from the TT Wunpol weights.

Definition:

MakeMapsPolarized(frame, pol_conv=spt3g.maps.MapPolConv.IAU)

spt3g.maps.map_modules.MakeMapsUnpolarized

Converts individual polarized maps to temperature-only versions of the same map.

Definition:

MakeMapsUnpolarized(frame)

spt3g.maps.map_modules.RemoveWeights

Remove weights from input maps. If zero_nans is True, empty pixels are skipped and pixels with zero weight are set to 0 instead of NaN. Operation is performed in place to minimize memory use.

Definition:

RemoveWeights(frame, zero_nans=False)

spt3g.maps.map_modules.ReplicateMaps

Clone the input map frame with Id input_map_id into new stub frames, one for each Id listed in output_map_ids.

Arguments:
input_map_idstring

ID of the map frame to replicate. The input frame is discarded after replication.

output_map_idslist of strings

List of IDs to assign to replicated map frames.

copy_weightsbool

If False, only the first output frame in the list includes a weights key (Wpol or Wunpol). If True, all output frames include a weights key.

Definition:

ReplicateMaps(frame, input_map_id, output_map_ids, copy_weights=False)

spt3g.maps.map_modules.ReprojectMaps

Reproject a map frame into a different projection. Original data are dropped and replaced by reprojected maps in the input frames. Maps can be changed between flat sky and healpix pixelizations, rotated between Equatorial and Galactic coordinates, and/or change polarization convention between COSMO and IAU, by setting the appropriate attributes of the input and stub maps. Attributes not defined in the stub map are assumed to be that of the input map. NB: coordinate rotation of polarized maps is not currently implemented.

Arguments:
map_stubG3SkyMap object

A stub (empty) sky map object to be used to construct the output maps. Can be a HealpixSkyMap or FlatSkyMap object. Setting the pol_conv and/or coord_ref attributes to values that differ from those of the input maps will result in output maps whose polarization convention and/or reference coordinate system have been changed.

rebinint

If supplied and >1, subdivide the output pixel by n x n with each sub-pixel taking on the input map values at pixel center (with interp or nearest neighbor). The output pixel takes on the average of the sub-pixel values. In the case that the input map has higher resolution than the output map (and that the input map is not low-pass filtered to remove information above the Nyquist freq. of the output map pixel), this reduces aliasing compared with direct sampling. But there would still be aliased power from the input map from freq above the ouput map pixel’s Nyquist.

interpbool

If True, use bilinear interpolation to extract values from the input map. Otherwise, the nearest-neighbor value is used.

Constructor:

ReprojectMaps(self, map_stub=None, rebin=1, interp=False)

spt3g.maps.map_modules.SetPolConv

Set or change the polarization convention of the input polarized map frame. If switching between IAU and COSMO conventions, flip the sign of the U map and the TU and QU weights. Otherwise, just set the polarization convention for all maps and weights.

Definition:

SetPolConv(frame, pol_conv=spt3g.maps.MapPolConv.IAU)

spt3g.maps.map_modules.ValidateMaps

Validate that the input map frame has all the necessary keys.

If ignore_missing_weights is False (default), a warning is issued when the frame contains weighted Stokes maps without a weights map. Set this option to True when feeding single bolometer map frames with common weights through a pipeline.

Definition:

ValidateMaps(frame, ignore_missing_weights=False)

spt3g.maps.quathelpers.AddTimingToPointingQuats

Transforms a quaternion pointing timestream <key> from a G3VectorQuat to a G3TimestreamQuat by adding timing information extracted from some other timestream-like object (scalar pointing timestreams, detector timestreams etc.) specified by <timing_ref>. Because this operation is backwards-compatible and involves no data loss, the transformation is done in-place – the previous data are deleted and replaced with the new one.

Definition:

AddTimingToPointingQuats(fr, key, timing_ref='RawBoresightAz')

Functions in spt3g.maps

spt3g.maps.apply_weights

apply_weights( (G3SkyMap)T, (G3SkyMap)Q, (G3SkyMap)U, (G3SkyMapWeights)W) -> None :

Apply weights to polarized maps.

spt3g.maps.apply_weights_t

apply_weights_t( (G3SkyMap)T, (G3SkyMapWeights)W) -> None :

Apply weights to unpolarized maps.

spt3g.maps.azel.convert_azel_to_radec

Convert timestreams of local azimuth and elevation to right ascension and declination.

Arguments:
az, elnp.ndarray or G3Timestream

Array of local coordinates. If inputs are G3Timestream objects, G3Timestreams are also returned.

locationastropy.coordinates.EarthLocation instance

The telescope location on Earth.

mjdnp.ndarray

An array of times for each az/el sample. If input az and el are not G3Timestreams, this argument is required.

Returns:

ra, dec : np.ndarray or G3Timestream

Definition:

convert_azel_to_radec(az, el, location=<EarthLocation (710.21505704, -701.59071905, -6359587.23641261) m>, mjd=None)

spt3g.maps.azel.convert_gal_to_radec

Convert timestreams of Galactic longitude and latitude to right ascension and declination.

Arguments:
glon, glatnp.ndarray or G3Timestream

Array of Galactic sky coordinates. If inputs are G3Timestream objects, G3Timestreams are also returned.

Returns:

ra, dec : np.ndarray or G3Timestream

Definition:

convert_gal_to_radec(glon, glat)

spt3g.maps.azel.convert_radec_to_azel

Convert timestreams of right ascension and declination to local azimuth and elevation.

Arguments:
ra, decnp.ndarray or G3Timestream

Array of Equatorial sky coordinates. If inputs are G3Timestream objects, G3Timestreams are also returned.

locationastropy.coordinates.EarthLocation instance

The telescope location on Earth.

mjdnp.ndarray

An array of times for each ra/dec sample. If input ra and dec are not G3Timestreams, this argument is required.

Returns:

az, el : np.ndarray or G3Timestream

Definition:

convert_radec_to_azel(ra, dec, location=<EarthLocation (710.21505704, -701.59071905, -6359587.23641261) m>, mjd=None)

spt3g.maps.azel.convert_radec_to_gal

Convert timestreams of right ascension and declination to Galactic longitude and latitude.

Arguments:
ra, decnp.ndarray or G3Timestream

Array of Equatorial sky coordinates. If inputs are G3Timestream objects, G3Timestreams are also returned.

Returns:

glon, glat : np.ndarray or G3Timestream

Definition:

convert_radec_to_gal(ra, dec)

spt3g.maps.convolve_map

convolve_map( (FlatSkyMap)map, (object)kernel) -> FlatSkyMap :

Convolve the input flat sky map with the given map-space kernel. The kernel must have odd dimensions and the same resolution as the map.

spt3g.maps.fitsio.load_skymap_fits

Load a fits file containing a sky map.

Arguments:
filenamestr

Path to fits file

hduint, optional

If supplied, the data are extract from the given HDU index.

keyslist of strings, optional

If supplied, return only these keys in the output dictionary. Options are: T, Q, U, W.

memmapbool, optional

Argument passed to astropy.io.fits.open. If True, the map is not read into memory, but only the required pixels are read when needed. Default: False.

apply_unitsbool, optional

If True, and input maps have known units, multiply by the appropriate conversion factor to return maps in G3Units.

Returns:

a dictionary of maps keyed with e.g. ‘T’, ‘Q’, ‘U’ and ‘W’.

Definition:

load_skymap_fits(filename, hdu=None, keys=None, memmap=False, apply_units=False)

spt3g.maps.fitsio.save_skymap_fits

Save G3 map objects to a fits file.

FlatSkyMap objects are stored in a series of (optionally compressed) ImageHDU entries, in which each HDU contains the projection information in its header in standard WCS format, along with the image data for a single map (one of the Stokes maps or a weight map component).

HealpixSkyMap objects are stored in a BinTableHDU extension, which contains the necessary header information for compatiblity with healpix map readers (e.g. healpix.read_map), and a single table with one column per Stokes map or weight map component. Sparse maps are stored as cut-sky pixel-indexed tables, while dense maps are stored with implicit indexing over all pixels. The former produces output that is equivalent to using healpy.write_map with the partial=True option.

Arguments:
filenamestr

Path to output file. Must not exist, unless overwrite is True.

T[, Q, U]FlatSkyMap or HealpixSkyMap

Maps to save

WG3SkyMapWeights

Weights to save with the maps

overwritebool

If True, any existing file with the same name will be ovewritten.

compressstr or bool

If defined, and if input maps are FlatSkyMap objects, store these in a series of compressed image HDUs, one per map. Otherwise, store input maps in a series of standard ImageHDUs, which are readable with older FITS readers (e.g. idlastro). If defined, the compression algorithm to be used by the Astropy class astropy.io.fits.CompImageHDU. Can be: ‘RICE_1’, ‘RICE_ONE’, ‘PLIO_1’, ‘GZIP_1’, ‘GZIP_2’ or ‘HCOMPRESS_1’. Only GZIP_1 and GZIP_2 are lossless, although only for integer data.

quantize_levelfloat

Floating point quantization level for compression. Higher values result in more accurate floating point representation, but worse compression ratio. See the astropy FITS image documention for details: https://docs.astropy.org/en/stable/io/fits/api/images.html

hdrdict

If defined, extra keywords to be appened to the FITS header. The dict can contain entries such as hdr['NEWKEY'] = 'New value' or hdr['NEWKEY'] = ('New value', "Comment for New value").

Definition:

save_skymap_fits(filename, T, Q=None, U=None, W=None, overwrite=False, compress=False, quantize_level=16.0, hdr=None)

spt3g.maps.flatten_pol

flatten_pol( (FlatSkyMap)Q, (FlatSkyMap)U [, (G3SkyMapWeights)W=None [, (float)h=0.001 [, (bool)invert=False]]]) -> None :

For maps defined on the sphere the direction of the polarization angle is is defined relative to the direction of North. When making maps we follow this definition.

For any flat sky estimators, the polarization angle is defined relative to the vertical axis. For some map projections the direction of north is not the same as the vertical axis. This function applies a rotation to the Q and U values to switch the curved sky Q/U definition to the flat sky Q/U definition.

If for whatever reason you want to reverse the process set the invert argument to True. Also applies the appropriate rotation to the Q and u elements of the associated weights.

spt3g.maps.get_boresight_rotator_timestream

get_boresight_rotator_timestream( (G3Timestream)arg1, (G3Timestream)arg2, (G3Timestream)arg3, (G3Timestream)arg4, (G3Timestream)arg5, (G3Timestream)arg6, (G3Timestream)arg7, (G3Timestream)arg8) -> G3TimestreamQuat :

Construct a transform quaternion timestream from timestreams of local and equatorial boresight pointing coordinates. Computes the transform from local (az_0, el_0) coordinates to equatorial (ra_0, dec_0), accounting for rotation about the boresight by including the second set of points.

spt3g.maps.get_fk5_j2000_to_gal_quat

get_fk5_j2000_to_gal_quat() -> quat :

Return the rotation quaternion to rotate from equatorial to galactic coordinates.

spt3g.maps.get_map_hist

get_map_hist( (object)map, (DoubleVector)bin_edges [, (G3SkyMapMask)mask=None [, (bool)ignore_zeros=False [, (bool)ignore_nans=False [, (bool)ignore_infs=False]]]]) -> DoubleVector :

Computes the histogram of the input map into bins defined by the array of bin edges, optionally ignoring zero, nan and/or inf values in the map. If a mask is supplied, then only the non-zero pixels in the mask are included.

spt3g.maps.get_map_moments

get_map_moments( (object)map [, (G3SkyMapMask)mask=None [, (int)order=2 [, (bool)ignore_zeros=False [, (bool)ignore_nans=False [, (bool)ignore_infs=False]]]]]) -> DoubleVector :

Computes moment statistics of the input map, optionally ignoring zero, nan and/or inf values in the map. If order = 1, only the mean is returned. If order = 2, 3 or 4 then the variance, skew and kurtosis are also included, respectively. If a mask is supplied, then only the non-zero pixels in the mask are included.

spt3g.maps.get_origin_rotator

get_origin_rotator( (float)alpha, (float)delta) -> quat :

Compute the transformation quaternion that would rotate the vector (1, 0, 0) to point in the given direction.

spt3g.maps.get_origin_rotator_timestream

get_origin_rotator_timestream( (G3Timestream)alpha, (G3Timestream)delta, (MapCoordReference)coord_sys) -> G3TimestreamQuat :

Construct a transform quaternion timestream from timestreams of sky coordinates. Equivalent to R_z(alpha) * R_y(delta).

spt3g.maps.get_ra_dec_map

get_ra_dec_map( (object)map_in) -> tuple :

Returns maps of the ra and dec angles for each pixel in the input map

spt3g.maps.get_ra_dec_mask

get_ra_dec_mask( (object)map_in, (float)ra_left, (float)ra_right, (float)dec_bottom, (float)dec_top) -> G3SkyMapMask :

Returns a mask that is nonzero for any pixels within the given ra and dec ranges

spt3g.maps.get_rot_ang

get_rot_ang( (quat)start_q, (quat)trans) -> float :

Computes the boresight rotation of the vector start_q when rotated by trans.

spt3g.maps.get_transform_quat

get_transform_quat( (float)arg1, (float)arg2, (float)arg3, (float)arg4, (float)arg5, (float)arg6, (float)arg7, (float)arg8) -> quat :

Computes a rotation that will take (as_0,ds_0) to (ae_0, de_0) and (as_1, ds_1) to (ae_1, de_1)

spt3g.maps.make_point_source_mask

make_point_source_mask( (object)map, (DoubleVector)ra, (DoubleVector)dec, (DoubleVector)radius) -> G3SkyMapMask :

Construct a mask from the input stub map with pixels within the given radius around each point source position set to 1.

spt3g.maps.map_modules.coadd_map_files

Coadd map files, optionally collating map Id’s into separate frames.

Arguments:
input_fileslist of str

List of input files to feed through the pipeline.

output_filestr

Output G3 filename. If not supplied, the output frames are returned without saving to disk.

coadderCoaddMaps instance

If set, use this instantiated module in the coadding pipeline. In this case, all other keyword arguments below are ignored.

map_idslist of str

A list of map Id’s to include in the coadd(s).

output_map_idstr

Id to use for the output map frame. If collate is True, this is the prefix applied to each output frame, with the input map Id appended to it.

collatebool

If True, coadd individual map Id’s into separate map frames. Otherwise, all map frames are coadded into one output frame.

weightedbool

If True, ensure that weights have been applied before coadding. Otherwise, the input maps are coadded as they are.

record_obs_idbool

If True, include source name and observation ID info in the output coadd frame InputMapIds key, along with the map ID for each input frame. If False, only the map frame ID is included.

Returns:
mapsG3Frame or dict of G3Frames

If collate is True, returns a dictionary of map frames keyed by Id. Otherwise, returns a single map frame.

Definition:

coadd_map_files(input_files, output_file=None, coadder=None, map_ids=None, output_map_id='Coadd', collate=False, weighted=True, record_obs_id=False)

spt3g.maps.maputils.flatsky_to_healpix

Re-pixelize a map to Healpix from one of the flat projections.

Parameters:
map_in: FlatSkyMap

The input map you want to reproject

map_stub[None]: HealpixSkyMap

Stub output map object to be used to construct the output map. If not supplied, one will be constructed using the remaining keyword arguments.

rebin[1]: int

If supplied and >1, subdivide the output pixel by n x n with each sub-pixel taking on the input map values at pixel center (with interp or nearest neighbor). The output pixel takes on the average of the sub-pixel values. In the case that the input map has higher resolution than the output map (and that the input map is not low-pass filtered to remove information above the Nyquist freq. of the output map pixel), this reduces aliasing compared with direct sampling. But there would still be aliased power from the input map from freq above the ouput map pixel’s Nyquist.

interp[false]: bool

If True, use bilinear interpolation to extract values from the input map. Otherwise, the nearest-neighbor value is used.

fullsky[false]: bool

If True a full-sky numpy array representation of the map is returned. Otherwise, a HealpixSkyMap instance is returned, containing only the pixels that overlap with the input map.

All additional keyword arguments are passed to HealpixSkyMap to construct the output map object. Required if map_stub is not supplied, otherwise ignored.

Maps can be rotated between Equatorial and Galactic coordinates, and/or change polarization convention between COSMO and IAU, by setting the appropriate attributes of the input and output maps. Attributes not defined in the output map are assumed to be that of the input map.

Returns:
output_map: numpy array or HealpixSkyMap

The array containing the healpix map. If fullsky is True, this is a numpy array, otherwise a HealpixSkyMap instance.

Definition:

flatsky_to_healpix(map_in, map_stub=None, rebin=1, interp=False, fullsky=False, **kwargs)

spt3g.maps.maputils.healpix_to_flatsky

Re-pixelize a map from Healpix to one of the flat sky projections.

Parameters:
map_in: numpy array or HealpixSkyMap

The array containing the input healpix map to reproject.

nest[False]: bool

Ordering of the healpix map, if the input is a numpy array. Ring ordering is assumed by default.

map_stub[None]: FlatSkyMap

Stub output map object to be used to construct the output map. If not supplied, one will be constructed using the remaining keyword arguments.

rebin[1]: int

If supplied and >1, subdivide the output pixel by n x n with each sub-pixel taking on the input map values at pixel center (with interp or nearest neighbor). The output pixel takes on the average of the sub-pixel values. In the case that the input map has higher resolution than the output map (and that the input map is not low-pass filtered to remove information above the Nyquist freq. of the output map pixel), this reduces aliasing compared with direct sampling. But there would still be aliased power from the input map from freq above the ouput map pixel’s Nyquist.

interp[false]: bool

If True, use bilinear interpolation to extract values from the 4 closest pixel centers of the healpix map. Otherwise, the nearest-neighbor value is used.

All additional keyword arguments are passed to FlatSkyMap to construct the output map object. Required if map_stub is not supplied, otherwise ignored.

Maps can be rotated between Equatorial and Galactic coordinates, and/or change polarization convention between COSMO and IAU, by setting the appropriate attributes of the input and output maps. Note that if the input map is a numpy array representation of a healpix map, the coordinate system and polarization convention are assumed to be that of the output map. Conversely, attributes not defined in the output map are assumed to be that of the input map.

Returns:
output_map: FlatSkyMap

The reprojected map

Definition:

healpix_to_flatsky(map_in, nest=False, map_stub=None, rebin=1, interp=False, **kwargs)

spt3g.maps.offsets_to_quat

offsets_to_quat( (float)x, (float)y) -> quat :

Returns the vector quaternion (0,1,0,0) rotated by the given x and y offsets. Equivalent to t * quat(0,1,0,0) / t, where t = get_origin_rotator(x, -y)

spt3g.maps.quathelpers.ang_to_quat

Convert a set of angles (or vector of them) specified as a (longitude, latitude) pair to a pointing quaternion (or vector of them). If start and stop are defined, the return value for vectors will be a G3TimestreamQuat with start and stop times set to the provided values.

Definition:

ang_to_quat(alpha, delta, start=None, stop=None)

spt3g.maps.quathelpers.quat_to_ang

Convert a pointing quaternion (or vector of them) to a set of angles (or vector of them) specified as a (longitude, latitude) pair.

Definition:

quat_to_ang(q)

spt3g.maps.remove_weights

remove_weights( (G3SkyMap)T, (G3SkyMap)Q, (G3SkyMap)U, (G3SkyMapWeights)W [, (bool)zero_nans=False]) -> None :

Remove weights from polarized maps. If zero_nans is true, empty pixels are skipped, and pixels with zero weight are set to 0 instead of nan.

spt3g.maps.remove_weights_t

remove_weights_t( (G3SkyMap)T, (G3SkyMapWeights)W [, (bool)zero_nans=False]) -> None :

Remove weights from unpolarized maps. If zero_nans is true, empty pixels are skipped, and pixels with zero weight are set to 0 instead of nan.

spt3g.maps.reproj_map

reproj_map( (object)in_map, (G3SkyMap)out_map [, (int)rebin=1 [, (bool)interp=False]]) -> None :

Reprojects the data from in_map onto out_map. out_map can have a different projection, size, resolution, etc. Optionally account for sub-pixel structure by setting rebin > 1 and/or enable bilinear interpolation of values from the input map by setting interp=True. Use the maps’ coord_ref attributes to rotate between Equatorial and Galactic coordinate systems. Use the maps’ pol_conv attributes to switch between COSMO and IAU polarization conventions. If output attributes are not set, they will be copied from the input map.