spherical_kde package

Submodules

spherical_kde.distributions module

Module containing the kernel function for the spherical KDE.

For more detail, see: https://en.wikipedia.org/wiki/Von_Mises-Fisher_distribution

spherical_kde.distributions.VonMisesFisher_distribution(phi, theta, phi0, theta0, sigma0)[source]

Von-Mises Fisher distribution function.

Parameters
phi, thetafloat or array_like

Spherical-polar coordinates to evaluate function at.

phi0, theta0float or array-like

Spherical-polar coordinates of the center of the distribution.

sigma0float

Width of the distribution.

Returns
float or array_like

log-probability of the vonmises fisher distribution.

Notes

Wikipedia:

https://en.wikipedia.org/wiki/Von_Mises-Fisher_distribution

spherical_kde.distributions.VonMisesFisher_sample(phi0, theta0, sigma0, size=None)[source]

Draw a sample from the Von-Mises Fisher distribution.

Parameters
phi0, theta0float or array-like

Spherical-polar coordinates of the center of the distribution.

sigma0float

Width of the distribution.

sizeint, tuple, array-like

number of samples to draw.

Returns
phi, thetafloat or array_like

Spherical-polar coordinates of sample from distribution.

spherical_kde.distributions.VonMises_mean(phi, theta)[source]

Von-Mises sample mean.

Parameters
phi, thetaarray-like

Spherical-polar coordinate samples to compute mean from.

Returns
float
..math::

sum_i^N x_i / || sum_i^N x_i ||

Notes

Wikipedia:

https://en.wikipedia.org/wiki/Von_Mises-Fisher_distribution#Estimation_of_parameters

spherical_kde.distributions.VonMises_std(phi, theta)[source]

Von-Mises sample standard deviation.

Parameters
phi, thetaarray-like

Spherical-polar coordinate samples to compute mean from.

Returns
solution for

..math:: 1/tanh(x) - 1/x = R,

where

..math:: R = || sum_i^N x_i || / N

Notes

Wikipedia:

https://en.wikipedia.org/wiki/Von_Mises-Fisher_distribution#Estimation_of_parameters but re-parameterised for sigma rather than kappa.

spherical_kde.utils module

Utilities

  • General stable functions

  • Transforming coordinates

  • Computing rotations

  • Performing spherical integrals

spherical_kde.utils.cartesian_from_polar(phi, theta)[source]

Embedded 3D unit vector from spherical polar coordinates.

Parameters
phi, thetafloat or numpy.array

azimuthal and polar angle in radians.

Returns
nhatnumpy.array

unit vector(s) in direction (phi, theta).

spherical_kde.utils.decra_from_polar(phi, theta)[source]

Convert from ra and dec to spherical polar coordinates.

Parameters
phi, thetafloat or numpy.array

azimuthal and polar angle in radians

Returns
ra, decfloat or numpy.array

Right ascension and declination in degrees.

spherical_kde.utils.logsinh(x)[source]

Compute log(sinh(x)), stably for large x.

Parameters
xfloat or numpy.array

argument to evaluate at, must be positive

Returns
float or numpy.array

log(sinh(x))

spherical_kde.utils.polar_from_cartesian(x)[source]

Embedded 3D unit vector from spherical polar coordinates.

Parameters
xarray_like

cartesian coordinates

Returns
phi, thetafloat or numpy.array

azimuthal and polar angle in radians.

spherical_kde.utils.polar_from_decra(ra, dec)[source]

Convert from spherical polar coordinates to ra and dec.

Parameters
ra, decfloat or numpy.array

Right ascension and declination in degrees.

Returns
phi, thetafloat or numpy.array

Spherical polar coordinates in radians

spherical_kde.utils.rotation_matrix(a, b)[source]

The rotation matrix that takes a onto b.

Parameters
a, bnumpy.array

Three dimensional vectors defining the rotation matrix

Returns
Mnumpy.array

Three by three rotation matrix

Notes

StackExchange post:

https://math.stackexchange.com/questions/180418/calculate-rotation-matrix-to-align-vector-a-to-vector-b-in-3d

spherical_kde.utils.spherical_integrate(f, log=False)[source]

Integrate an area density function over the sphere.

Parameters
fcallable

function to integrate (phi, theta) -> float

logbool

Should the function be exponentiated?

Returns
float

Spherical area integral

\[\int_0^{2\pi}d\phi\int_0^\pi d\theta f(\phi, \theta) \sin(\theta)\]
spherical_kde.utils.spherical_kullback_liebler(logp, logq)[source]

Compute the spherical Kullback-Liebler divergence.

Parameters
logp, logqcallable

log-probability distributions (phi, theta) -> float

Returns
float

Kullback-Liebler divergence

\[\int P(x)\log \frac{P(x)}{Q(x)} dx\]

Notes

Wikipedia post:

https://en.wikipedia.org/wiki/Kullback-Leibler_divergence

Module contents

The spherical kernel density estimator class.

class spherical_kde.SphericalKDE(phi_samples, theta_samples, weights=None, bandwidth=None, density=100)[source]

Bases: object

Spherical kernel density estimator

Parameters
phi_samples, theta_samplesarray_like

spherical-polar samples to construct the kde

weightsarray_like

Sample weighting default [1] * len(phi_samples))

bandwidthfloat

bandwidth of the KDE. Increasing bandwidth increases smoothness

densityint

number of grid points in theta and phi to draw contours.

Attributes
phi, thetanumpy.array

spherical polar samples

weightsnumpy.array

Sample weighting (normalised to sum to 1).

bandwidthfloat

Bandwidth of the kde. defaults to rule-of-thumb estimator https://en.wikipedia.org/wiki/Kernel_density_estimation Set to None to use this value

densityint

number of grid points in theta and phi to draw contours.

palefactorfloat

getdist-style colouration factor of sigma-contours.

Methods

__call__(phi, theta)

Log-probability density estimate

plot(ax[, colour])

Plot the KDE on an axis.

plot_samples(ax[, nsamples])

Plot equally weighted samples on an axis.

property bandwidth
plot(ax, colour='g', **kwargs)[source]

Plot the KDE on an axis.

Parameters
axmatplotlib.axes.Axes

matplotlib axis to plot on. This must be constructed with a cartopy.crs.projection:

>>> import cartopy
>>> import matplotlib.pyplot as plt
>>> fig = plt.subplots()
>>> ax = fig.add_subplot(111, projection=cartopy.crs.Mollweide())
color

Colour to plot the contours. arg can be an RGB or RGBA sequence or a string in any of several forms:

  1. a letter from the set ‘rgbcmykw’

  2. a hex color string, like ‘#00FFFF’

  3. a standard name, like ‘aqua’

  4. a string representation of a float, like ‘0.4’,

This is passed into matplotlib.colors.colorConverter.to_rgb

plot_samples(ax, nsamples=None, **kwargs)[source]

Plot equally weighted samples on an axis.

Parameters
axmatplotlib.axes.Axes

matplotlib axis to plot on. This must be constructed with a cartopy.crs.projection:

>>> import cartopy
>>> import matplotlib.pyplot as plt
>>> fig = plt.subplots()
>>> ax = fig.add_subplot(111, projection=cartopy.crs.Mollweide())
nsamplesint

Approximate number of samples to plot. Can only thin down to this number, not bulk up