API Reference

QPoint Class

The base class for computing boresight and offset astronomical pointing.

class qpoint.qpoint_class.QPoint(update_iers=False, **kwargs)[source]

Bases: object

__init__(update_iers=False, **kwargs)[source]

Initialize a QPoint memory instance for keeping track of pointing corrections over time.

If update_iers is True, attempt to call the method update_bulletin_a() to update the internal IERS-A database (requires astropy version 1.2 or newer).

Any remaining keyword arguments are passed to set() to update memory.

azel2bore(az, el, pitch, roll, lon, lat, ctime, q=None, **kwargs)[source]

Estimate the quaternion for the boresight orientation on the sky given the attitude (az/el/pitch/roll), location on the earth (lon/lat) and ctime. Input vectors must be numpy-array-like and broadcastable to the same shape.

Parameters:
  • az (array_like) – Boresight azimuth in degrees

  • el (array_like) – Boresight elevation in degrees

  • pitch (array_like) – Boresight pitch in degrees. If None, this term is ignored.

  • roll (array_like) – Boresight roll in degrees. If None, this term is ignored.

  • lon (array_like) – Observer longitude in degrees

  • lat (array_like) – Observer latitude in degrees

  • ctime (array_like) – Unix time in seconds UTC

  • q (array_like, optional) – Output quaternion array initialized by user. Supply this for in-place computation.

Returns:

q (array_like) – Nx4 numpy array of quaternions for each supplied timestamp.

Notes

Any keywords accepted by the set() method can also be passed here, and will be processed prior to calculation.

azel2radec(delta_az, delta_el, delta_psi, az, el, pitch, roll, lon, lat, ctime, hwp=None, sindec=False, return_pa=False, ra=None, dec=None, pa=None, sin2psi=None, cos2psi=None, **kwargs)[source]

Estimate the orientation on the sky for a detector offset from boresight, given the boresight attitude (az/el/pitch/roll), location on the earth (lon/lat) and UTC time. Input vectors must be numpy-array-like and broadcastable to the same shape. Detector offsets are defined assuming the boresight is pointed toward the horizon, and that the boresight polarization axis is along the horizontal.

Parameters:
  • delta_az (float) – Azimuthal offset of the detector in degrees

  • delta_el (float) – Elevation offset of the detector in degrees

  • delta_psi (float) – Polarization offset of the detector in degrees

  • az (array_like) – Boresight azimuth in degrees

  • el (array_like) – Boresight elevation in degrees

  • pitch (array_like) – Boresight pitch in degrees. If None, this term is ignored.

  • roll (array_like) – Boresight roll in degrees. If None, this term is ignored.

  • lon (array_like) – Observer longitude in degrees.

  • lat (array_like) – Observer latitude in degrees.

  • ctime (array_like) – Unix time in seconds UTC

  • hwp (array_like, optional) – HWP angles in degrees

  • sindec (bool, optional) – If True, return sin(dec) instead of dec in degrees (default False)

  • return_pa (bool, optional) – If True, return pa instead of sin2psi/cos2psi

Returns:

  • ra (array_like) – Detector right ascension in degrees

  • dec/sindec (array_like) – Detector declination in degrees

  • pa (array_like) – Detector position angle, if return_pa is True

  • sin2psi (array_like) – Detector polarization orientation, if return_pa is False

  • cos2psi (array_like) – Detector polarization orientation, if return_pa is False

Notes

Any keywords accepted by the set() method can also be passed here, and will be processed prior to calculation.

azelpsi2bore(az, el, psi, pitch, roll, lon, lat, ctime, q=None, **kwargs)[source]

Estimate the quaternion for the boresight orientation on the sky given the attitude (az/el/psi/pitch/roll), location on the earth (lon/lat) and ctime. Input vectors must be numpy-array-like and broadcastable to the same shape.

This is an augmented version of azel2bore to accept rotations of the focal plane about the optical axis.

Parameters:
  • az (array_like) – Boresight azimuth in degrees

  • el (array_like) – Boresight elevation in degrees

  • psi (array_like) – Boresight rotation angle in degrees

  • pitch (array_like) – Boresight pitch in degrees. If None, this term is ignored.

  • roll (array_like) – Boresight roll in degrees. If None, this term is ignored.

  • lon (array_like) – Observer longitude in degrees

  • lat (array_like) – Observer latitude in degrees

  • ctime (array_like) – Unix time in seconds UTC

  • q (array_like, optional) – Output quaternion array initialized by user. Supply this for in-place computation.

Returns:

q (array_like) – Nx4 numpy array of quaternions for each supplied timestamp.

Notes

Any keywords accepted by the set() method can also be passed here, and will be processed prior to calculation.

azelpsi2radec(delta_az, delta_el, delta_psi, az, el, psi, pitch, roll, lon, lat, ctime, hwp=None, sindec=False, return_pa=False, ra=None, dec=None, pa=None, sin2psi=None, cos2psi=None, **kwargs)[source]

Estimate the orientation on the sky for a detector offset from boresight, given the boresight attitude (az/el/pitch/roll), location on the earth (lon/lat) and UTC time. Input vectors must be numpy-array-like and broadcastable to the same shape. Detector offsets are defined assuming the boresight is pointed toward the horizon, and that the boresight polarization axis is along the horizontal.

This is agumented from azel2radec(), to accept rotations of the focal plane about the optical axis.

Parameters:
  • delta_az (float) – Azimuthal offset of the detector in degrees

  • delta_el (float) – Elevation offset of the detector in degrees

  • delta_psi (float) – Polarization offset of the detector in degrees

  • az (array_like) – Boresight azimuth in degrees

  • el (array_like) – Boresight elevation in degrees

  • psi (array_like) – Boresight rotation in degrees

  • pitch (array_like) – Boresight pitch in degrees. If None, this term is ignored.

  • roll (array_like) – Boresight roll in degrees. If None, this term is ignored.

  • lon (array_like) – Observer longitude in degrees.

  • lat (array_like) – Observer latitude in degrees.

  • ctime (array_like) – Unix time in seconds UTC

  • hwp (array_like, optional) – HWP angles in degrees

  • sindec (bool, optional) – If True, return sin(dec) instead of dec in degrees (default False)

  • return_pa (bool, optional) – If True, return pa instead of sin2psi/cos2psi

