View source: R/resolveTreeChar.R
resolveTreeChar | R Documentation |
This function resolves a set of given topology with less than fully-binary phylogenetic resolution so that
lineages are shifted and internal nodes added that minimize the number of independent character transitions needed to explain
an observed distribution of discrete character states for the taxa on such a tree, under various maximum-parsimony algorithms of
ancestral character reconstruction, powered ultimately by function ancestral.pars
in library phangorn
.
This function is mainly designed for use with poorly resolved trees which are being assessed with the function
minCharChange
.
resolveTreeChar( tree, trait, orderedChar = FALSE, stateBias = NULL, iterative = TRUE, cost = NULL, ambiguity = c(NA, "?"), dropAmbiguity = FALSE, polySymbol = "&", contrast = NULL )
tree |
A cladogram of type |
trait |
A vector of trait values for a discrete character, preferably named with taxon names identical to the tip labels on the input tree. |
orderedChar |
Is the character of interest given for |
stateBias |
This argument controls how |
iterative |
A logical argument which, if TRUE (the default), causes the function to repeat the polytomy-resolving functionality across the entire tree until the number of nodes stabilizes. If FALSE, polytomies are only passed a single time. |
cost |
A matrix of the cost (i.e. number of steps) necessary to
change between states of the input character trait.
If |
ambiguity |
A vector of values which indicate ambiguous
(i.e. missing or unknown) character state codings
in supplied |
dropAmbiguity |
A logical. If |
polySymbol |
A single symbol which separates alternative
states for polymorphic codings; the default symbol is
|
contrast |
A matrix of type integer with cells of 0
and 1, where each row is labeled with a string value
used for indicating character states in |
As shown in the example code below, this function offers a wide variety of options for manipulating the
maximum-parsimony algorithm used (i.e. MPR versus ACCTRAN), the ordering (or not) of character states,
and potential biasing of uncertainty character state reconstructions (when ordered characters are
assessed). This allows for a wide variety of possible resolutions for a given tree with polytomies
and a discrete character. In general, the author expects that use of this function will be optimal
when applied to ordered characters using one of the stateBias
options, perhaps
stateBias = "primitive"
(based on theoretical expectations for slow evolving characters). However,
anecdotal use of this function with various simulation datasets suggests that the results are quite
variable, and so the best option needs to be assessed based on the prior assumptions regarding the
data and the performance of the dataset with the various arguments of this function.
Returns the resulting tree, which may be fully resolved, partly more resolved or not more resolved at all
(i.e. have less polytomies) depending on what was possible, as constrained by ambiguities in character
reconstructions. Applying multi2di
is suggested as a post-step to obtain a fully-resolved
cladogram, if one is desired.
David W. Bapst
Hanazawa, M., H. Narushima, and N. Minaka. 1995. Generating most parsimonious reconstructions on a tree: A generalization of the Farris-Swofford-Maddison method. Discrete Applied Mathematics 56(2-3):245-265.
Narushima, H., and M. Hanazawa. 1997. A more efficient algorithm for MPR problems in phylogeny. Discrete Applied Mathematics 80(2-3):231-238.
Schliep, K. P. 2011. phangorn: phylogenetic analysis in R. Bioinformatics 27(4):592-593.
Swofford, D. L., and W. P. Maddison. 1987. Reconstructing ancestral character states under Wagner parsimony. Mathematical Biosciences 87(2):199-229.
ancPropStateMat
which is used internally by this function. This function was
intentionally designed for use with minCharChange
.
# let's write a quick&dirty ancestral trait plotting function quickAncPlot <- function(tree, trait, cex, orderedChar = FALSE, type = "MPR", cost = NULL){ ancData <- ancPropStateMat(tree = tree, trait = trait, orderedChar = orderedChar) ancCol <- (1:ncol(ancData))+1 plot(tree,show.tip.label = FALSE,no.margin = TRUE,direction = "upwards") tiplabels(pch = 16,pie = ancData[(1:Ntip(tree)),],cex = cex,piecol = ancCol, col = 0) nodelabels(pie = ancData[-(1:Ntip(tree)),],cex = cex,piecol = ancCol) } ########## # examples with simulated data set.seed(2) tree <- rtree(50) #simulate under a likelihood model trait <- rTraitDisc(tree,k = 3,rate = 0.7) tree <- degradeTree(tree,prop_collapse = 0.6) tree <- ladderize(tree,right = FALSE) #a bunch of type = MPR (default) examples treeUnord <- resolveTreeChar(tree,trait,orderedChar = FALSE) treeOrd <- resolveTreeChar(tree,trait,orderedChar = TRUE,stateBias = NULL) treeOrdPrim <- resolveTreeChar(tree,trait,orderedChar = TRUE,stateBias = "primitive") treeOrdDer <- resolveTreeChar(tree,trait,orderedChar = TRUE,stateBias = "derived") #compare number of nodes Nnode(tree) #original Nnode(treeUnord) #unordered, biasStates = NULL, MPR Nnode(treeOrd) #ordered, biasStates = NULL Nnode(treeOrdPrim) #ordered, biasStates = 'primitive' Nnode(treeOrdDer) #ordered, biasStates = 'derived' #let's compare original tree with unordered-resolved tree layout(1:2) quickAncPlot(tree,trait,orderedChar = FALSE,cex = 0.3) text(x = 43,y = 10,"Original",cex = 1.5) quickAncPlot(treeUnord,trait,orderedChar = FALSE,cex = 0.3) text(x = 43,y = 10,"orderedChar = FALSE",cex = 1.5) #some resolution gained #now let's compare the original and ordered, both biasStates = NULL layout(1:2) quickAncPlot(tree,trait,orderedChar = FALSE,cex = 0.3) text(x = 43,y = 10,"Original",cex = 1.5) quickAncPlot(treeOrd,trait,orderedChar = TRUE,cex = 0.3) text(x = 43,y = 10,"orderedChar = TRUE",cex = 1.5) #now let's compare the three ordered trees layout(1:3) quickAncPlot(treeOrd,trait,orderedChar = TRUE,cex = 0.3) text(x = 41,y = 8,"ordered, biasStates = NULL",cex = 1.5) quickAncPlot(treeOrdPrim,trait,orderedChar = TRUE,cex = 0.3) text(x = 41.5,y = 8,"ordered, biasStates = 'primitive'",cex = 1.5) quickAncPlot(treeOrdDer,trait,orderedChar = TRUE,cex = 0.3) text(x = 42,y = 8,"ordered, biasStates = 'derived'",cex = 1.5) #let's compare unordered with ordered, biasStates = 'primitive' layout(1:2) quickAncPlot(treeUnord,trait,orderedChar = FALSE,cex = 0.3) text(x = 41,y = 8,"orderedChar = FALSE",cex = 1.5) quickAncPlot(treeOrdPrim,trait,orderedChar = TRUE,cex = 0.3) text(x = 40,y = 11,"orderedChar = TRUE",cex = 1.5) text(x = 40,y = 4,"biasStates = 'primitive'",cex = 1.5) #these comparisons will differ greatly between datasets # need to try them on your own layout(1)
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.