speciesByRatesMatrix {BAMMtools}R Documentation

Compute species-specific rate through time trajectories

Description

Computes the mean of the marginal posterior density of speciation/extinction or phenotypic rates for equally spaced points along the root to tip path for each species.

Usage

speciesByRatesMatrix(ephy, nslices, index = NULL, spex = "s")

Arguments

ephy

An object of class bammdata.

nslices

An integer number of time slices. This determines the number of equally spaced points in time at which rates are computed for each species.

index

An integer or vector of mode integer indicating which posterior samples to use in the calculation. If NULL (default) all samples are used.

spex

A character string. "s" (default) calculates speciation rates; "e" calculates extinction rates; "netdiv" calculates diversification rates. Ignored if ephy$type = "trait".

Value

A list with two components:

Author(s)

Mike Grundler

References

http://bamm-project.org/

See Also

getRateThroughTimeMatrix

Examples

data(whales, events.whales)
ed <- getEventData(whales, events.whales, burnin=0.25, nsamples=500)
ratemat <- speciesByRatesMatrix(ed, nslices = 100)

dolphins <- extract.clade(whales, 140)$tip.label
plot.new()
plot.window(xlim = c(0, 35), ylim = c(0, 0.8))
for (i in 1:nrow(ratemat$rates)) {
	if (whales$tip.label[i] %in% dolphins) {
		lines(ratemat$times, ratemat$rates[i,], lwd = 2, col = 4)	
	} else {
		lines(ratemat$times, ratemat$rates[i,], lwd = 2, col = 8)
	}
}
axis(1, seq(-5, 35, 5))
axis(2, seq(-0.2, 0.8, 0.2), las = 1)
mtext("Time since root", 1, line = 2.5)
mtext("Speciation rate", 2, line = 2.5)


[Package BAMMtools version 2.1.11 Index]