exocartographer package

Submodules

Gaussian Process-based Map

Utilities to implement Gaussian process priors on healpix maps.

exocartographer.gp_map.draw_map(nside, mu, sigma, wn_rel_amp, lambda_spatial, nest=False)

Returns a map sampled from the Gaussian process with the given parameters.

Parameters:
  • mu – The mean of the GP (constant on the sphere).
  • sigma – The standard deviation of the GP at zero angular separation.
  • lambda_spatial – The angular correlation length of the process (radians).
  • nest – Healpix map ordering.
Returns:

Healpix map drawn from the GP with the given parameters.

exocartographer.gp_map.draw_map_cl(nside, mu, sigma, wn_rel_amp, lambda_spatial)

Returns a map statistically equivalent to draw_map(), but is much faster.

exocartographer.gp_map.exp_cov(nside, wn_rel_amp, lambda_angular, nest=False)

Returns a covariance function for a healpix map with nside.

\[C_{ij} = \left \langle p_i p_j \right\rangle = \exp\left( -\frac{1}{2} \left(\frac{\theta_{ij}}{\lambda}\right)^2 \right)\]

where \(\theta_{ij}\) is the great-circle angle between the \(i\) and \(j\) pixel centres.

Parameters:
  • nside – Healpix nside parameter.
  • wn_rel_amp – The relative amplitude of white-noise added to the covariance matrix. In other words, the pixel variance is always 1, but wn_rel_amp fraction of it comes from white noise, and 1-wn_rel_amp of it from the correlated noise with the above angular correlation function.
  • lambda_angular – The angular correlation length (radians).
  • nest – Ordering flag for the healpix map.
exocartographer.gp_map.exp_cov_cl(nside, wn_rel_amp, lambda_angular)

Returns the equivalent \(C_l\) angular power spectrum to the covariance returned by exp_cov(). The returned array is normalised by the Healpy convention.

exocartographer.gp_map.map_logprior(hpmap, mu, sigma, wn_rel_amp, lambda_angular, nest=False)

Returns the GP prior on the map with exponential covariance function.

Parameters:
  • hpmap – Healpix map on which the prior is to be evaluated.
  • mu – Mean of the GP
  • sigma – Standard deviation at zero angular separation.
  • lambda_angular – Angular correlation length.
  • nest – The ordering of the healpix map.
exocartographer.gp_map.map_logprior_cl(map, mu, sigma, wn_rel_amp, lambda_angular)

Returns the GP prior on maps, consistent with map_logprior(), but computed much more quickly. See map_logprior() for the description of the arguments.

The logic is the following: any Gaussian map whose correlation matrix depends only on angular separation between points will have a diagonal covariance matrix in Ylm space. (This is the spatial analog of stationarity in the time domain; stationary processes are completely described by their Fourier-space power spectrum, which is the diagonal of the covariance matrix in Fourier space.) Furthermore, any map that is statistically isotropic will have a Ylm covariance that is a funciton of l only (basically, rotational invariance demands that the variance of the different m’s at constant l be equal, since rotations mix these components together).

So, instead of computing the covariance matrix and using a multivariate Gaussian likelihood in pixel space, we compute the \(a_{lm}\) coefficients, and use a Gaussian likelihood in Spherical Harmonic space. Each \(a_{lm}\) is distributed like

\[a_{lm} \sim N\left(0, \sqrt{C_l}\right)\]

There is only one final wrinkle: the number of degrees of freedom in the \(a_{lm}\) coefficients is not equal to the number of pixels. So, we have a choice we can either

  • Under-constrain the pixels by using an \(l_\mathrm{max}\) that corresponds to fewer \(a_{lm}\) coefficients than pixels. Then there will be some linear combinations of pixels (corresponding to the \(Y_{lm}\) with \(l > l_\mathrm{max}\)) that are unconstrained by the prior.
  • Over-constrain the pixels by choosing an \(l_\mathrm{max}\) that corresponds to more \(a_{lm}\) coefficients than pixels. This means that the effective prior we are imposing in pixel space is not exactly the squared-exponential kernel.

We choose to do the latter, so the pixels are over-constrained. This makes sampling easier (no combinations of pixles that can run away to infinity).

exocartographer.gp_map.resolve(pix, new_nside, nest=False)

Illumination Map Posterior

class exocartographer.gp_illumination_map.IlluminationMapLikelihood(*args, **kwargs)

Bases: exocartographer.gp_illumination_map.IlluminationMapPosterior

class exocartographer.gp_illumination_map.IlluminationMapPosterior(times, reflectance, sigma_reflectance, nside=4, nside_illum=16, map_parameterization='pix')

Bases: object

A posterior class for mapping surfaces from a reflectance time series”

Parameters:
  • times – Array of times of photometric measurements.
  • reflectance – Array of photometric measurements at corresponding times.
  • sigma_reflectance – Single value or array of 1-sigma uncertainties on reflectance measurements.
  • nside – (optional) Resolution of HEALPix surface map. Has to be a power of 2. The number of pixels will be \(12 N_\mathrm{side}^2\). (default: 4)
  • nside_illum – (optional): Resolution of HEALPix illumination map, i.e., the “kernel” of illumination that is integrated against the pixel map. (default: 16)
  • map_parameterization – (optional): Parameterization of surface map to use. pix will parameterize the map with values of each pixel; alm will parameterize the map in spherical harmonic coefficients. (default: pix)
add_map_to_full_params(p, map)
analytic_visibility_illumination_matrix(p)
cos_inc(p)
cos_obl(p)
data_sigma_matrix(inv=False)
draw_map(p)
dtype
dtype_map
error_scale(p)
fix_params(params)

