R/MvtGLomax.R

Defines functions dmvglomax pmvglomax qmvglomax rmvglomax smvglomax

Documented in dmvglomax pmvglomax qmvglomax rmvglomax smvglomax

# ---------- Generalized Multivariate Lomax ---------------
#' Generalized Multivariate Lomax (Pareto Type II) Distribution
#' @name MvtGlomx
#' @description Calculation of density function, cumulative distribution function, equicoordinate quantile function and survival function, and random numbers generation for generalized multivariate Lomax distribution with a scalar parameter \code{parm1} and vectors of parameters \code{parm2} and \code{parm3}.
#' @param x vector or matrix of quantiles. If \eqn{x} is a matrix, each row vector constitutes a vector of quantiles for which the density \eqn{f(x)} is calculated (for \eqn{i}-th row \eqn{x_i}, \eqn{f(x_i)} is reported).
#' @param parm1 a scalar parameter, see parameter \eqn{a} in \strong{Details}.
#' @param parm2 a vector of parameters, see parameters \eqn{\theta_i} in \strong{Details}.
#' @param parm3 a vector of parameters, see parameters \eqn{l_i} in \strong{Details}.
#' @param k dimension of data or number of variates.
#' @param log logical; if TRUE, probability densities \eqn{f} are given as \eqn{log(f)}.
#'
#' @return \code{dmvglomax} gives the numerical values of the probability density.
#'
#' @examples
#' # Calculations for the generalized multivariate Lomax with parameters:
#' # a = 5, theta1 = 1, theta2 = 2, l1 = 4, l2 = 5
#' # Vector of quantiles: c(5, 6)
#'
#' dmvglomax(x = c(5, 6), parm1 = 5, parm2 = c(1, 2), parm3 = c(4, 5)) # Density
#'
#' @details
#' Generalized multivariate Lomax (Pareto type II) distribution was introduced by Nayak (1987) as a joint probability distribution of several skewed nonnegative random variables \eqn{X_1, X_2, \cdots, X_k}. Its probability density function is given by
#' \deqn{f(x_1, \cdots, x_k) = \frac{[ \prod_{i=1}^{k} \theta_i^{l_i}] \Gamma(\sum_{i=1}^{k} l_i + a ) \prod_{i=1}^{k} x_i^{l_i-1}}{\Gamma(a)[ \prod_{i=1}^{k} \Gamma(l_i)] (1+\sum_{i=1}^{k} \theta_i x_i )^{\sum_{i=1}^{k} l_i + a}},}
#' where \eqn{x_i>0, a,\theta_i, l_i>0, i=1,\cdots,k}.
#'
#' Cumulative distribution function \eqn{F(x_1, \dots, x_k)} is obtained by multiple integral
#' \deqn{F(x_1, \dots, x_k) = \int_{0}^{x_1} \cdots  \int_{0}^{x_k} f(y_1, \cdots, y_k) dy_k \cdots dy_1.}
#' This multiple integral is calculated by either adaptive multivariate integration using \code{\link{hcubature}} in package \strong{\link{cubature}} (Narasimhan et al., 2018) or via Monte Carlo method.
#'
#' Equicoordinate quantile is obtained by solving the following equation for \eqn{q} through the built-in one dimension root finding function \code{\link{uniroot}}:
#' \deqn{\int_{0}^{q} \cdots \int_{0}^{q} f(x_1, \cdots, x_k) dx_k \cdots dx_1 = p,}
#' where \eqn{p} is the given cumulative probability.
#'
#' The survival function \eqn{\bar{F}(x_1, \cdots, x_k)} is obtained either by the following formula related to cumulative distribution function \eqn{F(x_1, \dots, x_k)} (Joe, 1997)
#' \deqn{\bar{F}(x_1, \cdots, x_k) = 1 + \sum_{S \in \mathcal{S}} (-1)^{|S|} F_S(x_j, j \in S),}
#' or via Monte Carlo method.
#'
#' Random numbers from generalized multivariate Lomax distribution can be generated by simulating \eqn{k} independent gamma random variables having a common parameter following gamma distribution with shape parameter \eqn{a} and scale parameter \eqn{1}; see Nayak (1987).
#'
#' @seealso \code{\link{uniroot}} for one dimensional root (zero) finding.
#'
#' @references
#' Joe, H. (1997). \emph{Multivariate Models and Dependence Concepts}. London: Chapman & Hall.
#'
#' Narasimhan, B.,  Koller, M., Johnson, S. G., Hahn, T., Bouvier, A., Kiêu, K. and Gaure, S. (2018). cubature: Adaptive Multivariate Integration over Hypercubes. R package version 2.0.3.
#'
#' Nayak, T. K. (1987). Multivariate Lomax Distribution: Properties and Usefulness in Reliability Theory. \emph{Journal of Applied Probability}, Vol. 24, No. 1, 170-177.
#'
#' @importFrom stats rgamma
#' @export
dmvglomax <- function(x, parm1 = 1, parm2 = rep(1, k), parm3 = rep(1, k), log = FALSE) {
  dfun <- function(x) {
    if (any(x < 0)) {
      return(-Inf)
    }
    else {
      dlog <- sum(parm3 * log(parm2)) + lgamma(sum(parm3) + parm1) +
                  log(x) %*% (parm3 - 1)  - lgamma(parm1) - sum(lgamma(parm3)) -
                  (sum(parm3) + parm1) * log(1 + x %*% parm2)
      return(dlog)
    }
  }
  if (!is.vector(x, mode = "numeric") && !is.matrix(x)) {
    stop(sQuote("x"), " must be a vector or matrix of quantiles")
  }
  if (is.vector(x)) {
    x <- matrix(x, ncol = length(x))
  }
  k <- ncol(x)
  if (length(parm1) != 1) {
    stop(sQuote("parm1"), " must be a scalar")
  }
  if (parm1 <= 0) {
    stop(sQuote("parm1"), " must be a positive number")
  }
  if (any(parm2 <= 0) || !is.vector(parm2, mode = "numeric")) {
    stop(sQuote("parm2"), " must be a vector of positive numbers")
  }
  if (any(parm3 <= 0) || !is.vector(parm3, mode = "numeric")) {
    stop(sQuote("parm3"), " must be a vector of positive numbers")
  }
  if (k != length(parm2)) {
    stop("x and parm2 have non-conforming size")
  }
  if (k != length(parm3)) {
    stop("x and parm3 have non-conforming size")
  }

  logretval <- apply(x, 1, dfun)

  if (log)
    logretval
  else exp(logretval)
}

