TreeOfLife.jl
TreeOfLife — ModuleTreeOfLife.jl
TreeOfLife.jl defines data types for cladograms and chronograms in phylogenetics and methods analyzing these trees.
In development.
Alternative packages include:
- Phylo.jl ("in beta", and its outdated prototype Phylogenetics.jl)
 - PhyloTrees.jl
 - TreeTools.jl
 - Phylogenies.jl ("in development")
 
Examples
The followings are some simple examples illustrating the usage of this package. Please see the documentation for more details of TreeOfLife.jl.
julia> using TreeOfLife
julia> tree = fromnewick("(A:0.1,B:0.2,(C:0.3,D:0.4)E:0.5)F;")
ChronoTree(ChronoNode[ChronoNode("F", 0, 0, 2, 0.0, 0.0), ChronoNode("A", 1, 3, 0, 0.1, 0.1), ChronoNode("B", 1, 4, 0, 0.2, 0.2), ChronoNode("E", 1, 0, 5, 0.5, 0.5), ChronoNode("C", 4, 6, 0, 0.8, 0.3), ChronoNode("D", 4, 0, 0, 0.9, 0.4)])
julia> tonewick(tree)
"(A:0.1,B:0.2,(C:0.3,D:0.4)E:0.5)F;"
julia> tipnames = gettipnames(tree)
4-element Vector{String}:
 "A"
 "B"
 "C"
 "D"
julia> getmrca(tree, tipnames)
1
julia> phylodiv(tree, tipnames)
1.5
julia> ismonophyl(tree, ["A", "B"])
false
julia> ismonophyl(tree, ["C", "D"])
trueTypes
TreeOfLife.AbstractNode — TypeAbstractNodeAbstract supertype for ChronoNode and CladoNode.
TreeOfLife.ChronoNode — TypeChronoNode <: AbstractNodeType for nodes in a ChronoTree. Compare CladoNode.
TreeOfLife.CladoNode — TypeCladoNode <: AbstractNode
CladoNode(node::ChronoNode) :: CladoNodeType for nodes in a CladoTree. 
A ChronoNode can be converted to a CladoNode by removing all  information about time.
TreeOfLife.AbstractTree — TypeAbstractTreeAbstract supertype for ChronoTree and CladoTree.
TreeOfLife.ChronoTree — TypeChronoTree <: AbstractTree
ChronoTree(nodes::AbstractVector{ChronoNode}) :: ChronoTree
ChronoTree() :: ChronoTreeType for chronograms or dated phylogenetic trees, assumed to be rooted,  comprising ChronoNode instances. Compare CladoTree. 
When called with no arguments, the constructor returns an empty ChronoTree.
TreeOfLife.CladoTree — TypeCladoTree <: AbstractTree
CladoTree() :: CladoTree
CladoTree(nodes::AbstractVector{CladoNode}) :: CladoTree
CladoTree(tree::CladoTree) :: CladoTree
CladoTree(tree::ChronoTree) :: CladoTreeType for cladograms or undated phylogenetic trees, assumed to be rooted,  comprising CladoNode instances. 
When called with no arguments, the constructor returns an empty CladoTree.
A ChronoTree can be converted to a CladoTree by removing all  information about time.
Base.length — Functionlength(tree::AbstractTree) :: IntReturn the number of nodes in a phylogenetic tree.
Base.getindex — Functiongetindex(tree::AbstractTree, i) :: Any
getindex(tree::AbstractTree, i::Integer, trait::Symbol) :: Any
getindex(tree::AbstractTree, trait::Symbol) :: Dict{Int, <:Any}Get the i-th node of the tree, or its trait, or a Dict containing the  trait of all nodes of the tree. See also setindex!.
Base.setindex! — Functionsetindex!(tree::AbstractTree, node::AbstractNode, i::Integer) :: Vector
setindex!(tree::AbstractTree, tv, i::Integer, trait::Symbol) :: DictSet the i-th node of the tree or its trait. See also getindex.
Base.haskey — Functionhaskey(tree::AbstractTree, trait::Symbol) :: Bool
haskey(tree::AbstractTree, i::Integer, trait::Symbol) :: BoolTest if the tree or its i-th node has a trait or a trait value.
Base.empty — Functionempty(tree::AbstractTree) :: AbstractTreeConstruct a phylogenetic tree with only the root node of the same type.
Base.:== — Function==(node1::CladoNode, node2::CladoNode) :: Bool
==(node1::ChronoNode, node2::ChronoNode) :: BoolTest if two nodes are identical, in the sense that they have the same name, the same parent, sibling, and child, as well as they have approximate branch lengths.
==(tree1::AbstractTree, tree2::AbstractTree) :: BoolTest if two trees are identical, in the sense that they are of the same type,  isomorphic, and have the same node ordering; specifically, for dated trees or  ChronoTree instances, the node times are correspondingly equal. 
Identical trees are always isomorphic (can be tested by isisomorph).
Newick format
TreeOfLife.readnewick — Functionreadnewick(filename::AbstractString) :: AbstractTreeRead a phylogenetic tree from disk. Compare writenewick.
TreeOfLife.writenewick — Functionreadnewick(filename::AbstractString, tree::AbstractTree) :: AbstractTreeWrite a phylogenetic tree to disk. Compare readnewick.
TreeOfLife.fromnewick — Functionfromnewick(str::AbstractString; nocomments::Bool=false) :: AbstractTreeCreate a phylogenetic tree from a Newick-format tree string.
Its inverse function is tonewick.
The argument nocomments controls whether the comments (enclosed by square  brackets) are wiped out; by default it is set to false, i.e., all comments  are kept.
TreeOfLife.tonewick — Functiontonewick(tree::CladoTree) :: String
tonewick(tree::ChronoTree) :: StringExpress a phylogenetic tree as a Newick-format string.
Its inverse function is fromnewick.
Internal functions
TreeOfLife.ignore_comments — Functionignore_comments(str::AbstractString) :: StringRemove comments enclosed by square brackets (possibly nested) in a string.
TreeOfLife.from_newick — Functionfrom_newick(str::AbstractString) :: VectorParse a Newick-format tree string into segments of different types according  to their syntactical meanings. Used in fromnewick.
TreeOfLife.get_tree_type — Functionget_tree_type(elements::Vector) :: 
	Tuple{Type{<:AbstractTree}, Type{<:AbstractNode}}Test whether a vector of elements parsed from some Newick-syntax tree string  is a chronogram or a cladogram. Used in fromnewick.