Returns:

  • ra (array_like) – Detector right ascension in degrees

  • dec/sindec (array_like) – Detector declination in degrees

  • pa (array_like) – Detector position angle, if return_pa is True

  • sin2psi (array_like) – Detector polarization orientation, if return_pa is False

  • cos2psi (array_like) – Detector polarization orientation, if return_pa is False

Notes

Any keywords accepted by the set() method can also be passed here, and will be processed prior to calculation.

bore2dipole(q_off, ctime, q_bore, **kwargs)[source]

Calculate dipole timestream for given offset and boresight pointing.

Parameters:
  • q_off (quaternion) – Detector offset quaternion for a single detector, calculated using det_offset()

  • ctime (array_like) – Array of unix times in seconds UTC

  • q_bore (quaternion or array of quaternions) – Array of quaternions encoding the boresight orientation on the sky (as output by azel2radec() or similar). Broadcastable to the same length as ctime.

Returns:

dipole (array_like) – Dipole amplitude in K

Notes

Any keywords accepted by the set() method can also be passed here, and will be processed prior to calculation.

bore2pix(q_off, ctime, q_bore, q_hwp=None, nside=256, pol=True, return_pa=False, **kwargs)[source]

Calculate the orientation on the sky for a detector offset from the boresight. Detector offsets are defined assuming the boresight is pointed toward the horizon, and that the boresight polarization axis is along the vertical.

Parameters:
  • q_off (quaternion) – Detector offset quaternion for a single detector, calculated using det_offset().

  • ctime (array_like) – Unix times in seconds UTC, broadcastable to shape (N,), the long dimenions of q_bore.

  • q_bore (quaternion or array of quaternions) – Nx4 array of quaternions encoding the boresight orientation on the sky (as output by azel2radec() or equivalent)

  • q_hwp (quaternion or array of quaternions, optional) – HWP angle quaternions calculated using hwp_quat(). Must be broadcastable to the same shape as q_bore.

  • nside (int, optional) – HEALpix map dimension. Default: 256.

  • pol (bool, optional) – If False, return only the pixel timestream

  • return_pa (bool, optional) – If True, return pa instead of sin2psi / cos2psi

Returns:

  • pix (array_like) – Detector pixel number

  • pa/sin2psi (array_like) – Detector polarization orientation if return_pa is True, or sin(2*pa) if return_pa is False.

  • cos2psi (array_like) – detector polarization orientation cos(2*pa), if return_pa is False.

Notes

Any keywords accepted by the set() method can also be passed here, and will be processed prior to calculation.

bore2radec(q_off, ctime, q_bore, q_hwp=None, sindec=False, return_pa=False, ra=None, dec=None, pa=None, sin2psi=None, cos2psi=None, **kwargs)[source]

Calculate the orientation on the sky for a detector offset from the boresight. Detector offsets are defined assuming the boresight is pointed toward the horizon, and that the boresight polarization axis is along the vertical.

Parameters:
  • q_off (quaternion) – Detector offset quaternion for a single detector, calculated using det_offset().

  • ctime (array_like) – Unix time in seconds UTC, broadcastable to shape (N,), the long dimension of q_bore.

  • q_bore (quaternion or array of quaternions) – Nx4 array of quaternions encoding the boresight orientation on the sky (as output by azel2radec() or equivalent)

  • q_hwp (quaternion or array of quaternions, optional) – HWP angle quaternions calculated using hwp_quat(). Must be broadcastable to the same shape as q_bore.

  • sindec (bool, optional) – If True, return sin(dec) instead of dec in degrees (default False).

  • return_pa (bool, optional) – If True, return pa instead of sin2psi / cos2psi

Returns:

  • ra (array_like) – Detector right ascension in degrees

  • dec/sindec (array_like) – Detector declination in degrees or sin(dec) if sindec is True.

  • pa/sin2psi (array_like) – Detector polarization orientation if return_pa is True, or sin(2*pa) if return_pa is False.

  • cos2psi (array_like) – detector polarization orientation cos(2*pa), if return_pa is False.

Notes

Any keywords accepted by the set() method can also be passed here, and will be processed prior to calculation.

Pre-allocated output arguments can also be supplied as input keywords for in-place operation.

bore_offset(q_bore, ang1=None, ang2=None, ang3=None, post=False, inplace=False)[source]

Apply a fixed or variable offset to the boresight quaternion. Input arguments must be broadcastable to the same shape.

Parameters:
  • q_bore (array_like) – boresight pointing quaternion

  • ang1 (array_like, optional) – Azimuthal or ra offset in degrees

  • ang2 (array_like, optional) – Elevation or dec offset in degrees

  • ang3 (array_like, optional) – Position angle offset in degrees

  • post (bool, optional) – If False, apply offset as an az/el/pa pre-rotation If True, apply offset as an ra/dec/pa post-rotation

  • inplace (bool, optional) – If True, apply the rotation in-place in memory.

Returns:

q_bore (array_like) – Offset boresight quaternion

det_offset(delta_az, delta_el, delta_psi)[source]

Return quaternion corresponding to the requested detector centroid offset from boresight. Vectorized, input arguments must be broadcastable to the same shape.

Parameters:
  • delta_az (array_like) – Azimuthal centroid offset of the detector in degrees

  • delta_el (array_like) – Elevation centroid offset of the detector in degrees

  • delta_psi (array_like) – Polarization offset of the detector from vertical in degrees

Returns:

q (array_like) – Detector centroid offset quaternion for each detector

dipole(ctime, ra, dec, **kwargs)[source]

Return dipole amplitude in the given equatorial direction. Vectorized, input arguments must be broadcastable to the same shape.

Parameters:
  • ctime (array_like) – Unix time in seconds UTC

  • ra (array_like) – Right ascension on the sky, in degrees.

  • dec (array_like) – Declination on the sky, in degrees

Returns:

dipole (array_like) – Dipole amplitude in K

Notes

Any keywords accepted by the set() method can also be passed here, and will be processed prior to calculation.

gal2radec(ra, dec, pa=None, sin2psi=None, cos2psi=None, inplace=True, **kwargs)[source]

Rotate galactic coordinates to celestial (equatorial) coordinates. This is equivalent to calling rotate_coord() with coord=[‘G’, ‘C’]. Vectorized, input arguments must be broadcastable to the same shape.

