ApproxManifoldProducts.ManifoldKernelDensityType
struct ManifoldKernelDensity{M<:AbstractManifold, B<:BallTreeDensity, L, P<:AbstractArray}

On-manifold kernel density belief.

Notes

  • Allows partials as identified by list of coordinate dimensions e.g. ._partial = [1;3]
    • When building a partial belief, use full points with necessary information in the specified partial coords.

DevNotes

  • WIP AMP issue 41, use generic retractions during manifold products.
ApproxManifoldProducts.antimarginalMethod
antimarginal(newM, u0, mkd, newpartial)

If a marginal (statistics) of a probability reduce the dimensions (i.e. casts a shadow, or projects); then an antimarginal aims to increase dimension of the probability within reason.

Notes

  • Marginalization is a integration of for higher dimension to lower dimension, so antimarginal likely involve differentiation (anti-integral) instead.
  • In manifold language, this is an embedding into a higher dimensional space.
  • See structure from motion in machine vision, or stereo disparity for building depth clouds from 2D images.
  • Imagine combining three different partial embedding A=[1, nan, nan, 1.4], B=[nan,2.1,nan,4], C=[nan, nan, 3, nan] which should equal A+B+C = [notnan, notnan, notnan, notnan]
ApproxManifoldProducts.buildHybridManifoldCallbacksMethod
buildHybridManifoldCallbacks(manif)

Lots to do here, see RoME.jl #244 and standardized usage with Manifolds.jl.

Notes

  • diffop( test, reference ) <===> ΔX = inverse(test) * reference

DevNotes

  • FIXME replace with Manifolds.jl #41, RoME.jl #244
ApproxManifoldProducts.calcProductGaussiansMethod
calcProductGaussians(M, μ_, Σ_; dim, Λ_)

Calculate covariance weighted mean as product of incoming Gaussian points $μ_$ and coordinate covariances $Σ_$.

Notes

  • Return both weighted mean and new covariance (teh congruent product)
  • More efficient helper function allows passing keyword inverse covariances Λ_ instead.
  • Assume size(Σ_[1],1) == manifold_dimension(M).
  • calc lambdas first and use to calculate mean product second.
  • https://ccrma.stanford.edu/~jos/sasp/ProductTwoGaussian_PDFs.html
  • Pennec, X. Intrinsic Statistics on Riemannian Manifolds: Basic Tools for Geometric Measurements, HAL Archive, 2011, Inria, France.
ApproxManifoldProducts.evaluateManifoldNaive1D!Function
evaluateManifoldNaive1D!(ret, idx, pts, bw, x)
evaluateManifoldNaive1D!(ret, idx, pts, bw, x, loo)
evaluateManifoldNaive1D!(ret, idx, pts, bw, x, loo, diffop)

Evalute the KDE naively as equally weighted Gaussian kernels with common bandwidth. This function does, however, allow on-manifold evaluations.

ApproxManifoldProducts.getManifoldPartialFunction
getManifoldPartial(M, partial)
getManifoldPartial(M, partial, repr)
getManifoldPartial(M, partial, repr, offset; doError)

A so-called full dimension manifold can possibly be reduced to smaller partial manifolds over some of the dimensions, returning a new programatically generated <:AbstractManifold. This funciton can optionally also reduce a point representation for the desired partial dimensions too.

Example

using Manifolds
using ApproxManifoldProducts

# a familiar manifold of translation and rotation in 2D
M = SpecialEuclidean(2)


# make a new partial of only the translation components
M_x, _ = getManifoldPartial(M,[1;])
# returns a new TranslationGroup(1) corresponding to just x dimension

# representation is semidirect product of translation and rotation matrix
u0 = ArrayPartition([0.0;0],[1 0; 0 1.0])

# known coordinates are [x,y,θ], eg
#   vee(M,identity_element(M,u0),log(M,identity_element(M,u0),u0))
#   [0;0;0] in this example

# make another new partial of only the rotation component
M_rot, u_rot = getManifoldPartial(M,[3],u0)
# returns SpecialOrthogonal(2) information

# make another new partial of only  y and θ
M_yθ, u_yθ = getManifoldPartial(M,[2;3],u0)
# returns new manifold information as ProductArray(TranslationGroup(1),SpecialOrthogonal(2))

Notes

  • Partial dimensions of interest are defined via a AbstractVector{Int} of coordinate dimensions.
  • assumed to follow coordinates as transition step towards more general solution
  • assume ProductManifold is the only way to stitch multiple manifolds together
  • Still experimental.

DevNotes

  • FIXME any semidirect product or action information is lost and naive Product manifold assumptions are currently made.

Related

getManifold, [AbstractManifold], [ProductManifold], [GroupManifold], [ArrayPartition]

ApproxManifoldProducts.kerFunction
ker(M, p, q)
ker(M, p, q, sigma)
ker(M, p, q, sigma, distFnc)