TreeOfLife.to_newick! — Functionto_newick!(elements::AbstractVector, 
	tree::CladoTree, i::Integer) :: AbstractVector
to_newick!(elements::AbstractVector, 
	tree::ChronoTree, i::Integer) :: AbstractVectorConvert a phylogenetic tree to segments of a Newick-format string. Used in  tonewick.
Nexus format
TreeOfLife.readnexus — Functionreadnexus(filename::AbstractString; every=0) :: Vector{ChronoTree}Read dated phylogenetic trees from a Nexus-format file on disk.
Under development; correctness not guaranteed.
Methods involving one tree
TreeOfLife.gettips — Functiongettips(tree::AbstractTree) :: Vector{Int}Return the indices of all tip nodes of the tree.
TreeOfLife.gettipnames — Functiongettipnames(tree::AbstractTree) :: Vector{String}Return the names of all tip nodes of the tree.
TreeOfLife.alldistinct — Functionalldistinct(tree::AbstractTree) :: BoolTest if all tip nodes of the tree have distinct names.
TreeOfLife.getage — Functiongetage(tree::ChronoTree; 
	average=mean_, getrelerr::Bool=false, reltol=1e-6) :: Float64Return an average age (from the root node) of tip nodes of the tree.
The argument average defines the method for summarizing the ages to one; by  default it is set to mean_.
The argument getrelerr controls whether the relative standard deviation is  appended to the output ((mean, relstd)) or not (only mean); by default it  is set to false.
The argument reltol is a tolerance of relative error. By default it is set  to 1e-6. To suppress the judgment, set reltol=NaN.
TreeOfLife.getages — Functiongetages(tree::ChronoTree; average=mean_, reltol=1e-6) :: Vector{Float64}Return ages (from the root node) of all nodes of the tree.
The argument average defines the method for summarizing the ages to one; by  default it is set to mean_.
The argument reltol is a tolerance of relative error. By default it is set  to 1e-6. To suppress the judgment, set reltol=NaN.
TreeOfLife.getparent — Functiongetparent(tree::AbstractTree, i::Integer) :: IntFind the parent or direct ancestor of tree[i] and return its index.
TreeOfLife.getchildren — Functiongetchildren(tree::AbstractTree, i::Integer) :: Vector{Int}Find all children or direct descendents of tree[i] and return their indices.
TreeOfLife.preorder — Functionpreorder(tree::AbstractTree, i=1) :: Vector{Int}Return the pre-order traversal sequence of the whole tree, or its subtree  with root node tree[i].
TreeOfLife.postorder — Functionpostorder(tree::AbstractTree, i=1) :: Vector{Int}Return the post-order traversal sequence of the whole tree, or its subtree  with root node tree[i].
TreeOfLife.isroot — Functionisroot(tree::AbstractTree, i::Integer) :: BoolTest if the i-th node of the tree is the root.
TreeOfLife.istip — Functionistip(tree::AbstractTree, i::Integer) :: BoolTest if the i-th node of the tree is a tip or leaf node.
TreeOfLife.getname — Functiongetname(node::AbstractNode) :: String
getname(tree::AbstractTree, i::Integer) :: StringExtract the name of the given node as a string.
TreeOfLife.hassibling — Functionhassibling(tree::AbstractTree, i::Integer) :: BoolTest if the i-th node of the tree has following sibling(s).
TreeOfLife.rename — Functionrename(oldtree::AbstractTree, 
	oldtonew::Dict{<:AbstractString,<:AbstractString}Create a new tree whose nodes are respectively renamed from the oldtree by  a mapping from old names to new names. Specifically, nodes with empty names  remain.
TreeOfLife.rename! — Functionrename!(tree::Tree, 
	oldtonew::Dict{<:AbstractString,<:AbstractString}Rename nodes of the tree in place by a mapping from old names to new names. Specifically, nodes with empty names remain.
TreeOfLife.getsubtree — Functiongetsubtree(oldtree::AbstractTree, tipset; 
	simplify::Bool=true, keeproot::Bool=false) :: AbstractTreeExtract the subtree generated from a given set of tips of the tree.
The argument simplify controls whether internal node with only one child  needs to be reduced, i.e., connecting directly its child and its parent; by  default it is set to true. 
The argument keeproot controls whether the original root node needs to be  contained in the subtree; by default it is set to false, in other words,  yielding a truly minimum spanning tree (MST). 
When simplify is set to false, the value of keeproot has no effect.
TreeOfLife.getmrca — Functiongetmrca(tree::AbstractTree, tipset) :: IntFind the index of the most recent common ancestor node for a set of nodes.
TreeOfLife.ismonophyl — Functionismonophyl(tree::AbstractTree, tipset) :: BoolTest if a given set of tip nodes are monophyletic based on the tree.
TreeOfLife.getphylodiv — Functiongetphylodiv(tree::ChronoTree, tipset; keeproot::Bool=false)Compute the phylogenetic diversity (PD) of a given set of tips of the tree, i.e., the sum of branch lengths of the subtree generated from the set.
The argument keeproot controls whether the original root node needs to be  contained in the subtree; by default it is set to false.
TreeOfLife.cutfromroot — Functioncutfromroot(tree::ChronoTree, dist::Real; keep::Symbol=:both)
	:: Union{Vector{NTuple{2,Int}}, Vector{Int}}Find the temporal section by dist time units after the root is born.  The argument keep can be set among four options, i.e., :both (tuples  containing parents and childs), :parent, :child, or :closer. 
TreeOfLife.cutfromtips — Functioncutfromtips(tree::ChronoTree, dist::Real; 
	keep::Symbol=:both, average=minimum, reltol=1e-6) 
	:: Union{Vector{NTuple{2,Int}}, Vector{Int}}Find the temporal section by dist time units before the root is born. The argument keep can be set among four options, i.e., :both (tuples  containing parents and childs), :parent, :child, or :closer. 
TreeOfLife.isbinary — Functionisbinary(tree::AbstractTree) :: BoolTest if the tree is strictly binary or dichotonous, i.e., all non-tip nodes have exactly two descendents.
TreeOfLife.isisomorph — Functionisisomorph(tree1::CladoTree, tree2::CladoTree) :: Bool
isisomorph(tree1::ChronoTree, tree2::ChronoTree) :: BoolTest if two trees are isomorphic.
When both phylogenetic tree are dated, the isomorphism implies that branch lengths are correspondingly equal; otherwise, only the tree topology are compared.
TreeOfLife.getdescs — Functiongetdescs(tree::AbstractTree, mrca::Int) :: Vector{Int}Find the indices of all tip descendents from a common ancestor node.
TreeOfLife.getdescnames — Functiongetdescnames(tree::AbstractTree, mrca::Int) :: Vector{String}
getdescnames(tree::AbstractTree) :: Vector{Vector{String}}Find the names of all tip descendents from a common ancestor node, or such descendent name lists for all nodes of the tree.
Internal functions
TreeOfLife.calibrate_t_root! — Functioncalibrate_t_root!(tree::ChronoTree) :: ChronoTree
calibrate_t_root!(tree::AbstractTree) :: AbstractTreeCalculate all t_root values according to t_branch values. Used in  fromnewick.
TreeOfLife.calibrate_t_branch! — Functioncalibrate_t_branch!(tree::ChronoTree) :: ChronoTree
calibrate_t_branch!(tree::AbstractTree) :: AbstractTreeCalibrate all t_root values of nodes of the tree so that the root's t_root  is zero, and then recalculate all t_branch values according to the new  t_root values. Used in getsubtree.
TreeOfLife.mean_ — Functionmean_(a::Vector{Float64}) :: Float64Compute the arithmetic mean of a vector of 64-bit float numbers.
TreeOfLife.preorder! — Functionpreorder!(sequence, tree::AbstractTree, i=1) :: NothingAppend the pre-order traversal sequence of the whole tree, or its subtree  with root node tree[i].
TreeOfLife.postorder! — Functionpostorder!(sequence, tree::AbstractTree, i=1) :: NothingAppend the post-order traversal sequence of the whole tree, or its subtree  with root node tree[i].
TreeOfLife.get_counts — Functionget_selected(oldtree::AbstractTree, tipset; 
	simplify::Bool=true, keeproot::Bool=false) :: Vector{Int}Select nodes of a subtree generated from a given set of tips of the tree. Used  in getsubtree and isbinary. 
Arguments simplify and keeproot have same meanings as in getsubtree. 
TreeOfLife.sum_t_branch — Functionsum_t_branch(tree::ChronoTree)Compute the sum of branch lengths of the tree. Used in getphylodiv. 
TreeOfLife.tree_hash — Functiontree_hash(tree::CladoTree, h::UInt=zero(UInt)) :: UInt
tree_hash(tree::ChronoTree, h::UInt=zero(UInt)) :: UIntCompute a hash value for a phylogenetic tree so that isomorphic trees  necessarily have the same hash value (tested by isisomorph).
Methods involving multiple trees
TreeOfLife.getconsensus — Functiongetconsensus(trees::Vector{<:AbstractTree}; 
	threshold::Float64=0.5, every::Int=0) :: CladoTreeSummarize phylogenetic trees to one such that clades with support rate no less than a threshold remain in the consensus tree.
The argument threshold should be a number between 0.5 and 1.0; the default  value is 0.5.
The argument every, if set to a positive integer, makes the (probably long)  analyzing process print a log per every trees; the default value is 0.
Internal functions
TreeOfLife.count_clade_ages — Functioncount_clade_ages(trees::Vector{ChronoTree}; every::Int=0) 
	:: Dict{Set{String}, Vector{Float64}}Count the age for every possible tip node combination as long as occured at least once in some of the trees.
The argument every, if set to a positive integer, makes the (probably long)  analyzing process print a log per every trees; the default value is 0.
TreeOfLife.count_clades — Functioncount_clades(trees::Vector{<:AbstractTree}; every::Int=0)
	:: Dict{Set{String}, Int}Count every possible tip node combination as long as occured at least once in  some of the trees. Used in getconsensus.
The argument every, if set to a positive integer, makes the (probably long)  analyzing process print a log per every trees; the default value is 0.
TreeOfLife.construct_tree — Functionconstruct_tree(parents::Vector{Int}) :: Tuple{CladoTree, Int}Construct a cladogram from parent indices. Used in getconsensus.
Index
TreeOfLifeTreeOfLife.AbstractNodeTreeOfLife.AbstractTreeTreeOfLife.ChronoNodeTreeOfLife.ChronoTreeTreeOfLife.CladoNodeTreeOfLife.CladoTreeBase.:==Base.emptyBase.getindexBase.haskeyBase.lengthBase.setindex!TreeOfLife.alldistinctTreeOfLife.calibrate_t_branch!TreeOfLife.calibrate_t_root!TreeOfLife.construct_treeTreeOfLife.count_clade_agesTreeOfLife.count_cladesTreeOfLife.cutfromrootTreeOfLife.cutfromtipsTreeOfLife.from_newickTreeOfLife.fromnewickTreeOfLife.get_countsTreeOfLife.get_tree_typeTreeOfLife.getageTreeOfLife.getagesTreeOfLife.getchildrenTreeOfLife.getconsensusTreeOfLife.getdescnamesTreeOfLife.getdescsTreeOfLife.getmrcaTreeOfLife.getnameTreeOfLife.getparentTreeOfLife.getphylodivTreeOfLife.getsubtreeTreeOfLife.gettipnamesTreeOfLife.gettipsTreeOfLife.hassiblingTreeOfLife.ignore_commentsTreeOfLife.isbinaryTreeOfLife.isisomorphTreeOfLife.ismonophylTreeOfLife.isrootTreeOfLife.istipTreeOfLife.mean_TreeOfLife.postorderTreeOfLife.postorder!TreeOfLife.preorderTreeOfLife.preorder!TreeOfLife.readnewickTreeOfLife.readnexusTreeOfLife.renameTreeOfLife.rename!TreeOfLife.sum_t_branchTreeOfLife.to_newick!TreeOfLife.tonewickTreeOfLife.tree_hashTreeOfLife.writenewick