Parameters:
  • l (array_like) – – or –

  • b (array_like) – – or –

  • pa (array_like) – – or –

  • l – arrays of coordinates, of shape (n,). If none of pa, sin2psi or cos2psi are supplied, pa = 0 is assumed.

  • b – arrays of coordinates, of shape (n,). If none of pa, sin2psi or cos2psi are supplied, pa = 0 is assumed.

  • sin2psi (array_like) – arrays of coordinates, of shape (n,). If none of pa, sin2psi or cos2psi are supplied, pa = 0 is assumed.

  • cos2psi (array_like) – arrays of coordinates, of shape (n,). If none of pa, sin2psi or cos2psi are supplied, pa = 0 is assumed.

  • inplace (bool, optional) – If True, apply the rotation in-place on the input quaternion. Otherwise, return a copy of the input array. Default: True.

Returns:

  • ra, dec, pa (array_like) – – or –

  • ra, dec, sin2psi, cos2psi (array_like) – rotated coordinate arrays, same form as input

Notes

Any keywords accepted by the set() method can also be passed here, and will be processed prior to calculation.

get(*args)[source]

Returns a dictionary of the requested state parameters. If no parameters are supplied, then all are returned. If a single parameter is supplied, then just that value is returned. See set() for a list of parameter names.

Can also select ‘options’, ‘rates’, ‘weather’, or ‘params’ to return all of that subset of parameters.

get_bulletin_a(mjd)[source]

Return dut1/x/y for given mjd. Numpy-vectorized.

get_interp_val(map_in, ra, dec, nest=False)[source]

Interpolate map pixels to these coordinates. Uses a C implementation of the bilinear interpolation method get_interpol() as implemented in the equivalent healpix_cxx / healpy function.

Parameters:
  • map_in (array_like) – A single healpix map or list of maps which to interpolate from.

  • ra (array_like) – Timestreams of coordinates to interpolate to, in degrees, broadcastable to a common shape (nsample,)

  • dec (array_like) – Timestreams of coordinates to interpolate to, in degrees, broadcastable to a common shape (nsample,)

  • nest (bool, optional) – If True, input map is in the nested pixel ordering scheme. Otherwise, ring ordering is assumed. Default: False.

Returns:

values (array_like) – Array of interpolated map values, of shape (nmap, nsample).

gmst(ctime, **kwargs)[source]

Return Greenwich mean sidereal time for given ctimes. Vectorized.

Parameters:

ctime (array_like) – Unix time in seconds UTC

Returns:

gmst (array_like) – Greenwich mean sidereal time of the observer

Notes

Any keywords accepted by the set() method can also be passed here, and will be processed prior to calculation.

hwp_quat(theta)[source]

Return quaternion corresponding to rotation by 2 * theta, where theta is the physical waveplate angle. Vectorized.

Parameters:

theta (array_like) – HWP physical angle in degrees

Returns:

q (array_like) – Quaternion for each hwp angle

lmst(ctime, lon, **kwargs)[source]

Return local mean sidereal time for given ctimes and longitudes. Vectorized, input arguments must be broadcastable to the same shape.

Parameters:
  • ctime (array_like) – Unix time in seconds UTC

  • lon (array_like) – Observer longitude (degrees)

Returns:

lmst (array_like) – Local mean sidereal time of the observer

Notes

Any keywords accepted by the set() method can also be passed here, and will be processed prior to calculation.

load_bulletin_a(filename, columns=['mjd', 'dut1', 'x', 'y'], **kwargs)[source]

Load IERS Bulletin A from file and store in memory. The file must be readable using numpy.loadtxt with unpack=True, and is assumed to be sorted by mjd.

Parameters:
  • filename (string) – Name of the text file containing IERS Bulletin A parameters.

  • columns (list of strings) – list of columns as they appear in the file. A KeyError is raise if the list does not contain each of [‘mjd’, ‘dut1’, ‘x’, ‘y’].

  • function (Any other keyword arguments are passed to the numpy.loadtxt)

Returns:

  • mjd (array_like) – Modified Julian date

  • dut1 (array_like) – UT1-UTC time correction

  • x (array_like)

  • y (array_like) – Polar motion corrections

print_memory()[source]

Print current memory state in C.

quat2pix(quat, nside=256, pol=True, **kwargs)[source]

Calculate HEALpix pixel number and optional polarization angle for a given orientation.

Parameters:
  • quat (quaternion or array of quaternions) – Pointing orientation(s)

  • nside (int, optional) – HEALpix resolution parameter

  • pol (bool, optional) – If True, return sin2psi and cos2psi along with the pixel number(s)

Returns:

  • pix (array_like) – Pixel number(s) for the given input quaternion(s)

  • sin2psi (array_like)

  • cos2psi (array_like) – Polarization coefficients, if pol is True.

Notes

Any keywords accepted by the set() method can also be passed here, and will be processed prior to calculation.

quat2pixpa(quat, nside=256, **kwargs)[source]

Calculate ra/dec/pa for input quaternion(s).

Parameters:

q (array_like) – Pointing quaternion

Returns:

  • pix (array_like) – Pixel number(s) corresponding to the input positions(s).

  • pa (array_like) – Position angle

Notes

Any keywords accepted by the set() method can also be passed here, and will be processed prior to calculation.

quat2radecpa(quat, **kwargs)[source]

Calculate ra/dec/pa for input quaternion(s).

Parameters:

q (array_like) – Pointing quaternion

Returns:

  • ra (array_like) – Right ascension angle

  • dec (array_like) – Declination angle

  • pa (array_like) – Position angle

Notes

Any keywords accepted by the set() method can also be passed here, and will be processed prior to calculation.

radec2azel(ra, dec, pa, lon, lat, ctime, az=None, el=None, hpa=None, **kwargs)[source]

Estimate the horizon coordinates for a given set of equatorial coordinates (ra/dec/psi), location on the earth (lon/lat) and UTC time. Input vectors must be numpy-array-like and broadcastable to the same shape.

Parameters:
  • ra (array_like) – Right ascension angle

  • dec (array_like) – Declination angle

  • pa (array_like) – Position angle in equatorial coordinates

  • lon (array_like) – Observer longitude in degrees.

  • lat (array_like) – Observer latitude in degrees.

  • ctime (array_like) – Unix time in seconds UTC

Returns:

  • az (array_like) – Azimuth in degrees

  • el (array_like) – Elevation in degrees

  • hpa (array_like) – Position angle in horizon coordinates

Notes

