Title: | Base Classes and Functions for Phylogenetic Tree Input and Output |
---|---|
Description: | 'treeio' is an R package to make it easier to import and store phylogenetic tree with associated data; and to link external data from different sources to phylogeny. It also supports exporting phylogenetic tree with heterogeneous associated data to a single tree file and can be served as a platform for merging tree with associated data and converting file formats. |
Authors: | Guangchuang Yu [aut, cre] , Tommy Tsan-Yuk Lam [ctb, ths], Shuangbin Xu [ctb] , Bradley Jones [ctb], Casey Dunn [ctb], Tyler Bradley [ctb], Konstantinos Geles [ctb] |
Maintainer: | Guangchuang Yu <[email protected]> |
License: | Artistic-2.0 |
Version: | 1.31.0 |
Built: | 2024-10-30 09:17:31 UTC |
Source: | https://github.com/YuLab-SMU/treeio |
convert phylo to treedata
## S3 method for class 'phylo' as.treedata(tree, boot = NULL, ...)
## S3 method for class 'phylo' as.treedata(tree, boot = NULL, ...)
tree |
input tree, a |
boot |
optional, can be bootstrap value from ape::boot.phylo |
... |
additional parameters |
converting phylo object to treedata object
Guangchuang Yu
find the hierarchical cluster analysis among the nodes of graph based on the length of all the shortest paths in the graph.
find.hclust( x, graph.mst = FALSE, weights = NULL, hclust.method = "average", ... )
find.hclust( x, graph.mst = FALSE, weights = NULL, hclust.method = "average", ... )
x |
a igraph object |
graph.mst |
logical whether obtain the minimum spanning tree first then find.hclust, default is FALSE. |
weights |
a numeric vector giving edge weights or a character.
If this is |
hclust.method |
the agglomeration method to be used, This should be (an
unambiguous abbreviation of) one of |
... |
additional parameters |
hclust object
library(igraph) set.seed(123) g <- igraph::sample_gnp(100, .1) %>% set_edge_attr(name='weight', value=abs(rnorm(E(.),3))) tr1 <- find.hclust(g, weights = NA) tr2 <- find.hclust(g) tr3 <- find.hclust(g, graph.mst = TRUE)
library(igraph) set.seed(123) g <- igraph::sample_gnp(100, .1) %>% set_edge_attr(name='weight', value=abs(rnorm(E(.),3))) tr1 <- find.hclust(g, weights = NA) tr2 <- find.hclust(g) tr3 <- find.hclust(g, graph.mst = TRUE)
access placement information
get.placements(tree, ...) ## S3 method for class 'jplace' get.placements(tree, by = "best", ...)
get.placements(tree, ...) ## S3 method for class 'jplace' get.placements(tree, by = "best", ...)
tree |
tree object |
... |
additional parameters |
by |
one of 'best' and 'all' |
placement tibble
access phylo slot
get.tree(x, ...)
get.tree(x, ...)
x |
tree object |
... |
additional parameters |
phylo object
Guangchuang Yu
access tree text (newick text) from tree object
get.treetext(object, ...)
get.treetext(object, ...)
object |
treedata object |
... |
additional parameter |
phylo object
calculate total number of nodes
getNodeNum(tree) Nnode2(tree)
getNodeNum(tree) Nnode2(tree)
tree |
tree object |
number
Guangchuang Yu
getNodeNum(rtree(30)) Nnode2(rtree(30))
getNodeNum(rtree(30)) Nnode2(rtree(30))
test whether input object is produced by ggtree function
is.ggtree(x)
is.ggtree(x)
x |
object |
TRUE or FALSE
Guangchuang Yu
Class "jplace" This class stores phylogenetic placements
phylo
phylo object for tree structure
treetext
newick tree string
data
associated data
extraInfo
extra information, reserve for merge_tree
file
tree file
placements
reserve for jplace file to store placement information
info
extra information, e.g. metadata, software version etc.
Guangchuang Yu https://guangchuangyu.github.io
label branch for PAML to infer selection pressure using branch model
label_branch_paml(tree, node, label)
label_branch_paml(tree, node, label)
tree |
phylo object |
node |
node number |
label |
label of branch, e.g. #1 |
updated phylo object
Guangchuang Yu
site mask
mask(tree_object, field, site, mask_site = FALSE)
mask(tree_object, field, site, mask_site = FALSE)
tree_object |
tree object |
field |
selected field |
site |
site |
mask_site |
if TRUE, site will be masked. if FALSE, selected site will not be masked, while other sites will be masked. |
updated tree object
Guangchuang Yu
merge two tree object
merge_tree(obj1, obj2)
merge_tree(obj1, obj2)
obj1 |
tree object 1 |
obj2 |
tree object 2 |
tree object
Guangchuang Yu
print information of a list of treedata objects
## S3 method for class 'treedataList' print(x, ...)
## S3 method for class 'treedataList' print(x, ...)
x |
a list of treedata objects |
... |
no used |
message
convert raxml bootstrap tree to newick format
raxml2nwk(infile, outfile = "raxml.tree")
raxml2nwk(infile, outfile = "raxml.tree")
infile |
input file |
outfile |
output file |
newick file
Guangchuang Yu
parse ASTRAL output newick text
read.astral(file)
read.astral(file)
file |
ASTRAL Newick file |
treedata object
Guangchuang Yu
tt <- paste0( "((species1,(species2,species3)'[pp1=0.75;pp2=0.24;pp3=0.01]':", "1.2003685744180805)'[pp1=0.98;pp2=0.02;pp3=0]':0.9679599282730038,", "((species4,species5)'[pp1=0.88;pp2=0.11;pp3=0.01]':1.2454851536484994))" ) read.astral(textConnection(tt))
tt <- paste0( "((species1,(species2,species3)'[pp1=0.75;pp2=0.24;pp3=0.01]':", "1.2003685744180805)'[pp1=0.98;pp2=0.02;pp3=0]':0.9679599282730038,", "((species4,species5)'[pp1=0.88;pp2=0.11;pp3=0.01]':1.2454851536484994))" ) read.astral(textConnection(tt))
read beast/mrbayes/mega Nexus output
read beast/mrbayes/mega newick file format
read.beast(file, threads = 1, verbose = FALSE) read.mrbayes(file, threads = 1, verbose = FALSE) read.beast.newick(file, threads = 1, verbose = FALSE) read.mega(file, threads = 1, verbose = FALSE)
read.beast(file, threads = 1, verbose = FALSE) read.mrbayes(file, threads = 1, verbose = FALSE) read.beast.newick(file, threads = 1, verbose = FALSE) read.mega(file, threads = 1, verbose = FALSE)
file |
newick file |
threads |
number of threads for multithreading (default: 1) |
verbose |
set TRUE to log progress (default: FALSE) |
treedata object
treedata object
Guangchuang Yu https://guangchuangyu.github.io
Bradley R Jones
file <- system.file("extdata/BEAST", "beast_mcc.tree", package="treeio") read.beast(file) file <- system.file("extdata/MrBayes", "Gq_nxs.tre", package="treeio") read.mrbayes(file) tr <- read.beast.newick( textConnection( '(a[&rate=1]:2,(b[&rate=1.1]:1,c[&rate=0.9]:1)[&rate=1]:1);' ) )
file <- system.file("extdata/BEAST", "beast_mcc.tree", package="treeio") read.beast(file) file <- system.file("extdata/MrBayes", "Gq_nxs.tre", package="treeio") read.mrbayes(file) tr <- read.beast.newick( textConnection( '(a[&rate=1]:2,(b[&rate=1.1]:1,c[&rate=0.9]:1)[&rate=1]:1);' ) )
read baseml output
read.codeml(rstfile, mlcfile, tree = "mlc", type = "Joint")
read.codeml(rstfile, mlcfile, tree = "mlc", type = "Joint")
rstfile |
rst file |
mlcfile |
mlc file |
tree |
one of 'mlc' or 'rst' |
type |
one of 'Marginal' or 'Joint' |
A treedata
object
Guangchuang Yu
rstfile <- system.file("extdata/PAML_Codeml", "rst", package="treeio") mlcfile <- system.file("extdata/PAML_Codeml", "mlc", package="treeio") read.codeml(rstfile, mlcfile)
rstfile <- system.file("extdata/PAML_Codeml", "rst", package="treeio") mlcfile <- system.file("extdata/PAML_Codeml", "mlc", package="treeio") read.codeml(rstfile, mlcfile)
read mlc file of codeml output
read.codeml_mlc(mlcfile)
read.codeml_mlc(mlcfile)
mlcfile |
mlc file |
A codeml_mlc
object
Guangchuang Yu
mlcfile <- system.file("extdata/PAML_Codeml", "mlc", package="treeio") read.codeml_mlc(mlcfile)
mlcfile <- system.file("extdata/PAML_Codeml", "mlc", package="treeio") read.codeml_mlc(mlcfile)
read FASTA file
read.fasta(fasta, type = "auto")
read.fasta(fasta, type = "auto")
fasta |
fasta file |
type |
sequence type of the input file, one of 'NT' or 'AA'. Default is 'auto' and guess the sequence type automatically |
This function supports both DNA or AA sequences
DNAbin or AAbin object
Guangchuang Yu
read HYPHY output
read.hyphy(nwk, ancseq, tip.fasfile = NULL)
read.hyphy(nwk, ancseq, tip.fasfile = NULL)
nwk |
tree file in nwk format, one of hyphy output |
ancseq |
ancestral sequence file in nexus format, one of hyphy output |
tip.fasfile |
tip sequence file |
A hyphy object
Guangchuang Yu https://guangchuangyu.github.io
nwk <- system.file("extdata/HYPHY", "labelledtree.tree", package="treeio") ancseq <- system.file("extdata/HYPHY", "ancseq.nex", package="treeio") read.hyphy(nwk, ancseq)
nwk <- system.file("extdata/HYPHY", "labelledtree.tree", package="treeio") ancseq <- system.file("extdata/HYPHY", "ancseq.nex", package="treeio") read.hyphy(nwk, ancseq)
parse sequences from hyphy output
read.hyphy.seq(file)
read.hyphy.seq(file)
file |
output of hyphy ancestral sequence inference; nexus format |
DNAbin object
Guangchuang Yu
ancseq <- system.file("extdata/HYPHY", "ancseq.nex", package="treeio") read.hyphy.seq(ancseq)
ancseq <- system.file("extdata/HYPHY", "ancseq.nex", package="treeio") read.hyphy.seq(ancseq)
parse IQ-TREE output
read.iqtree(file)
read.iqtree(file)
file |
IQ-TREE Newick text |
treedata object
Guangchuang Yu
read jplace file
read.jplace(file)
read.jplace(file)
file |
jplace file |
jplace
instance
Guangchuang Yu
jp <- system.file("extdata", "sample.jplace", package="treeio") read.jplace(jp)
jp <- system.file("extdata", "sample.jplace", package="treeio") read.jplace(jp)
Import tree data from jtree file, which is JSON-based text and probably output by write.jtree
read.jtree(file)
read.jtree(file)
file |
tree file |
treedata object
Guangchuang Yu
read MCMCTree output Tree
read.mcmctree(file, force.ultrametric = FALSE)
read.mcmctree(file, force.ultrametric = FALSE)
file |
the output tree file of MCMCTree |
force.ultrametric |
logical whether convert the tree to be ultrametric, if it is not ultrametric, default is FALSE. When the tree is ultrametric, branch times will be calculated automatically. |
treedata object
file <- system.file("extdata/MCMCTree", "mcmctree_output.tree", package="treeio") tr <- read.mcmctree(file) tr
file <- system.file("extdata/MCMCTree", "mcmctree_output.tree", package="treeio") tr <- read.mcmctree(file) tr
parse tabular output of MEGA
read.mega_tabular(file)
read.mega_tabular(file)
file |
MEGA tabular file |
treedata object
Guangchuang Yu
read newick tree
read.newick(file, node.label = "label", ...)
read.newick(file, node.label = "label", ...)
file |
newick file |
node.label |
parse node label as 'label' or 'support' value |
... |
additional parameter, passed to 'read.tree' |
phylo or treedata object
Guangchuang Yu
read.nextstrain.json
read.nextstrain.json(x)
read.nextstrain.json(x)
x |
the json tree file of auspice from nextstrain. |
treedata object
Shuangbin Xu
file1 <- system.file("extdata/nextstrain.json", "minimal_v2.json", package="treeio") tr <- read.nextstrain.json(file1) tr
file1 <- system.file("extdata/nextstrain.json", "minimal_v2.json", package="treeio") tr <- read.nextstrain.json(file1) tr
read nhx tree file
read.nhx(file)
read.nhx(file)
file |
nhx file |
nhx object
Guangchuang Yu https://guangchuangyu.github.io
nhxfile <- system.file("extdata/NHX", "ADH.nhx", package="treeio") read.nhx(nhxfile)
nhxfile <- system.file("extdata/NHX", "ADH.nhx", package="treeio") read.nhx(nhxfile)
read rst file from paml (both baseml and codeml) output
read.paml_rst(rstfile, type = "Joint")
read.paml_rst(rstfile, type = "Joint")
rstfile |
rst file |
type |
one of 'Marginal' or 'Joint' |
A treedata
object
Guangchuang Yu https://guangchuangyu.github.io
rstfile <- system.file("extdata/PAML_Baseml", "rst", package="treeio") read.paml_rst(rstfile)
rstfile <- system.file("extdata/PAML_Baseml", "rst", package="treeio") read.paml_rst(rstfile)
parsing phylip tree format
read.phylip(file)
read.phylip(file)
file |
phylip file |
an instance of 'phylip'
Guangchuang Yu
phyfile <- system.file("extdata", "sample.phy", package="treeio") read.phylip(phyfile)
phyfile <- system.file("extdata", "sample.phy", package="treeio") read.phylip(phyfile)
read aligned sequences from phylip format
read.phylip.seq(file)
read.phylip.seq(file)
file |
phylip file, currently only sequential format is supported |
DNAbin object
Guangchuang Yu
http://evolution.genetics.washington.edu/phylip/doc/sequence.html
parse tree from phylip file
read.phylip.tree(file)
read.phylip.tree(file)
file |
phylip file |
phylo or multiPhylo object
Guangchuang Yu
read.phyloxml
read.phyloxml(file)
read.phyloxml(file)
file |
phyloxml file |
treedata class or treedataList class
xmlfile1 <- system.file("extdata/phyloxml", "test_x2.xml", package="treeio") px1 <- read.phyloxml(xmlfile1) px1 xmlfile2 <- system.file("extdata/phyloxml", "phyloxml_examples.xml", package="treeio") px2 <- read.phyloxml(xmlfile2) px2
xmlfile1 <- system.file("extdata/phyloxml", "test_x2.xml", package="treeio") px1 <- read.phyloxml(xmlfile1) px1 xmlfile2 <- system.file("extdata/phyloxml", "phyloxml_examples.xml", package="treeio") px2 <- read.phyloxml(xmlfile2) px2
parse output from r8s
read.r8s(file)
read.r8s(file)
file |
r8s output log file |
multiPhylo object
Guangchuang Yu
read.r8s(system.file("extdata/r8s", "H3_r8s_output.log", package="treeio"))
read.r8s(system.file("extdata/r8s", "H3_r8s_output.log", package="treeio"))
parse RAxML bootstrapping analysis output
read.raxml(file)
read.raxml(file)
file |
RAxML bootstrapping analysis output |
treedata object
Guangchuang Yu
raxml_file <- system.file("extdata/RAxML", "RAxML_bipartitionsBranchLabels.H3", package="treeio") read.raxml(raxml_file)
raxml_file <- system.file("extdata/RAxML", "RAxML_bipartitionsBranchLabels.H3", package="treeio") read.raxml(raxml_file)
read.treeqza
read.treeqza(treeqza, node.label = "label", ...)
read.treeqza(treeqza, node.label = "label", ...)
treeqza |
the qiime2 output file contained tree file. |
node.label |
parse node label as 'label' or 'support' value. |
... |
additional parameter, passed to 'read.tree'. |
phylo tree object or treedata object when node.label was parsed 'support'.
qzafile1 <- system.file("extdata/qiime2treeqza", "fasttree-tree.qza", package="treeio") qzafile2 <- system.file("extdata/qiime2treeqza", "iqt-tree.qza", package="treeio") qzafile3 <- system.file("extdata/qiime2treeqza", "raxml-cat-tree.qza", package="treeio") tr1 <- read.treeqza(qzafile1) tr1 tr2 <- read.treeqza(qzafile2) tr2 tr3 <- read.treeqza(qzafile3) tr3 # parse node label as 'support' value. qzafile4 <- system.file("extdata/qiime2treeqza", "raxml-cat-bootstrap-tree.qza", package="treeio") tr4 <- read.treeqza(qzafile4, node.label="support") tr4
qzafile1 <- system.file("extdata/qiime2treeqza", "fasttree-tree.qza", package="treeio") qzafile2 <- system.file("extdata/qiime2treeqza", "iqt-tree.qza", package="treeio") qzafile3 <- system.file("extdata/qiime2treeqza", "raxml-cat-tree.qza", package="treeio") tr1 <- read.treeqza(qzafile1) tr1 tr2 <- read.treeqza(qzafile2) tr2 tr3 <- read.treeqza(qzafile3) tr3 # parse node label as 'support' value. qzafile4 <- system.file("extdata/qiime2treeqza", "raxml-cat-bootstrap-tree.qza", package="treeio") tr4 <- read.treeqza(qzafile4, node.label="support") tr4
read timetree output
read.treetime(file) read.timetree(file)
read.treetime(file) read.timetree(file)
file |
the output tree file of timetree |
treedata object
rename tip label of phylogenetic tree
rename_taxa(tree, data, key = 1, value = 2)
rename_taxa(tree, data, key = 1, value = 2)
tree |
tree object, either treedata or phylo |
data |
data frame |
key |
column in data that match tip label (use 1st column by default) |
value |
column in data for rename tip label (use 2nd column by default) |
tree object
Guangchuang Yu
tree <- rtree(3) d <- data.frame(old = paste0('t', 1:3), new = LETTERS[1:3]) rename_taxa(tree, d) rename_taxa(tree, d, old, new)
tree <- rtree(3) d <- data.frame(old = paste0('t', 1:3), new = LETTERS[1:3]) rename_taxa(tree, d) rename_taxa(tree, d, old, new)
rescale branch length of tree object
rescale_tree(tree_object, branch.length)
rescale_tree(tree_object, branch.length)
tree_object |
tree object |
branch.length |
numerical features (e.g. dN/dS) |
update tree object
Guangchuang Yu
spt method
spt(x, from, to, weights = NULL, ...)
spt(x, from, to, weights = NULL, ...)
x |
a igraph object |
from |
a specific node of network. |
to |
other nodes of the network, length of it must be larger than 2. |
weights |
a numeric vector giving edge weights or a character.
If this is |
... |
additional parameters |
phylo object
library(igraph) set.seed(123) g <- igraph::sample_gnp(100, .1) %>% set_edge_attr(name='weight', value=abs(rnorm(E(.),3))) tr1 <- spt(g, from = 6, to=V(g), weights = 'weight') tr1 tr2 <- spt(g, from = 6, to = V(g), weights = NA) tr2
library(igraph) set.seed(123) g <- igraph::sample_gnp(100, .1) %>% set_edge_attr(name='weight', value=abs(rnorm(E(.),3))) tr1 <- spt(g, from = 6, to=V(g), weights = 'weight') tr1 tr2 <- spt(g, from = 6, to = V(g), weights = NA) tr2
Export treedata
object to BEAST NEXUS file. This function was adopted and modified from ape::write.nexus
write.beast(treedata, file = "", translate = TRUE, tree.name = NULL)
write.beast(treedata, file = "", translate = TRUE, tree.name = NULL)
treedata |
|
file |
output file. If file = "", print the output content on screen |
translate |
whether to translate taxa labels |
tree.name |
names of the trees, NULL to use existing tree names |
output file or file content on screen
Guangchuang Yu
nhxfile <- system.file("extdata/NHX", "phyldog.nhx", package="treeio") nhx <- read.nhx(nhxfile) write.beast(nhx)
nhxfile <- system.file("extdata/NHX", "phyldog.nhx", package="treeio") nhx <- read.nhx(nhxfile) write.beast(nhx)
Export treedata
object to BEAST Newick file. This is useful for making BEAST starting trees with metadata
write.beast.newick( treedata, file = "", append = FALSE, digits = 10, tree.prefix = "" )
write.beast.newick( treedata, file = "", append = FALSE, digits = 10, tree.prefix = "" )
treedata |
|
file |
output file. If file = "", print the output content on screen |
append |
logical. Only used if the argument 'file' is the name of file (and not a connection or "|cmd"). If 'TRUE' output will be appended to 'file'; otherwise, it will overwrite the contents of file. |
digits |
integer, the indicating the number of decimal places, default is 10. |
tree.prefix |
character the tree prefix, default is "". |
output file or file content on screen
Guangchuang Yu
nhxfile <- system.file("extdata/NHX", "phyldog.nhx", package="treeio") nhx <- read.nhx(nhxfile) write.beast.newick(nhx)
nhxfile <- system.file("extdata/NHX", "phyldog.nhx", package="treeio") nhx <- read.nhx(nhxfile) write.beast.newick(nhx)
Export jplace
object to jplace file.
write.jplace(x, outfile)
write.jplace(x, outfile)
x |
a jplace object. |
outfile |
the output file name |
jp <- system.file("extdata", "sample.jplace", package="treeio") tr1 <- read.jplace(jp) outfile <- tempfile() write.jplace(tr1, outfile) tr2 <- read.jplace(outfile) tr2
jp <- system.file("extdata", "sample.jplace", package="treeio") tr1 <- read.jplace(jp) outfile <- tempfile() write.jplace(tr1, outfile) tr2 <- read.jplace(outfile) tr2
Export treedata
object to json tree file
write.jtree(treedata, file = "")
write.jtree(treedata, file = "")
treedata |
|
file |
output file. If file = "", print the output content on screen |
output file or file content on screen
Guangchuang Yu