msc2time | R Documentation |
Calibrate a BPP A00 MCMC sample from a multi-species coalescent phylogeny to absolute time using a fossil calibration or a prior on the molecular rate.
msc2time.t(mcmc, node.name, calf, ...)
msc2time.r(mcmc, u.mean, u.sd, g.mean, g.sd)
mcmc |
A data frame containing the MCMC output of a BPP A00 analysis |
node.name |
A character vector of length one with the name of the node to which the calibration will be applied |
calf |
A calibration function to generate random numbers |
... |
Further parameters passed to |
u.mean |
Numeric vector of length one with the mean for the per-generation molecular rate calibration |
u.sd |
Numeric vector of length one with the SD for for the per-generation molecular rate calibration |
g.mean |
Numeric vector of length one with the mean for the generation time calibration |
g.sd |
Numeric vector of length one with the SD for the generation time calibration |
msc2time.t
will calibarte a BPP A00 analysis to geological
time using a fossil calibration and msc2time.r
will do the same but
using a prior on the rate.
msc2time.t
will obtain a sample of times from the random
distribution in calf
. Suitable choices for calf
are
runif
, rgamma
, and rslnorm
. The sampled times are then
used to calculate the molecular rate, and then re-scale the relative times
(tau's) for the other nodes in mcmc
to geological time.
u.mean
and u.sd
, and g.mean
and g.sd
, are used
to construct gamma density calibrations for the per-generation molecular
rate and generation time respectively. The gamma density with mean m
and s.d. s
has shape a = (m / s)^2
and rate b = m / s^2
.
In msc2time.r
, the gamma densities are used to obtain random samples
of the per-generation rate and generation time. From these the molecular
rate per absolute time unit is calculated, and then used to convert the
relative times (tau's) to absolute divergence times. The relative
population sizes (theta's) are converted to effective population sizes in
number of individuals.
Angelis and dos Reis (2015) give the random sampling procedure used in these functions. Yoder et al. (2016) gives an example of calibrating a mouse lemur phylogeny using a prior on the rate. The BPP A00 analysis is described in Yang (2015).
Note that setting g.mean = 1
and g.sd = 0
allows the use of a
rate calibration in geological time: In this case set u.mean
to be
the desire rate in calendar time units, such as substitutions per site per
year or per millions of years, with u.sd
set to give the desired
standard deviation in the calendar time unit.
A data frame with a posterior sample of the calibrated times and
molecular rate, and additionally, in the case of msc2time.r
, the
population sizes.
Mario dos Reis
K. Angelis and M. dos Reis (2015) The impact of ancestral population size and incomplete lineage sorting on Bayesian estimation of species divergence times. Curr. Zool., 61: 874–885.
C. R. Campbell, G. P. Tiley, J. W. Poelstra, K. E. Hunnicutt, P. A. Larsen, H.-J. Lee, J. L. Thorne, M. dos Reis, and A. D. Yoder. (2021) Pedigree-based and phylogenetic methods support surprising patterns of mutation rate and spectrum in the gray mouse lemur Heredity, 127: 233–244.
Z. Yang (2015) The BPP program for species tree estimation and species delimitation. Curr. Zool., 61: 854–865.
A. D. Yoder, C. R. Campbell, M. B. Blanco, M. dos Reis, J. U. Ganzhorn, S. M. Goodman, K. E. Hunnicutt, P. A. Larsen, P. M. Kappeler, R. M. Rasoloarison, J. M. Ralison, D. L. Swofford, and D. W. Weisrock. (2016) Geogenetic patterns in mouse lemurs (genus Microcebus) reveal the ghosts of Madagascar's forests past. Proc. Nat. Acad. Sci. USA., 113: 8049–8056.
data(hominids)
# Calibrate the hominid phylogeny with a uniform fossil calibration of
# between 6.5 to 10 Ma for the human-chimp divergence.
calmsc <- msc2time.t(mcmc=hominids$mcmc, node="7humanchimp", calf=runif,
min=6.5, max=10)
# posterior age of human-chimp (this is the same as the calibration)
plot(density(calmsc$t_7humanchimp, adj=.1), xlab="Time (Ma)",
main="Human-chimp age")
rug(calmsc$t_7humanchimp)
## Not run:
# calculate posterior summary (requires CODA package)
mcmc.summary(calmsc)
## End(Not run)
# Calibrate the hominid phylogeny using a shifted-lognormal fossil calibration
# with minimum at 6.5 Ma for the human-chimp divergence.
calmsc <- msc2time.t(mcmc=hominids$mcmc, node="7humanchimp", calf=rslnorm,
shift=6.5, sdlog=.5)
plot(density(calmsc$t_7humanchimp, adj=.1), xlab="Time (Ma)",
main="Human-chimp age")
rug(calmsc$t_7humanchimp)
data(microcebus)
# Calibrate the Microcebus phylogeny to absoluate divergence times using a
# prior on the per-generation rate and generation time
calmsc <- msc2time.r(mcmc=microcebus$mcmc, u.mean=8.7e-9, u.sd=1.65e-9,
g.mean=3.75, g.sd=0.375)
# posterior age of the phylogeny's root (in thousans of years)
plot(density(calmsc$t_7OLMXRB / 1e3, adj=.1), xlab="Time (Ka)",
main="Root age (Microcebus)")
rug(calmsc$t_7OLMXRB / 1e3)
# Posterior of the ancestral effective population at the root (in
# thousands of individuals)
plot(density(calmsc$Ne_7OLMXRB / 1e3, adj=.1),
xlab="Ne (x 10^3 individuals)", main = "Ne at root (Microcebus)")
rug(calmsc$Ne_7OLMXRB / 1e3)
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.