]> git.donarmstrong.com Git - ape.git/blob - R/as.phylo.R
some changes for ape 3.0-6
[ape.git] / R / as.phylo.R
1 ## as.phylo.R (2011-03-25)
2
3 ##     Conversion Among Tree Objects
4
5 ## Copyright 2005-2011 Emmanuel Paradis
6
7 ## This file is part of the R-package `ape'.
8 ## See the file ../COPYING for licensing issues.
9
10 old2new.phylo <- function(phy)
11 {
12     mode(phy$edge) <- "numeric"
13     phy$Nnode <- -min(phy$edge)
14     n <- length(phy$tip.label)
15     NODES <- phy$edge < 0
16     phy$edge[NODES] <- n - phy$edge[NODES]
17     phy
18 }
19
20 new2old.phylo <- function(phy)
21 {
22     NTIP <- length(phy$tip.label)
23     NODES <- phy$edge > NTIP
24     phy$edge[NODES] <- NTIP - phy$edge[NODES]
25     mode(phy$edge) <- "character"
26     phy$Nnode <- NULL
27     phy
28 }
29
30 as.phylo <- function (x, ...)
31 {
32     if (length(class(x)) == 1 && class(x) == "phylo")
33         return(x)
34     UseMethod("as.phylo")
35 }
36
37 as.phylo.hclust <- function(x, ...)
38 {
39     N <- dim(x$merge)[1]
40     edge <- matrix(0L, 2*N, 2)
41     edge.length <- numeric(2*N)
42     ## `node' gives the number of the node for the i-th row of x$merge
43     node <- integer(N)
44     node[N] <- N + 2L
45     cur.nod <- N + 3L
46     j <- 1L
47     for (i in N:1) {
48         edge[j:(j + 1), 1] <- node[i]
49         for (l in 1:2) {
50             k <- j + l - 1L
51             y <- x$merge[i, l]
52             if (y > 0) {
53                 edge[k, 2] <- node[y] <- cur.nod
54                 cur.nod <- cur.nod + 1L
55                 edge.length[k] <- x$height[i] - x$height[y]
56             } else {
57                 edge[k, 2] <- -y
58                 edge.length[k] <- x$height[i]
59             }
60         }
61         j <- j + 2L
62     }
63     if (is.null(x$labels))
64         x$labels <- as.character(1:(N + 1))
65     obj <- list(edge = edge, edge.length = edge.length / 2,
66                 tip.label = x$labels, Nnode = N)
67     class(obj) <- "phylo"
68     reorder(obj)
69 }
70
71 as.phylo.phylog <- function(x, ...)
72 {
73     tr <- read.tree(text = x$tre)
74     n <- length(tr$tip.label)
75     edge.length <- numeric(dim(tr$edge)[1])
76     term  <- which(tr$edge[, 2] <= n)
77     inte  <- which(tr$edge[, 2] > n)
78     edge.length[term] <- x$leaves[tr$tip.label]
79     edge.length[inte] <- x$nodes[tr$node.label][-1]
80     tr$edge.length <- edge.length
81     if (x$nodes["Root"] != 0) {
82         tr$edge.root <- x$nodes["Root"]
83         names(tr$edge.root) <- NULL
84     }
85     tr
86 }
87
88 as.hclust.phylo <- function(x, ...)
89 {
90     if (!is.ultrametric(x)) stop("the tree is not ultrametric")
91     if (!is.binary.tree(x)) stop("the tree is not binary")
92     if (!is.rooted(x)) stop("the tree is not rooted")
93     n <- length(x$tip.label)
94     x$node.label <- NULL # by Jinlong Zhang (2010-12-15)
95     bt <- sort(branching.times(x))
96     inode <- as.numeric(names(bt))
97     N <- n - 1L
98     nm <- numeric(N + n) # hash table
99     nm[inode] <- 1:N
100     merge <- matrix(NA, N, 2)
101     for (i in 1:N) {
102         ind <- which(x$edge[, 1] == inode[i])
103         for (k in 1:2) {
104             tmp <- x$edge[ind[k], 2]
105             merge[i, k] <- if (tmp <= n) -tmp else nm[tmp]
106         }
107     }
108     names(bt) <- NULL
109     obj <- list(merge = merge, height = bt, order = 1:n, labels = x$tip.label,
110                 call = match.call(), method = "unknown")
111     class(obj) <- "hclust"
112     obj
113 }
114
115 as.network.phylo <- function(x, directed = is.rooted(x), ...)
116 {
117     if (is.null(x$node.label)) x <- makeNodeLabel(x)
118     res <- network(x$edge, directed = directed, ...)
119     network.vertex.names(res) <- c(x$tip.label, x$node.label)
120     res
121 }
122
123 as.igraph <- function(x, ...) UseMethod("as.igraph")
124
125 as.igraph.phylo <- function(x, directed = is.rooted(x), use.labels = TRUE, ...)
126 {
127     ## local copy because x will be changed before evaluating is.rooted(x):
128     directed <- directed
129     if (use.labels) {
130         if (is.null(x$node.label)) x <- makeNodeLabel(x)
131         x$edge <- matrix(c(x$tip.label, x$node.label)[x$edge], ncol = 2)
132     } else x$edge <- x$edge - 1L
133     graph.edgelist(x$edge, directed = directed, ...)
134 }