Any keywords accepted by the set() method can also be passed here, and will be processed prior to calculation.

radec2gal(ra, dec, pa=None, sin2psi=None, cos2psi=None, inplace=True, **kwargs)[source]

Rotate celestial (equatorial) coordinates to galactic coordinates. This is equivalent to calling rotate_coord() with coord=[‘C’, ‘G’]. Vectorized, input arguments must be broadcastable to the same shape.

Parameters:
  • ra (array_like) – – or –

  • dec (array_like) – – or –

  • pa (array_like) – – or –

  • ra – arrays of coordinates, of shape (n,). If none of pa, sin2psi or cos2psi are supplied, pa = 0 is assumed.

  • dec – arrays of coordinates, of shape (n,). If none of pa, sin2psi or cos2psi are supplied, pa = 0 is assumed.

  • sin2psi (array_like) – arrays of coordinates, of shape (n,). If none of pa, sin2psi or cos2psi are supplied, pa = 0 is assumed.

  • cos2psi (array_like) – arrays of coordinates, of shape (n,). If none of pa, sin2psi or cos2psi are supplied, pa = 0 is assumed.

  • inplace (bool, optional) – If True, apply the rotation in-place on the input quaternion. Otherwise, return a copy of the input array. Default: True.

Returns:

  • l, b, pa (array_like) – – or –

  • l, b, sin2psi, cos2psi (array_like) – rotated coordinate arrays, same form as input

Notes

Any keywords accepted by the set() method can also be passed here, and will be processed prior to calculation.

radec2pix(ra, dec, nside=256, **kwargs)[source]

Calculate HEALpix pixel number for given ra/dec and nside. Vectorized, input arguments must be broadcastable to the same shape.

Parameters:
  • ra (array_like) – Right ascension angle

  • dec (array_like) – Declination angle

  • nside (int) – HEALpix resolution parameter

Returns:

pix (array_like) – Pixel number(s) corresponding to the input positions(s).

Notes

Any keywords accepted by the set() method can also be passed here, and will be processed prior to calculation.

radecpa2quat(ra, dec, pa, **kwargs)[source]

Calculate quaternion for input ra/dec/pa. Vectorized, input arguments must be broadcastable to the same shape.

Parameters:
  • ra (array_like) – Right ascension angle

  • dec (array_like) – Declination angle

  • pa (array_like) – Position angle

Returns:

q (array_like) – Quaternion constructed from the input angles.

Notes

Any keywords accepted by the set() method can also be passed here, and will be processed prior to calculation.

refraction(q, **kwargs)[source]
refraction(delta) None

Update refraction parameters

Parameters:
  • q (quaternion or array of quaternions) – Observer orientation in horizon coordinates

  • temperature (float) – Ambient temperature, Celcius

  • pressure (float) – Ambient pressure, mbar

  • humidity (float) – Ambient relative humidity, fraction

  • frequency (float) – Observing frequency, GHz

  • delta (float) – The refraction correction itself, in degrees

Returns:

delta (array_like) – Refraction correction computed at each input orientation

Notes

If q is given, then the refraction correction in degrees is calculated, stored and returned after updating any other given parameters. Otherwise, the correction is returned w/out recalculating.

Alternatively, if a single numerical argument, or the delta keyword argument is given, then the correction is stored with this value instead of being recalculated.

Numpy-vectorized for the q argument. Note that this is not an efficient vectorization, and only the last calculated value is stored for use in the coordinate conversion functions.

reset_inv_rates()[source]

Reset update counters for each inverse state. Useful to force an updated correction term at the beginning of each chunk.

reset_rates()[source]

Reset update counters for each state. Useful to force an updated correction term at the beginning of each chunk.

rotate_coord(ra, dec, pa=None, sin2psi=None, cos2psi=None, coord=['C', 'G'], inplace=True, **kwargs)[source]

Rotate coordinates from one coordinate system to another. Vectorized, input arguments must be broadcastable to the same shape. Supported coordinates:

C: celestial (equatorial) coordinates G: galactic coordinates

Parameters:
  • ra (array_like) – – or –

  • dec (array_like) – – or –

  • pa (array_like) – – or –

  • ra – arrays of coordinates, of shape (n,) If none of pa, sin2psi or cos2psi are supplied, pa = 0 is assumed.

  • dec – arrays of coordinates, of shape (n,) If none of pa, sin2psi or cos2psi are supplied, pa = 0 is assumed.

  • sin2psi (array_like) – arrays of coordinates, of shape (n,) If none of pa, sin2psi or cos2psi are supplied, pa = 0 is assumed.

  • cos2psi (array_like) – arrays of coordinates, of shape (n,) If none of pa, sin2psi or cos2psi are supplied, pa = 0 is assumed.

  • coord (list, optional) – 2-element list of input and output coordinates. Supported systems: C, G.

  • inplace (bool, optional) – If True, apply the rotation in-place on the input coordinates. Otherwise, return a copy of the input array. Default: True.

Returns:

  • ra, dec, pa (array_like) – – or –

  • ra, dec, sin2psi, cos2psi (array_like) – rotated coordinate arrays, same form as input

Notes

Any keywords accepted by the set() method can also be passed here, and will be processed prior to calculation.

rotate_map(map_in, coord=['C', 'G'], map_out=None, interp_pix=True, **kwargs)[source]

Rotate a polarized 3-x-npix map from one coordinate system to another. Supported coordinates:

C = celestial (equatorial J2000) G = galactic

Parameters:
  • map_in (array_like) – Input map, of shape (3, N)

  • coord (list, optional) – 2-element list of input and output coordinates. Supported systems: C, G.

  • map_out (array_like, optional) – Rotated output map, for inplace operation. Same shape as map_in.

  • interp_pix (bool, optional) – If True, interpolate the rotated map.

Returns:

map_out (array_like) – Rotated output map.

Notes

Any keywords accepted by the set() method can also be passed here, and will be processed prior to calculation.

Only full-sky maps are currently supported.

rotate_quat(quat, coord=['C', 'G'], inplace=True, **kwargs)[source]

Rotate a quaternion from one coordinate system to another. Supported coordinates:

C: celestial (equatorial) coordinates G: galactic coordinates

Parameters:
  • quat (array_like) – array of quaternions, of shape (n, 4)

  • coord (list, optional) – 2-element list of input and output coordinates

  • inplace (bool, optional) – If True, apply the rotation in-place on the input quaternion. Otherwise, return a copy of the input array. Default: True.

