]> git.donarmstrong.com Git - ape.git/blob - R/cophenetic.phylo.R
bug fixed in write.tree
[ape.git] / R / cophenetic.phylo.R
1 ## cophenetic.phylo.R (2007-01-23)
2
3 ##   Pairwise Distances from a Phylogenetic Tree
4
5 ## Copyright 2006-2007 Emmanuel Paradis
6
7 ## This file is part of the R-package `ape'.
8 ## See the file ../COPYING for licensing issues.
9
10 dist.nodes <- function(x)
11 {
12     if (is.null(x$edge.length))
13       stop("your tree has no branch lengths")
14
15     if (!is.binary.tree(x) || !is.rooted(x))
16       x <- multi2di(x, random = FALSE)
17     n <- length(x$tip.label)
18     n.node <- x$Nnode
19     N <- n + n.node
20     x <- reorder(x, order = "pruningwise")
21
22     res <- matrix(NA, N, N)
23     res[cbind(1:N, 1:N)] <- 0 # implicit mode conversion
24
25     ## I like the simplicity of this one:
26     res[x$edge] <- res[x$edge[, 2:1]] <- x$edge.length
27
28     ## compute the distances ...
29     for (i in seq(from = 1, by = 2, length.out = n.node)) {
30         j <- i + 1
31         anc <- x$edge[i, 1]
32         des1 <- x$edge[i, 2]
33         des2 <- x$edge[j, 2]
34
35         ## If `des1' is a node, we look for the nodes and tips for
36         ## which the distance up to `des1' has already been
37         ## computed, including `des1' itself. For all these, we can
38         ## compute the distance up to `anc' and all node(s) and
39         ## tip(s) in `des2'.
40         if (des1 > n) des1 <- which(!is.na(res[des1, ]))
41
42         ## id. for `des2'
43         if (des2 > n) des2 <- which(!is.na(res[des2, ]))
44
45         ## The following expression is vectorized only on `des2' and
46         ## not on `des1' because they may be of different lengths.
47         for (y in des1)
48           res[y, des2] <- res[des2, y] <- res[anc, y] + res[anc, des2]
49         ## compute the distances between the tip(s) and node(s)
50         ## in `des2' and the ancestor of `anc'; id. for `des2'
51         ## (only if it is not the root)
52         if (anc != n + 1) {
53             ind <- which(x$edge[, 2] == anc)
54             nod <- x$edge[ind, 1] # the ancestor of `anc'
55             l <- x$edge.length[ind]
56             res[des2, nod] <- res[nod, des2] <- res[anc, des2] + l
57             res[des1, nod] <- res[nod, des1] <- res[anc, des1] + l
58         }
59     }
60     dimnames(res)[1:2] <- list(1:N)
61     res
62 }
63
64 cophenetic.phylo <- function(x)
65 {
66     n <- length(x$tip.label)
67     ans <- dist.nodes(x)[1:n, 1:n]
68     dimnames(ans)[1:2] <- list(x$tip.label)
69     ans
70 }