Title: | Elastic Full Procrustes Means for Sparse and Irregular Planar Curves |
---|---|
Description: | Provides functions for the computation of functional elastic shape means over sets of open planar curves. The package is particularly suitable for settings where these curves are only sparsely and irregularly observed. It uses a novel approach for elastic shape mean estimation, where planar curves are treated as complex functions and a full Procrustes mean is estimated from the corresponding smoothed Hermitian covariance surface. This is combined with the methods for elastic mean estimation proposed in Steyer, Stöcker, Greven (2022) <doi:10.1111/biom.13706>. See Stöcker et. al. (2022) <arXiv:2203.10522> for details. |
Authors: | Manuel Pfeuffer [aut, cre], Lisa Steyer [aut], Almond Stoecker [aut] |
Maintainer: | Manuel Pfeuffer <[email protected]> |
License: | GPL (>= 3) |
Version: | 0.1.7.9000 |
Built: | 2025-02-18 05:43:13 UTC |
Source: | https://github.com/mpff/elastes |
Computes an elastic full Procrustes mean for curves stored in data_curves
.
Constructor function for class elastic_shape_mean
.
compute_elastic_shape_mean( data_curves, knots = seq(0, 1, len = 13), type = c("smooth", "polygon"), penalty = 2, var_type = c("smooth", "constant", "zero"), pfit_method = c("smooth", "polygon"), smooth_warp = function(i) 0, eps = 0.05, max_iter = 50, verbose = FALSE, cluster = NULL )
compute_elastic_shape_mean( data_curves, knots = seq(0, 1, len = 13), type = c("smooth", "polygon"), penalty = 2, var_type = c("smooth", "constant", "zero"), pfit_method = c("smooth", "polygon"), smooth_warp = function(i) 0, eps = 0.05, max_iter = 50, verbose = FALSE, cluster = NULL )
data_curves |
list of |
knots |
set of knots for the mean spline curve |
type |
if "smooth" linear srv-splines are used which results in a differentiable mean curve if "polygon" the mean will be piecewise linear. |
penalty |
the penalty to use in the covariance smoothing step. use '-1' for no penalty. |
var_type |
(experimental) assume "smooth", "constant" or "zero" measurement-error variance along t |
pfit_method |
(experimental) "smooth" or "polygon" |
smooth_warp |
(experimental) controls the weighting of original and smoothed observations over the iterations, if pfit_method == "smooth". |
eps |
the algorithm stops if L2 norm of coefficients changes by less than |
max_iter |
maximal number of iterations |
verbose |
print iterations |
cluster |
(experimental) use the parallel package for faster computation |
an object of class elastic_shape_mean
, which is a list
with entries
type |
"smooth" if mean was modeled using linear srv-splines, "polygon" if constant srv-splines |
coefs |
spline coefficients |
knots |
spline knots |
variance |
sample elastic shape variance |
data_curves |
list of |
fit |
see |
curve <- function(t){ rbind(t*cos(13*t), t*sin(13*t)) } set.seed(18) data_curves <- lapply(1:4, function(i){ m <- sample(10:15, 1) delta <- abs(rnorm(m, mean = 1, sd = 0.05)) t <- cumsum(delta)/sum(delta) data.frame(t(curve(t)) + 0.07*t*matrix(cumsum(rnorm(2*length(delta))), ncol = 2)) }) #randomly rotate and scale curves rand_scale <- function(curve){ ( 0.5 + runif(1) ) * curve } rand_rotate <- function(curve){ names <- colnames(curve) theta <- 2*pi*runif(1) mat <- matrix(c(cos(theta), sin(theta), -sin(theta), cos(theta)), nrow = 2, ncol = 2) curve.rot <- as.matrix(curve) %*% t(mat) curve.rot <- as.data.frame(curve.rot) colnames(curve.rot) <- names return(curve.rot) } data_curves <- lapply(data_curves, rand_scale) data_curves <- lapply(data_curves, rand_rotate) #compute smooth procrustes mean with 2 order penalty knots <- seq(0,1, length = 11) elastic_shape_mean <- compute_elastic_shape_mean( data_curves, knots = knots, type = "smooth", penalty = 2 ) plot(elastic_shape_mean)
curve <- function(t){ rbind(t*cos(13*t), t*sin(13*t)) } set.seed(18) data_curves <- lapply(1:4, function(i){ m <- sample(10:15, 1) delta <- abs(rnorm(m, mean = 1, sd = 0.05)) t <- cumsum(delta)/sum(delta) data.frame(t(curve(t)) + 0.07*t*matrix(cumsum(rnorm(2*length(delta))), ncol = 2)) }) #randomly rotate and scale curves rand_scale <- function(curve){ ( 0.5 + runif(1) ) * curve } rand_rotate <- function(curve){ names <- colnames(curve) theta <- 2*pi*runif(1) mat <- matrix(c(cos(theta), sin(theta), -sin(theta), cos(theta)), nrow = 2, ncol = 2) curve.rot <- as.matrix(curve) %*% t(mat) curve.rot <- as.data.frame(curve.rot) colnames(curve.rot) <- names return(curve.rot) } data_curves <- lapply(data_curves, rand_scale) data_curves <- lapply(data_curves, rand_rotate) #compute smooth procrustes mean with 2 order penalty knots <- seq(0,1, length = 11) elastic_shape_mean <- compute_elastic_shape_mean( data_curves, knots = knots, type = "smooth", penalty = 2 ) plot(elastic_shape_mean)
Finds optimal rotation and scaling alignment for a discrete open srv curve to a smooth curve
fit_alignment_proc2d( q, type, knots, var_type, coefs.compl, method, cov_fit, pca, L )
fit_alignment_proc2d( q, type, knots, var_type, coefs.compl, method, cov_fit, pca, L )
q |
complex srv curve with parametrization, needs to be vectorized.
The result of a call to |
type |
spline degree |
knots |
basis knots |
var_type |
either "smooth" or "constant" measurement error in cov_fit object |
coefs.compl |
complex coefficients of smooth curve |
method |
temp |
cov_fit |
temp |
pca |
temp |
L |
temp |
optimal rotation G and scaling b
Fits an elastic full Procrustes mean for open, planar curves.
Is usually called from compute_elastic_shape_mean
.
fit_mean( srv_data_curves, knots, penalty, var_type, pfit_method, max_iter, type, eps, cluster, verbose, smooth_warp )
fit_mean( srv_data_curves, knots, penalty, var_type, pfit_method, max_iter, type, eps, cluster, verbose, smooth_warp )
srv_data_curves |
list of |
knots |
set of knots for the mean spline curve |
penalty |
the penalty to use in the covariance smoothing step. use '-1' for no penalty. |
var_type |
(experimental) assume "smooth", "constant" or "zero" measurement-error variance along t |
pfit_method |
(experimental) "smooth" or "polygon" |
max_iter |
maximal number of iterations |
type |
if "smooth" linear srv-splines are used which results in a differentiable mean curve if "polygon" the mean will be piecewise linear. |
eps |
the algorithm stops if L2 norm of coefficients changes less |
cluster |
a cluster object for use in the |
verbose |
print iterations |
smooth_warp |
(experimental) controls the weighting of original and smoothed observations over the iterations, if pfit_method == "smooth". |
a list
with entries
type |
"smooth" or "polygon" |
coefs |
|
knots |
spline knots |
penalty |
penalty used in the covariance estimation |
distances |
distances to mean |
fit |
a |
Calculate the center of a curve
get_center(curve)
get_center(curve)
curve |
a |
The average of observed points in curve
.
Finds the distance of a discrete open srv curve to a smooth curve
get_distance(srv_curve, s, q, eps = 10 * .Machine$double.eps)
get_distance(srv_curve, s, q, eps = 10 * .Machine$double.eps)
srv_curve |
srv transformation of the smooth curve, needs to be vectorized |
s |
time points for q, first has to be 0, last has to be 1 |
q |
square root velocity vectors, one less than time points in s |
eps |
convergence tolerance |
distance between srv_curve and q
Evaluate a curve on a grid
get_evals(curve, t_grid = NULL, ...) ## S3 method for class 'data.frame' get_evals(curve, t_grid = NULL, ...) ## S3 method for class 'elastic_shape_mean' get_evals(curve, t_grid = NULL, centering = TRUE, srv = FALSE, ...)
get_evals(curve, t_grid = NULL, ...) ## S3 method for class 'data.frame' get_evals(curve, t_grid = NULL, ...) ## S3 method for class 'elastic_shape_mean' get_evals(curve, t_grid = NULL, centering = TRUE, srv = FALSE, ...)
curve |
a one parameter function which is to be evaluated on a grid |
t_grid |
the curve is evaluated at the values in t_grid, first value needs to be 0, last value needs to be 1. If t_grid = NULL, a default regular grid with grid length 0.01 is chosen |
... |
other arguments |
centering |
TRUE if curves shall be centered |
srv |
TRUE if SRV curve shall be evaluated |
a data.frame
with evaluations of the curve
at the values in t_grid
in its rows.
See get_evals
for the original code.
curve <- function(t){c(t*sin(10*t), t*cos(10*t))} plot(get_evals(curve), type = "b")
curve <- function(t){c(t*sin(10*t), t*cos(10*t))} plot(get_evals(curve), type = "b")
Finds optimal aligned time points for srv curve q to srv curve p using coordinate wise optimization.
get_optimal_t(srv_procrustes_curves, coefs, t_optims, type, knots, eps, i)
get_optimal_t(srv_procrustes_curves, coefs, t_optims, type, knots, eps, i)
srv_procrustes_curves |
scaling and rotation aligned srv curves |
coefs |
mean coefficients |
t_optims |
current optimal parametrization |
type |
"smooth" or "polygon" |
knots |
mean basis knots |
eps |
convergence tolerance |
i |
current iteration |
optimal time points for srv_data_curves, without first value 0 and last value 1 optimal time points have the distance of the observation to the srv_curve as an attribute
Calculate the polygon length of a curve
get_polygon_length(curve)
get_polygon_length(curve)
curve |
a |
The length of curve
, treating it as a polygon.
Compute the Procrustes aligned data curve...
get_procrustes_fit(data_curve)
get_procrustes_fit(data_curve)
data_curve |
A |
Aligned data_curve as a data.frame
.
Compute the Procrustes fit given optimal rotation, scaling and translation.
get_procrustes_fit_from_param( data_curve, rot, scale, plength, trans, norm_factor )
get_procrustes_fit_from_param( data_curve, rot, scale, plength, trans, norm_factor )
data_curve |
A |
rot |
The rotation (in radian). |
scale |
The scaling. |
plength |
The polygon length of the original curve. |
trans |
The translation. |
norm_factor |
The normalization factor from the smooth curve estimate. |
Plots objects of class elastic_shape_mean
.
## S3 method for class 'elastic_shape_mean' plot(x, srv = FALSE, centering = TRUE, asp = 1, col = "red", ...)
## S3 method for class 'elastic_shape_mean' plot(x, srv = FALSE, centering = TRUE, asp = 1, col = "red", ...)
x |
object of class |
srv |
TRUE if the SRV curve should be plotted |
centering |
TRUE if mean and pfits should be centered |
asp |
numeric, giving the aspect ratio of the two coordinates,
see |
col |
color of the mean curve. |
... |
further plotting parameters. |
No return value, called for side effects.
For examples see documentation of compute_elastic_shape_mean
.
See plot.elastic_mean
for the original code.