Fix parameters to the specified values and remove them from the Posterior’s dtype.

Parameters:params – A dictionary of parameters to fix, and the values to fix them to.
fixed_params
full_dtype
full_dtype_map
gp_sigma_matrix(p)
hpmap(p)
inc(p)
lightcurve(p)
lmax
log_map_prior(p)
log_period_prior(logP)
log_prior(p)
loglikelihood(p)
map_parameterization
mmax
nalms
nmap_params
nparams
nparams_full
nparams_full_map
nparams_map
npix
npix_illum
nside
nside_illum
ntimes
obl(p)
obl_orientation(p)
orbital_period(p)
params_map_to_params(pm, include_fixed_params=False)
phi_orb(p)
reflectance
resolved_visibility_illumination_matrix(p)
rotation_period(p)
set_params(p, dict)
sigma(p)
sigma_reflectance
spatial_scale(p)
spatial_scale_high
spatial_scale_low
sub_observer_latlong(p, times)
times
to_params(p)

Return a typed version of ndarray p.

unfix_params(params=None)

Let fixed parameters vary.

Parameters:params – A list of parameters to unfix. If None, all parameters will be allowed to vary. (default: None)
visibility_illumination_matrix(p)

Produce the “kernel” of illumination that is integrated against the pixel map.

Parameters:p – Array of values for unfixed parameters.

The kernel is composed of a product of cosines: the illumination is proportional to \(\vec{n} \cdot \vec{n_s}\), where \(\vec{n}\) is the pixel normal and \(\vec{n_s}\) is the vector to the star. The visibility is proportional to \(\vec{n} \cdot \vec{n_o}\), where \(\vec{n_o}\) is the vector to the observer. The contribution of any pixel of value p to the lightcurve is therefore \(p (\vec{n} \cdot \vec{n_s}) (\vec{n} \cdot \vec{n_o})\) if both \(\vec{n} \cdot \vec{n_s}\) and \(\vec{n} \cdot \vec{n_o}\) are > 0, and zero otherwise. So, we need to evaluate these dot products.

Fix a coordinate system in which to evaluate these dot products as follows:

The orbit is in the x-y plane (so \(\hat{L} \parallel \hat{z}\)), with the x-axis pointing along superior conjunction. The observer is therefore in the x-z plane, and has an inclination angle \(\iota\) in \([0, \pi/2]\) to the orbital plane. So \(\vec{n_o} = (-\sin(\iota), 0, \cos(\iota))\).

The planet star vector, \(\vec{n_s}\), is given by \(\vec{n_s} = (-\sin(x_i), -\cos(x_i), 0)\), where \(x_i\) is the orbital phase. If the orbit has phase \(x_{i0}\) at t = 0, then \(n_s = R_z(2\pi/P_\mathrm{orb}t + x_{i0}) \cdot (-1,0,0)\).

For the normal vector to the planet, we must describe the series of rotations that maps the orbital coordinate system into the planet-centred coordinate system. Imagine that the planet spin axis is at first aligned with the z-axis. We apply \(R_y(\mathrm{obl})\), where obl is the obliquity angle in \([0, \pi/2]\), and then \(R_z(\phi_\mathrm{rot})\), where \(\phi_\mathrm{rot}\) is the azimuthal angle of the planet’s spin axis in \([0, 2\pi]\). Now the planet’s spin axis points to \(S = (\cos(\phi_\mathrm{rot})*\sin(\mathrm{obl}), \sin(\phi_\mathrm{rot})*\sin(\mathrm{obl}), \cos(\mathrm{obl}))\). At time \(t\), the planet’s normals are given by \(n(t) = R_S(2\pi/P_\mathrm{rot} t) n(0)\), where \(n(0) = R_z(\phi_\mathrm{rot}) R_y(\mathrm{obl}) n\), with n the body-centred normals to the pixels. We can now evaluate dot products in the fixed orbital frame.

In principle, we are done, but there is an efficiency consideration. For each time at which we have an observation, we need to perform rotations on the various vectors in the problem. There are, in general, a large number of n vectors (we have many pixels covering the planet’s surface), but only one n_o and n_s vector. It will be more efficient to apply rotations only to the n_o and n_s vectors instead of rotating the many, many, n vectors. (This is equivalent to performing the dot products in the planet-centred frame; we could have done that from the beginning, but it is easier to describe the geometry in the orbit-centred frame.) We can accomplish this via a trick: when vectors in a dot product are rotated, the rotations can be moved from one vector to the other:

\[(R_A a) (R_B b) = \langle R_A a | R_B b \rangle = \langle R_B^{-1} R_A a | b \rangle = \langle a | R_A^{-1} R_B b \rangle\]

So, instead of rotating the normal vectors to the planet pixels by

\[n(t) = R_S(\omega_\mathrm{rot}(t)) R_z(\phi_\mathrm{rot}) R_y(\mathrm{obl}) n\]

We can just rotate the vectors that are inner-producted with these vectors by

\[R_\mathrm{inv} = R_y(-\mathrm{obl}) R_z(-\phi_\mathrm{rot}) R_S(-\omega_\mathrm{rot})\]
wn_high
wn_low
wn_rel_amp(p)
class exocartographer.gp_illumination_map.IlluminationMapPrior(*args, **kwargs)

Bases: exocartographer.gp_illumination_map.IlluminationMapPosterior

exocartographer.gp_illumination_map.quaternion_multiply(qa, qb)
exocartographer.gp_illumination_map.rotate_vector(rqs, v)
exocartographer.gp_illumination_map.rotation_quaternions(axis, angles)