phy.build {pez} | R Documentation |
Build a novel phylogeny from existing data
Description
Build a novel phylogeny from existing data
congeneric.impute
sequentially add species to a phylogeny to
form an _imputed_ bifurcating tree. Makes use of a result from
Steel & Mooers (2010) that gives the expected branch-length under a
Yule model whether the rate of diversification has been
estimated. The intention of this is to approximate the method by
which phylogenetic structure is sampled from the prior in BEAST;
i.e., to approximate the standard Kuhn et al. (2011) method for
imputing a phylogeny. When using congeneric.impute
you
should (1) repeat your analyses across many (if in doubt,
thousands) of separate runs (see Kuhn et al. 2011) and (2) check
for yourself that your trees are unbiased for your purpose - I make
no guarantee this is appropriate, and in many cases I think it
would not be. See also 'notes' below.
Usage
congeneric.merge(tree, species, split = "_", ...)
bind.replace(backbone, donor, replacing.tip.label, donor.length = NA)
congeneric.impute(tree, species, split = "_", max.iter = 1000, ...)
Arguments
tree |
|
species |
vector of species names to be bound into the tree if missing from it |
split |
the character that splits genus and species names in
your phylogeny. Default is |
... |
ignored |
backbone |
backbone phylogeny ( |
donor |
phylogeny ( |
replacing.tip.label |
the species in the donor phylogeny that's being replaced by the donor phylogeny |
donor.length |
how deep the donor phylogeny should be cut into the backbone phylogeny. If NA (default), then the bladj algorithm is followed (or, in plain English, it's put half-way along the branch) |
max.iter |
Sometimes the random draw for the new branch length
to be added will be too large to allow it to be added to the
tree. In such cases, |
Details
congeneric.merge
Binds missing species into a
phylogeny by replacing all members of the clade it belongs to with
a polytomy. Assumes the tip.labels
represent Latin
binomials, split by the split
argument. This code was
originally shipped with phyloGenerator - this is the merge
method in that program.
bind.replace
Binds a phylogeny (donor) into a
bigger phylogeny ('backbone'); useful if you're building a
phylogeny a la Phylomatic. A version of this R code was shipped
with phyloGenerator (Pearse & Purvis 2013). This is really an
internal function for congeneric.merge
, but hopefully it's
of some use to you!
Value
phylo
phylogeny
phylogeny (phylo
)
Note
Thank you to Josep Padulles Cubino, who found that the genus
name splitting in a previous version of this function could
cause incorrect placement in oddly named cases. As with all
phylogenetic construction tools, the output from
congeneric.merge
should be checked for accuracy before
use.
Caveats for congeneric.impute
: something I noticed is
that BEAST randomly picks an edge to break when adding species
(starting from a null tree), and this is the behaviour I have
(attempted to) replicate here. It is not clear to me that this
is unbiased, since a clade undergoing rapid diversification
will have many edges but these will be short (and so cannot
have an edge inserted into them following the method below). My
understanding is this is a known problem, and I simply cannot
think of a better way of doing this that doesn't incorporate
what I consider to be worse pathology. Thus this method, even
if it works (which I can't guarantee), it should tend to break
long branches.
Author(s)
Will Pearse
Will Pearse
References
Pearse W.D. & Purvis A. phyloGenerator: an automated phylogeny generation tool for ecologists. Methods in Ecology and Evolution 4(7): 692–698.
Steel, M., & Mooers, A. (2010). The expected length of pendant and interior edges of a Yule tree. Applied Mathematics Letters, 23(11), 1315-1319.
Kuhn, T. S., Mooers, A. O., & Thomas, G. H. (2011). A simple polytomy resolver for dated phylogenies. Methods in Ecology and Evolution, 2(5), 427-436.
Examples
tree <- read.tree(text="((a_a:1,b_b:1):1, c_c:2):1;")
tree <- congeneric.merge(tree, c("a_nother", "a_gain", "b_sharp"))
tree <- read.tree(text="((a_a:1,b_b:1):1, c_c:2):1;")
tree <- congeneric.impute(tree, c("a_nother", "a_gain", "b_sharp"))