| Title: | Distributions Compatible with Automatic Differentiation by 'RTMB' |
|---|---|
| Description: | Extends the functionality of the 'RTMB' <https://kaskr.r-universe.dev/RTMB> package by providing a collection of non-standard probability distributions compatible with automatic differentiation (AD). While 'RTMB' enables flexible and efficient modelling, including random effects, its built-in support is limited to standard distributions. The package adds additional AD-compatible distributions, broadening the range of models that can be implemented and estimated using 'RTMB'. Automatic differentiation and Laplace approximation are described in Kristensen et al. (2016) <doi:10.18637/jss.v070.i05>. |
| Authors: | Jan-Ole Fischer [aut, cre] (ORCID: <https://orcid.org/0009-0004-1556-9053>) |
| Maintainer: | Jan-Ole Fischer <[email protected]> |
| License: | MIT + file LICENSE |
| Version: | 1.0.4 |
| Built: | 2026-05-28 10:59:38 UTC |
| Source: | https://github.com/janolefi/rtmbdist |
Smooth approximation to the absolute value function
abs_smooth(x, epsilon = 1e-06)abs_smooth(x, epsilon = 1e-06)
x |
vector of evaluation points |
epsilon |
smoothing constant |
We approximate the absolute value here as
Smooth absolute value of x.
abs(0) abs_smooth(0, 1e-4)abs(0) abs_smooth(0, 1e-4)
Density, distribution function, quantile function, and random generation for the Box–Cox Cole and Green distribution.
dbccg(x, mu = 1, sigma = 0.1, nu = 1, log = FALSE) pbccg(q, mu = 1, sigma = 0.1, nu = 1, lower.tail = TRUE, log.p = FALSE) qbccg(p, mu = 1, sigma = 0.1, nu = 1, lower.tail = TRUE, log.p = FALSE) rbccg(n, mu = 1, sigma = 0.1, nu = 1)dbccg(x, mu = 1, sigma = 0.1, nu = 1, log = FALSE) pbccg(q, mu = 1, sigma = 0.1, nu = 1, lower.tail = TRUE, log.p = FALSE) qbccg(p, mu = 1, sigma = 0.1, nu = 1, lower.tail = TRUE, log.p = FALSE) rbccg(n, mu = 1, sigma = 0.1, nu = 1)
x, q
|
vector of quantiles |
mu |
location parameter, must be positive. |
sigma |
scale parameter, must be positive. |
nu |
skewness parameter (real). |
log, log.p
|
logical; if |
lower.tail |
logical; if |
p |
vector of probabilities |
n |
number of random values to return |
This implementation of dbccg and pbccg allows for automatic differentiation with RTMB while the other functions are imported from gamlss.dist package.
See gamlss.dist::BCCG for more details.
The density is
where for and for , and is the standard normal CDF.
dbccg gives the density, pbccg gives the distribution function, qbccg gives the quantile function, and rbccg generates random deviates.
Rigby, R. A., Stasinopoulos, D. M., Heller, G. Z., and De Bastiani, F. (2019) Distributions for modeling location, scale, and shape: Using GAMLSS in R, Chapman and Hall/CRC, doi:10.1201/9780429298547. An older version can be found in https://www.gamlss.com/.
x <- rbccg(5, mu = 10, sigma = 0.2, nu = 0.5) d <- dbccg(x, mu = 10, sigma = 0.2, nu = 0.5) p <- pbccg(x, mu = 10, sigma = 0.2, nu = 0.5) q <- qbccg(p, mu = 10, sigma = 0.2, nu = 0.5)x <- rbccg(5, mu = 10, sigma = 0.2, nu = 0.5) d <- dbccg(x, mu = 10, sigma = 0.2, nu = 0.5) p <- pbccg(x, mu = 10, sigma = 0.2, nu = 0.5) q <- qbccg(p, mu = 10, sigma = 0.2, nu = 0.5)
Density, distribution function, quantile function, and random generation for the Box-Cox Power Exponential distribution.
dbcpe(x, mu = 5, sigma = 0.1, nu = 1, tau = 2, log = FALSE) pbcpe(q, mu = 5, sigma = 0.1, nu = 1, tau = 2, lower.tail = TRUE, log.p = FALSE) qbcpe(p, mu = 5, sigma = 0.1, nu = 1, tau = 2, lower.tail = TRUE, log.p = FALSE) rbcpe(n, mu = 5, sigma = 0.1, nu = 1, tau = 2)dbcpe(x, mu = 5, sigma = 0.1, nu = 1, tau = 2, log = FALSE) pbcpe(q, mu = 5, sigma = 0.1, nu = 1, tau = 2, lower.tail = TRUE, log.p = FALSE) qbcpe(p, mu = 5, sigma = 0.1, nu = 1, tau = 2, lower.tail = TRUE, log.p = FALSE) rbcpe(n, mu = 5, sigma = 0.1, nu = 1, tau = 2)
x, q
|
vector of quantiles |
mu |
location parameter, must be positive. |
sigma |
scale parameter, must be positive. |
nu |
vector of |
tau |
vector of |
log, log.p
|
logical; if |
lower.tail |
logical; if |
p |
vector of probabilities |
n |
number of random values to return |
This implementation of dbcpe and pbcpe allows for automatic differentiation with RTMB while the other functions are imported from gamlss.dist package.
See gamlss.dist::BCPE for more details.
The density is
where for and for , and and are the PDF and CDF of the power exponential (PE) distribution with shape .
dbcpe gives the density, pbcpe gives the distribution function, qbcpe gives the quantile function, and rbcpe generates random deviates.
Rigby, R. A., Stasinopoulos, D. M., Heller, G. Z., and De Bastiani, F. (2019) Distributions for modeling location, scale, and shape: Using GAMLSS in R, Chapman and Hall/CRC, doi:10.1201/9780429298547. An older version can be found in https://www.gamlss.com/.
x <- rbcpe(1, mu = 5, sigma = 0.1, nu = 1, tau = 1) d <- dbcpe(x, mu = 5, sigma = 0.1, nu = 1, tau = 1) p <- pbcpe(x, mu = 5, sigma = 0.1, nu = 1, tau = 1) q <- qbcpe(p, mu = 5, sigma = 0.1, nu = 1, tau = 1)x <- rbcpe(1, mu = 5, sigma = 0.1, nu = 1, tau = 1) d <- dbcpe(x, mu = 5, sigma = 0.1, nu = 1, tau = 1) p <- pbcpe(x, mu = 5, sigma = 0.1, nu = 1, tau = 1) q <- qbcpe(p, mu = 5, sigma = 0.1, nu = 1, tau = 1)
Density, distribution function, quantile function, and random generation for the Box–Cox t distribution.
dbct(x, mu = 5, sigma = 0.1, nu = 1, tau = 2, log = FALSE) pbct(q, mu = 5, sigma = 0.1, nu = 1, tau = 2, lower.tail = TRUE, log.p = FALSE) qbct(p, mu = 5, sigma = 0.1, nu = 1, tau = 2, lower.tail = TRUE, log.p = FALSE) rbct(n, mu = 5, sigma = 0.1, nu = 1, tau = 2)dbct(x, mu = 5, sigma = 0.1, nu = 1, tau = 2, log = FALSE) pbct(q, mu = 5, sigma = 0.1, nu = 1, tau = 2, lower.tail = TRUE, log.p = FALSE) qbct(p, mu = 5, sigma = 0.1, nu = 1, tau = 2, lower.tail = TRUE, log.p = FALSE) rbct(n, mu = 5, sigma = 0.1, nu = 1, tau = 2)
x, q
|
vector of quantiles |
mu |
location parameter, must be positive. |
sigma |
scale parameter, must be positive. |
nu |
skewness parameter (real). |
tau |
degrees of freedom, must be positive. |
log, log.p
|
logical; if |
lower.tail |
logical; if |
p |
vector of probabilities |
n |
number of random values to return |
This implementation of dbct and pbct allows for automatic differentiation with RTMB while the other functions are imported from gamlss.dist package.
See gamlss.dist::BCT for more details.
The density is
where for and for , and and are the PDF and CDF of Student's distribution with degrees of freedom.
dbct gives the density, pbct gives the distribution function, qbct gives the quantile function, and rbct generates random deviates.
Rigby, R. A., Stasinopoulos, D. M., Heller, G. Z., and De Bastiani, F. (2019) Distributions for modeling location, scale, and shape: Using GAMLSS in R, Chapman and Hall/CRC, doi:10.1201/9780429298547. An older version can be found in https://www.gamlss.com/.
x <- rbct(1, mu = 10, sigma = 0.2, nu = 0.5, tau = 4) d <- dbct(x, mu = 10, sigma = 0.2, nu = 0.5, tau = 4) p <- pbct(x, mu = 10, sigma = 0.2, nu = 0.5, tau = 4) q <- qbct(p, mu = 10, sigma = 0.2, nu = 0.5, tau = 4)x <- rbct(1, mu = 10, sigma = 0.2, nu = 0.5, tau = 4) d <- dbct(x, mu = 10, sigma = 0.2, nu = 0.5, tau = 4) p <- pbct(x, mu = 10, sigma = 0.2, nu = 0.5, tau = 4) q <- qbct(p, mu = 10, sigma = 0.2, nu = 0.5, tau = 4)
Density, distribution function, quantile function, and random generation for the beta distribution reparameterised in terms of mean and concentration.
dbeta(x, shape1, shape2, log = FALSE, eps = 0) dbeta2(x, mu, phi, log = FALSE, eps = 0) pbeta2(q, mu, phi, lower.tail = TRUE, log.p = FALSE) qbeta2(p, mu, phi, lower.tail = TRUE, log.p = FALSE) rbeta2(n, mu, phi)dbeta(x, shape1, shape2, log = FALSE, eps = 0) dbeta2(x, mu, phi, log = FALSE, eps = 0) pbeta2(q, mu, phi, lower.tail = TRUE, log.p = FALSE) qbeta2(p, mu, phi, lower.tail = TRUE, log.p = FALSE) rbeta2(n, mu, phi)
x, q
|
vector of quantiles |
shape1, shape2
|
non-negative parameters |
log, log.p
|
logical; if |
eps |
for internal use only, don't change. |
mu |
mean parameter, must be in the interval from 0 to 1. |
phi |
concentration parameter, must be positive. |
lower.tail |
logical; if |
p |
vector of probabilities |
n |
number of random values to return. |
This implementation allows for automatic differentiation with RTMB.
Currently, dbeta masks RTMB::dbeta because the latter has a numerically unstable gradient.
dbeta and dbeta2 compute the beta density. dbeta2 reparameterises by
mean and concentration :
dbeta2 gives the density, pbeta2 gives the distribution function, qbeta2 gives the quantile function, and rbeta2 generates random deviates.
set.seed(123) x <- rbeta2(1, 0.5, 1) d <- dbeta2(x, 0.5, 1) p <- pbeta2(x, 0.5, 1) q <- qbeta2(p, 0.5, 1)set.seed(123) x <- rbeta2(1, 0.5, 1) d <- dbeta2(x, 0.5, 1) p <- pbeta2(x, 0.5, 1) q <- qbeta2(p, 0.5, 1)
Density and random generation for the beta-binomial distribution.
dbetabinom(x, size, shape1, shape2, log = FALSE) rbetabinom(n, size, shape1, shape2)dbetabinom(x, size, shape1, shape2, log = FALSE) rbetabinom(n, size, shape1, shape2)
x |
vector of non-negative counts. |
size |
vector of total counts (number of trials). Needs to be >= |
shape1 |
positive shape parameter 1 of the Beta prior. |
shape2 |
positive shape parameter 2 of the Beta prior. |
log |
logical; if |
n |
number of random values to return (for |
This implementation of dbetabinom allows for automatic differentiation with RTMB.
dbetabinom gives the density and rbetabinom generates random samples.
set.seed(123) x <- rbetabinom(1, 10, 2, 5) d <- dbetabinom(x, 10, 2, 5)set.seed(123) x <- rbetabinom(1, 10, 2, 5) d <- dbetabinom(x, 10, 2, 5)
Density, distribution function, quantile function, and random generation for the Beta prime distribution.
dbetaprime(x, shape1, shape2, log = FALSE) pbetaprime(q, shape1, shape2, lower.tail = TRUE, log.p = FALSE) qbetaprime(p, shape1, shape2, lower.tail = TRUE, log.p = FALSE) rbetaprime(n, shape1, shape2)dbetaprime(x, shape1, shape2, log = FALSE) pbetaprime(q, shape1, shape2, lower.tail = TRUE, log.p = FALSE) qbetaprime(p, shape1, shape2, lower.tail = TRUE, log.p = FALSE) rbetaprime(n, shape1, shape2)
x, q
|
vector of quantiles |
shape1, shape2
|
non-negative shape parameters of the corresponding Beta distribution |
log, log.p
|
logical; if |
lower.tail |
logical; if |
p |
vector of probabilities |
n |
number of random values to return. |
This implementation allows for automatic differentiation with RTMB.
If , then
dbetaprime gives the density, pbetaprime gives the distribution function, qbetaprime gives the quantile function, and rbetaprime generates random deviates.
x <- rbetaprime(1, 2, 1) d <- dbetaprime(x, 2, 1) p <- pbetaprime(x, 2, 1) q <- qbetaprime(p, 2, 1)x <- rbetaprime(1, 2, 1) d <- dbetaprime(x, 2, 1) p <- pbetaprime(x, 2, 1) q <- qbetaprime(p, 2, 1)
Construct a function that computes the log density or CDF of the bivariate Clayton copula,
intended to be used with dcopula.
cclayton(theta) Cclayton(theta)cclayton(theta) Cclayton(theta)
theta |
Positive dependence parameter ( |
The Clayton copula density is
A function of two arguments (u,v) returning log copula density (cclayton) or copula CDF (Cclayton).
cgaussian(), cgumbel(), cfrank()
x <- c(0.5, 1); y <- c(0.2, 0.8) d1 <- dnorm(x, 1, log = TRUE); d2 <- dbeta(y, 2, 1, log = TRUE) p1 <- pnorm(x, 1); p2 <- pbeta(y, 2, 1) dcopula(d1, d2, p1, p2, copula = cclayton(2), log = TRUE) # CDF version (for discrete copulas) Cclayton(1.5)(0.5, 0.4)x <- c(0.5, 1); y <- c(0.2, 0.8) d1 <- dnorm(x, 1, log = TRUE); d2 <- dbeta(y, 2, 1, log = TRUE) p1 <- pnorm(x, 1); p2 <- pbeta(y, 2, 1) dcopula(d1, d2, p1, p2, copula = cclayton(2), log = TRUE) # CDF version (for discrete copulas) Cclayton(1.5)(0.5, 0.4)
Returns a function computing the log density of the bivariate Frank copula,
intended to be used with dcopula.
cfrank(theta) Cfrank(theta)cfrank(theta) Cfrank(theta)
theta |
Dependence parameter ( |
The Frank copula density is
A function of two arguments (u, v) returning either
the log copula density (cfrank) or the copula CDF (Cfrank).
cgaussian(), cclayton(), cgumbel()
x <- c(0.5, 1); y <- c(1, 2) d1 <- dnorm(x, 1, log = TRUE); d2 <- dexp(y, 2, log = TRUE) p1 <- pnorm(x, 1); p2 <- pexp(y, 2) dcopula(d1, d2, p1, p2, copula = cfrank(2), log = TRUE)x <- c(0.5, 1); y <- c(1, 2) d1 <- dnorm(x, 1, log = TRUE); d2 <- dexp(y, 2, log = TRUE) p1 <- pnorm(x, 1); p2 <- pexp(y, 2) dcopula(d1, d2, p1, p2, copula = cfrank(2), log = TRUE)
Returns a function computing the log density of the bivariate Gaussian copula,
intended to be used with dcopula.
cgaussian(rho = 0)cgaussian(rho = 0)
rho |
Correlation parameter ( |
Function of two arguments (u,v) returning log copula density.
The Gaussian copula density is
where , , and .
cclayton(), cgumbel(), cfrank()
x <- c(0.5, 1); y <- c(1, 2) d1 <- dnorm(x, 1, log = TRUE); d2 <- dexp(y, 2, log = TRUE) p1 <- pnorm(x, 1); p2 <- pexp(y, 2) dcopula(d1, d2, p1, p2, copula = cgaussian(0.5), log = TRUE)x <- c(0.5, 1); y <- c(1, 2) d1 <- dnorm(x, 1, log = TRUE); d2 <- dexp(y, 2, log = TRUE) p1 <- pnorm(x, 1); p2 <- pexp(y, 2) dcopula(d1, d2, p1, p2, copula = cgaussian(0.5), log = TRUE)
Returns a function computing the log density of the multivariate Gaussian copula, parameterised by the inverse correlation matrix.
cgmrf(Q)cgmrf(Q)
Q |
Inverse of a positive definite correlation matrix with unit diagonal. Can either be sparse or dense matrix. |
Caution: Parameterising the inverse correlation directly is difficult, as inverting it needs to yield a positive definite matrix with unit diagonal.
Hence we still advise parameterising the correaltion matrix R and computing its inverse.
This function is useful when you need access to the precision (i.e. inverse correlation) in your likelihood function.
Function with matrix argument U returning log copula density.
x <- c(0.5, 1); y <- c(1, 2); z <- c(0.2, 0.8) d1 <- dnorm(x, 1, log = TRUE); d2 <- dexp(y, 2, log = TRUE); d3 <- dbeta(z, 2, 1, log = TRUE) p1 <- pnorm(x, 1); p2 <- pexp(y, 2); p3 <- pbeta(z, 2, 1) R <- matrix(c(1,0.5,0.3,0.5,1,0.4,0.3,0.4,1), nrow = 3) ## Based on correlation matrix dmvcopula(cbind(d1, d2, d3), cbind(p1, p2, p3), copula = cmvgauss(R), log = TRUE) ## Based on precision matrix Q <- solve(R) dmvcopula(cbind(d1, d2, d3), cbind(p1, p2, p3), copula = cgmrf(Q), log = TRUE) ## Parameterisation inside a model # using RTMB::unstructured to get a valid correlation matrix library(RTMB) d <- 5 # dimension cor_func <- unstructured(d) npar <- length(cor_func$parms()) R <- cor_func$corr(rep(0.1, npar))x <- c(0.5, 1); y <- c(1, 2); z <- c(0.2, 0.8) d1 <- dnorm(x, 1, log = TRUE); d2 <- dexp(y, 2, log = TRUE); d3 <- dbeta(z, 2, 1, log = TRUE) p1 <- pnorm(x, 1); p2 <- pexp(y, 2); p3 <- pbeta(z, 2, 1) R <- matrix(c(1,0.5,0.3,0.5,1,0.4,0.3,0.4,1), nrow = 3) ## Based on correlation matrix dmvcopula(cbind(d1, d2, d3), cbind(p1, p2, p3), copula = cmvgauss(R), log = TRUE) ## Based on precision matrix Q <- solve(R) dmvcopula(cbind(d1, d2, d3), cbind(p1, p2, p3), copula = cgmrf(Q), log = TRUE) ## Parameterisation inside a model # using RTMB::unstructured to get a valid correlation matrix library(RTMB) d <- 5 # dimension cor_func <- unstructured(d) npar <- length(cor_func$parms()) R <- cor_func$corr(rep(0.1, npar))
Construct functions that compute either the log density or the CDF
of the bivariate Gumbel copula, intended for use with dcopula.
cgumbel(theta) Cgumbel(theta)cgumbel(theta) Cgumbel(theta)
theta |
Dependence parameter ( |
The Gumbel copula density
where contains the derivative terms ensuring the function is a density.
A function of two arguments (u, v) returning either
the log copula density (cgumbel) or the copula CDF (Cgumbel).
cgaussian(), cclayton(), cfrank()
x <- c(0.5, 1); y <- c(0.2, 0.4) d1 <- dnorm(x, 1, log = TRUE); d2 <- dbeta(y, 2, 1, log = TRUE) p1 <- pnorm(x, 1); p2 <- pbeta(y, 2, 1) dcopula(d1, d2, p1, p2, copula = cgumbel(1.5), log = TRUE) # CDF version (for discrete copulas) Cgumbel(1.5)(0.5, 0.4)x <- c(0.5, 1); y <- c(0.2, 0.4) d1 <- dnorm(x, 1, log = TRUE); d2 <- dbeta(y, 2, 1, log = TRUE) p1 <- pnorm(x, 1); p2 <- pbeta(y, 2, 1) dcopula(d1, d2, p1, p2, copula = cgumbel(1.5), log = TRUE) # CDF version (for discrete copulas) Cgumbel(1.5)(0.5, 0.4)
Returns a function computing the log density of the multivariate Gaussian copula,
intended to be used with dmvcopula.
cmvgauss(R)cmvgauss(R)
R |
Positive definite correlation matrix (unit diagonal) |
Function with matrix argument U returning log copula density.
x <- c(0.5, 1); y <- c(1, 2); z <- c(0.2, 0.8) d1 <- dnorm(x, 1, log = TRUE); d2 <- dexp(y, 2, log = TRUE); d3 <- dbeta(z, 2, 1, log = TRUE) p1 <- pnorm(x, 1); p2 <- pexp(y, 2); p3 <- pbeta(z, 2, 1) R <- matrix(c(1,0.5,0.3,0.5,1,0.4,0.3,0.4,1), nrow = 3) ## Based on correlation matrix dmvcopula(cbind(d1, d2, d3), cbind(p1, p2, p3), copula = cmvgauss(R), log = TRUE) ## Based on precision matrix Q <- solve(R) dmvcopula(cbind(d1, d2, d3), cbind(p1, p2, p3), copula = cgmrf(Q), log = TRUE) ## Parameterisation inside a model # using RTMB::unstructured to get a valid correlation matrix library(RTMB) d <- 5 # dimension cor_func <- unstructured(d) npar <- length(cor_func$parms()) R <- cor_func$corr(rep(0.1, npar))x <- c(0.5, 1); y <- c(1, 2); z <- c(0.2, 0.8) d1 <- dnorm(x, 1, log = TRUE); d2 <- dexp(y, 2, log = TRUE); d3 <- dbeta(z, 2, 1, log = TRUE) p1 <- pnorm(x, 1); p2 <- pexp(y, 2); p3 <- pbeta(z, 2, 1) R <- matrix(c(1,0.5,0.3,0.5,1,0.4,0.3,0.4,1), nrow = 3) ## Based on correlation matrix dmvcopula(cbind(d1, d2, d3), cbind(p1, p2, p3), copula = cmvgauss(R), log = TRUE) ## Based on precision matrix Q <- solve(R) dmvcopula(cbind(d1, d2, d3), cbind(p1, p2, p3), copula = cgmrf(Q), log = TRUE) ## Parameterisation inside a model # using RTMB::unstructured to get a valid correlation matrix library(RTMB) d <- 5 # dimension cor_func <- unstructured(d) npar <- length(cor_func$parms()) R <- cor_func$corr(rep(0.1, npar))
Computes the joint density (or log-density) of a bivariate distribution constructed from two arbitrary margins combined with a specified copula.
dcopula(d1, d2, p1, p2, copula = cgaussian(0), log = FALSE)dcopula(d1, d2, p1, p2, copula = cgaussian(0), log = FALSE)
d1, d2
|
Marginal density values. If |
p1, p2
|
Marginal CDF values. Need not be supplied on log scale. |
copula |
A function of two arguments ( |
log |
Logical; if |
The joint density is
where are the marginal CDFs, are the marginal densities,
and is the copula density.
The marginal densities d1, d2 and CDFs p1, p2
must be differentiable for automatic differentiation (AD) to work.
Available copula constructors are:
Joint density (or log-density) under the bivariate copula.
# Normal + Exponential margins with Gaussian copula x <- c(0.5, 1); y <- c(1, 2) d1 <- dnorm(x, 1, log = TRUE); d2 <- dexp(y, 2, log = TRUE) p1 <- pnorm(x, 1); p2 <- pexp(y, 2) dcopula(d1, d2, p1, p2, copula = cgaussian(0.5), log = TRUE) # Normal + Beta margins with Clayton copula x <- c(0.5, 1); y <- c(0.2, 0.8) d1 <- dnorm(x, 1, log = TRUE); d2 <- dbeta(y, 2, 1, log = TRUE) p1 <- pnorm(x, 1); p2 <- pbeta(y, 2, 1) dcopula(d1, d2, p1, p2, copula = cclayton(2), log = TRUE) # Normal + Beta margins with Gumbel copula x <- c(0.5, 1); y <- c(0.2, 0.4) d1 <- dnorm(x, 1, log = TRUE); d2 <- dbeta(y, 2, 1, log = TRUE) p1 <- pnorm(x, 1); p2 <- pbeta(y, 2, 1) dcopula(d1, d2, p1, p2, copula = cgumbel(1.5), log = TRUE) # Normal + Exponential margins with Frank copula x <- c(0.5, 1); y <- c(1, 2) d1 <- dnorm(x, 1, log = TRUE); d2 <- dexp(y, 2, log = TRUE) p1 <- pnorm(x, 1); p2 <- pexp(y, 2) dcopula(d1, d2, p1, p2, copula = cfrank(2), log = TRUE)# Normal + Exponential margins with Gaussian copula x <- c(0.5, 1); y <- c(1, 2) d1 <- dnorm(x, 1, log = TRUE); d2 <- dexp(y, 2, log = TRUE) p1 <- pnorm(x, 1); p2 <- pexp(y, 2) dcopula(d1, d2, p1, p2, copula = cgaussian(0.5), log = TRUE) # Normal + Beta margins with Clayton copula x <- c(0.5, 1); y <- c(0.2, 0.8) d1 <- dnorm(x, 1, log = TRUE); d2 <- dbeta(y, 2, 1, log = TRUE) p1 <- pnorm(x, 1); p2 <- pbeta(y, 2, 1) dcopula(d1, d2, p1, p2, copula = cclayton(2), log = TRUE) # Normal + Beta margins with Gumbel copula x <- c(0.5, 1); y <- c(0.2, 0.4) d1 <- dnorm(x, 1, log = TRUE); d2 <- dbeta(y, 2, 1, log = TRUE) p1 <- pnorm(x, 1); p2 <- pbeta(y, 2, 1) dcopula(d1, d2, p1, p2, copula = cgumbel(1.5), log = TRUE) # Normal + Exponential margins with Frank copula x <- c(0.5, 1); y <- c(1, 2) d1 <- dnorm(x, 1, log = TRUE); d2 <- dexp(y, 2, log = TRUE) p1 <- pnorm(x, 1); p2 <- pexp(y, 2) dcopula(d1, d2, p1, p2, copula = cfrank(2), log = TRUE)
Computes the joint probability mass function of two discrete margins combined with a copula CDF.
ddcopula(d1, d2, p1, p2, Copula, log = FALSE)ddcopula(d1, d2, p1, p2, Copula, log = FALSE)
d1, d2
|
Marginal p.m.f. values at the observed points, not on log-scale. |
p1, p2
|
Marginal CDF values at the observed points. |
Copula |
A function of two arguments returning the copula CDF. |
log |
Logical; if |
The joint probability mass function for two discrete margins is
where are the marginal CDFs, and is the copula CDF.
Available copula CDF constructors are:
Joint probability (or log-probability) under chosen copula
x <- c(3,5); y <- c(2,4) d1 <- dpois(x, 4); d2 <- dpois(y, 3) p1 <- ppois(x, 4); p2 <- ppois(y, 3) ddcopula(d1, d2, p1, p2, Copula = Cclayton(2), log = FALSE)x <- c(3,5); y <- c(2,4) d1 <- dpois(x, 4); d2 <- dpois(y, 3) p1 <- ppois(x, 4); p2 <- ppois(y, 3) ddcopula(d1, d2, p1, p2, Copula = Cclayton(2), log = FALSE)
Density and and random generation for the Dirichlet distribution.
ddirichlet(x, alpha, log = FALSE) rdirichlet(n, alpha)ddirichlet(x, alpha, log = FALSE) rdirichlet(n, alpha)
x |
vector or matrix of quantiles. If |
alpha |
vector or matrix of positive shape parameters |
log |
logical; if |
n |
number of random values to return. |
This implementation of ddirichlet allows for automatic differentiation with RTMB.
where lies on the unit simplex and .
ddirichlet gives the density, rdirichlet generates random deviates.
# single alpha alpha <- c(1,2,3) x <- rdirichlet(1, alpha) d <- ddirichlet(x, alpha) # vectorised over alpha alpha <- rbind(alpha, 2*alpha) x <- rdirichlet(2, alpha)# single alpha alpha <- c(1,2,3) x <- rdirichlet(1, alpha) d <- ddirichlet(x, alpha) # vectorised over alpha alpha <- rbind(alpha, 2*alpha) x <- rdirichlet(2, alpha)
Density and and random generation for the Dirichlet-multinomial distribution.
ddirmult(x, size, alpha, log = FALSE) rdirmult(n, size, alpha)ddirmult(x, size, alpha, log = FALSE) rdirmult(n, size, alpha)
x |
vector or matrix of non-negative counts, where rows are observations and columns are categories. |
size |
vector of total counts for each observation. Needs to match the row sums of |
alpha |
vector or matrix of positive shape parameters |
log |
logical; if |
n |
number of random values to return. |
This implementation of ddirmult allows for automatic differentiation with RTMB.
where and .
ddirmult gives the density and rdirmult generates random samples.
# single alpha alpha <- c(1,2,3) size <- 10 x <- rdirmult(1, size, alpha) d <- ddirmult(x, size, alpha) # vectorised over alpha and size alpha <- rbind(alpha, 2*alpha) size <- c(size, 3*size) x <- rdirmult(2, size, alpha)# single alpha alpha <- c(1,2,3) size <- 10 x <- rdirmult(1, size, alpha) d <- ddirmult(x, size, alpha) # vectorised over alpha and size alpha <- rbind(alpha, 2*alpha) size <- c(size, 3*size) x <- rdirmult(2, size, alpha)
Computes the joint density (or log-density) of a distribution constructed from any number of arbitrary margins combined with a specified copula.
dmvcopula(D, P, copula = cmvgauss(diag(ncol(D))), log = FALSE)dmvcopula(D, P, copula = cmvgauss(diag(ncol(D))), log = FALSE)
D |
Matrix of marginal density values of with rows corresponding to observations and columns corresponding to dimensions.
If |
P |
Matrix of marginal CDF values of the same dimension as |
copula |
A function of a matrix argument |
log |
Logical; if |
The joint density is
where are the marginal CDFs, are the marginal densities,
and is the copula density.
The marginal densities d_1, ..., d_d and CDFs p_1, ..., p_d
must be differentiable for automatic differentiation (AD) to work.
Available multivariate copula constructors are:
cmvgauss (Multivariate Gaussian copula)
cgmrf (Multivariate Gaussian copula parameterised by precision (inverse correlation) matrix)
Joint density (or log-density) under the chosen copula.
x <- c(0.5, 1); y <- c(1, 2); z <- c(0.2, 0.8) d1 <- dnorm(x, 1, log = TRUE); d2 <- dexp(y, 2, log = TRUE); d3 <- dbeta(z, 2, 1, log = TRUE) p1 <- pnorm(x, 1); p2 <- pexp(y, 2); p3 <- pbeta(z, 2, 1) R <- matrix(c(1,0.5,0.3,0.5,1,0.4,0.3,0.4,1), nrow = 3) ### Multivariate Gaussian copula ## Based on correlation matrix dmvcopula(cbind(d1, d2, d3), cbind(p1, p2, p3), copula = cmvgauss(R), log = TRUE) ## Based on precision matrix Q <- solve(R) dmvcopula(cbind(d1, d2, d3), cbind(p1, p2, p3), copula = cgmrf(Q), log = TRUE) ## Parameterisation inside a model # using RTMB::unstructured to get a valid correlation matrix library(RTMB) d <- 5 # dimension cor_func <- unstructured(d) npar <- length(cor_func$parms()) R <- cor_func$corr(rep(0.1, npar))x <- c(0.5, 1); y <- c(1, 2); z <- c(0.2, 0.8) d1 <- dnorm(x, 1, log = TRUE); d2 <- dexp(y, 2, log = TRUE); d3 <- dbeta(z, 2, 1, log = TRUE) p1 <- pnorm(x, 1); p2 <- pexp(y, 2); p3 <- pbeta(z, 2, 1) R <- matrix(c(1,0.5,0.3,0.5,1,0.4,0.3,0.4,1), nrow = 3) ### Multivariate Gaussian copula ## Based on correlation matrix dmvcopula(cbind(d1, d2, d3), cbind(p1, p2, p3), copula = cmvgauss(R), log = TRUE) ## Based on precision matrix Q <- solve(R) dmvcopula(cbind(d1, d2, d3), cbind(p1, p2, p3), copula = cgmrf(Q), log = TRUE) ## Parameterisation inside a model # using RTMB::unstructured to get a valid correlation matrix library(RTMB) d <- 5 # dimension cor_func <- unstructured(d) npar <- length(cor_func$parms()) R <- cor_func$corr(rep(0.1, npar))
AD-compatible error function and complementary error function
erf(x) erfc(x)erf(x) erfc(x)
x |
vector of evaluation points |
erf(x) returns the error function and erfc(x) returns the complementary error function.
erf(1) erfc(1)erf(1) erfc(1)
Density, distribution function, quantile function, and random generation for the exponentially modified Gaussian distribution.
dexgauss(x, mu = 0, sigma = 1, lambda = 1, log = FALSE) pexgauss(q, mu = 0, sigma = 1, lambda = 1, lower.tail = TRUE, log.p = FALSE) qexgauss(p, mu = 0, sigma = 1, lambda = 1, lower.tail = TRUE, log.p = FALSE) rexgauss(n, mu = 0, sigma = 1, lambda = 1)dexgauss(x, mu = 0, sigma = 1, lambda = 1, log = FALSE) pexgauss(q, mu = 0, sigma = 1, lambda = 1, lower.tail = TRUE, log.p = FALSE) qexgauss(p, mu = 0, sigma = 1, lambda = 1, lower.tail = TRUE, log.p = FALSE) rexgauss(n, mu = 0, sigma = 1, lambda = 1)
x, q
|
vector of quantiles |
mu |
mean parameter of the Gaussian part |
sigma |
standard deviation parameter of the Gaussian part, must be positive. |
lambda |
rate parameter of the exponential part, must be positive. |
log, log.p
|
logical; if |
lower.tail |
logical; if |
p |
vector of probabilities |
n |
number of random values to return |
This implementation of dexgauss and pexgauss allows for automatic differentiation with RTMB.
qexgauss and rexgauss are reparameterised imports from gamlss.dist::exGAUS.
If and , then
follows the exponentially modified Gaussian distribution with parameters , , and .
The density is
where is the standard normal CDF.
dexgauss gives the density, pexgauss gives the distribution function, qexgauss gives the quantile function, and rexgauss generates random deviates.
Rigby, R. A., Stasinopoulos, D. M., Heller, G. Z., and De Bastiani, F. (2019) Distributions for modeling location, scale, and shape: Using GAMLSS in R, Chapman and Hall/CRC, doi:10.1201/9780429298547. An older version can be found in https://www.gamlss.com/.
x <- rexgauss(1, 1, 2, 2) d <- dexgauss(x, 1, 2, 2) p <- pexgauss(x, 1, 2, 2) q <- qexgauss(p, 1, 2, 2)x <- rexgauss(1, 1, 2, 2) d <- dexgauss(x, 1, 2, 2) p <- pexgauss(x, 1, 2, 2) q <- qexgauss(p, 1, 2, 2)
Density, distribution function, and random generation for the folded normal distribution.
dfoldnorm(x, mu = 0, sigma = 1, log = FALSE) pfoldnorm(q, mu = 0, sigma = 1, lower.tail = TRUE, log.p = FALSE) rfoldnorm(n, mu = 0, sigma = 1)dfoldnorm(x, mu = 0, sigma = 1, log = FALSE) pfoldnorm(q, mu = 0, sigma = 1, lower.tail = TRUE, log.p = FALSE) rfoldnorm(n, mu = 0, sigma = 1)
x, q
|
vector of quantiles |
mu |
location parameter |
sigma |
scale parameter, must be positive. |
log, log.p
|
logical; if |
lower.tail |
logical; if |
n |
number of random values to return |
p |
vector of probabilities |
This implementation of dfoldnorm allows for automatic differentiation with RTMB.
dfoldnorm gives the density, pfoldnorm gives the distribution function, and rfoldnorm generates random deviates.
x <- rfoldnorm(1, 1, 2) d <- dfoldnorm(x, 1, 2) p <- pfoldnorm(x, 1, 2)x <- rfoldnorm(1, 1, 2) d <- dfoldnorm(x, 1, 2) p <- pfoldnorm(x, 1, 2)
Density, distribution function, quantile function, and random generation for the gamma distribution reparameterised in terms of mean and standard deviation.
dgamma2(x, mean = 1, sd = 1, log = FALSE) pgamma2(q, mean = 1, sd = 1, lower.tail = TRUE, log.p = FALSE) qgamma2(p, mean = 1, sd = 1, lower.tail = TRUE, log.p = FALSE) rgamma2(n, mean = 1, sd = 1)dgamma2(x, mean = 1, sd = 1, log = FALSE) pgamma2(q, mean = 1, sd = 1, lower.tail = TRUE, log.p = FALSE) qgamma2(p, mean = 1, sd = 1, lower.tail = TRUE, log.p = FALSE) rgamma2(n, mean = 1, sd = 1)
x, q
|
vector of quantiles |
mean |
mean parameter, must be positive. |
sd |
standard deviation parameter, must be positive. |
log, log.p
|
logical; if |
lower.tail |
logical; if |
p |
vector of probabilities |
n |
number of random values to return. |
This implementation allows for automatic differentiation with RTMB.
Reparameterises the gamma distribution via and :
dgamma2 gives the density, pgamma2 gives the distribution function, qgamma2 gives the quantile function, and rgamma2 generates random deviates.
x <- rgamma2(1) d <- dgamma2(x) p <- pgamma2(x) q <- qgamma2(p)x <- rgamma2(1) d <- dgamma2(x) p <- pgamma2(x) q <- qgamma2(p)
Density, distribution function, quantile function, and random generation for the generalised Gamma distribution.
dgengamma(x, mu = 1, sigma = 0.5, nu = 1, log = FALSE) pgengamma(q, mu = 1, sigma = 0.5, nu = 1, lower.tail = TRUE, log.p = FALSE) qgengamma(p, mu = 1, sigma = 0.5, nu = 1, lower.tail = TRUE, log.p = FALSE) rgengamma(n, mu = 1, sigma = 0.5, nu = 1)dgengamma(x, mu = 1, sigma = 0.5, nu = 1, log = FALSE) pgengamma(q, mu = 1, sigma = 0.5, nu = 1, lower.tail = TRUE, log.p = FALSE) qgengamma(p, mu = 1, sigma = 0.5, nu = 1, lower.tail = TRUE, log.p = FALSE) rgengamma(n, mu = 1, sigma = 0.5, nu = 1)
x, q
|
vector of quantiles |
mu |
location parameter, must be positive. |
sigma |
scale parameter, must be positive. |
nu |
skewness parameter (real). |
log, log.p
|
logical; if |
lower.tail |
logical; if |
p |
vector of probabilities |
n |
number of random values to return |
This implementation of dgengamma, pgengamma, and qgengamma allows for automatic differentiation with RTMB.
For the density is that of a generalised gamma: setting ,
For the distribution reduces to a log-normal with log-mean and log-standard-deviation .
dgengamma gives the density, pgengamma gives the distribution function, qgengamma gives the quantile function, and rgengamma generates random deviates.
Rigby, R. A., Stasinopoulos, D. M., Heller, G. Z., and De Bastiani, F. (2019) Distributions for modeling location, scale, and shape: Using GAMLSS in R, Chapman and Hall/CRC, doi:10.1201/9780429298547. An older version can be found in https://www.gamlss.com/.
x <- rgengamma(5, mu = 4, sigma = 0.5, nu = 0.5) d <- dgengamma(x, mu = 4, sigma = 0.5, nu = 0.5) p <- pgengamma(x, mu = 4, sigma = 0.5, nu = 0.5) q <- qgengamma(p, mu = 4, sigma = 0.5, nu = 0.5)x <- rgengamma(5, mu = 4, sigma = 0.5, nu = 0.5) d <- dgengamma(x, mu = 4, sigma = 0.5, nu = 0.5) p <- pgengamma(x, mu = 4, sigma = 0.5, nu = 0.5) q <- qgengamma(p, mu = 4, sigma = 0.5, nu = 0.5)
Probability mass function, distribution function, and random generation for the generalised Poisson distribution.
dgenpois(x, lambda = 1, phi = 1, log = FALSE) pgenpois(q, lambda = 1, phi = 1, lower.tail = TRUE, log.p = FALSE) qgenpois(p, lambda = 1, phi = 1, lower.tail = TRUE, log.p = FALSE, max.value = 10000) rgenpois(n, lambda = 1, phi = 1, max.value = 10000)dgenpois(x, lambda = 1, phi = 1, log = FALSE) pgenpois(q, lambda = 1, phi = 1, lower.tail = TRUE, log.p = FALSE) qgenpois(p, lambda = 1, phi = 1, lower.tail = TRUE, log.p = FALSE, max.value = 10000) rgenpois(n, lambda = 1, phi = 1, max.value = 10000)
x, q
|
integer vector of counts |
lambda |
vector of positive means |
phi |
vector of non-negative dispersion parameters |
log, log.p
|
logical; return log-density if TRUE |
lower.tail |
logical; if |
p |
vector of probabilities |
max.value |
a constant, set to the default value of 10000 for how far the algorithm should look for |
n |
number of random values to return. |
This implementation of dgenpois allows for automatic differentiation with RTMB.
The other functions are imported from gamlss.dist::GPO.
The distribution has mean and variance .
For it reduces to the Poisson distribution, however must be strictly positive here.
dgenpois gives the probability mass function, pgenpois gives the distribution function, qgenpois gives the quantile function, and rgenpois generates random deviates.
set.seed(123) x <- rgenpois(1, 2, 3) d <- dgenpois(x, 2, 3) p <- pgenpois(x, 2, 3) q <- qgenpois(p, 2, 3)set.seed(123) x <- rgenpois(1, 2, 3) d <- dgenpois(x, 2, 3) p <- pgenpois(x, 2, 3) q <- qgenpois(p, 2, 3)
Density, distribution function, quantile function, and random generation for the Gompertz distribution.
dgompertz(x, eta = 1, b = 1, log = FALSE) pgompertz(q, eta = 1, b = 1, lower.tail = TRUE, log.p = FALSE) qgompertz(p, eta = 1, b = 1, lower.tail = TRUE, log.p = FALSE) rgompertz(n, eta = 1, b = 1)dgompertz(x, eta = 1, b = 1, log = FALSE) pgompertz(q, eta = 1, b = 1, lower.tail = TRUE, log.p = FALSE) qgompertz(p, eta = 1, b = 1, lower.tail = TRUE, log.p = FALSE) rgompertz(n, eta = 1, b = 1)
x, q
|
vector of quantiles (non-negative) |
eta |
shape parameter, must be positive |
b |
rate parameter, must be positive |
log, log.p
|
logical; if |
lower.tail |
logical; if |
p |
vector of probabilities |
n |
number of random values to return |
The Gompertz distribution with shape and rate has density
with CDF and quantile function
.
dgompertz gives the density, pgompertz gives the distribution function,
qgompertz gives the quantile function, and rgompertz generates random deviates.
https://en.wikipedia.org/wiki/Gompertz_distribution
x <- rgompertz(1, eta = 1, b = 1) d <- dgompertz(x, eta = 1, b = 1) p <- pgompertz(x, eta = 1, b = 1) q <- qgompertz(p, eta = 1, b = 1)x <- rgompertz(1, eta = 1, b = 1) d <- dgompertz(x, eta = 1, b = 1) p <- pgompertz(x, eta = 1, b = 1) q <- qgompertz(p, eta = 1, b = 1)
Density, distribution function, quantile function, and random generation for the Gumbel distribution.
dgumbel(x, location = 0, scale = 1, log = FALSE) pgumbel(q, location = 0, scale = 1, lower.tail = TRUE, log.p = FALSE) qgumbel(p, location = 0, scale = 1, lower.tail = TRUE, log.p = FALSE) rgumbel(n, location = 0, scale = 1)dgumbel(x, location = 0, scale = 1, log = FALSE) pgumbel(q, location = 0, scale = 1, lower.tail = TRUE, log.p = FALSE) qgumbel(p, location = 0, scale = 1, lower.tail = TRUE, log.p = FALSE) rgumbel(n, location = 0, scale = 1)
x, q
|
vector of quantiles |
location |
location parameter |
scale |
scale parameter, must be positive. |
log, log.p
|
logical; if |
lower.tail |
logical; if |
p |
vector of probabilities |
n |
number of random values to return |
This implementation of dgumbel allows for automatic differentiation with RTMB.
dgumbel gives the density, pgumbel gives the distribution function, qgumbel gives the quantile function, and rgumbel generates random deviates.
x <- rgumbel(1, 0.5, 2) d <- dgumbel(x, 0.5, 2) p <- pgumbel(x, 0.5, 2) q <- qgumbel(p, 0.5, 2)x <- rgumbel(1, 0.5, 2) d <- dgumbel(x, 0.5, 2) p <- pgumbel(x, 0.5, 2) q <- qgumbel(p, 0.5, 2)
Density, distribution function, quantile function, and random generation for the inverse Chi-squared distribution.
dinvchisq(x, df, scale = 1/df, log = FALSE) pinvchisq(q, df, scale = 1/df, lower.tail = TRUE, log.p = FALSE) qinvchisq(p, df, scale = 1/df, lower.tail = TRUE, log.p = FALSE) rinvchisq(n, df, scale = 1/df)dinvchisq(x, df, scale = 1/df, log = FALSE) pinvchisq(q, df, scale = 1/df, lower.tail = TRUE, log.p = FALSE) qinvchisq(p, df, scale = 1/df, lower.tail = TRUE, log.p = FALSE) rinvchisq(n, df, scale = 1/df)
x, q
|
vector of quantiles, must be positive. |
df |
degrees of freedom ( |
scale |
optional positive scale parameter. Default value of |
log, log.p
|
logical; if |
lower.tail |
logical; if |
p |
vector of probabilities |
n |
number of random values to return |
If , then .
The inverse Chi-squared distribution with degrees of freedom has
density
This implementation of dinvchisq, pinvchisq, and qinvchisq allows for automatic differentiation with RTMB.
dinvchisq gives the density, pinvchisq gives the distribution function,
qinvchisq gives the quantile function, and rinvchisq generates random deviates.
x <- rinvchisq(1, df = 5) d <- dinvchisq(x, df = 5) p <- pinvchisq(x, df = 5) q <- qinvchisq(p, df = 5)x <- rinvchisq(1, df = 5) d <- dinvchisq(x, df = 5) p <- pinvchisq(x, df = 5) q <- qinvchisq(p, df = 5)
Density, distribution function, and random generation for the inverse Gamma distribution.
dinvgamma(x, shape, rate, scale = 1/rate, log = FALSE) pinvgamma(q, shape, rate, scale = 1/rate, lower.tail = TRUE, log.p = FALSE) qinvgamma(p, shape, rate, scale = 1/rate, lower.tail = TRUE, log.p = FALSE) rinvgamma(n, shape, rate, scale = 1/rate)dinvgamma(x, shape, rate, scale = 1/rate, log = FALSE) pinvgamma(q, shape, rate, scale = 1/rate, lower.tail = TRUE, log.p = FALSE) qinvgamma(p, shape, rate, scale = 1/rate, lower.tail = TRUE, log.p = FALSE) rinvgamma(n, shape, rate, scale = 1/rate)
x, q
|
vector of quantiles, must be positive. |
shape, rate, scale
|
positive parameters of corresponding gamma distribution |
log, log.p
|
logical; if |
lower.tail |
logical; if |
p |
vector of probabilities |
n |
number of random values to return |
This implementation of dinvgamma, pinvgamma, and qinvgamma allows for automatic differentiation with RTMB.
If , then .
where .
dinvgamma gives the density, pinvgamma gives the distribution function, qinvgamma gives the quantile function, and rinvgamma generates random deviates.
x <- rinvgamma(1, 1, 0.5) d <- dinvgamma(x, 1, 0.5) p <- pinvgamma(x, 1, 0.5) q <- qinvgamma(p, 1, 0.5)x <- rinvgamma(1, 1, 0.5) d <- dinvgamma(x, 1, 0.5) p <- pinvgamma(x, 1, 0.5) q <- qinvgamma(p, 1, 0.5)
Density, distribution function, and random generation for the inverse Gaussian distribution.
dinvgauss(x, mean = 1, shape = 1, log = FALSE) pinvgauss(q, mean = 1, shape = 1, lower.tail = TRUE, log.p = FALSE) qinvgauss(p, mean = 1, shape = 1, lower.tail = TRUE, log.p = FALSE, ...) rinvgauss(n, mean = 1, shape = 1)dinvgauss(x, mean = 1, shape = 1, log = FALSE) pinvgauss(q, mean = 1, shape = 1, lower.tail = TRUE, log.p = FALSE) qinvgauss(p, mean = 1, shape = 1, lower.tail = TRUE, log.p = FALSE, ...) rinvgauss(n, mean = 1, shape = 1)
x, q
|
vector of quantiles, must be positive. |
mean |
location parameter |
shape |
shape parameter, must be positive. |
log, log.p
|
logical; if |
lower.tail |
logical; if |
p |
vector of probabilities |
... |
additional parameter passed to |
n |
number of random values to return |
This implementation of dinvgauss allows for automatic differentiation with RTMB.
qinvgauss and rinvgauss are imported from the statmod package.
dinvgauss gives the density, pinvgauss gives the distribution function, qinvgauss gives the quantile function, and rinvgauss generates random deviates.
x <- rinvgauss(1, 1, 0.5) d <- dinvgauss(x, 1, 0.5) p <- pinvgauss(x, 1, 0.5) q <- qinvgauss(p, 1, 0.5)x <- rinvgauss(1, 1, 0.5) d <- dinvgauss(x, 1, 0.5) p <- pinvgauss(x, 1, 0.5) q <- qinvgauss(p, 1, 0.5)
Density, distribution function, quantile function, and random generation for the Kumaraswamy distribution.
dkumar(x, a, b, log = FALSE) pkumar(q, a, b, lower.tail = TRUE, log.p = FALSE) qkumar(p, a, b, lower.tail = TRUE, log.p = FALSE) rkumar(n, a, b)dkumar(x, a, b, log = FALSE) pkumar(q, a, b, lower.tail = TRUE, log.p = FALSE) qkumar(p, a, b, lower.tail = TRUE, log.p = FALSE) rkumar(n, a, b)
x, q
|
vector of quantiles in |
a, b
|
positive shape parameters |
log, log.p
|
logical; if |
lower.tail |
logical; if |
p |
vector of probabilities |
n |
number of random values to return. |
dkumar gives the density, pkumar gives the distribution function, qkumar gives the quantile function, and rkumar generates random deviates.
x <- rkumar(1, a = 1, b = 2) d <- dkumar(x, a = 1, b = 2) p <- pkumar(x, a = 1, b = 2) q <- qkumar(p, a = 1, b = 2)x <- rkumar(1, a = 1, b = 2) d <- dkumar(x, a = 1, b = 2) p <- pkumar(x, a = 1, b = 2) q <- qkumar(p, a = 1, b = 2)
Density, distribution function, quantile function, and random generation for the Laplace distribution.
dlaplace(x, mu = 0, b = 1, eps = NULL, log = FALSE) plaplace(q, mu = 0, b = 1, lower.tail = TRUE, log.p = FALSE) qlaplace(p, mu = 0, b = 1, lower.tail = TRUE, log.p = FALSE) rlaplace(n, mu = 0, b = 1)dlaplace(x, mu = 0, b = 1, eps = NULL, log = FALSE) plaplace(q, mu = 0, b = 1, lower.tail = TRUE, log.p = FALSE) qlaplace(p, mu = 0, b = 1, lower.tail = TRUE, log.p = FALSE) rlaplace(n, mu = 0, b = 1)
x, q
|
vector of quantiles |
mu |
location parameter |
b |
scale parameter, must be positive. |
eps |
optional smoothing parameter for |
log, log.p
|
logical; if |
lower.tail |
logical; if |
p |
vector of probabilities |
n |
number of random values to return |
This implementation of dlaplace allows for automatic differentiation with RTMB.
dlaplace gives the density, plaplace gives the distribution function, qlaplace gives the quantile function, and rlaplace generates random deviates.
x <- rlaplace(1, 1, 1) d <- dlaplace(x, 1, 1) p <- plaplace(x, 1, 1) q <- qlaplace(p, 1, 1)x <- rlaplace(1, 1, 1) d <- dlaplace(x, 1, 1) p <- plaplace(x, 1, 1) q <- qlaplace(p, 1, 1)
Density, distribution function, quantile function, and random generation for the log-logistic distribution.
dllogis(x, alpha = 1, beta = 1, log = FALSE) pllogis(q, alpha = 1, beta = 1, lower.tail = TRUE, log.p = FALSE) qllogis(p, alpha = 1, beta = 1, lower.tail = TRUE, log.p = FALSE) rllogis(n, alpha = 1, beta = 1)dllogis(x, alpha = 1, beta = 1, log = FALSE) pllogis(q, alpha = 1, beta = 1, lower.tail = TRUE, log.p = FALSE) qllogis(p, alpha = 1, beta = 1, lower.tail = TRUE, log.p = FALSE) rllogis(n, alpha = 1, beta = 1)
x, q
|
vector of quantiles ( |
alpha |
scale parameter ( |
beta |
shape parameter ( |
log |
logical; if |
lower.tail |
logical; if |
log.p |
logical; if |
p |
vector of probabilities. |
n |
number of random values to return. |
The log-logistic distribution has density
where is the scale parameter and is the shape parameter.
The scale parameter equals the median. Larger gives lighter tails; the distribution
has a finite mean only when and a finite variance only when .
The log-logistic arises naturally as the distribution of where
, which yields the
numerically convenient log-density
The CDF is
and the quantile function is
dllogis allows for automatic differentiation with RTMB.
dllogis gives the density, pllogis gives the distribution function,
qllogis gives the quantile function, and rllogis generates random deviates.
x <- rllogis(5, alpha = 1, beta = 2) dllogis(x, alpha = 1, beta = 2) pllogis(c(0.5, 1, 2), alpha = 1, beta = 2) qllogis(c(0.25, 0.5, 0.75), alpha = 1, beta = 2)x <- rllogis(5, alpha = 1, beta = 2) dllogis(x, alpha = 1, beta = 2) pllogis(c(0.5, 1, 2), alpha = 1, beta = 2) qllogis(c(0.25, 0.5, 0.75), alpha = 1, beta = 2)
Efficient Monte Carlo sampling of parameters (and REPORTed quantities) from the approximate posterior of an (R)TMB model.
See sdreport for details on posterior variance-covariance in random effects models.
mcreport( obj, nSamples = 1000, include_random_pars = TRUE, report = FALSE, Q = NULL, ... )mcreport( obj, nSamples = 1000, include_random_pars = TRUE, report = FALSE, Q = NULL, ... )
obj |
Optimised |
nSamples |
Number of samples to draw |
include_random_pars |
Logical; Should random parameters be included in the output? |
report |
Logical; Should |
Q |
Optional precalculated sparse precision matrix returned by |
... |
For internal use only |
Caution: This has nothing to do with Bayesian posterior sampling. It is simple Monte Carlo sampling from a Gaussian distribution with mean at the MLE and covariance given by the inverse of the Hessian (for fixed effects models) or the joint precision matrix (for random effects models). This is a common approach to get an approximate idea of parameter uncertainty around the MLE, but it relies on large-sample asymptotics.
Sampling random effects from their posterior as compared to calculating marginal standard devitions like sdreport is particularly useful for multidimensional random effects (e.g. for locations and )
where pointwise confidence intervals (e.g. along a path) based on standard deviations are not possible.
A list structured like the original parameter list used in the MakeADFun call (potentially including additional REPORTed quantities). Each entry is a list with nSamples entries.
set.seed(123) ### Example 1: Without random effects ### # simulate data x <- rgumbel(100, location = 5, scale = 2) # negative log-likelihood function nll <- function(par) { getAll(par) x <- OBS(x) # mark x as the response scale <- exp(log_scale) ADREPORT(scale); REPORT(scale) -sum(dgumbel(x, loc, scale, log = TRUE)) } # RTMB AD object par <- list(loc = 5, log_scale = log(2)) obj <- MakeADFun(nll, par, silent = TRUE) # model fitting using AD gradient opt <- nlminb(obj$par, obj$fn, obj$gr) # sample from asymptotic distribution of the MLE samples <- mcreport(obj, report = TRUE) # parameters mean(unlist(samples$loc)); sd(unlist(samples$loc)) # reported transformed parameters mean(unlist(samples$scale)); sd(unlist(samples$scale)) # compare to delta method sd from sdreport sdr <- sdreport(obj) summary(sdr) ### Example 2: With random effects ### # simulate random effects id <- rep(1:10, each = 100) # 10 groups with 100 observations each true_gamma <- rnorm(10, 0, 2) # random coefficients true_probs <- plogis(1 + true_gamma) # linear predictor # simulate Bernoulli data conditional on random coefficients x <- rbinom(1000, 1, true_probs[id]) # negative log-likelihood function nll <- function(par) { getAll(par) x <- OBS(x) # mark x as the response # random effect likelihood sd_gamma <- exp(log_sd_gamma) ADREPORT(sd_gamma); REPORT(sd_gamma) nll <- -sum(dnorm(gamma, 0, sd_gamma, log = TRUE)) # data likelihood probs <- plogis(beta0 + gamma) # linear predictor ADREPORT(probs); REPORT(probs) nll <- nll - sum(dbinom(x, 1, probs[id], log = TRUE)) return(nll) } # RTMB Laplace approximation par <- list(beta0 = 1, log_sd_gamma = log(2), gamma = rep(0,10)) obj <- MakeADFun(nll, par, random = "gamma", silent = TRUE) # model fitting using AD gradient opt <- nlminb(obj$par, obj$fn, obj$gr) # draw samples samples <- mcreport(obj, report = TRUE) # run sdreport sdr <- sdreport(obj) # standard deviation using delta method as.list(sdr, "Std", report = TRUE)$probs # standard deviation from samples apply(simplify2array(samples$probs), 1, sd)set.seed(123) ### Example 1: Without random effects ### # simulate data x <- rgumbel(100, location = 5, scale = 2) # negative log-likelihood function nll <- function(par) { getAll(par) x <- OBS(x) # mark x as the response scale <- exp(log_scale) ADREPORT(scale); REPORT(scale) -sum(dgumbel(x, loc, scale, log = TRUE)) } # RTMB AD object par <- list(loc = 5, log_scale = log(2)) obj <- MakeADFun(nll, par, silent = TRUE) # model fitting using AD gradient opt <- nlminb(obj$par, obj$fn, obj$gr) # sample from asymptotic distribution of the MLE samples <- mcreport(obj, report = TRUE) # parameters mean(unlist(samples$loc)); sd(unlist(samples$loc)) # reported transformed parameters mean(unlist(samples$scale)); sd(unlist(samples$scale)) # compare to delta method sd from sdreport sdr <- sdreport(obj) summary(sdr) ### Example 2: With random effects ### # simulate random effects id <- rep(1:10, each = 100) # 10 groups with 100 observations each true_gamma <- rnorm(10, 0, 2) # random coefficients true_probs <- plogis(1 + true_gamma) # linear predictor # simulate Bernoulli data conditional on random coefficients x <- rbinom(1000, 1, true_probs[id]) # negative log-likelihood function nll <- function(par) { getAll(par) x <- OBS(x) # mark x as the response # random effect likelihood sd_gamma <- exp(log_sd_gamma) ADREPORT(sd_gamma); REPORT(sd_gamma) nll <- -sum(dnorm(gamma, 0, sd_gamma, log = TRUE)) # data likelihood probs <- plogis(beta0 + gamma) # linear predictor ADREPORT(probs); REPORT(probs) nll <- nll - sum(dbinom(x, 1, probs[id], log = TRUE)) return(nll) } # RTMB Laplace approximation par <- list(beta0 = 1, log_sd_gamma = log(2), gamma = rep(0,10)) obj <- MakeADFun(nll, par, random = "gamma", silent = TRUE) # model fitting using AD gradient opt <- nlminb(obj$par, obj$fn, obj$gr) # draw samples samples <- mcreport(obj, report = TRUE) # run sdreport sdr <- sdreport(obj) # standard deviation using delta method as.list(sdr, "Std", report = TRUE)$probs # standard deviation from samples apply(simplify2array(samples$probs), 1, sd)
Density and and random generation for the multivariate t distribution
dmvt(x, mu, Sigma, df, log = FALSE) rmvt(n, mu, Sigma, df)dmvt(x, mu, Sigma, df, log = FALSE) rmvt(n, mu, Sigma, df)
x |
vector or matrix of quantiles |
mu |
vector or matrix of location parameters (mean if |
Sigma |
positive definite scale matrix (proportional to the covariance matrix if |
df |
degrees of freedom; must be positive |
log |
logical; if |
n |
number of random values to return. |
This implementation of dmvt allows for automatic differentiation with RTMB.
Note: for df the mean is undefined, and for df the covariance is infinite.
For df > 2, the covariance is df/(df-2) * Sigma.
where is the dimension.
dmvt gives the density, rmvt generates random deviates.
# single mu mu <- c(1,2,3) Sigma <- diag(c(1,1,1)) df <- 5 x <- rmvt(2, mu, Sigma, df) d <- dmvt(x, mu, Sigma, df) # vectorised over mu mu <- rbind(c(1,2,3), c(0, 0.5, 1)) x <- rmvt(2, mu, Sigma, df) d <- dmvt(x, mu, Sigma, df)# single mu mu <- c(1,2,3) Sigma <- diag(c(1,1,1)) df <- 5 x <- rmvt(2, mu, Sigma, df) d <- dmvt(x, mu, Sigma, df) # vectorised over mu mu <- rbind(c(1,2,3), c(0, 0.5, 1)) x <- rmvt(2, mu, Sigma, df) d <- dmvt(x, mu, Sigma, df)
Probability mass function, distribution function, quantile function, and random generation for the negative binomial distribution reparameterised in terms of mean and size.
dnbinom2(x, mu, size, log = FALSE) pnbinom2(q, mu, size, lower.tail = TRUE, log.p = FALSE) qnbinom2(p, mu, size, lower.tail = TRUE, log.p = FALSE) rnbinom2(n, mu, size)dnbinom2(x, mu, size, log = FALSE) pnbinom2(q, mu, size, lower.tail = TRUE, log.p = FALSE) qnbinom2(p, mu, size, lower.tail = TRUE, log.p = FALSE) rnbinom2(n, mu, size)
x, q
|
vector of quantiles |
mu |
mean parameter, must be positive. |
size |
size parameter, must be positive. |
log, log.p
|
logical; if |
lower.tail |
logical; if |
p |
vector of probabilities |
n |
number of random values to return. |
This implementation allows for automatic differentiation with RTMB.
pnbinom is an AD-compatible implementation of the standard parameterisation of the CDF, missing from RTMB.
Reparameterises the negative binomial by mean via :
dnbinom2 gives the density, pnbinom2 gives the distribution function, qnbinom2 gives the quantile function, and rnbinom2 generates random deviates.
set.seed(123) x <- rnbinom2(1, 1, 2) d <- dnbinom2(x, 1, 2) p <- pnbinom2(x, 1, 2) q <- qnbinom2(p, 1, 2)set.seed(123) x <- rnbinom2(1, 1, 2) d <- dnbinom2(x, 1, 2) p <- pnbinom2(x, 1, 2) q <- qnbinom2(p, 1, 2)
Density, distribution function, and random generation for the one-inflated beta distribution.
doibeta(x, shape1, shape2, oneprob = 0, log = FALSE) poibeta(q, shape1, shape2, oneprob = 0, lower.tail = TRUE, log.p = FALSE) roibeta(n, shape1, shape2, oneprob = 0)doibeta(x, shape1, shape2, oneprob = 0, log = FALSE) poibeta(q, shape1, shape2, oneprob = 0, lower.tail = TRUE, log.p = FALSE) roibeta(n, shape1, shape2, oneprob = 0)
x, q
|
vector of quantiles |
shape1, shape2
|
non-negative shape parameters of the beta distribution |
oneprob |
one-inflation probability between 0 and 1. |
log, log.p
|
logical; if |
lower.tail |
logical; if |
n |
number of random values to return. |
This implementation allows for automatic differentiation with RTMB.
where is the beta density and is oneprob.
doibeta gives the density, poibeta gives the distribution function, and roibeta generates random deviates.
set.seed(123) x <- roibeta(1, 2, 2, 0.5) d <- doibeta(x, 2, 2, 0.5) p <- poibeta(x, 2, 2, 0.5)set.seed(123) x <- roibeta(1, 2, 2, 0.5) d <- doibeta(x, 2, 2, 0.5) p <- poibeta(x, 2, 2, 0.5)
Density, distribution function, and random generation for the one-inflated beta distribution reparameterised in terms of mean and concentration.
doibeta2(x, mu, phi, oneprob = 0, log = FALSE) poibeta2(q, mu, phi, oneprob = 0, lower.tail = TRUE, log.p = FALSE) roibeta2(n, mu, phi, oneprob = 0)doibeta2(x, mu, phi, oneprob = 0, log = FALSE) poibeta2(q, mu, phi, oneprob = 0, lower.tail = TRUE, log.p = FALSE) roibeta2(n, mu, phi, oneprob = 0)
x, q
|
vector of quantiles |
mu |
mean parameter, must be in the interval from 0 to 1. |
phi |
concentration parameter, must be positive. |
oneprob |
one-inflation probability between 0 and 1. |
log, log.p
|
logical; if |
lower.tail |
logical; if |
n |
number of random values to return. |
This implementation allows for automatic differentiation with RTMB.
Uses the same density as oibeta with and :
doibeta2 gives the density, poibeta2 gives the distribution function, and roibeta2 generates random deviates.
set.seed(123) x <- roibeta2(1, 0.6, 2, 0.5) d <- doibeta2(x, 0.6, 2, 0.5) p <- poibeta2(x, 0.6, 2, 0.5)set.seed(123) x <- roibeta2(1, 0.6, 2, 0.5) d <- doibeta2(x, 0.6, 2, 0.5) p <- poibeta2(x, 0.6, 2, 0.5)
Density, distribution function, quantile function, and random generation for the pareto distribution.
dpareto(x, mu = 1, log = FALSE) ppareto(q, mu = 1, lower.tail = TRUE, log.p = FALSE) qpareto(p, mu = 1, lower.tail = TRUE, log.p = FALSE) rpareto(n, mu = 1)dpareto(x, mu = 1, log = FALSE) ppareto(q, mu = 1, lower.tail = TRUE, log.p = FALSE) qpareto(p, mu = 1, lower.tail = TRUE, log.p = FALSE) rpareto(n, mu = 1)
x, q
|
vector of quantiles |
mu |
location parameter, must be positive. |
log, log.p
|
logical; if |
lower.tail |
logical; if |
p |
vector of probabilities |
n |
number of random values to return |
This implementation of dpareto and ppareto allows for automatic differentiation with RTMB while the other functions are imported from gamlss.dist package.
See gamlss.dist::PARETO for more details.
dpareto gives the density, ppareto gives the distribution function, qpareto gives the quantile function, and rpareto generates random deviates.
Rigby, R. A., Stasinopoulos, D. M., Heller, G. Z., and De Bastiani, F. (2019) Distributions for modeling location, scale, and shape: Using GAMLSS in R, Chapman and Hall/CRC, doi:10.1201/9780429298547. An older version can be found in https://www.gamlss.com/.
set.seed(123) x <- rpareto(1, mu = 5) d <- dpareto(x, mu = 5) p <- ppareto(x, mu = 5) q <- qpareto(p, mu = 5)set.seed(123) x <- rpareto(1, mu = 5) d <- dpareto(x, mu = 5) p <- ppareto(x, mu = 5) q <- qpareto(p, mu = 5)
Survival, hazard, cumulative distribution,
density, quantile and sampling function for the power generalized
Weibull (PgW) distribution with parameters scale, shape and powershape.
spgweibull(q, scale = 1, shape = 1, powershape = 1, log = FALSE) hpgweibull(x, scale = 1, shape = 1, powershape = 1, log = FALSE) ppgweibull(q, scale = 1, shape = 1, powershape = 1, lower.tail = TRUE, log.p = FALSE) dpgweibull(x, scale = 1, shape = 1, powershape = 1, log = FALSE) qpgweibull( p, scale = 1, shape = 1, powershape = 1, lower.tail = TRUE, log.p = FALSE ) rpgweibull(n, scale = 1, shape = 1, powershape = 1)spgweibull(q, scale = 1, shape = 1, powershape = 1, log = FALSE) hpgweibull(x, scale = 1, shape = 1, powershape = 1, log = FALSE) ppgweibull(q, scale = 1, shape = 1, powershape = 1, lower.tail = TRUE, log.p = FALSE) dpgweibull(x, scale = 1, shape = 1, powershape = 1, log = FALSE) qpgweibull( p, scale = 1, shape = 1, powershape = 1, lower.tail = TRUE, log.p = FALSE ) rpgweibull(n, scale = 1, shape = 1, powershape = 1)
scale |
positive scale parameter |
shape |
positive shape parameter |
powershape |
positive power shape parameter |
log, log.p
|
logical; if |
x, q
|
vector of quantiles |
lower.tail |
logical; if |
p |
vector of probabilities |
n |
number of observations |
The survival function of the PgW distribution is:
The hazard function is
The cumulative distribution function is then and the density function
is .
If both shape parameters equal 1, the PgW distribution reduces to the exponential distribution
(see dexp) with
If the power shape parameter equals 1, the PgW distribution simplifies to the Weibull distribution
(see dweibull) with the same parametrization.
dpgweibull gives the density, ppgweibull gives the distribution function, qpgweibull gives the quantile function, and rpgweibull generates random deviates.
spgweibull gives the survival function and hpgweibull gives the hazard function.
x <- rpgweibull(1, 2, 2, 3) d <- dpgweibull(x, 2, 2, 3) p <- ppgweibull(x, 2, 2, 3) q <- qpgweibull(p, 2, 2, 3) s <- spgweibull(x, 2, 2, 3) h <- hpgweibull(x, 2, 2, 3)x <- rpgweibull(1, 2, 2, 3) d <- dpgweibull(x, 2, 2, 3) p <- ppgweibull(x, 2, 2, 3) q <- qpgweibull(p, 2, 2, 3) s <- spgweibull(x, 2, 2, 3) h <- hpgweibull(x, 2, 2, 3)
Density, distribution function, quantile function, and random generation for the Power Exponential distribution (two versions).
dpowerexp(x, mu = 0, sigma = 1, nu = 2, log = FALSE) ppowerexp(q, mu = 0, sigma = 1, nu = 2, lower.tail = TRUE, log.p = FALSE) qpowerexp(p, mu = 0, sigma = 1, nu = 2, lower.tail = TRUE, log.p = FALSE) rpowerexp(n, mu = 0, sigma = 1, nu = 2) dpowerexp2(x, mu = 0, sigma = 1, nu = 2, log = FALSE) ppowerexp2(q, mu = 0, sigma = 1, nu = 2, lower.tail = TRUE, log.p = FALSE) qpowerexp2(p, mu = 0, sigma = 1, nu = 2, lower.tail = TRUE, log.p = FALSE) rpowerexp2(n, mu = 0, sigma = 1, nu = 2)dpowerexp(x, mu = 0, sigma = 1, nu = 2, log = FALSE) ppowerexp(q, mu = 0, sigma = 1, nu = 2, lower.tail = TRUE, log.p = FALSE) qpowerexp(p, mu = 0, sigma = 1, nu = 2, lower.tail = TRUE, log.p = FALSE) rpowerexp(n, mu = 0, sigma = 1, nu = 2) dpowerexp2(x, mu = 0, sigma = 1, nu = 2, log = FALSE) ppowerexp2(q, mu = 0, sigma = 1, nu = 2, lower.tail = TRUE, log.p = FALSE) qpowerexp2(p, mu = 0, sigma = 1, nu = 2, lower.tail = TRUE, log.p = FALSE) rpowerexp2(n, mu = 0, sigma = 1, nu = 2)
x, q
|
vector of quantiles |
mu |
location parameter |
sigma |
scale parameter, must be positive |
nu |
shape parameter (real) |
log, log.p
|
logical; if |
lower.tail |
logical; if |
p |
vector of probabilities |
n |
number of random values to return |
This implementation of the densities and distribution functions allow for automatic differentiation with RTMB while the other functions are imported from gamlss.dist package.
For powerexp, mu is the mean and sigma is the standard deviation while this does not hold for powerexp2.
See gamlss.dist::PE for more details.
For powerexp (PE), is the standard deviation; the density is
where .
For powerexp2 (PE2), is a scale parameter; the density is
dpowerexp gives the density, ppowerexp gives the distribution function, qpowerexp gives the quantile function, and rpowerexp generates random deviates.
Rigby, R. A., Stasinopoulos, D. M., Heller, G. Z., and De Bastiani, F. (2019) Distributions for modeling location, scale, and shape: Using GAMLSS in R, Chapman and Hall/CRC, doi:10.1201/9780429298547. An older version can be found in https://www.gamlss.com/.
# PE x <- rpowerexp(1, mu = 0, sigma = 1, nu = 2) d <- dpowerexp(x, mu = 0, sigma = 1, nu = 2) p <- ppowerexp(x, mu = 0, sigma = 1, nu = 2) q <- qpowerexp(p, mu = 0, sigma = 1, nu = 2) # PE2 x <- rpowerexp2(1, mu = 0, sigma = 1, nu = 2) d <- dpowerexp2(x, mu = 0, sigma = 1, nu = 2) p <- ppowerexp2(x, mu = 0, sigma = 1, nu = 2) q <- qpowerexp2(p, mu = 0, sigma = 1, nu = 2)# PE x <- rpowerexp(1, mu = 0, sigma = 1, nu = 2) d <- dpowerexp(x, mu = 0, sigma = 1, nu = 2) p <- ppowerexp(x, mu = 0, sigma = 1, nu = 2) q <- qpowerexp(p, mu = 0, sigma = 1, nu = 2) # PE2 x <- rpowerexp2(1, mu = 0, sigma = 1, nu = 2) d <- dpowerexp2(x, mu = 0, sigma = 1, nu = 2) p <- ppowerexp2(x, mu = 0, sigma = 1, nu = 2) q <- qpowerexp2(p, mu = 0, sigma = 1, nu = 2)
Draw samples from a multivariate Gaussian distribution specified by a sparse precision matrix. This is numerically efficient for high-dimensional but sparse systems.
rgmrf(n, mean = 0, Q)rgmrf(n, mean = 0, Q)
n |
Number of samples to draw. |
mean |
Mean vector (or scalar, which will be recycled to match the dimension of |
Q |
Sparse precision matrix ( |
A matrix of samples with rows corresponding to samples and columns to dimensions.
rgmrf(3, mean = c(1, 2, 3), Q = Matrix::Diagonal(3))rgmrf(3, mean = c(1, 2, 3), Q = Matrix::Diagonal(3))
Probability mass function, distribution function, and random generation for the Skellam distribution.
dskellam(x, mu1, mu2, log = FALSE) rskellam(n, mu1, mu2)dskellam(x, mu1, mu2, log = FALSE) rskellam(n, mu1, mu2)
x |
integer vector of counts |
mu1, mu2
|
Poisson means |
log |
logical; return log-density if TRUE |
n |
number of random values to return. |
The Skellam distribution is the distribution of the difference of two Poisson random variables. Specifically, if and , then .
This implementation of dskellam allows for automatic differentiation with RTMB.
dskellam gives the probability mass function and rskellam generates random deviates.
x <- rskellam(1, 2, 3) d <- dskellam(x, 2, 3)x <- rskellam(1, 2, 3) d <- dskellam(x, 2, 3)
Density, distribution function, quantile function, and random generation for the skew normal distribution.
dskewnorm(x, xi = 0, omega = 1, alpha = 0, log = FALSE) pskewnorm(q, xi = 0, omega = 1, alpha = 0, lower.tail = TRUE, log.p = FALSE) qskewnorm(p, xi = 0, omega = 1, alpha = 0, lower.tail = TRUE, log.p = FALSE) rskewnorm(n, xi = 0, omega = 1, alpha = 0)dskewnorm(x, xi = 0, omega = 1, alpha = 0, log = FALSE) pskewnorm(q, xi = 0, omega = 1, alpha = 0, lower.tail = TRUE, log.p = FALSE) qskewnorm(p, xi = 0, omega = 1, alpha = 0, lower.tail = TRUE, log.p = FALSE) rskewnorm(n, xi = 0, omega = 1, alpha = 0)
x, q
|
vector of quantiles |
xi |
location parameter |
omega |
scale parameter, must be positive. |
alpha |
skewness parameter, +/- |
log, log.p
|
logical; if |
lower.tail |
logical; if |
p |
vector of probabilities |
n |
number of random values to return |
This implementation of dskewnorm allows for automatic differentiation with RTMB while the other functions are imported from the sn package.
See sn::dsn for more details.
where and are the standard normal PDF and CDF.
dskewnorm gives the density, pskewnorm gives the distribution function, qskewnorm gives the quantile function, and rskewnorm generates random deviates.
# alpha is skew parameter x <- rskewnorm(1, alpha = 1) d <- dskewnorm(x, alpha = 1) p <- pskewnorm(x, alpha = 1) q <- qskewnorm(p, alpha = 1)# alpha is skew parameter x <- rskewnorm(1, alpha = 1) d <- dskewnorm(x, alpha = 1) p <- pskewnorm(x, alpha = 1) q <- qskewnorm(p, alpha = 1)
Density, distribution function, quantile function and random generation for the skew normal distribution reparameterised in terms of mean, standard deviation and skew magnitude
dskewnorm2(x, mean = 0, sd = 1, alpha = 0, log = FALSE) pskewnorm2(q, mean = 0, sd = 1, alpha = 0, lower.tail = TRUE, log.p = FALSE) qskewnorm2(p, mean = 0, sd = 1, alpha = 0, lower.tail = TRUE, log.p = FALSE) rskewnorm2(n, mean = 0, sd = 1, alpha = 0)dskewnorm2(x, mean = 0, sd = 1, alpha = 0, log = FALSE) pskewnorm2(q, mean = 0, sd = 1, alpha = 0, lower.tail = TRUE, log.p = FALSE) qskewnorm2(p, mean = 0, sd = 1, alpha = 0, lower.tail = TRUE, log.p = FALSE) rskewnorm2(n, mean = 0, sd = 1, alpha = 0)
x, q
|
vector of quantiles |
mean |
mean parameter |
sd |
standard deviation, must be positive. |
alpha |
skewness parameter, +/- |
log, log.p
|
logical; if |
lower.tail |
logical; if |
p |
vector of probabilities |
n |
number of random values to return |
This implementation of dskewnorm2 allows for automatic differentiation with RTMB while the other functions are imported from the sn package.
Uses the same density as skewnorm with location , scale , and shape
reparameterised from the mean , standard deviation , and skewness parameter :
dskewnorm2 gives the density, pskewnorm2 gives the distribution function, qskewnorm2 gives the quantile function, and rskewnorm2 generates random deviates.
# alpha is skew parameter x <- rskewnorm2(1, alpha = 1) d <- dskewnorm2(x, alpha = 1) p <- pskewnorm2(x, alpha = 1) q <- qskewnorm2(p, alpha = 1)# alpha is skew parameter x <- rskewnorm2(1, alpha = 1) d <- dskewnorm2(x, alpha = 1) p <- pskewnorm2(x, alpha = 1) q <- qskewnorm2(p, alpha = 1)
Density, distribution function, quantile function, and random generation for the skew t distribution (type 2).
dskewt(x, mu = 0, sigma = 1, skew = 0, df = 100, log = FALSE) pskewt(q, mu = 0, sigma = 1, skew = 0, df = 100, method = 0, lower.tail = TRUE, log.p = FALSE) qskewt(p, mu = 0, sigma = 1, skew = 0, df = 100, tol = 1e-8, method = 0, lower.tail = TRUE, log.p = FALSE) rskewt(n, mu = 0, sigma = 1, skew = 0, df = 100)dskewt(x, mu = 0, sigma = 1, skew = 0, df = 100, log = FALSE) pskewt(q, mu = 0, sigma = 1, skew = 0, df = 100, method = 0, lower.tail = TRUE, log.p = FALSE) qskewt(p, mu = 0, sigma = 1, skew = 0, df = 100, tol = 1e-8, method = 0, lower.tail = TRUE, log.p = FALSE) rskewt(n, mu = 0, sigma = 1, skew = 0, df = 100)
x, q
|
vector of quantiles |
mu |
location parameter |
sigma |
scale parameter, must be positive. |
skew |
skewness parameter, can be positive or negative. |
df |
degrees of freedom, must be positive. |
log, log.p
|
logical; if |
method |
an integer value between 0 and 5 which selects the computing method; see ‘Details’ in the |
lower.tail |
logical; if |
p |
vector of probabilities |
tol |
a scalar value which regulates the accuracy of the result of qsn, measured on the probability scale. |
n |
number of random values to return. |
This corresponds to the skew t type 2 distribution in GAMLSS (ST2), see pp. 411-412 of Rigby et al. (2019) and the version implemented in the sn package.
This implementation of dskewt allows for automatic differentiation with RTMB while the other functions are imported from the sn package.
See sn::dst for more details.
Caution: In a numerial optimisation, the skew parameter should NEVER be initialised with exactly zero.
This will cause the initial and all subsequent derivatives to be exactly zero and hence the parameter will remain at its initial value.
where , is the Student- PDF, and is its CDF.
dskewt gives the density, pskewt gives the distribution function, qskewt gives the quantile function, and rskewt generates random deviates.
x <- rskewt(1, 1, 2, 5, 2) d <- dskewt(x, 1, 2, 5, 2) p <- pskewt(x, 1, 2, 5, 2) q <- qskewt(p, 1, 2, 5, 2)x <- rskewt(1, 1, 2, 5, 2) d <- dskewt(x, 1, 2, 5, 2) p <- pskewt(x, 1, 2, 5, 2) q <- qskewt(p, 1, 2, 5, 2)
Density, distribution function, quantile function, and random generation
for the skew t distribution reparameterised so that mean and sd
correspond to the true mean and standard deviation.
dskewt2(x, mean = 0, sd = 1, skew = 0, df = 100, log = FALSE) pskewt2(q, mean = 0, sd = 1, skew = 0, df = 100, method = 0, lower.tail = TRUE, log.p = FALSE) qskewt2( p, mean = 0, sd = 1, skew = 0, df = 100, tol = 1e-08, method = 0, lower.tail = TRUE, log.p = FALSE ) rskewt2(n, mean = 0, sd = 1, skew = 0, df = 100)dskewt2(x, mean = 0, sd = 1, skew = 0, df = 100, log = FALSE) pskewt2(q, mean = 0, sd = 1, skew = 0, df = 100, method = 0, lower.tail = TRUE, log.p = FALSE) qskewt2( p, mean = 0, sd = 1, skew = 0, df = 100, tol = 1e-08, method = 0, lower.tail = TRUE, log.p = FALSE ) rskewt2(n, mean = 0, sd = 1, skew = 0, df = 100)
x, q
|
vector of quantiles |
mean |
mean parameter |
sd |
standard deviation parameter, must be positive. |
skew |
skewness parameter, can be positive or negative. |
df |
degrees of freedom, must be positive. |
log, log.p
|
logical; if |
method |
an integer value between 0 and 5 which selects the computing method; see ‘Details’ in the |
lower.tail |
logical; if |
p |
vector of probabilities |
tol |
a scalar value which regulates the accuracy of the result of qsn, measured on the probability scale. |
n |
number of random values to return. |
This corresponds to the skew t type 2 distribution in GAMLSS (ST2), see pp. 411-412 of Rigby et al. (2019) and the version implemented in the sn package.
However, it is reparameterised in terms of a standard deviation parameter sd rather than just a scale parameter sigma. Details of this reparameterisation are given below.
This implementation of dskewt allows for automatic differentiation with RTMB while the other functions are imported from the sn package.
See sn::dst for more details.
Caution: In a numerial optimisation, the skew parameter should NEVER be initialised with exactly zero.
This will cause the initial and all subsequent derivatives to be exactly zero and hence the parameter will remain at its initial value.
For given skew and df = , define
then
dskewt2 gives the density, pskewt2 gives the distribution function, qskewt2 gives the quantile function, and rskewt2 generates random deviates.
x <- rskewt2(1, 1, 2, 5, 5) d <- dskewt2(x, 1, 2, 5, 5) p <- pskewt2(x, 1, 2, 5, 5) q <- qskewt2(p, 1, 2, 5, 5)x <- rskewt2(1, 1, 2, 5, 5) d <- dskewt2(x, 1, 2, 5, 5) p <- pskewt2(x, 1, 2, 5, 5) q <- qskewt2(p, 1, 2, 5, 5)
Density, distribution function, quantile function, and random generation for the t distribution with location and scale parameters.
dt2(x, mu, sigma, df, log = FALSE) pt2(q, mu, sigma, df, lower.tail = TRUE, log.p = FALSE) rt2(n, mu, sigma, df) qt2(p, mu, sigma, df, lower.tail = TRUE, log.p = FALSE) pt(q, df)dt2(x, mu, sigma, df, log = FALSE) pt2(q, mu, sigma, df, lower.tail = TRUE, log.p = FALSE) rt2(n, mu, sigma, df) qt2(p, mu, sigma, df, lower.tail = TRUE, log.p = FALSE) pt(q, df)
x, q
|
vector of quantiles |
mu |
location parameter |
sigma |
scale parameter, must be positive. |
df |
degrees of freedom, must be positive. |
log, log.p
|
logical; if |
lower.tail |
logical; if |
n |
number of random values to return. |
p |
vector of probabilities |
This implementation of dt2 allows for automatic differentiation with RTMB.
where is the Student- PDF with degrees of freedom.
dt2 gives the density, pt2 gives the distribution function, qt2 gives the quantile function, and rt2 generates random deviates.
x <- rt2(1, 1, 2, 5) d <- dt2(x, 1, 2, 5) p <- pt2(x, 1, 2, 5) q <- qt2(p, 1, 2, 5)x <- rt2(1, 1, 2, 5) d <- dt2(x, 1, 2, 5) p <- pt2(x, 1, 2, 5) q <- qt2(p, 1, 2, 5)
Density, distribution function, quantile function, and random generation for the truncated normal distribution.
dtruncnorm(x, mean = 0, sd = 1, min = -Inf, max = Inf, log = FALSE) ptruncnorm(q, mean = 0, sd = 1, min = -Inf, max = Inf, lower.tail = TRUE, log.p = FALSE) qtruncnorm(p, mean = 0, sd = 1, min = -Inf, max = Inf, lower.tail = TRUE, log.p = FALSE) rtruncnorm(n, mean = 0, sd = 1, min = -Inf, max = Inf)dtruncnorm(x, mean = 0, sd = 1, min = -Inf, max = Inf, log = FALSE) ptruncnorm(q, mean = 0, sd = 1, min = -Inf, max = Inf, lower.tail = TRUE, log.p = FALSE) qtruncnorm(p, mean = 0, sd = 1, min = -Inf, max = Inf, lower.tail = TRUE, log.p = FALSE) rtruncnorm(n, mean = 0, sd = 1, min = -Inf, max = Inf)
x, q
|
vector of quantiles |
mean |
mean parameter, must be positive. |
sd |
standard deviation parameter, must be positive. |
min, max
|
truncation bounds. |
log, log.p
|
logical; if |
lower.tail |
logical; if |
p |
vector of probabilities |
n |
number of random values to return. |
This implementation of dtruncnorm allows for automatic differentiation with RTMB.
where and are the standard normal PDF and CDF.
dtruncnorm gives the density, ptruncnorm gives the distribution function, qtruncnorm gives the quantile function, and rtruncnorm generates random deviates.
x <- rtruncnorm(1, mean = 2, sd = 2, min = -1, max = 5) d <- dtruncnorm(x, mean = 2, sd = 2, min = -1, max = 5) p <- ptruncnorm(x, mean = 2, sd = 2, min = -1, max = 5) q <- qtruncnorm(p, mean = 2, sd = 2, min = -1, max = 5)x <- rtruncnorm(1, mean = 2, sd = 2, min = -1, max = 5) d <- dtruncnorm(x, mean = 2, sd = 2, min = -1, max = 5) p <- ptruncnorm(x, mean = 2, sd = 2, min = -1, max = 5) q <- qtruncnorm(p, mean = 2, sd = 2, min = -1, max = 5)
Density, distribution function, quantile function, and random generation for the truncated t distribution.
dtrunct(x, df, min = -Inf, max = Inf, log = FALSE) ptrunct(q, df, min = -Inf, max = Inf, lower.tail = TRUE, log.p = FALSE) qtrunct(p, df, min = -Inf, max = Inf, lower.tail = TRUE, log.p = FALSE) rtrunct(n, df, min = -Inf, max = Inf)dtrunct(x, df, min = -Inf, max = Inf, log = FALSE) ptrunct(q, df, min = -Inf, max = Inf, lower.tail = TRUE, log.p = FALSE) qtrunct(p, df, min = -Inf, max = Inf, lower.tail = TRUE, log.p = FALSE) rtrunct(n, df, min = -Inf, max = Inf)
x, q
|
vector of quantiles |
df |
degrees of freedom parameter, must be positive. |
min, max
|
truncation bounds. |
log, log.p
|
logical; if |
lower.tail |
logical; if |
p |
vector of probabilities |
n |
number of random values to return. |
This implementation of dtrunct allows for automatic differentiation with RTMB.
where and are the Student- PDF and CDF.
dtrunct gives the density, ptrunct gives the distribution function,
qtrunct gives the quantile function, and rtrunct generates random deviates.
x <- rtrunct(1, df = 5, min = -1, max = 5) d <- dtrunct(x, df = 5, min = -1, max = 5) p <- ptrunct(x, df = 5, min = -1, max = 5) q <- qtrunct(p, df = 5, min = -1, max = 5)x <- rtrunct(1, df = 5, min = -1, max = 5) d <- dtrunct(x, df = 5, min = -1, max = 5) p <- ptrunct(x, df = 5, min = -1, max = 5) q <- qtrunct(p, df = 5, min = -1, max = 5)
Density, distribution function, quantile function, and random generation for
the truncated t distribution with location mu and scale sigma.
dtrunct2(x, df, mu = 0, sigma = 1, min = -Inf, max = Inf, log = FALSE) ptrunct2(q, df, mu = 0, sigma = 1, min = -Inf, max = Inf, lower.tail = TRUE, log.p = FALSE) qtrunct2(p, df, mu = 0, sigma = 1, min = -Inf, max = Inf, lower.tail = TRUE, log.p = FALSE) rtrunct2(n, df, mu = 0, sigma = 1, min = -Inf, max = Inf)dtrunct2(x, df, mu = 0, sigma = 1, min = -Inf, max = Inf, log = FALSE) ptrunct2(q, df, mu = 0, sigma = 1, min = -Inf, max = Inf, lower.tail = TRUE, log.p = FALSE) qtrunct2(p, df, mu = 0, sigma = 1, min = -Inf, max = Inf, lower.tail = TRUE, log.p = FALSE) rtrunct2(n, df, mu = 0, sigma = 1, min = -Inf, max = Inf)
x, q
|
vector of quantiles |
df |
degrees of freedom parameter, must be positive. |
mu |
location parameter. |
sigma |
scale parameter, must be positive. |
min, max
|
truncation bounds. |
log, log.p
|
logical; if |
lower.tail |
logical; if |
p |
vector of probabilities |
n |
number of random values to return. |
This implementation of dtrunct2 allows for automatic differentiation with RTMB.
dtrunct2 gives the density, ptrunct2 gives the distribution function,
qtrunct2 gives the quantile function, and rtrunct2 generates random deviates.
x <- rtrunct2(1, df = 5, mu = 2, sigma = 3, min = -1, max = 5) d <- dtrunct2(x, df = 5, mu = 2, sigma = 3, min = -1, max = 5) p <- ptrunct2(x, df = 5, mu = 2, sigma = 3, min = -1, max = 5) q <- qtrunct2(p, df = 5, mu = 2, sigma = 3, min = -1, max = 5)x <- rtrunct2(1, df = 5, mu = 2, sigma = 3, min = -1, max = 5) d <- dtrunct2(x, df = 5, mu = 2, sigma = 3, min = -1, max = 5) p <- ptrunct2(x, df = 5, mu = 2, sigma = 3, min = -1, max = 5) q <- qtrunct2(p, df = 5, mu = 2, sigma = 3, min = -1, max = 5)
Density, distribution function, and random generation for the von Mises distribution.
dvm(x, mu = 0, kappa = 1, log = FALSE) pvm( q, mu = 0, kappa = 1, from = NULL, tol = 1e-20, lower.tail = TRUE, log.p = FALSE ) rvm(n, mu = 0, kappa = 1, wrap = TRUE)dvm(x, mu = 0, kappa = 1, log = FALSE) pvm( q, mu = 0, kappa = 1, from = NULL, tol = 1e-20, lower.tail = TRUE, log.p = FALSE ) rvm(n, mu = 0, kappa = 1, wrap = TRUE)
x, q
|
vector of angles measured in radians at which to evaluate the density function. |
mu |
mean direction of the distribution measured in radians. |
kappa |
non-negative numeric value for the concentration parameter of the distribution. |
log |
logical; if |
from |
value from which the integration for CDF starts. If |
tol |
the precision in evaluating the distribution function |
lower.tail |
logical; if |
log.p |
logical; if |
n |
number of random values to return. |
wrap |
logical; if |
This implementation of dvm allows for automatic differentiation with RTMB.
rvm and pvm are simply wrappers of the corresponding functions from circular.
where is the modified Bessel function of the first kind of order 0.
dvm gives the density, pvm gives the distribution function, and rvm generates random deviates.
set.seed(1) x <- rvm(10, 0, 1) d <- dvm(x, 0, 1) p <- pvm(x, 0, 1)set.seed(1) x <- rvm(10, 0, 1) d <- dvm(x, 0, 1) p <- pvm(x, 0, 1)
Density, distribution function, and random generation for the von Mises-Fisher distribution.
dvmf(x, mu, kappa, log = FALSE) rvmf(n, mu, kappa)dvmf(x, mu, kappa, log = FALSE) rvmf(n, mu, kappa)
x |
unit vector or matrix (with each row being a unit vector) of evaluation points |
mu |
unit mean vector |
kappa |
non-negative numeric value for the concentration parameter of the distribution. |
log |
logical; if |
n |
number of random values to return. |
This implementation of dvmf allows for automatic differentiation with RTMB. rvmf is a reparameterised import from movMF::rmovMF.
where is the normalising constant and is the modified Bessel function of the first kind of order .
dvmf gives the density and rvm generates random deviates.
set.seed(123) # single parameter set mu <- rep(1, 3) / sqrt(3) kappa <- 4 x <- rvmf(1, mu, kappa) d <- dvmf(x, mu, kappa) # vectorised over parameters mu <- matrix(mu, nrow = 1) mu <- mu[rep(1,10), ] kappa <- rep(kappa, 10) x <- rvmf(10, mu, kappa) d <- dvmf(x, mu, kappa)set.seed(123) # single parameter set mu <- rep(1, 3) / sqrt(3) kappa <- 4 x <- rvmf(1, mu, kappa) d <- dvmf(x, mu, kappa) # vectorised over parameters mu <- matrix(mu, nrow = 1) mu <- mu[rep(1,10), ] kappa <- rep(kappa, 10) x <- rvmf(10, mu, kappa) d <- dvmf(x, mu, kappa)
Density, distribution function, and random generation for the von Mises-Fisher distribution.
dvmf2(x, theta, log = FALSE) rvmf2(n, theta)dvmf2(x, theta, log = FALSE) rvmf2(n, theta)
x |
unit vector or matrix (with each row being a unit vector) of evaluation points |
theta |
direction and concentration vector. The direction of |
log |
logical; if |
n |
number of random values to return. |
In this parameterisation, , where is a unit vector and is the concentration parameter.
dvmf2 allows for automatic differentiation with RTMB. rvmf2 is imported from movMF::rmovMF.
dvmf gives the density and rvm generates random deviates.
set.seed(123) # single parameter set theta <- c(1,2,3) x <- rvmf2(1, theta) d <- dvmf2(x, theta) # vectorised over parameters theta <- matrix(theta, nrow = 1) theta <- theta[rep(1,10), ] x <- rvmf2(10, theta) d <- dvmf2(x, theta)set.seed(123) # single parameter set theta <- c(1,2,3) x <- rvmf2(1, theta) d <- dvmf2(x, theta) # vectorised over parameters theta <- matrix(theta, nrow = 1) theta <- theta[rep(1,10), ] x <- rvmf2(10, theta) d <- dvmf2(x, theta)
Density and random generation for the wishart distribution
dwishart(x, nu, Sigma, log = FALSE) rwishart(n, nu, Sigma)dwishart(x, nu, Sigma, log = FALSE) rwishart(n, nu, Sigma)
x |
positive definite |
nu |
degrees of freedom, needs to be greater than |
Sigma |
scale matrix, needs to be positive definite and match the dimension of |
log |
logical; if |
n |
number of random deviates to return |
where is the dimension and is the multivariate gamma function.
dwishart gives the density,
rwishart generates random deviates (matrix for n = 1, array with n slices for n > 1)
# single input: matrix-valued x <- rwishart(1, nu = 5, Sigma = diag(3)) d <- dwishart(x, nu = 5, Sigma = diag(3)) # multiple inputs: array of matrices x <- rwishart(4, nu = 5, Sigma = diag(3)) d <- dwishart(x, nu = 5, Sigma = diag(3)) # multiple inputs for x, nu and Sigma nu <- c(7,5,8,9) Sigma <- array(dim = c(3,3,4)) for(i in 1:4) Sigma[,,i] <- (i + 3) * diag(3) x <- rwishart(4, nu, Sigma) d <- dwishart(x, nu, Sigma)# single input: matrix-valued x <- rwishart(1, nu = 5, Sigma = diag(3)) d <- dwishart(x, nu = 5, Sigma = diag(3)) # multiple inputs: array of matrices x <- rwishart(4, nu = 5, Sigma = diag(3)) d <- dwishart(x, nu = 5, Sigma = diag(3)) # multiple inputs for x, nu and Sigma nu <- c(7,5,8,9) Sigma <- array(dim = c(3,3,4)) for(i in 1:4) Sigma[,,i] <- (i + 3) * diag(3) x <- rwishart(4, nu, Sigma) d <- dwishart(x, nu, Sigma)
Density and random generation for the wrapped Cauchy distribution.
dwrpcauchy(x, mu = 0, rho, log = FALSE) rwrpcauchy(n, mu = 0, rho, wrap = TRUE)dwrpcauchy(x, mu = 0, rho, log = FALSE) rwrpcauchy(n, mu = 0, rho, wrap = TRUE)
x |
vector of angles measured in radians at which to evaluate the density function. |
mu |
mean direction of the distribution measured in radians. |
rho |
concentration parameter of the distribution, must be in the interval from 0 to 1. |
log |
logical; if |
n |
number of random values to return. |
wrap |
logical; if |
This implementation of dwrpcauchy allows for automatic differentiation with RTMB.
rwrpcauchy is simply a wrapper for rwrappedcauchyimported from circular.
dwrpcauchy gives the density and rwrpcauchy generates random deviates.
set.seed(1) x <- rwrpcauchy(10, 0, 0.5) d <- dwrpcauchy(x, 0, 0.5)set.seed(1) x <- rwrpcauchy(10, 0, 0.5) d <- dwrpcauchy(x, 0, 0.5)
Constructs a zero-inflated density function from a given probability density function
zero_inflate(dist, discrete = NULL)zero_inflate(dist, discrete = NULL)
dist |
either a probability density function or a probability mass function |
discrete |
logical; if |
The definition of zero-inflation is different for discrete and continuous distributions.
For discrete distributions with p.m.f. and zero-inflation probability , we have
and
For continuous distributions with p.d.f. , we have
where is the Dirac delta function at zero.
zero-inflated density function with first argument x, second argument zeroprob, and additional arguments ... that will be passed to dist.
# Zero-inflated normal distribution dzinorm <- zero_inflate(dnorm) dzinorm(c(NA, 0, 2), 0.5, mean = 1, sd = 1) # Zero-inflated Poisson distribution zipois <- zero_inflate(dpois) zipois(c(NA, 0, 1), 0.5, 1) # Non-standard case: Zero-inflated reparametrised beta distribution dzibeta2 <- zero_inflate(dbeta2, discrete = FALSE)# Zero-inflated normal distribution dzinorm <- zero_inflate(dnorm) dzinorm(c(NA, 0, 2), 0.5, mean = 1, sd = 1) # Zero-inflated Poisson distribution zipois <- zero_inflate(dpois) zipois(c(NA, 0, 1), 0.5, 1) # Non-standard case: Zero-inflated reparametrised beta distribution dzibeta2 <- zero_inflate(dbeta2, discrete = FALSE)
Density, distribution function, and random generation for the zero-inflated beta distribution.
dzibeta(x, shape1, shape2, zeroprob = 0, log = FALSE) pzibeta(q, shape1, shape2, zeroprob = 0, lower.tail = TRUE, log.p = FALSE) rzibeta(n, shape1, shape2, zeroprob = 0)dzibeta(x, shape1, shape2, zeroprob = 0, log = FALSE) pzibeta(q, shape1, shape2, zeroprob = 0, lower.tail = TRUE, log.p = FALSE) rzibeta(n, shape1, shape2, zeroprob = 0)
x, q
|
vector of quantiles |
shape1, shape2
|
non-negative shape parameters of the beta distribution |
zeroprob |
zero-inflation probability between 0 and 1. |
log, log.p
|
logical; if |
lower.tail |
logical; if |
n |
number of random values to return. |
This implementation allows for automatic differentiation with RTMB.
where is zeroprob.
dzibeta gives the density, pzibeta gives the distribution function, and rzibeta generates random deviates.
set.seed(123) x <- rzibeta(1, 2, 2, 0.5) d <- dzibeta(x, 2, 2, 0.5) p <- pzibeta(x, 2, 2, 0.5)set.seed(123) x <- rzibeta(1, 2, 2, 0.5) d <- dzibeta(x, 2, 2, 0.5) p <- pzibeta(x, 2, 2, 0.5)
Density, distribution function, and random generation for the zero-inflated beta distribution reparameterised in terms of mean and concentration.
dzibeta2(x, mu, phi, zeroprob = 0, log = FALSE) pzibeta2(q, mu, phi, zeroprob = 0, lower.tail = TRUE, log.p = FALSE) rzibeta2(n, mu, phi, zeroprob = 0)dzibeta2(x, mu, phi, zeroprob = 0, log = FALSE) pzibeta2(q, mu, phi, zeroprob = 0, lower.tail = TRUE, log.p = FALSE) rzibeta2(n, mu, phi, zeroprob = 0)
x, q
|
vector of quantiles |
mu |
mean parameter, must be in the interval from 0 to 1. |
phi |
concentration parameter, must be positive. |
zeroprob |
zero-inflation probability between 0 and 1. |
log, log.p
|
logical; if |
lower.tail |
logical; if |
n |
number of random values to return. |
p |
vector of probabilities |
This implementation allows for automatic differentiation with RTMB.
Uses the same density as zibeta with and :
dzibeta2 gives the density, pzibeta2 gives the distribution function, and rzibeta2 generates random deviates.
set.seed(123) x <- rzibeta2(1, 0.5, 1, 0.5) d <- dzibeta2(x, 0.5, 1, 0.5) p <- pzibeta2(x, 0.5, 1, 0.5)set.seed(123) x <- rzibeta2(1, 0.5, 1, 0.5) d <- dzibeta2(x, 0.5, 1, 0.5) p <- pzibeta2(x, 0.5, 1, 0.5)
Probability mass function, distribution function, and random generation for the zero-inflated binomial distribution.
dzibinom(x, size, prob, zeroprob = 0, log = FALSE) pzibinom(q, size, prob, zeroprob = 0, lower.tail = TRUE, log.p = FALSE) rzibinom(n, size, prob, zeroprob = 0)dzibinom(x, size, prob, zeroprob = 0, log = FALSE) pzibinom(q, size, prob, zeroprob = 0, lower.tail = TRUE, log.p = FALSE) rzibinom(n, size, prob, zeroprob = 0)
x, q
|
vector of quantiles |
size |
number of trials (zero or more). |
prob |
probability of success on each trial. |
zeroprob |
zero-inflation probability between 0 and 1 |
log, log.p
|
logical; return log-density if TRUE |
lower.tail |
logical; if |
n |
number of random values to return. |
This implementation allows for automatic differentiation with RTMB.
where is zeroprob.
dzibinom gives the probability mass function, pzibinom gives the distribution function, and rzibinom generates random deviates.
set.seed(123) x <- rzibinom(1, size = 10, prob = 0.5, zeroprob = 0.5) d <- dzibinom(x, size = 10, prob = 0.5, zeroprob = 0.5) p <- pzibinom(x, size = 10, prob = 0.5, zeroprob = 0.5)set.seed(123) x <- rzibinom(1, size = 10, prob = 0.5, zeroprob = 0.5) d <- dzibinom(x, size = 10, prob = 0.5, zeroprob = 0.5) p <- pzibinom(x, size = 10, prob = 0.5, zeroprob = 0.5)
Density, distribution function, and random generation for the zero-inflated gamma distribution.
dzigamma(x, shape, scale, zeroprob = 0, log = FALSE) pzigamma(q, shape, scale, zeroprob = 0, lower.tail = TRUE, log.p = FALSE) rzigamma(n, shape, scale, zeroprob = 0)dzigamma(x, shape, scale, zeroprob = 0, log = FALSE) pzigamma(q, shape, scale, zeroprob = 0, lower.tail = TRUE, log.p = FALSE) rzigamma(n, shape, scale, zeroprob = 0)
x, q
|
vector of quantiles |
shape |
positive shape parameter |
scale |
positive scale parameter |
zeroprob |
zero-inflation probability between 0 and 1. |
log, log.p
|
logical; if |
lower.tail |
logical; if |
n |
number of random values to return |
This implementation allows for automatic differentiation with RTMB.
where is zeroprob.
dzigamma gives the density, pzigamma gives the distribution function, and rzigamma generates random deviates.
x <- rzigamma(1, 1, 1, 0.5) d <- dzigamma(x, 1, 1, 0.5) p <- pzigamma(x, 1, 1, 0.5)x <- rzigamma(1, 1, 1, 0.5) d <- dzigamma(x, 1, 1, 0.5) p <- pzigamma(x, 1, 1, 0.5)
Density, distribution function, and random generation for the zero-inflated gamma distribution reparameterised in terms of mean and standard deviation.
dzigamma2(x, mean = 1, sd = 1, zeroprob = 0, log = FALSE) pzigamma2(q, mean = 1, sd = 1, zeroprob = 0, lower.tail = TRUE, log.p = FALSE) rzigamma2(n, mean = 1, sd = 1, zeroprob = 0)dzigamma2(x, mean = 1, sd = 1, zeroprob = 0, log = FALSE) pzigamma2(q, mean = 1, sd = 1, zeroprob = 0, lower.tail = TRUE, log.p = FALSE) rzigamma2(n, mean = 1, sd = 1, zeroprob = 0)
x, q
|
vector of quantiles |
mean |
mean parameter, must be positive. |
sd |
standard deviation parameter, must be positive. |
zeroprob |
zero-inflation probability between 0 and 1. |
log, log.p
|
logical; if |
lower.tail |
logical; if |
n |
number of random values to return |
This implementation allows for automatic differentiation with RTMB.
Uses the same density as zigamma with and :
dzigamma2 gives the density, pzigamma2 gives the distribution function, and rzigamma2 generates random deviates.
x <- rzigamma2(1, 2, 1, 0.5) d <- dzigamma2(x, 2, 1, 0.5) p <- pzigamma2(x, 2, 1, 0.5)x <- rzigamma2(1, 2, 1, 0.5) d <- dzigamma2(x, 2, 1, 0.5) p <- pzigamma2(x, 2, 1, 0.5)
Density, distribution function, and random generation for the zero-inflated inverse Gaussian distribution.
dziinvgauss(x, mean = 1, shape = 1, zeroprob = 0, log = FALSE) pziinvgauss(q, mean = 1, shape = 1, zeroprob = 0, lower.tail = TRUE, log.p = FALSE) rziinvgauss(n, mean = 1, shape = 1, zeroprob = 0)dziinvgauss(x, mean = 1, shape = 1, zeroprob = 0, log = FALSE) pziinvgauss(q, mean = 1, shape = 1, zeroprob = 0, lower.tail = TRUE, log.p = FALSE) rziinvgauss(n, mean = 1, shape = 1, zeroprob = 0)
x, q
|
vector of quantiles |
mean |
location parameter |
shape |
shape parameter, must be positive. |
zeroprob |
zero-probability, must be in |
log, log.p
|
logical; if |
lower.tail |
logical; if |
n |
number of random values to return |
This implementation of zidinvgauss allows for automatic differentiation with RTMB.
where is zeroprob and is the inverse Gaussian density.
dziinvgauss gives the density, pziinvgauss gives the distribution function, and rziinvgauss generates random deviates.
x <- rziinvgauss(1, 1, 2, 0.5) d <- dziinvgauss(x, 1, 2, 0.5) p <- pziinvgauss(x, 1, 2, 0.5)x <- rziinvgauss(1, 1, 2, 0.5) d <- dziinvgauss(x, 1, 2, 0.5) p <- pziinvgauss(x, 1, 2, 0.5)
Density, distribution function, and random generation for the zero-inflated log normal distribution.
dzilnorm(x, meanlog = 0, sdlog = 1, zeroprob = 0, log = FALSE) pzilnorm(q, meanlog = 0, sdlog = 1, zeroprob = 0, lower.tail = TRUE, log.p = FALSE) rzilnorm(n, meanlog = 0, sdlog = 1, zeroprob = 0) plnorm(q, meanlog = 0, sdlog = 1, lower.tail = TRUE, log.p = FALSE)dzilnorm(x, meanlog = 0, sdlog = 1, zeroprob = 0, log = FALSE) pzilnorm(q, meanlog = 0, sdlog = 1, zeroprob = 0, lower.tail = TRUE, log.p = FALSE) rzilnorm(n, meanlog = 0, sdlog = 1, zeroprob = 0) plnorm(q, meanlog = 0, sdlog = 1, lower.tail = TRUE, log.p = FALSE)
x, q
|
vector of quantiles |
meanlog, sdlog
|
mean and standard deviation of the distribution on the log scale with default values of 0 and 1 respectively. |
zeroprob |
zero-inflation probability between 0 and 1. |
log, log.p
|
logical; if |
lower.tail |
logical; if |
n |
number of random values to return |
This implementation allows for automatic differentiation with RTMB.
where is zeroprob, = meanlog, = sdlog, and is the log-normal density.
dzilnorm gives the density, pzilnorm gives the distribution function, and rzilnorm generates random deviates.
x <- rzilnorm(1, 1, 1, 0.5) d <- dzilnorm(x, 1, 1, 0.5) p <- pzilnorm(x, 1, 1, 0.5)x <- rzilnorm(1, 1, 1, 0.5) d <- dzilnorm(x, 1, 1, 0.5) p <- pzilnorm(x, 1, 1, 0.5)
Probability mass function, distribution function, quantile function, and random generation for the zero-inflated negative binomial distribution.
dzinbinom(x, size, prob, zeroprob = 0, log = FALSE) pzinbinom(q, size, prob, zeroprob = 0, lower.tail = TRUE, log.p = FALSE) rzinbinom(n, size, prob, zeroprob = 0)dzinbinom(x, size, prob, zeroprob = 0, log = FALSE) pzinbinom(q, size, prob, zeroprob = 0, lower.tail = TRUE, log.p = FALSE) rzinbinom(n, size, prob, zeroprob = 0)
x, q
|
vector of (non-negative integer) quantiles |
size |
size parameter, must be positive. |
prob |
mean parameter, must be positive. |
zeroprob |
zero-inflation probability between 0 and 1. |
log, log.p
|
logical; if |
lower.tail |
logical; if |
n |
number of random values to return. |
p |
vector of probabilities |
This implementation allows for automatic differentiation with RTMB.
where is zeroprob.
dzinbinom gives the density, pzinbinom gives the distribution function, and rzinbinom generates random deviates.
set.seed(123) x <- rzinbinom(1, size = 2, prob = 0.5, zeroprob = 0.5) d <- dzinbinom(x, size = 2, prob = 0.5, zeroprob = 0.5) p <- pzinbinom(x, size = 2, prob = 0.5, zeroprob = 0.5)set.seed(123) x <- rzinbinom(1, size = 2, prob = 0.5, zeroprob = 0.5) d <- dzinbinom(x, size = 2, prob = 0.5, zeroprob = 0.5) p <- pzinbinom(x, size = 2, prob = 0.5, zeroprob = 0.5)
Probability mass function, distribution function, quantile function and random generation for the zero-inflated negative binomial distribution reparameterised in terms of mean and size.
dzinbinom2(x, mu, size, zeroprob = 0, log = FALSE) pzinbinom2(q, mu, size, zeroprob = 0, lower.tail = TRUE, log.p = FALSE) rzinbinom2(n, mu, size, zeroprob = 0)dzinbinom2(x, mu, size, zeroprob = 0, log = FALSE) pzinbinom2(q, mu, size, zeroprob = 0, lower.tail = TRUE, log.p = FALSE) rzinbinom2(n, mu, size, zeroprob = 0)
x, q
|
vector of (non-negative integer) quantiles |
mu |
mean parameter, must be positive. |
size |
size parameter, must be positive. |
zeroprob |
zero-inflation probability between 0 and 1. |
log, log.p
|
logical; if |
lower.tail |
logical; if |
n |
number of random values to return. |
p |
vector of probabilities |
This implementation allows for automatic differentiation with RTMB.
Uses the same density as zinbinom with :
dzinbinom2 gives the density, pzinbinom2 gives the distribution function, and rzinbinom2 generates random deviates.
set.seed(123) x <- rzinbinom2(1, 2, 1, zeroprob = 0.5) d <- dzinbinom2(x, 2, 1, zeroprob = 0.5) p <- pzinbinom2(x, 2, 1, zeroprob = 0.5)set.seed(123) x <- rzinbinom2(1, 2, 1, zeroprob = 0.5) d <- dzinbinom2(x, 2, 1, zeroprob = 0.5) p <- pzinbinom2(x, 2, 1, zeroprob = 0.5)
Probability mass function, distribution function, and random generation for the zero-inflated Poisson distribution.
dzipois(x, lambda, zeroprob = 0, log = FALSE) pzipois(q, lambda, zeroprob = 0, lower.tail = TRUE, log.p = FALSE) rzipois(n, lambda, zeroprob = 0)dzipois(x, lambda, zeroprob = 0, log = FALSE) pzipois(q, lambda, zeroprob = 0, lower.tail = TRUE, log.p = FALSE) rzipois(n, lambda, zeroprob = 0)
x, q
|
integer vector of counts |
lambda |
vector of (non-negative) means |
zeroprob |
zero-inflation probability between 0 and 1 |
log, log.p
|
logical; return log-density if TRUE |
lower.tail |
logical; if |
n |
number of random values to return. |
This implementation allows for automatic differentiation with RTMB.
where is zeroprob.
dzipois gives the probability mass function, pzipois gives the distribution function, and rzipois generates random deviates.
set.seed(123) x <- rzipois(1, 0.5, 1) d <- dzipois(x, 0.5, 1) p <- pzipois(x, 0.5, 1)set.seed(123) x <- rzipois(1, 0.5, 1) d <- dzipois(x, 0.5, 1) p <- pzipois(x, 0.5, 1)
Density, distribution function, and random generation for the zero-inflated Weibull distribution.
dziweibull(x, shape, scale, zeroprob = 0, log = FALSE) pziweibull(q, shape, scale, zeroprob = 0, lower.tail = TRUE, log.p = FALSE) rziweibull(n, shape, scale, zeroprob = 0)dziweibull(x, shape, scale, zeroprob = 0, log = FALSE) pziweibull(q, shape, scale, zeroprob = 0, lower.tail = TRUE, log.p = FALSE) rziweibull(n, shape, scale, zeroprob = 0)
x, q
|
vector of quantiles |
shape |
positive shape parameter |
scale |
positive scale parameter |
zeroprob |
zero-inflation probability between 0 and 1. |
log, log.p
|
logical; if |
lower.tail |
logical; if |
n |
number of random values to return |
This implementation allows for automatic differentiation with RTMB.
where is zeroprob, = shape, = scale.
dziweibull gives the density, pziweibull gives the distribution function, and rziweibull generates random deviates.
x <- rziweibull(1, 1, 1, 0.5) d <- dziweibull(x, 1, 1, 0.5) p <- pziweibull(x, 1, 1, 0.5)x <- rziweibull(1, 1, 1, 0.5) d <- dziweibull(x, 1, 1, 0.5) p <- pziweibull(x, 1, 1, 0.5)
Density, distribution function, and random generation for the zero-one-inflated beta distribution.
dzoibeta(x, shape1, shape2, zeroprob = 0, oneprob = 0, log = FALSE) pzoibeta(q, shape1, shape2, zeroprob = 0, oneprob = 0, lower.tail = TRUE, log.p = FALSE) rzoibeta(n, shape1, shape2, zeroprob = 0, oneprob = 0)dzoibeta(x, shape1, shape2, zeroprob = 0, oneprob = 0, log = FALSE) pzoibeta(q, shape1, shape2, zeroprob = 0, oneprob = 0, lower.tail = TRUE, log.p = FALSE) rzoibeta(n, shape1, shape2, zeroprob = 0, oneprob = 0)
x, q
|
vector of quantiles |
shape1, shape2
|
non-negative shape parameters of the beta distribution |
zeroprob |
zero-inflation probability between 0 and 1. |
oneprob |
one-inflation probability between 0 and 1. |
log, log.p
|
logical; if |
lower.tail |
logical; if |
n |
number of random values to return. |
This implementation allows for automatic differentiation with RTMB.
where = zeroprob and = oneprob.
dzoibeta gives the density, pzoibeta gives the distribution function, and rzoibeta generates random deviates.
set.seed(123) x <- rzoibeta(1, 2, 2, 0.2, 0.3) d <- dzoibeta(x, 2, 2, 0.2, 0.3) p <- pzoibeta(x, 2, 2, 0.2, 0.3)set.seed(123) x <- rzoibeta(1, 2, 2, 0.2, 0.3) d <- dzoibeta(x, 2, 2, 0.2, 0.3) p <- pzoibeta(x, 2, 2, 0.2, 0.3)
Density, distribution function, and random generation for the zero-one-inflated beta distribution reparameterised in terms of mean and concentration.
dzoibeta2(x, mu, phi, zeroprob = 0, oneprob = 0, log = FALSE) pzoibeta2(q, mu, phi, zeroprob = 0, oneprob = 0, lower.tail = TRUE, log.p = FALSE) rzoibeta2(n, mu, phi, zeroprob = 0, oneprob = 0)dzoibeta2(x, mu, phi, zeroprob = 0, oneprob = 0, log = FALSE) pzoibeta2(q, mu, phi, zeroprob = 0, oneprob = 0, lower.tail = TRUE, log.p = FALSE) rzoibeta2(n, mu, phi, zeroprob = 0, oneprob = 0)
x, q
|
vector of quantiles |
mu |
mean parameter, must be in the interval from 0 to 1. |
phi |
concentration parameter, must be positive. |
zeroprob |
zero-inflation probability between 0 and 1. |
oneprob |
one-inflation probability between 0 and 1. |
log, log.p
|
logical; if |
lower.tail |
logical; if |
n |
number of random values to return. |
This implementation allows for automatic differentiation with RTMB.
Uses the same density as zoibeta with and .
dzoibeta2 gives the density, pzoibeta2 gives the distribution function, and rzoibeta2 generates random deviates.
set.seed(123) x <- rzoibeta2(1, 0.6, 2, 0.2, 0.3) d <- dzoibeta2(x, 0.6, 2, 0.2, 0.3) p <- pzoibeta2(x, 0.6, 2, 0.2, 0.3)set.seed(123) x <- rzoibeta2(1, 0.6, 2, 0.2, 0.3) d <- dzoibeta2(x, 0.6, 2, 0.2, 0.3) p <- pzoibeta2(x, 0.6, 2, 0.2, 0.3)
Probability mass function, distribution function, and random generation for the zero-truncated Binomial distribution.
dztbinom(x, size, prob, log = FALSE) pztbinom(q, size, prob, lower.tail = TRUE, log.p = FALSE) rztbinom(n, size, prob)dztbinom(x, size, prob, log = FALSE) pztbinom(q, size, prob, lower.tail = TRUE, log.p = FALSE) rztbinom(n, size, prob)
x, q
|
integer vector of counts |
size |
number of trials |
prob |
success probability in each trial |
log, log.p
|
logical; return log-density if TRUE |
lower.tail |
logical; if |
n |
number of random values to return. |
This implementation allows for automatic differentiation with RTMB.
By definition, this distribution only has support on the positive integers (1, ..., size). Any zero-truncated distribution is defined as
where is the probability mass function of the corresponding untruncated distribution.
dztbinom gives the probability mass function, pztbinom gives the distribution function, and rztbinom generates random deviates.
set.seed(123) x <- rztbinom(1, size = 10, prob = 0.3) d <- dztbinom(x, size = 10, prob = 0.3) p <- pztbinom(x, size = 10, prob = 0.3)set.seed(123) x <- rztbinom(1, size = 10, prob = 0.3) d <- dztbinom(x, size = 10, prob = 0.3) p <- pztbinom(x, size = 10, prob = 0.3)
Probability mass function, distribution function, and random generation for the zero-truncated Negative Binomial distribution.
dztnbinom(x, size, prob, log = FALSE) pztnbinom(q, size, prob, lower.tail = TRUE, log.p = FALSE) rztnbinom(n, size, prob)dztnbinom(x, size, prob, log = FALSE) pztnbinom(q, size, prob, lower.tail = TRUE, log.p = FALSE) rztnbinom(n, size, prob)
x, q
|
integer vector of counts |
size |
target for number of successful trials, or dispersion parameter (the shape parameter of the gamma mixing distribution). Must be strictly positive, need not be integer. |
prob |
probability of success in each trial. 0 < prob <= 1. |
log, log.p
|
logical; return log-density if TRUE |
lower.tail |
logical; if |
n |
number of random values to return. |
This implementation allows for automatic differentiation with RTMB.
By definition, this distribution only has support on the positive integers (1, 2, ...). Any zero-truncated distribution is defined as
where is the probability mass function of the corresponding untruncated distribution.
dztnbinom gives the probability mass function, pztnbinom gives the distribution function, and rztnbinom generates random deviates.
set.seed(123) x <- rztnbinom(1, size = 2, prob = 0.5) d <- dztnbinom(x, size = 2, prob = 0.5) p <- pztnbinom(x, size = 2, prob = 0.5)set.seed(123) x <- rztnbinom(1, size = 2, prob = 0.5) d <- dztnbinom(x, size = 2, prob = 0.5) p <- pztnbinom(x, size = 2, prob = 0.5)
Probability mass function, distribution function, quantile function, and random generation for the zero-truncated negative binomial distribution reparameterised in terms of mean and size.
dztnbinom2(x, mu, size, log = FALSE) pztnbinom2(q, mu, size, lower.tail = TRUE, log.p = FALSE) rztnbinom2(n, mu, size)dztnbinom2(x, mu, size, log = FALSE) pztnbinom2(q, mu, size, lower.tail = TRUE, log.p = FALSE) rztnbinom2(n, mu, size)
x, q
|
integer vector of counts |
mu |
mean parameter, must be positive |
size |
size/dispersion parameter, must be positive |
log, log.p
|
logical; return log-density if TRUE |
lower.tail |
logical; if |
n |
number of random values to return. |
This implementation allows for automatic differentiation with RTMB.
By definition, this distribution only has support on the positive integers (1, 2, ...). Any zero-truncated distribution is defined as
where is the probability mass function of the corresponding untruncated distribution.
dztnbinom2 gives the probability mass function, pztnbinom2 gives the distribution function, and rztnbinom2 generates random deviates.
set.seed(123) x <- rztnbinom2(1, mu = 2, size = 1) d <- dztnbinom2(x, mu = 2, size = 1) p <- pztnbinom2(x, mu = 2, size = 1)set.seed(123) x <- rztnbinom2(1, mu = 2, size = 1) d <- dztnbinom2(x, mu = 2, size = 1) p <- pztnbinom2(x, mu = 2, size = 1)
Probability mass function, distribution function, and random generation for the zero-truncated Poisson distribution.
dztpois(x, lambda, log = FALSE) pztpois(q, lambda, lower.tail = TRUE, log.p = FALSE) rztpois(n, lambda)dztpois(x, lambda, log = FALSE) pztpois(q, lambda, lower.tail = TRUE, log.p = FALSE) rztpois(n, lambda)
x, q
|
integer vector of counts |
lambda |
vector of (non-negative) means |
log, log.p
|
logical; return log-density if TRUE |
lower.tail |
logical; if |
n |
number of random values to return. |
This implementation allows for automatic differentiation with RTMB.
By definition, this distribution only has support on the positive integers (1, 2, ...). Any zero-truncated distribution is defined as
where is the probability mass function of the corresponding untruncated distribution.
dztpois gives the probability mass function, pztpois gives the distribution function, and rztpois generates random deviates.
set.seed(123) x <- rztpois(1, 0.5) d <- dztpois(x, 0.5) p <- pztpois(x, 0.5)set.seed(123) x <- rztpois(1, 0.5) d <- dztpois(x, 0.5) p <- pztpois(x, 0.5)