#' @rdname MvtGlomx
#' @param q a vector of quantiles.
#' @param algorithm method to be used for calculating cumulative probability. Two options are provided as (i) \code{numerical} using adaptive multivariate integral and (ii) \code{MC} using Monte Carlo method. Recommend algorithm \code{numerical} for \eqn{(k <= 4)} dimension and \code{MC} for \eqn{(k > 4)} dimension based on running time consumption. Default option is set as \code{numerical}.
#' @param nsim number of simulations used in algorithm \code{MC}.
#' @return \code{pmvglomax} gives a list of two items:
#'
#' \eqn{\quad} \code{value} cumulative probability
#'
#' \eqn{\quad} \code{error} the estimated relative error by \code{algorithm = "numerical"} or the estimated standard error by \code{algorithm = "MC"}
#' @importFrom cubature hcubature
#' @importFrom stats sd
#' @examples
#'
#' # Cumulative Probability using adaptive multivariate integral
#' pmvglomax(q = c(5, 6), parm1 = 5, parm2 = c(1, 2), parm3 = c(4, 5))
#'
#' \donttest{
#' # Cumulative Probability using Monte Carlo method
#' pmvglomax(q = c(5, 6), parm1 = 5, parm2 = c(1, 2), parm3 = c(4, 5), algorithm = "MC")}
#'
#' @export
pmvglomax <- function(q, parm1 = 1, parm2 = rep(1, k), parm3 = rep(1, k), algorithm = c("numerical", "MC"), nsim = 1e7) {
  if (any(q < 0) ||  !is.vector(q, mode = "numeric")) {
    stop(sQuote("q"), " must be a vector of non-negative quantiles")
  }
  k <- length(q)
  if (length(parm1) != 1) {
    stop(sQuote("parm1"), " must be a scalar")
  }
  if (parm1 <= 0) {
    stop(sQuote("parm1"), " must be a positive number")
  }
  if (any(parm2 <= 0) || !is.vector(parm2, mode = "numeric")) {
    stop(sQuote("parm2"), " must be a vector of positive numbers")
  }
  if (any(parm3 <= 0) || !is.vector(parm3, mode = "numeric")) {
    stop(sQuote("parm3"), " must be a vector of positive numbers")
  }
  if (k != length(parm2)) {
    stop("q and parm2 have non-conforming size")
  }
  if (k != length(parm3)) {
    stop("q and parm3 have non-conforming size")
  }

  method <- match.arg(algorithm, c("numerical", "MC"))

  if (method == "numerical") {
    CDF <- hcubature(dmvglomax, lowerLimit=rep(0, k), upperLimit = rep(1, k), parm1 = parm1, parm2 = q * parm2, parm3 = parm3)
    return(list(value = as.vector(CDF$integral), error=as.vector(CDF$error)))
  }
  if (method == "MC") {
	if (nsim %% 1 != 0 || nsim <= 0){
		stop(sQuote("nsim"), " must be a positive integer")
	}
    mcSample <- rmvglomax(n = nsim, parm1 = parm1, parm2 = parm2, parm3 = parm3)
    trueTable <- sweep(mcSample, 2, q, FUN = "<")
    numTable <- 1 * trueTable
    MC.vec <- rep(1, nsim)

    for (i in 1:k) {
      MC.vec <- MC.vec * as.vector(numTable[,i])
    }

    prob = mean(MC.vec)
  }

  return(list(value = prob, error = sd(MC.vec)/sqrt(nsim)))
}