Normal kernel used for Hilbert space embeddings.

ApproxManifoldProducts.makePointFromCoordsFunction
makePointFromCoords(M, coords)
makePointFromCoords(M, coords, u0)
makePointFromCoords(M, coords, u0, ϵ)
makePointFromCoords(M, coords, u0, ϵ, retraction_method)

Helper function to convert coordinates to a desired on-manifold point.

DevNotes

  • FIXME need much better consolidation or even removal of this function entirely.
    • This function makes too strong an assumption of groups

Notes

  • u0 is used to identify the data type for a point
  • Pass in a different exp if needed.
ApproxManifoldProducts.manifoldLooCrossValidationMethod
manifoldLooCrossValidation(pts, bw; own, diffop)

Calculate negative entropy with leave one out (j'th element) cross validation.

Background

From: Silverman, B.: Density Estimation for Statistics and Data Analysis, 1986, p.52

Probability density function p(x), as estimated by kernels

\[hatp_{-j}(x) = 1/(N-1) Σ_{i != j}^N frac{1}{sqrt{2pi}σ } exp{ -frac{(x-μ)^2}{2 σ^2} }\]

and has Cross Validation number as the average log evaluations of leave one out hatp_{-j}(x):

\[CV(p) = 1/N Σ_i^N log hat{p}_{-j}(x_i)\]

This quantity CV is related to an entropy H(p) estimate via:

\[H(p) = -CV(p)\]

ApproxManifoldProducts.manifoldProductMethod
manifoldProduct(ff)
manifoldProduct(
    ff,
    mani;
    makeCopy,
    Niter,
    ndims,
    N,
    u0,
    oldPoints,
    addEntropy,
    recordLabels,
    selectedLabels,
    _randU,
    _randN,
    logger
)

Approximate the pointwise the product of functionals on manifolds using KernelDensityEstimate.

Notes:

  • Always pass full beliefs, for partials use for e.g. partialDimsWorkaround=[1;3;6]
  • Can also multiply different partials together

Example

# setup
M = TranslationGroup(3)
N = 75
p = manikde!(M, [randn(3) for _ in 1:N])
q = manikde!(M, [randn(3) .+ 1 for _ in 1:N])

# approximate the product between hybrid manifold densities
pq = manifoldProduct([p;q])

# direct histogram plot
using Gadfly
plot( x=getPoints(pq)[1,:], y=getPoints(pq)[2,:], Geom.histogram2d )

# TODO, convenient plotting (work in progress...)
ApproxManifoldProducts.mmdFunction
mmd(MF, a, b)
mmd(MF, a, b, N)
mmd(MF, a, b, N, M)
mmd(MF, a, b, N, M, threads; bw)

MMD disparity (i.e. 'distance') measure based on Kernel Hilbert Embeddings between two beliefs.

Notes:

  • This is a wrapper to the in-place mmd! function.

Related

mmd!, ker

ApproxManifoldProducts.mmd!Function
mmd!(MF, val, a, b)
mmd!(MF, val, a, b, N)
mmd!(MF, val, a, b, N, M)
mmd!(MF, val, a, b, N, M, threads; bw)

MMD disparity (i.e. 'distance') measure based on Kernel Hilbert Embeddings between two beliefs.

Notes:

  • This is the in-place version (well more in-place than mmd)

DevNotes:

  • TODO dont assume equally weighted particles
  • TODO profile SIMD vs SLEEF
  • TODO optimize memory
  • TODO make multithreaded

See also: mmd, ker

ApproxManifoldProducts.normDistAccAt!Method
normDistAccAt!(ret, idx, x, sigma)
normDistAccAt!(ret, idx, x, sigma, w)

Probability density function p(x), as estimated by kernels

\[hatp_{-j}(x) = 1/(N-1) Σ_{i != j}^N frac{1}{sqrt{2pi}σ } exp{ -frac{(x-μ)^2}{2 σ^2} }\]

ApproxManifoldProducts.productbeliefMethod
productbelief(
    denspts,
    manifold,
    dens,
    N;
    asPartial,
    dbg,
    logger
)

Take product of dens (including optional partials beliefs) as proposals to be multiplied together.

Notes

  • Return points of full dimension, even if only partial dimensions in proposals.
    • 'Other' dimensions left unchanged from incoming denspts
  • d dimensional product approximation
  • Incorporate ApproxManifoldProducts to process variables in individual batches.

DevNotes

KernelDensityEstimate.getPointsMethod
getPoints(x)
getPoints(x, )

Return underlying points used to construct the ManifoldKernelDensity.

Notes

  • Return type is ::Vector{P} where P represents a Manifold point type (e.g. group element or coordinates).
  • Second argument controls whether partial dimensions only should be returned (=true default).

DevNotes

  • Currently converts down to manifold from matrix of coordinates (legacy), to be deprecated TODO