| getEventData {BAMMtools} | R Documentation | 
Create bammdata object from MCMC output
Description
getEventData Reads shift configuration data (the
"event data" output) from a BAMM analysis and creates a
bammdata object. The bammdata object is fundamental
for extracting information about macroevolutionary rate variation
through time and among lineages.
Usage
getEventData(
  phy,
  eventdata,
  burnin = 0,
  nsamples = NULL,
  verbose = FALSE,
  type = "diversification"
)
Arguments
phy | 
 An object of class   | 
eventdata | 
 A character string specifying the path to a   | 
burnin | 
 A numeric indicating the fraction of posterior samples to discard as burn-in.  | 
nsamples | 
 An integer indicating the number of posterior samples to
include in the   | 
verbose | 
 A logical. If   | 
type | 
 A character string. Either "diversification" or "trait"
depending on your   | 
Details
In the BAMM framework, an "event" defines a
macroevolutionary process of diversification or trait evolution. Every
sample from the posterior includes at least one process, defined by
such an "event". If a given sample includes just a single event, then
the dynamics of diversification or trait evolution can be described
entirely by a single time-constant or time-varying process that begins
at the root of the tree. Any sample from the posterior distribution
may include a complex mixture of distinct processes. To represent
temporal heterogeneity in macroevolutionary rates, BAMM models
a rate R, e.g. speciation, as a function that changes
exponentially with time:
R(t) = R(0)*exp(b*t).
Here R(0) is the initial rate and b is a parameter
determining how quickly that rate grows or decays with time. 
The eventdata file (or data frame) is a record of events and
associated parameters that were sampled with BAMM during
simulation of the posterior with reversible jump MCMC. This complex,
information-rich file is processed into a bammdata object,
which serves as the core data object for numerous downstream analyses.
From a bammdata object, you can summarize rate variation
through time, among clades, extract locations of rate shifts,
summarize clade-specific rates of speciation and extinction, and more.
In general, the user does not need to be concerned with the details of
a bammdata object. The object is used as input by a number of
BAMMtools functions. 
The parameter nsamples can be used to reduce the total amount
of data included in the raw eventdata output from a BAMM run.
The final bammdata object will consist of all data for
nsamples from the posterior. These nsamples are equally
spaced after discarding some burnin fraction as "burn-in". If
nsamples is set to NULL, the bammdata object will
include all samples in the posterior after discarding the
burnin fraction.
Value
A list with many components:
edge: See documentation for class
phyloin package ape.Nnode: See documentation for class
phyloin package ape.tip.label: See documentation for class
phyloin package ape.edge.length: See documentation for class
phyloin package ape.begin: The beginning time of each branch in absolute time (the root is set to time zero)
end: The ending time of each branch in absolute time.
numberEvents: An integer vector with the number of events contained in
phyfor each posterior sample. The length of this vector is equal to the number of posterior samples in thebammdataobject.eventData: A list of dataframes. Each element is a single posterior sample. Each row in a dataframe holds the data for a single event. Data associated with an event are:
node- a node number. This identifies the branch where the event originates.time- this is the absolute time on that branch where the event originates (with the root at time 0).lam1- an initial rate of speciation or trait evolution.lam2- a decay/growth parameter.mu1- an initial rate of extinction.mu2- a decay/growth parameter.index- a unique integer associated with the event. See 'Details'.eventVectors: A list of integer vectors. Each element is a single posterior sample. For each branch in
phythe index of the event that occurs along that branch. Branches are ordered increasing here and elsewhere.eventBranchSegs: A list of matrices. Each element is a single posterior sample. Each matrix has four columns:
Column 1identifies a node inphy.Column 2identifies the beginning time of the branch or segment of the branch that subtends the node inColumn 1.Column 3identifies the ending time of the branch or segment of the branch that subtends the node inColumn 1.Column 4identifies the index of the event that occurs along the branch or segment of the branch that subtends the node inColumn 1.tipStates: A list of integer vectors. Each element is a single posterior sample. For each tip the index of the event that occurs along the branch subtending the tip. Tips are ordered increasing here and elsewhere.
tipLambda: A list of numeric vectors. Each element is a single posterior sample. For each tip the rate of speciation or trait evolution at the end of the terminal branch subtending that tip.
tipMu: A list of numeric vectors. Each element is a single posterior sample. For each tip the rate of extinction at the end of the terminal branch subtending that tip. Meaningless if working with
BAMMtrait results.meanTipLambda: For each tip the mean of the marginal posterior density of the rate of speciation or trait evolution at the end of the terminal branch subtending that tip.
meanTipMu: For each tip the mean of the marginal posterior density of the rate of extinction at the end of the terminal branch subtending that tip. Meaningless if working with
BAMMtrait results.type: A character string. Either "diversification" or "trait" depending on your
BAMManalysis.downseq: An integer vector holding the nodes of
phy. The order corresponds to the order in which nodes are visited by a pre-order tree traversal.lastvisit: An integer vector giving the index of the last node visited by the node in the corresponding position in
downseq.downseqandlastvisitcan be used to quickly retrieve the descendants of any node. e.g. the descendants of node 89 can be found bydownseq[which(downseq==89):which(downseq==lastvisit[89]).
Note
Currently the function does not check for duplicate tip labels in
phy, which may cause the function to choke.
Author(s)
Dan Rabosky, Mike Grundler
References
See Also
summary.bammdata, plot.bammdata,
dtRates.
Examples
data(primates, events.primates)
xx <- getEventData(primates, events.primates, burnin=0.25, nsamples=500,
                   type = 'trait')
# compute mean phenotypic rate for primate body size evolution:
brates <- getCladeRates(xx)
mean(brates$beta)
# Plot rates:
plot(xx)