#' @rdname MvtGlomx
#' @param p a scalar value corresponding to probability.
#' @param interval a vector containing the end-points of the interval to be searched. Default value is set as \code{c(1e-8, 1e8)}.
#' @return \code{qmvglomax} gives the equicoordinate quantile. \code{NaN} is returned for no solution found in the given interval. The result is seed dependent if Monte Carlo algorithm is chosen  (\code{algorithm = "MC"}).
#' @examples
#' \donttest{
#' # Equicoordinate quantile of cumulative probability 0.5
#' qmvglomax(p = 0.5, parm1 = 5, parm2 = c(1, 2), parm3 = c(4, 5))}
#'
#' @importFrom stats uniroot
#' @export
qmvglomax <- function(p, parm1 = 1, parm2 = rep(1, k), parm3 = rep(1, k), interval = c(1e-8, 1e8), algorithm = c("numerical", "MC"), nsim = 1e6) {

  if (length(p) != 1 || p <= 0 || p >= 1) {
    stop(sQuote("p"), " is not a double between zero and one")
  }
  if (length(parm1) != 1) {
    stop(sQuote("parm1"), " must be a scalar")
  }
  if (parm1 <= 0) {
    stop(sQuote("parm1"), " must be a positive number")
  }
  if (missing(parm2)) {
    stop(sQuote("parm2"), " is a required argument")
  }
  if (any(parm2 <= 0) || !is.vector(parm2, mode = "numeric")) {
    stop(sQuote("parm2"), " must be a vector of positive numbers")
  }
  k <- length(parm2)
  if (any(parm3 <= 0) || !is.vector(parm3, mode = "numeric")) {
    stop(sQuote("parm3"), " must be a vector of positive numbers")
  }
  if (k != length(parm3)) {
    stop("parm2 and parm3 have non-conforming size")
  }
  if (!is.vector(interval) || length(interval) != 2 || interval[1] > interval[2]) {
    stop(sQuote("interval"), " must be a vector of valid range")
  }

  method <- match.arg(algorithm, c("numerical", "MC"))

  rootfun <- function(x) {
    equi.x <- rep(x, k)
    val <- pmvglomax(equi.x, parm1 = parm1, parm2 = parm2, parm3 = parm3,
                     algorithm = method, nsim = nsim)$value - p
  }

  find.root <- tryCatch(uniroot(rootfun, interval = interval), error = function(e) e)

  if (inherits(find.root, "error")) {
    if (find.root$message == "f() values at end points not of opposite sign")
      print("No quantile was found in provided interval")
      return(NaN)
  }

  return(find.root$root)

}

