This is the high-level entry point to the POUMM method. The POUMM function fits the POUMM method to a tree and observed trait-values at its tips and returns an object of class "POUMM".

POUMM(
  z,
  tree,
  se = 0,
  zName = "z",
  treeName = "tree",
  parDigits = 6,
  usempfr = 0,
  useCpp = TRUE,
  ...,
  spec = NULL,
  doMCMC = TRUE,
  likPOUMM_lowLevelFun = likPOUMMGivenTreeVTipsC,
  verbose = FALSE,
  debug = FALSE
)

Arguments

z

Either a numeric vector containing the phenotypic values at the tips of tree or a named list containing named elements z - a numeric vector and tree - a phylo object (it is possible to specify different element names using the arguments zName and treeName).

tree

A phylo object or NULL in case z is a list.

se

A non-negative numerical vector (or single number) indicating known measurement standard error (defaults to 0). Note the elements of this vector are assumed to describe the measurement error at individual nodes independent of the environmental contribution (described by the parameter sigmae). The total error standard deviation is thus sqrt(sigmae2+se^2).

zName, treeName

Character strings used when the parameter z is a list; indicate the names in the list of the values-vector and the tree. Default: 'z' and 'tree'.

parDigits

Integer specifying rounding to be done on the parameter vector before likelihood calculation. Defaults to 6 decimal digits. This can be useful during maximum likelihood optimization to prevent likelihood calculation on very small but positive values of alpha, but should be used with caution since specifying a small number of digits, i.e. 2 or 3 can result in an infinite loop during optim. Specify a negative number to disable rounding.

usempfr

integer indicating if and how mpfr should be used for small parameter values (`any(c(alpha, sigma, sigmae) < 0.01)`). Using the mpfr package can be forced by specifying an integer greater or equal to 2. Setting usempfr=0 (default) causes high precision likelihood calculation to be done on each encounter of parameters with at least 1 bigger log-likelihood value than any of the currently found maximum log-likelihood or the previously calculated log-likelihood value Requires the Rmpfr package. Note that using mpfr may increase the time for one likelihood calculation more than 100-fold. Set usempfr to -1 or less to completely disable Rmpfr functionality.

useCpp

Logical indicating whether C++ likelihood calculation should be used for faster vector operations. Defaults to TRUE. Since the C++ likelihood implementation does not support mpfr, useCpp gets disabled when usempfr is bigger than 0.

...

additional arguments passed to the `likPOUMMGivenTreeVTips()` function (`?dVGivenTreeOU` for details).

spec

A named list specifying how the ML and MCMC fit should be done. See `?specifyPOUMM`.

doMCMC

Deprecated - replaced by specifying nSamplesMCMC as a member of spec instead (see `?specifyPOUMM`). logical: should a MCMC fit be performed. An MCMC fit provides a sample from the posterior distribution of the parameters given a prior distribution and the data. Unlike the ML-fit, it allows to estimate confidence intervals for the estimated parameters. This argument is TRUE by default. The current implementation uses the adaptive Metropolis sampler from the package `adaptMCMC` written by Andreas Scheidegger. To obtain meaningful estimates MCMC may need to run for several millions of iterations (parameter nSamplesMCMC set to 1e5 by default). See parameters ending at MCMC in `?specifyPOUMM` for details.

likPOUMM_lowLevelFun

the low-level function used for POUMM - likelihood calculation. Default value is likPOUMMGivenTreeVTipsC.

verbose, debug

Logical flags indicating whether to print informative and/or debug information on the standard output (both are set to to FALSE by default).

Value

An object of S3 class 'POUMM'. This object can be analyzed using S3 generic functions: summary, plot, AIC, BIC, coef, logLik, fitted.

References

Mitov, V., and Stadler, T. (2017). Fast and Robust Inference of Phylogenetic Ornstein-Uhlenbeck Models Using Parallel Likelihood Calculation. bioRxiv, 115089. https://doi.org/10.1101/115089 Mitov, V., & Stadler, T. (2017). Fast Bayesian Inference of Phylogenetic Models Using Parallel Likelihood Calculation and Adaptive Metropolis Sampling. Systematic Biology, 235739. http://doi.org/10.1101/235739 Vihola, M. (2012). Robust adaptive Metropolis algorithm with coerced acceptance rate. Statistics and Computing, 22(5), 997-1008. http://doi.org/10.1007/s11222-011-9269-5

Scheidegger, A. (2012). adaptMCMC: Implementation of a generic adaptive Monte Carlo Markov Chain sampler. http://CRAN.R-project.org/package=adaptMCMC

See also

specifyPOUMM for parametrizations and custom settings of the POUMM fit.

Examples

if (FALSE) { # Please, read the package vignette for more detailed examples. N <- 500 tr <- ape::rtree(N) z <- rVNodesGivenTreePOUMM(tr, 0, 2, 3, 1, 1)[1:N] fit <- POUMM(z, tr, spec = specifyPOUMM(nSamplesMCMC = 5e4)) plot(fit) summary(fit) AIC(fit) BIC(fit) coef(fit) logLik(fit) fitted(fit) plot(resid(fit)) abline(h=0) # fit PMM to the same data and do a likelihood ratio test fitPMM <- POUMM(z, tr, spec = specifyPMM(nSamplesMCMC = 5e4)) lmtest::lrtest(fitPMM, fit) }