]> git.donarmstrong.com Git - ape.git/blob - R/dist.topo.R
f2a10ae26352e16325763cea2cfda19fa5bcad46
[ape.git] / R / dist.topo.R
1 ## dist.topo.R (2010-01-22)
2
3 ##      Topological Distances, Tree Bipartitions,
4 ##   Consensus Trees, and Bootstrapping Phylogenies
5
6 ## Copyright 2005-2010 Emmanuel Paradis
7
8 ## This file is part of the R-package `ape'.
9 ## See the file ../COPYING for licensing issues.
10
11 dist.topo <- function(x, y, method = "PH85")
12 {
13     if (method == "score" && (is.null(x$edge.length) || is.null(y$edge.length)))
14         stop("trees must have branch lengths for Billera et al.'s distance.")
15     nx <- length(x$tip.label)
16     bp1 <- .Call("bipartition", x$edge, nx, x$Nnode, PACKAGE = "ape")
17     bp1 <- lapply(bp1, function(xx) sort(x$tip.label[xx]))
18     ny <- length(y$tip.label) # fix by Otto Cordero
19     ## fix by Tim Wallstrom:
20     bp2.tmp <- .Call("bipartition", y$edge, ny, y$Nnode, PACKAGE = "ape")
21     bp2 <- lapply(bp2.tmp, function(xx) sort(y$tip.label[xx]))
22     bp2.comp <- lapply(bp2.tmp, function(xx) setdiff(1:ny, xx))
23     bp2.comp <- lapply(bp2.comp, function(xx) sort(y$tip.label[xx]))
24     ## End
25     q1 <- length(bp1)
26     q2 <- length(bp2)
27     if (method == "PH85") {
28         p <- 0
29         for (i in 1:q1) {
30             for (j in 1:q2) {
31                 if (identical(bp1[[i]], bp2[[j]]) | identical(bp1[[i]], bp2.comp[[j]])) {
32                     p <- p + 1
33                     break
34                 }
35             }
36         }
37         dT <- q1 + q2 - 2 * p # same than:
38         ##dT <- if (q1 == q2) 2*(q1 - p) else 2*(min(q1, q2) - p) + abs(q1 - q2)
39     }
40     if (method == "score") {
41         dT <- 0
42         found1 <- FALSE
43         found2 <- logical(q2)
44         found2[1] <- TRUE
45         for (i in 2:q1) {
46             for (j in 2:q2) {
47                 if (identical(bp1[[i]], bp2[[j]]) | identical(bp1[[i]], bp2.comp[[j]])) {
48                     if (i == 19) browser()
49                     dT <- dT + (x$edge.length[which(x$edge[, 2] == nx + i)] -
50                                 y$edge.length[which(y$edge[, 2] == ny + j)])^2
51                     found1 <- found2[j] <- TRUE
52                     break
53                 }
54             }
55             if (found1) found1 <- FALSE
56             else dT <- dT + (x$edge.length[which(x$edge[, 2] == nx + i)])^2
57         }
58         if (!all(found2))
59             dT <- dT + sum((y$edge.length[y$edge[, 2] %in% (ny + which(!found2))])^2)
60         dT <- sqrt(dT)
61     }
62     dT
63 }
64
65 .compressTipLabel <- function(x)
66 {
67     ## 'x' is a list of objects of class "phylo" possibly with no class
68     if (!is.null(attr(x, "TipLabel"))) return(x)
69     ref <- x[[1]]$tip.label
70     if (any(table(ref) != 1))
71         stop("some tip labels are duplicated in tree no. 1")
72     n <- length(ref)
73     for (i in 2:length(x)) {
74         if (identical(x[[i]]$tip.label, ref)) next
75         ilab <- match(x[[i]]$tip.label, ref)
76         ## can use tabulate here because 'ilab' contains integers
77         if (any(tabulate(ilab) > 1))
78             stop(paste("some tip labels are duplicated in tree no.", i))
79         if (any(is.na(ilab)))
80             stop(paste("tree no.", i, "has different tip labels"))
81         ie <- match(1:n, x[[i]]$edge[, 2])
82         x[[i]]$edge[ie, 2] <- ilab
83     }
84     for (i in 1:length(x)) x[[i]]$tip.label <- NULL
85     attr(x, "TipLabel") <- ref
86     x
87 }
88
89 prop.part <- function(..., check.labels = TRUE)
90 {
91     obj <- list(...)
92     if (length(obj) == 1 && class(obj[[1]]) != "phylo")
93         obj <- obj[[1]]
94     ## <FIXME>
95     ## class(obj) <- NULL # needed?
96     ## </FIXME>
97     ntree <- length(obj)
98     if (ntree == 1) check.labels <- FALSE
99     if (check.labels) obj <- .compressTipLabel(obj)
100     for (i in 1:ntree) storage.mode(obj[[i]]$Nnode) <- "integer"
101     ## <FIXME>
102     ## The 1st must have tip labels
103     ## Maybe simply pass the number of tips to the C code??
104     if (!is.null(attr(obj, "TipLabel")))
105         for (i in 1:ntree) obj[[i]]$tip.label <- attr(obj, "TipLabel")
106     ## </FIXME>
107     clades <- .Call("prop_part", obj, ntree, TRUE, PACKAGE = "ape")
108     attr(clades, "number") <- attr(clades, "number")[1:length(clades)]
109     attr(clades, "labels") <- obj[[1]]$tip.label
110     class(clades) <- "prop.part"
111     clades
112 }
113
114 print.prop.part <- function(x, ...)
115 {
116     if (is.null(attr(x, "labels"))) {
117         for (i in 1:length(x)) {
118             cat("==>", attr(x, "number")[i], "time(s):")
119             print(x[[i]], quote = FALSE)
120         }
121     } else {
122         for (i in 1:length(attr(x, "labels")))
123           cat(i, ": ", attr(x, "labels")[i], "\n", sep = "")
124         cat("\n")
125         for (i in 1:length(x)) {
126             cat("==>", attr(x, "number")[i], "time(s):")
127             print(x[[i]], quote = FALSE)
128         }
129     }
130 }
131
132 summary.prop.part <- function(object, ...) attr(object, "number")
133
134 plot.prop.part <- function(x, barcol = "blue", leftmar = 4, ...)
135 {
136     if (is.null(attr(x, "labels")))
137       stop("cannot plot this partition object; see ?prop.part for details.")
138     L <- length(x)
139     n <- length(attr(x, "labels"))
140     layout(matrix(1:2, 2, 1), heights = c(1, 3))
141     par(mar = c(0.1, leftmar, 0.1, 0.1))
142     plot(1:L, attr(x, "number"), type = "h", col = barcol, xlim = c(1, L),
143          xlab = "", ylab = "Frequency", xaxt = "n", bty = "n")
144     plot(0, type = "n", xlim = c(1, L), ylim = c(1, n),
145          xlab = "", ylab = "", xaxt = "n", yaxt = "n")
146     for (i in 1:L) points(rep(i, length(x[[i]])), x[[i]], ...)
147     mtext(attr(x, "labels"), side = 2, at = 1:n, las = 1)
148 }
149
150 prop.clades <- function(phy, ..., part = NULL)
151 {
152     if (is.null(part)) {
153         obj <- list(...)
154         if (length(obj) == 1 && class(obj[[1]]) != "phylo")
155           obj <- unlist(obj, recursive = FALSE)
156         part <- prop.part(obj, check.labels = TRUE)
157     }
158     bp <- .Call("bipartition", phy$edge, length(phy$tip.label),
159                 phy$Nnode, PACKAGE = "ape")
160     if (!is.null(attr(part, "labels")))
161       for (i in 1:length(part))
162         part[[i]] <- sort(attr(part, "labels")[part[[i]]])
163     bp <- lapply(bp, function(xx) sort(phy$tip.label[xx]))
164     n <- numeric(phy$Nnode)
165     for (i in 1:phy$Nnode) {
166         for (j in 1:length(part)) {
167             if (identical(all.equal(bp[[i]], part[[j]]), TRUE)) {
168                 n[i] <- attr(part, "number")[j]
169                 done <-  TRUE
170                 break
171             }
172         }
173     }
174     n
175 }
176
177 boot.phylo <- function(phy, x, FUN, B = 100, block = 1, trees = FALSE)
178 {
179     if (is.list(x) && !is.data.frame(x)) {
180         if (inherits(x, "DNAbin")) x <- as.matrix(x)
181         else {
182             nm <- names(x)
183             n <- length(x)
184             x <- unlist(x)
185             nL <- length(x)
186             x <- matrix(x, n, nL/n, byrow = TRUE)
187             rownames(x) <- nm
188         }
189     }
190     boot.tree <- vector("list", B)
191     for (i in 1:B) {
192         if (block > 1) {
193             y <- seq(block, ncol(x), block)
194             boot.i <- sample(y, replace = TRUE)
195             boot.samp <- numeric(ncol(x))
196             boot.samp[y] <- boot.i
197             for (j in 1:(block - 1))
198               boot.samp[y - j] <- boot.i - j
199         } else boot.samp <- sample(ncol(x), replace = TRUE)
200         boot.tree[[i]] <- FUN(x[, boot.samp])
201     }
202     for (i in 1:B) storage.mode(boot.tree[[i]]$Nnode) <- "integer"
203     storage.mode(phy$Nnode) <- "integer"
204     ans <- attr(.Call("prop_part", c(list(phy), boot.tree),
205                       B + 1, FALSE, PACKAGE = "ape"), "number") - 1
206     if (trees) ans <- list(BP = ans, trees = boot.tree)
207     ans
208 }
209
210 consensus <- function(..., p = 1, check.labels = TRUE)
211 {
212     foo <- function(ic, node) {
213         ## ic: index of 'pp'
214         ## node: node number in the final tree
215         pool <- pp[[ic]]
216         if (ic < m) {
217             for (j in (ic + 1):m) {
218                 wh <- match(pp[[j]], pool)
219                 if (!any(is.na(wh))) {
220                     edge[pos, 1] <<- node
221                     pool <- pool[-wh]
222                     edge[pos, 2] <<- nextnode <<- nextnode + 1L
223                     pos <<- pos + 1L
224                     foo(j, nextnode)
225                 }
226             }
227         }
228         size <- length(pool)
229         if (size) {
230             ind <- pos:(pos + size - 1)
231             edge[ind, 1] <<- node
232             edge[ind, 2] <<- pool
233             pos <<- pos + size
234         }
235     }
236     obj <- list(...)
237     if (length(obj) == 1) {
238         ## better than unlist(obj, recursive = FALSE)
239         ## because "[[" keeps the class of 'obj':
240         obj <- obj[[1]]
241         if (class(obj) == "phylo") return(obj)
242     }
243     if (!is.null(attr(obj, "TipLabel")))
244         labels <- attr(obj, "TipLabel")
245     else {
246         labels <- obj[[1]]$tip.label
247         if (check.labels) obj <- .compressTipLabel(obj)
248     }
249     ntree <- length(obj)
250     ## Get all observed partitions and their frequencies:
251     pp <- prop.part(obj, check.labels = FALSE)
252     ## Drop the partitions whose frequency is less than 'p':
253     pp <- pp[attr(pp, "number") >= p * ntree]
254     ## Get the order of the remaining partitions by decreasing size:
255     ind <- sort(unlist(lapply(pp, length)), decreasing = TRUE,
256                 index.return = TRUE)$ix
257     pp <- pp[ind]
258     n <- length(labels)
259     m <- length(pp)
260     edge <- matrix(0L, n + m - 1, 2)
261     if (m == 1) {
262         edge[, 1] <- n + 1L
263         edge[, 2] <- 1:n
264     } else {
265         nextnode <- n + 1L
266         pos <- 1L
267         foo(1, nextnode)
268     }
269     structure(list(edge = edge, tip.label = labels,
270               Nnode = m), class = "phylo")
271 }