rbase.mean {RiemBase}R Documentation

Fréchet Mean of Manifold-valued Data

Description

For manifold-valued data, Fréchet mean is the solution of following cost function,

\textrm{min}_x \sum_{i=1}^n \rho^2 (x, x_i),\quad x\in\mathcal{M}

for a given data \{x_i\}_{i=1}^n and \rho(x,y) is the geodesic distance between two points on manifold \mathcal{M}. It uses a gradient descent method with a backtracking search rule for updating.

Usage

rbase.mean(input, maxiter = 496, eps = 1e-06, parallel = FALSE)

Arguments

input

a S3 object of riemdata class. See riemfactory for more details.

maxiter

maximum number of iterations for gradient descent algorithm.

eps

stopping criterion for the norm of gradient.

parallel

a flag for enabling parallel computation.

Value

a named list containing

x

an estimate Fréchet mean.

iteration

number of iterations until convergence.

Author(s)

Kisung You

References

Karcher H (1977). “Riemannian center of mass and mollifier smoothing.” Communications on Pure and Applied Mathematics, 30(5), 509–541. ISSN 00103640, 10970312.

Kendall WS (1990). “Probability, Convexity, and Harmonic Maps with Small Image I: Uniqueness and Fine Existence.” Proceedings of the London Mathematical Society, s3-61(2), 371–406. ISSN 00246115.

Afsari B, Tron R, Vidal R (2013). “On the Convergence of Gradient Descent for Finding the Riemannian Center of Mass.” SIAM Journal on Control and Optimization, 51(3), 2230–2260. ISSN 0363-0129, 1095-7138.

Examples


### Generate 100 data points on Sphere S^2 near (0,0,1).
ndata = 100
theta = seq(from=-0.99,to=0.99,length.out=ndata)*pi
tmpx  = cos(theta) + rnorm(ndata,sd=0.1)
tmpy  = sin(theta) + rnorm(ndata,sd=0.1)

### Wrap it as 'riemdata' class
data  = list()
for (i in 1:ndata){
  tgt = c(tmpx[i],tmpy[i],1)
  data[[i]] = tgt/sqrt(sum(tgt^2)) # project onto Sphere
}
data = riemfactory(data, name="sphere")

### Compute Fréchet Mean
out1 = rbase.mean(data)
out2 = rbase.mean(data,parallel=TRUE) # test parallel implementation



[Package RiemBase version 0.2.5 Index]