bvm {ridgetorus} | R Documentation |
Density evaluation, sampling, and parameter estimation of the bivariate sine von Mises distribution
Description
Computation of the density and normalizing constant
of the bivariate sine von Mises
Simulation of samples from a bivariate sine von Mises.
Maximum likelihood and method of moments estimation of the
parameters .
Usage
d_bvm(x, mu, kappa, log_const = NULL)
const_bvm(kappa, M = 25, MC = 10000)
r_bvm(n, mu, kappa)
fit_bvm_mm(
x,
lower = c(0, 0, -30),
upper = c(30, 30, 30),
start = NULL,
M = 25,
hom = FALSE,
indep = FALSE,
...
)
fit_bvm_mle(
x,
start = NULL,
M = 25,
lower = c(-pi, -pi, 0, 0, -30),
upper = c(pi, pi, 30, 30, 30),
hom = FALSE,
indep = FALSE,
...
)
Arguments
x |
matrix of size |
mu |
circular means of the density, a vector of length |
kappa |
vector of length |
log_const |
logarithm of the normalizing constant. Computed internally
if |
M |
truncation of the series expansion for computing the normalizing
constant. Defaults to |
MC |
Monte Carlo replicates for computing the normalizing
constant when there is no series expansion. Defaults to |
n |
sample size. |
lower , upper |
vectors of length |
start |
a vector of length |
hom |
assume a homogeneous distribution with equal marginal
concentrations? Defaults to |
indep |
set the dependence parameter to zero? Defaults to |
... |
further parameters passed to
|
Value
-
d_bvm
: a vector of lengthnx
with the density evaluated atx
. -
const_bvm
: the value of the normalizing constant.
-
r_bvm
: a matrix of sizec(n, 2)
with the random sample. -
fit_mme_bvm, fit_mle_bvm
: a list with the estimated parametersand the object
opt
containing the optimization summary.
References
Mardia, K. V., Hughes, G., Taylor, C. C., and Singh, H. (2008). A multivariate von Mises with applications to bioinformatics. Canadian Journal of Statistics, 36(1):99–109. doi:10.1002/cjs.5550360110
Singh, H., Hnizdo, V., and Demchuk, E. (2002). Probabilistic model for two dependent circular variables. Biometrika, 89(3):719–723. doi:10.1093/biomet/89.3.719
Examples
## Density evaluation
mu <- c(0, 0)
kappa <- 3:1
nth <- 50
th <- seq(-pi, pi, l = nth)
x <- as.matrix(expand.grid(th, th))
const <- const_bvm(kappa = kappa)
d <- d_bvm(x = x, mu = mu, kappa = kappa, log_const = log(const))
filled.contour(th, th, matrix(d, nth, nth), col = viridisLite::viridis(31),
levels = seq(0, max(d), l = 30))
## Sampling and estimation
n <- 100
samp <- r_bvm(n = n, mu = mu, kappa = kappa)
(param_mm <- fit_bvm_mm(samp)$par)
(param_mle <- fit_bvm_mle(samp)$par)