#' @rdname MvtGlomx
#' @param n number of observations.
#' @return \code{rmvglomax} generates random numbers.
#' @examples
#' # Random numbers generation with sample size 100
#' rmvglomax(n = 100, parm1 = 5, parm2 = c(1, 2), parm3 = c(4, 5))
#'
#' @export
rmvglomax <- function(n, parm1 = 1, parm2 = rep(1, k), parm3 = rep(1, k)) {
  if (!is.numeric(n) || n <= 0 || n %% 1 != 0 || length(n) != 1 ) {
    stop("sample size n must be a positive integer")
  }
  if (length(parm1) != 1) {
    stop(sQuote("parm1"), " must be a scalar")
  }
  if (parm1 <= 0) {
    stop(sQuote("parm1"), " must be a positive number")
  }
  if (missing(parm2)) {
    stop(sQuote("parm2"), " is a required argument")
  }
  if (any(parm2 <= 0) || !is.vector(parm2, mode = "numeric")) {
    stop(sQuote("parm2"), " must be a vector of positive numbers")
  }
  k <- length(parm2)
  if (any(parm3 <= 0) || !is.vector(parm3, mode = "numeric")) {
    stop(sQuote("parm3"), " must be a vector of positive numbers")
  }
  if (k != length(parm3)) {
    stop("parm2 and parm3 have non-conforming size")
  }

  retval <- matrix(0, nrow = n, ncol = k)
  eta <- rgamma(n, shape = parm1, rate = 1)
  for (i in 1:k) {
    retval[, i] <- rgamma(n, shape = parm3[i], rate = eta * parm2[i])
  }
  return(retval)
}

#' @rdname MvtGlomx
#' @return \code{smvglomax} gives a list of two items:
#'
#' \eqn{\quad} \code{value} the value of survial function
#'
#' \eqn{\quad} \code{error} the estimated relative error by \code{algorithm = "numerical"} or the estimated standard error by \code{algorithm = "MC"}
#' @examples
#' smvglomax(q = c(5, 6), parm1 = 5, parm2 = c(1, 2), parm3 = c(4, 5)) # Survival function
#'
#' @importFrom stats sd
#'
#' @export
smvglomax <- function(q, parm1 = 1, parm2 = rep(1, k), parm3 = rep(1, k), algorithm = c("numerical", "MC"), nsim = 1e7) {
  if (any(q < 0) ||  !is.vector(q, mode = "numeric")) {
    stop(sQuote("q"), " must be a vector of non-negative quantiles")
  }
  k <- length(q)
  if (length(parm1) != 1) {
    stop(sQuote("parm1"), " must be a scalar")
  }
  if (parm1 <= 0) {
    stop(sQuote("parm1"), " must be a positive number")
  }
  if (any(parm2 <= 0) || !is.vector(parm2, mode = "numeric")) {
    stop(sQuote("parm2"), " must be a vector of positive numbers")
  }
  if (any(parm3 <= 0) || !is.vector(parm3, mode = "numeric")) {
    stop(sQuote("parm3"), " must be a vector of positive numbers")
  }
  if (k != length(parm2)) {
    stop("q and parm2 have non-conforming size")
  }
  if (k != length(parm3)) {
    stop("q and parm3 have non-conforming size")
  }

  method <- match.arg(algorithm, c("numerical", "MC"))

  if (method == "numerical") {
    SF <- 1
    SF.error <- 0
    Index <- seq(1, k)
    for (i in 1:k) {
      combn.Index <- as.matrix(combn(Index, i), length(combn(Index, i)))
      card <- ncol(combn.Index)
      coef.One <- rep((-1)^i, card)
      sub.cdf <- rep(0, card)
      sub.error <- rep(0, card)
      for (j in 1:card) {
        sub.index <- combn.Index[,j]
        sub.q <- as.vector(q[sub.index])
        sub.parm2 <- parm2[sub.index]
        sub.parm3 <- parm3[sub.index]
        sub.calc <- pmvglomax(sub.q, parm1 = parm1, parm2 = sub.parm2, parm3 = sub.parm3, algorithm = method)
        sub.cdf[j] <- sub.calc$value
        sub.error[j] <- sub.calc$error
      }
      SF <- SF + sum(coef.One * sub.cdf)
      SF.error <- SF.error + sum(sub.error)
    }
    return(list(value = SF, error = SF.error))
  }

  if (method == "MC") {
    if (nsim %% 1 != 0 || nsim <= 0){
      stop(sQuote("nsim"), " must be a positive integer")
    }
    mcSample <- rmvglomax(n = nsim, parm1 = parm1, parm2 = parm2, parm3 = parm3)
    trueTable <- sweep(mcSample, 2, q, FUN = ">")
    numTable <- 1 * trueTable
    MC.vec <- rep(1, nsim)

    for (i in 1:k) {
      MC.vec <- MC.vec * as.vector(numTable[,i])
    }

    SF = mean(MC.vec)
    return(list(value = SF, error = sd(MC.vec)/sqrt(nsim)))
  }
}

Try the NonNorMvtDist package in your browser

Any scripts or data that you put into this service are public.

NonNorMvtDist documentation built on March 23, 2020, 5:06 p.m.