Returns:

quat (array_like) – rotated quaternion array

Notes

Any keywords accepted by the set() method can also be passed here, and will be processed prior to calculation.

set(**kwargs)[source]

Set computation options.

Parameters:
  • rate_daber ({'never', 'once', 'always'}, or float) – Rate at which the diurnal aberration correction is applied in seconds (NB: this can only be applied always or never)

  • rate_lonlat ({'never', 'once', 'always'}, or float) – Rate at which observer’s lon and lat are updated

  • rate_wobble ({'never', 'once', 'always'}, or float) – Rate at which the polar motion correction is updated (NB: these are not estimated for dates beyond a year from now)

  • rate_dut1 ({'never', 'once', 'always'}, or float) – Rate at which the ut1-utc correction is updated (NB: this is not estimated for dates beyond a year from now)

  • rate_erot ({'never', 'once', 'always'}, or float) – Rate at which the earth’s rotation angle is updated

  • rate_npb ({'never', 'once', 'always'}, or float) – Rate at which the nutation/precession/frame-bias terms are updated

  • rate_aaber ({'never', 'once', 'always'}, or float) – Rate at which the annual aberration correction (due to the earth’s orbital velocity) is updated

  • rate_ref ({'never', 'once', 'always'}, or float) – Rate at which the refaction correction is updated (NB: this correction can also be updated manually – see refraction())

  • accuracy ('low' or 'high') – If ‘low’, use a truncated form (2000b) for the NPB correction, which is much faster but less accurate. If ‘high’ (default), use the full 2006/2000a form.

  • mean_aber (bool) – If True, apply the aberration correction as an average for the entire field of view. This is gives a 1-2 arcsec deviation at the edges of the SPIDER field of view.

  • fast_aber (bool) – If True, apply the aberration correction using the exact angle and quaternion rotation. Otherwise, use a small-angle approximation for the aberration quaternion (default).

  • fast_math (bool) – If True, use polynomial approximations for trig functions

  • polconv ('cosmo' or 'iau') – Specify the ‘cosmo’ or ‘iau’ polarization convention

  • pix_order ('nest' or 'ring') – HEALPix pixel ordering

  • interp_pix (bool) – If True, interpolate between pixels in scanning the source map.

  • fast_pix (bool) – If True, use vec2pix to get pixel number directly from the quaternion instead of ang2pix from ra/dec.

  • error_missing (bool) – If True, raise an error if reading/writing missing pixels.

  • nan_missing (bool) – If True, fill samples from missing pixels with NaN. Only used if error_missing is False.

  • interp_missing (bool) – If True and interp_pix is True, drop missing neighbors and reweight remaining neighbors. Overrides nan_missing.

  • num_threads (bool) – Number of openMP threads to use for mapmaking.

  • temperature (float) – Ambient temperature, Celcius. For computing refraction corrections.

  • pressure (float) – Ambient pressure, mbar. For computing refraction corrections.

  • humidity (float) – Relative humidity, fraction. For computing refraction corrections.

  • frequency (float) – Observer frequency, GHz. For computing refraction corrections.

  • dut1 (float) – UT1 correction

  • ref_delta (float) – Refraction correction

update_bulletin_a(start_year=2000)[source]

Update the IERS Bulletin A database using astropy, and return the stored entries. Issues an ImportWarning if astropy version 1.2 or newer is not found.

Parameters:

start_year (int, optional) – Oldest year for which data should be stored.

Returns:

  • mjd (array_like) – Modified Julian date

  • dut1 (array_like) – UT1-UTC time correction

  • x (array_like)

  • y (array_like) – Polar motion (wobble) corrections

QMap Class

The class for projecting time-ordered data onto maps and vice-versa.

class qpoint.qmap_class.QMap(nside=None, pol=True, vpol=False, source_map=None, source_pol=True, source_vpol=False, q_bore=None, ctime=None, q_hwp=None, **kwargs)[source]

Bases: QPoint

Quaternion-based mapmaker that generates per-channel pointing on-the-fly.

__init__(nside=None, pol=True, vpol=False, source_map=None, source_pol=True, source_vpol=False, q_bore=None, ctime=None, q_hwp=None, **kwargs)[source]

Initialize the internal structures and data depo for mapmaking and/or timestream generation.

Parameters:
  • nside (int, optional) – Resolution of the output maps. Default: None If None, dest is not initialized.

  • pol (bool, optional) – If True, output maps are polarized.

  • vpol (bool, optional) – If True, output maps contain V polarization.

  • source_map (array_like, optional) – If supplied, passed to init_source() to initialize the source map structure.

  • source_pol (bool, optional) – If True, source_map is polarized. See init_source() for details.

  • source_vpol (bool, optional) – If True, source_map contains V polarization.

  • q_bore (array_like, optional) – Boresight pointing data. See init_point() for details. If not supplied, the pointing structure is left uninitialized.

  • ctime (array_like, optional) – Boresight pointing data. See init_point() for details. If not supplied, the pointing structure is left uninitialized.

  • q_hwp (array_like, optional) – Boresight pointing data. See init_point() for details. If not supplied, the pointing structure is left uninitialized.

Notes

Remaining keyword arguments are passed to the set() method.

An internal depo dictionary attribute stores the source and output maps, timestreams, and pointing data for retrieval by the user. Only pointers to these arrays in memory are passed to the C library. To ensure that extraneous copies of data are not made, supply these methods with C-contiguous arrays of the correct shape.

depo

Dictionary of source and output maps, timetreams and pointing data. Pointers to these arrays in memory are passed to the C library.

dest_is_init()[source]

Return True if the dest map is initialized, otherwise False.

dest_is_pol()[source]

Return True if the destination map is polarized, otherwise False. Raise an error if destination map is not initialized.

dest_is_vpol()[source]

Return True if the destination map contains V polarization, otherwise False. Raise an error if destination map is not initialized.

from_tod(q_off, tod=None, count_hits=True, weight=None, gain=None, mueller=None, flag=None, weights=None, do_diff=False, **kwargs)[source]

Calculate signal and hits maps for given detectors.

