getRateThroughTimeMatrix {BAMMtools}R Documentation

Generate rate-through-time matrix from bammdata object

Description

Computes a matrix of macroevolutionary rates at specified timepoints from a bammdata object. These rates can be used for plotting speciation rates (and other rates) through time.

Usage

getRateThroughTimeMatrix(
  ephy,
  start.time = NULL,
  end.time = NULL,
  nslices = 100,
  node = NULL,
  nodetype = "include"
)

Arguments

ephy

An object of class bammdata.

start.time

The start time (in units before present) for the time. sequence over which rates should be computed. If NULL, starts at the root.

end.time

The end time (in units before present) for the time sequence over which rates should be computed. If NULL, ends in the present.

nslices

The number of time points at which to compute rate estimates (between start.time and end.time).

node

Allows user to extract rate-through-time information for the subtree descended from a specific node. Alternatively, a specified subtree can be excluded from the rate matrix calculations.

nodetype

Two options: "include" and "exclude". If "include", computes rate matrix only for the descendants of subtree defined by node. If "exclude", computes rate matrix for all background lineages in tree after excluding the subtree defined by node. Ignored if node = NULL.

Details

Computes evolutionary rates for each sample in the posterior included as part of the bammdata object. Rates are computed by draping an imaginary grid over the phylogeny, where the grid begins at start.time and ends at end.time, with nslices vertical lines through the phylogeny. The mean rate at each point in time (for a given sample from the posterior) is simply the mean rate at that time for all branches that are intersected by the grid (see the grid plot in the examples section).

This function is used by plotRateThroughTime, but the user can work directly with the bamm-ratematrix object for greater control in plotting rate-through-time trajectories for individual clades. See examples for an example of how this can be used to plot confidence intervals on a rate trajectory using shaded polygons.

The node options are particularly useful. If you have run BAMM on a large phylogeny, you can easily generate the rate-through-time data for a particular subtree by specifying the node number along with nodetype = "include". Likewise, if you want to look at just the background rate - excluding some particular lineage - just specify nodetype = "exclude".

Value

An object of class bamm-ratematrix with the following components:

lambda

A nsamples x nslices matrix of speciation rates, where nsamples is the number of posterior samples in the bammdata object.

mu

A nsamples x nslices matrix of extinction rates.

beta

A nsamples x nslices matrix of phenotypic rates (if applicable).

times

A vector of timepoints where rates were computed.

times

A vector of timepoints where rates were computed (see Examples).

type

Either "diversification" or "trait", depending on the input data.

Author(s)

Dan Rabosky

See Also

plotRateThroughTime

Examples

## Not run: 
# Plot a rate-through-time curve with 
# confidence intervals for the whale dataset:

data(whales, events.whales)
ed <- getEventData(whales, events.whales)

rmat <- getRateThroughTimeMatrix(ed)

plot.new()
plot.window(xlim=c(0, 36), ylim=c(0, .7))

## Speciation quantiles: plot 90% CIs
qq <- apply(rmat$lambda, 2, quantile, c(0.05, 0.5, 0.95))

xv <- c(rmat$times, rev(rmat$times))
yv <- c(qq[1,], rev(qq[3,]))

## Add the confidence polygon on rate distributions:
polygon(xv, yv, col='gray80', border=FALSE)

## Add the median rate line:
lines(rmat$times, qq[2,], lwd=3, col='red')

## Add axes
axis(1, at=seq(-5, 35, by=5))
axis(2, at=seq(-0.2, 1, by=0.2), las=1)

####### Now we will show the actual grid used for rate calculations:

plot(whales, show.tip.label=FALSE)
axisPhylo()

mbt <- max(branching.times(whales))
tvec <- mbt - rmat$times;
tvec <- rmat$times;

for (i in 1:length(tvec)){
    lines(x=c(tvec[i], tvec[i]), y=c(0, 90), lwd=0.7, col='gray70')
}

## This shows the grid of time slices over the phylogeny
## End(Not run)

[Package BAMMtools version 2.1.11 Index]