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