Parameters:
  • q_off (array_like) – quaternion offset array, of shape (ndet, 4)

  • tod (array_like, optional) – output array for timestreams, of shape (ndet, nsamp) if not supplied, only the projection map is populated.

  • count_hits (bool, optional) – if True (default), populate projection map.

  • weight (array_like, optional) – array of channel weights, of shape (ndet,). Defaults to 1 if not supplied.

  • gain (array_like, optional) – Per-channel gains, of shape (ndet,) or a constant. Default : 1.

  • mueller (array_like, optional) – array of Mueller matrix A/B/C elements, of shape (ndet,3). Defaults to [1, 1, 0] per channel if not supplied.

  • flag (array_like, optional) – array of flag timestreams for each channel, of shape (ndet, nsamp).

  • weights (array_like, optional) – array of weight timestreams for each channel, of shape (ndet, nsamp).

  • do_diff (do timestream differencing. Assumes first half of tods are) – one pair and the second half are the other.

Returns:

  • vec (array_like, optional) – binned signal map, if tod is supplied

  • proj (array_like, optional) – binned projection matrix map, if count_hits is True

Notes

The remaining keyword arguments are passed to the set() method.

init_dest(nside=None, pol=True, vec=None, proj=None, pixels=None, vpol=False, copy=False, reset=False, update=False)[source]

Initialize the destination map structure. Timestreams are binned and projection matrices accumulated into this structure.

Parameters:
  • nside (int, optional) – map dimension. If pixels is supplied, this argument is required. Otherwise, the default is 256.

  • pol (bool, optional) – If True, a polarized map will be created.

  • vec (array_like or bool, optional, shape (N, npix)) – If supplied, nside and pol are determined from this map, and the vector (binned signal) map is initialized from this. If False, accumulation of this map from timestreams is disabled.

  • proj (array_like or bool, optional, shape (N*(N+1)/2, npix)) – Array of upper-triangular elements of the projection matrix for each pixel. If not supplied, a blank map of the appropriate shape is created. If False, accumulation of the projection matrix is disabled.

  • pixels (1D array_like, optional) – Array of pixel numbers for each map index, if vec and proj are partial maps.

  • vpol (bool, optional) – If True, a polarized map including V polarization will be created.

  • copy (bool, optional) – If True and vec/proj are supplied, make copies of these inputs to avoid in-place operations.

  • reset (bool, optional) – If True, and if the structure has already been initialized, it is reset and re-initialized with the new map. If False, a RuntimeError is raised if the structure has already been initialized.

  • update (bool, optional) – If True, and if the structure has already been initialized, the supplied vec and proj are replaced in the existing dest structure rather than reinitializing from scratch.

init_detarr(q_off, weight=None, gain=None, mueller=None, tod=None, flag=None, weights=None, do_diff=False, write=False)[source]

Initialize the detector listing structure. Detector properties and timestreams are passed to and from the mapmaker through this structure.

Parameters:
  • q_off (array_like) – Array of offset quaternions, of shape (ndet, 4).

  • weight (array_like, optional) – Per-channel mapmaking weights, of shape (ndet,) or a constant. Default : 1.

  • gain (array_like, optional) – Per-channel gains, of shape (ndet,) or a constant. Default : 1.

  • mueller (array_like, optional) – Per-channel polarization efficiencies, of shape(ndet, 4). Default : [1., 1., 0., 1.] per det.

  • tod (array_like, optional) – Timestream array, of shape (ndet, nsamp). nsamp must match that of the pointing structure. If not supplied and write is True, then a zero-filled timestream array is initialized.

  • flag (array_like, optional) – Flag array, of shape (ndet, nsamp), for excluding data from mapmaking. nsamp must match that of the pointing structure. If not supplied, a zero-filled array is initialized (i.e. no flagged samples).

  • weights (array_like, optional) – Weight array, of shape (ndet, nsamp), for weighting each sample of data. nsamp must match that of the pointing structure. If not supplied, this option is not used.

  • do_diff (bool, optional) – If True, initialize pairs of arrays for pair-differenced mapmaking.

  • write (bool, optional) – If True, the timestreams are ensured writable and created if necessary.

init_point(q_bore=None, ctime=None, q_hwp=None)[source]

Initialize or update the boresight pointing data structure.

Parameters:
  • q_bore (array_like, optional) – Boresight pointing quaternion, of shape (nsamp, 4). If supplied, the pointing structure is reset if already initialized.

  • ctime (array_like, optional) – time since the UTC epoch. If not None, the time array is updated to this. Shape must be (nsamp,)

  • q_hwp (array_like, optional) – Waveplate quaternion. If not None, the quaternion is updated to this. Shape must be (nsamp, 4)

init_source(source_map, pol=True, pixels=None, nside=None, vpol=False, reset=False, update=False)[source]

Initialize the source map structure. Timestreams are produced by scanning this map.

Parameters:
  • source_map (array_like) – Input map. Must be of shape (N, npix), where N can be 1, 3, 6, 9, or 18.

  • pol (bool, optional) – If True, and the map shape is (3, npix), then input is a polarized map (and not T + first derivatives).

  • pixels (1D array_like, optional) – Array of pixel numbers for each map index, if source_map is a partial map.

  • nside (int, optional) – map dimension. If pixels is supplied, this argument is required. Otherwise, the nside is determined from the input map.

  • vpol (bool, optional) – If True, and the input map shape is (4, npix), then input is a polarized map that includes V polarization.

  • reset (bool, optional) – If True, and if the structure has already been initialized, it is reset and re-initialized with the new map. If False, a RuntimeError is raised if the structure has already been initialized.

  • update (bool, optional) – If True, and if the structure has already been initialized, the supplied source_map is replaced in the existing source structure rather than reinitializing from scratch.

Notes

This method will automatically determine the type of map given its shape. Note that for N=3, the pol keyword argument should be used to disambiguate the two map types. By default, a polarized map with (T,Q,U) components is assumed.

  • A map of shape (1, npix) or (npix,) contains only a T map.

  • A map of shape (3, npix) contains (T, Q, U) if pol is True, or (T, dTdt, dTdp) if pol is False.

  • A map of shape (4, npix) contains (T, Q, U, V).

  • A map of shape (6, npix) contains (T, dTdt, dTdp, dT2dt2, dT2dpdt, dT2dp2).

  • A map of shape (9, npix) contains (T, Q, U, dTdt, dQdt, dUdt, dTdp, dQdp, dUdp).

  • A map of shape (18, npix) contains all the columns of the 9-column map, followed by (dT2dt2, dQ2dt2, dU2dt2, dT2dpdt, dQ2dpdt, dU2dpdt, dT2dp2, dQ2dp2, dU2dp2).

point_is_init()[source]

Return True if the point map is initialized, otherwise False.

proj_cond(proj=None, mode=None, partial=False)[source]

Hits-normalized projection matrix condition number for each pixel.

Parameters:
  • proj (array_like) – Projection matrix, of shape (N*(N+1)/2, npix). If None, obtained from the depo.

  • mode ({None, 1, -1, 2, -2, inf, -inf, 'fro'}, optional) – condition number order. See numpy.linalg.cond. Default: None (2-norm from SVD)

  • partial (bool, optional) – If True, the map is not checked to ensure a proper healpix nside.

Returns:

cond (array_like) – Condition number of each pixel.

reset()[source]

Reset the internal data structures, and clear the data depo.

reset_dest()[source]

Reset the destination map structure. Must be reinitialized to continue mapmaking.

reset_detarr()[source]

Reset the detector array structure.

reset_point()[source]

Reset the pointing data structure.

reset_source()[source]

Reset the source map structure. Must be reinitialized to produce more timestreams.

solve_map(vec=None, proj=None, mask=None, copy=True, return_proj=False, return_mask=False, partial=None, fill=0, cond=None, cond_thresh=1000000.0, method='exact')[source]

Solve for a map, given the binned map and the projection matrix for each pixel.

Parameters:
  • vec (array_like, optional) – A map or list of N maps. Default to the “vec” entry in depo.

  • proj (array_like, optional) – An array of upper-triangular projection matrices for each pixel, of shape (N*(N+1)/2, npix). Default to the “proj” entry in depo.

  • mask (array_like, optional) – A mask of shape (npix,), evaluates to True where pixels are valid. The input mask in converted to a boolean array if supplied.

  • copy (bool, optional) – if False, do the computation in-place so that the input maps are modified. Otherwise, a copy is created prior to solving. Default: False.

  • return_proj (bool, optional) – if True, return the Cholesky-decomposed projection matrix. if False, and inplace is True, the input projection matrix is not modified.

  • return_mask (bool, optional) – if True, return the mask array, updated with any pixels that could not be solved.

  • partial (bool, optional) – If True, the map is not checked to ensure a proper healpix nside.

  • fill (scalar, optional) – Fill the solved map where proj == 0 with this value. Default: 0.

  • cond (array_like, optional) – A map of condition number per pixel. If not supplied, this will be calculated using proj_cond()

  • cond_thresh (scalar, optional) – A threshold to place on the condition number to exclude pixels prior to solving. Reduce this to avoid LinAlgError due to singular matrices.

  • method (string, optional) – Map inversion method. If “exact”, invert the pointing matrix directly If “cho”, use Cholesky decomposition to solve. Default: “exact”.

Returns:

  • map (array_like) – A solved map or set of maps, in shape (N, npix).

  • proj_out (array_like) – The upper triangular elements of the decomposed projection matrix, (if method is ‘cho’) or of the matrix inverse (if method is ‘exact’), if requested, in shape (N*(N+1)/2, npix).

  • mask (array_like) – 1-d array, True for valid pixels, if return_mask is True

solve_map_cho(*args, **kwargs)[source]

Solve for a map, given the binned map and the projection matrix for each pixel, using Cholesky decomposition. This method uses the scipy.linalg.cho_factor and scipy.linalg.cho_solve functions internally.

Parameters:
  • vec (array_like, optional) – A map or list of N maps. Default to the “vec” entry in depo.

  • proj (array_like, optional) – An array of upper-triangular projection matrices for each pixel, of shape (N*(N+1)/2, npix). Default to the “proj” entry in depo.

  • mask (array_like, optional) – A mask of shape (npix,), evaluates to True where pixels are valid. The input mask in converted to a boolean array if supplied.

  • copy (bool, optional) – if False, do the computation in-place so that the input maps are modified. Otherwise, a copy is created prior to solving. Default: False.

  • return_proj (bool, optional) – if True, return the Cholesky-decomposed projection matrix. if False, and inplace is True, the input projection matrix is not modified.

  • return_mask (bool, optional) – if True, return the mask array, updated with any pixels that could not be solved.

  • partial (bool, optional) – If True, the map is not checked to ensure a proper healpix nside.

  • fill (scalar, optional) – Fill the solved map where proj == 0 with this value. Default: 0.

  • cond (array_like, optional) – A map of condition number per pixel. If not supplied, this will be calculated using proj_cond()

  • cond_thresh (scalar, optional) – A threshold to place on the condition number to exclude pixels prior to solving. Reduce this to avoid LinAlgError due to singular matrices.

Returns:

  • map (array_like) – A solved map or set of maps, in shape (N, npix).

  • proj_out (array_like) – The upper triangular elements of the decomposed projection matrix, if requested, in shape (N*(N+1)/2, npix).

  • mask (array_like) – 1-d array, True for valid pixels, if return_mask is True

source_is_init()[source]

Return True if the source map is initialized, otherwise False.

source_is_pol()[source]

Return True if the source map is polarized, otherwise False. Raise an error if source map is not initialized.

source_is_vpol()[source]

Return True if the source map contains V polarization, otherwise False. Raise an error if source map is not initialized.

to_tod(q_off, gain=None, mueller=None, tod=None, flag=None, **kwargs)[source]

Calculate signal TOD from source map for multiple channels.

Parameters:
  • q_off (array_like) – quaternion offset array, of shape (ndet, 4)

  • gain (array_like, optional) – Per-channel gains, of shape (ndet,) or a constant. Default : 1.

  • mueller (array_like, optional) – array of Mueller matrix A/B/C/D elements, of shape (ndet, 4). Defaults t [1, 1, 0, 1] per channel if not supplied.

  • tod (array_like, optional) – output array for timestreams, of shape (ndet, nsamp) use this keyword argument for in-place computation.

Returns:

tod (array_like) – A timestream sampled from the input map for each requested detector. The output array shape is (ndet, nsamp).

Notes

The remaining keyword arguments are passed to the set() method.

unsolve_map(map_in, proj=None, mask=None, copy=True, return_proj=False, return_mask=False, partial=None, fill=0)[source]

Invert the solved map to recover the binned vec array.

Parameters:
  • map_in (array_like) – A map or list of N maps.

  • proj (array_like, optional) – An array of upper-triangular projection matrices for each pixel, of shape (N*(N+1)/2, npix). Default to the “proj” entry in depo.

  • mask (array_like, optional) – A mask of shape (npix,), evaluates to True where pixels are valid. The input mask in converted to a boolean array if supplied.

  • copy (bool, optional) – if False, do the computation in-place so that the input maps are modified. Otherwise, a copy is created prior to solving. Default: False.

  • return_proj (bool, optional) – if True, return the Cholesky-decomposed projection matrix. if False, and inplace is True, the input projection matrix is not modified.

  • return_mask (bool, optional) – if True, return the mask array, updated with any pixels that could not be solved.

  • partial (bool, optional) – If True, the map is not checked to ensure a proper healpix nside.

  • fill (scalar, optional) – Fill the solved map where proj == 0 with this value. Default: 0.

Returns:

  • vec (array_like) – A binned map or set of maps, in shape (N, npix).

  • proj_out (array_like) – The upper triangular elements of the projection matrix, if requested, in shape (N*(N+1)/2, npix).

  • mask (array_like) – 1-d array, True for valid pixels, if return_mask is True

Tools

Other tools used for testing and argument checks.

qpoint.version()[source]

Print version string

qpoint.qpoint_class.check_input(name, arg, shape=None, quat=False, dtype=<class 'numpy.float64'>, inplace=True, fill=0, allow_transpose=True, allow_tuple=True, output=False)[source]

Ensure input argument is an aligned array of the right type and shape.

Parameters:
  • name (string) – Name of the argument

  • arg (array_like) – The argument itself

  • shape (tuple, optional) – If supplied, ensure arg has this shape.

  • quat (bool, optional) – If True, ensure that the last dimension of the array has length 4.

  • dtype (numpy.dtype, optional) – Ensure arg is of this dtype. Default: numpy.double.

  • inplace (bool, optional) – If False make sure that a copy of the input arg is made prior to returning. Otherwise, the output arg may share memory with the input.

  • fill (scalar, optional) – If the input is missing or empty, fill with this value. If None, leave the array empty.

  • allow_transpose (bool, optional) – If True, transpose the input array if the transposed shape matches shape.

  • allow_tuple (bool, optional) – If True, numpy.vstack the input arg if it is a tuple.

  • output (bool, optional) – If True, ensure that the output arg is a writeable array. Otherwise, the array is only ensured to be aligned and C-contiguous.

Returns:

arg – Aligned, contiguous and properly typed and shaped array for passing to the C library.

qpoint.qpoint_class.check_inputs(*args, **kwargs)[source]

Ensure that a group of input arguments have the same shape by broadcasting.

Parameters:
  • args – A list of broadcastable array_like arguments.

  • kwargs – Dictionary of arguments to pass to check_input.

Returns:

args – A list of broadcast and properly aligned/shaped/typed arrays for passing to the C library as a set of timestreams.

qpoint.qpoint_class.check_output(name, arg=None, shape=None, quat=False, dtype=<class 'numpy.float64'>, inplace=True, fill=None, allow_transpose=True, allow_tuple=True, **kwargs)[source]

Ensure that the output argument is properly aligned/shaped/typed. Check input kwargs to see if a pointer to the output array has been passed in.

Parameters:
  • name (string) – Name of the argument

  • arg (array_like, optional) – The argument itself. If not supplied, kwargs is checked.

  • shape (tuple, optional) – If supplied, ensure arg has this shape. If arg is None, this argument is required.

  • dtype (numpy.dtype, optional) – Ensure arg is of this dtype. Default: numpy.double.

  • inplace (bool, optional) – If False make sure that a copy of the input arg is made prior to returning. Otherwise, the output arg may share memory with the input.

  • fill (scalar, optional) – If the input is missing or empty, fill with this value. If None, leave the array empty.

  • allow_transpose (bool, optional) – If True, transpose the input array if the transposed shape matches shape.

  • allow_tuple (bool, optional) – If True, numpy.vstack the input arg if it is a tuple.

  • kwargs – Any remaining input arguments. If arg is None, kwargs is searched for the name key. If not found, a an empty array is created of the appropriate shape.

Returns:

arg – Aligned, contiguous, writeable and properly typed and shaped array for passing to the C library.

qpoint.qmap_class.check_map(map_in, copy=False, partial=False, dtype=<class 'numpy.float64'>)[source]

Return a properly transposed and memory-aligned map and its nside.

Parameters:
  • map_in (map or list of maps) – Input map(s)

  • copy (bool, optional) – If True, ensure that output map does not share memory with the input map. Use if you do not want in-place operations to modify the map contents.

  • partial (bool, optional) – If True, the map is not checked to ensure a proper healpix nside, and the number of pixels is returned instead.

Returns:

  • map_out (numpy.ndarray) – Properly shaped and memory-aligned map, copied from the input if necessary.

  • nside or npix (int) – If partial is False, the map nside. Otherwise, the number of pixels.

qpoint.qmap_class.check_proj(proj_in, copy=False, partial=False)[source]

Return a properly transposed and memory-aligned projection map, its nside, and the map dimension.

Parameters:
  • proj_in (map or list of maps) – Input projection matrix map

  • copy (bool, optional) – If True, ensure that output map does not share memory with the input map. Use if you do not want in-place operations to modify the map contents.

  • partial (bool, optional) – If True, the map is not checked to ensure a proper healpix nside, and the number of pixels is returned instead.

Returns:

  • map_out (numpy.ndarray) – Properly shaped and memory-aligned map, copied from the input if necessary.

  • nside or npix (int) – If partial is False, the map nside. Otherwise, the number of pixels.

  • nmap (int) – The map size this projection matrix is intended to invert, i.e. the solution to len(proj) = nmap * (nmap + 1) / 2. Raises an error if an integer solution is not found.

qpoint.tools.refraction(el, temp, press, hum, freq=150.0)[source]

Standalone function for calculating the refraction correction without storing any parameters. Useful for testing, numpy-vectorized.

Parameters:
  • el (array_like) – Observer elevation angle, degrees

  • temperature (array_like) – Ambient temperature, Celcius

  • pressure (array_like) – Ambient pressure, mbar

  • humidity (array_like) – Relative humidity, fraction

  • frequency (array_like) – Observing frequency, GHz

Returns:

delta (array_like) – Refraction correction, in degrees