]> git.donarmstrong.com Git - ape.git/blob - R/drop.tip.R
fixing drop.tip and bionj
[ape.git] / R / drop.tip.R
1 ## drop.tip.R (2011-05-19)
2
3 ##   Remove Tips in a Phylogenetic Tree
4
5 ## Copyright 2003-2011 Emmanuel Paradis
6
7 ## This file is part of the R-package `ape'.
8 ## See the file ../COPYING for licensing issues.
9
10 extract.clade <- function(phy, node, root.edge = 0, interactive = FALSE)
11 {
12     Ntip <- length(phy$tip.label)
13     ROOT <- Ntip + 1
14     Nedge <- dim(phy$edge)[1]
15     wbl <- !is.null(phy$edge.length)
16     if (interactive) node <- identify(phy)$nodes else {
17         if (length(node) > 1) {
18             node <- node[1]
19             warning("only the first value of 'node' has been considered")
20         }
21         if (is.character(node)) {
22             if (is.null(phy$node.label))
23                 stop("the tree has no node labels")
24             node <- which(phy$node.label %in% node) + Ntip
25         }
26         if (node <= Ntip)
27             stop("node number must be greater than the number of tips")
28     }
29     if (node == ROOT) return(phy)
30     phy <- reorder(phy) # insure it is in cladewise order
31     root.node <- which(phy$edge[, 2] == node)
32     start <- root.node + 1 # start of the clade looked for
33     anc <- phy$edge[root.node, 1] # the ancestor of 'node'
34     next.anc <- which(phy$edge[-(1:start), 1] <= anc) # find the next occurence of 'anc' or an 'older' node
35
36     keep <- if (length(next.anc)) start + 0:(next.anc[1] - 1) else start:Nedge
37
38     if (root.edge) {
39         NewRootEdge <- phy$edge.length[root.node]
40         root.edge <- root.edge - 1
41         while (root.edge) {
42             if (anc == ROOT) break
43             i <- which(phy$edge[, 2] ==  anc)
44             NewRootEdge <- NewRootEdge + phy$edge.length[i]
45             root.edge <- root.edge - 1
46             anc <- phy$edge[i, 1]
47         }
48         if (root.edge && !is.null(phy$root.edge))
49             NewRootEdge <- NewRootEdge + phy$root.edge
50         phy$root.edge <- NewRootEdge
51     }
52
53     phy$edge <- phy$edge[keep, ]
54     if (wbl) phy$edge.length <- phy$edge.length[keep]
55     TIPS <- phy$edge[, 2] <= Ntip
56     tip <- phy$edge[TIPS, 2]
57     phy$tip.label <- phy$tip.label[sort(tip)] # <- added sort to avoid shuffling of tip labels (2010-07-21)
58     ## keep the ordering so no need to reorder tip.label:
59     phy$edge[TIPS, 2] <- order(tip)
60     if (!is.null(phy$node.label))
61         phy$node.label <- phy$node.label[sort(unique(phy$edge[, 1])) - Ntip]
62     Ntip <- length(phy$tip.label)
63     phy$Nnode <- dim(phy$edge)[1] - Ntip + 1L
64     ## The block below renumbers the nodes so that they conform
65     ## to the "phylo" format -- same as in root()
66     newNb <- integer(Ntip + phy$Nnode)
67     newNb[node] <- Ntip + 1L
68     sndcol <- phy$edge[, 2] > Ntip
69     ## executed from right to left, so newNb is modified before phy$edge:
70     phy$edge[sndcol, 2] <- newNb[phy$edge[sndcol, 2]] <-
71         (Ntip + 2):(Ntip + phy$Nnode)
72     phy$edge[, 1] <- newNb[phy$edge[, 1]]
73     phy
74 }
75
76 drop.tip <-
77     function(phy, tip, trim.internal = TRUE, subtree = FALSE,
78              root.edge = 0, rooted = is.rooted(phy), interactive = FALSE)
79 {
80     if (!inherits(phy, "phylo"))
81         stop('object "phy" is not of class "phylo"')
82     if (!length(tip)) return(phy)
83
84     Ntip <- length(phy$tip.label)
85     ## find the tips to drop:
86     if (interactive) {
87         cat("Left-click close to the tips you want to drop; right-click when finished...\n")
88         xy <- locator()
89         nToDrop <- length(xy$x)
90         tip <- integer(nToDrop)
91         lastPP <- get("last_plot.phylo", envir = .PlotPhyloEnv)
92         for (i in 1:nToDrop) {
93             d <- sqrt((xy$x[i] - lastPP$xx)^2 + (xy$y[i] - lastPP$yy)^2)
94             tip[i] <- which.min(d)
95         }
96     } else {
97         if (is.character(tip))
98             tip <- which(phy$tip.label %in% tip)
99     }
100     if (any(tip > Ntip))
101         warning("some tip numbers were higher than the number of tips")
102
103     if (!rooted && subtree) {
104         phy <- root(phy, (1:Ntip)[-tip][1])
105         root.edge <- 0
106     }
107
108     phy <- reorder(phy)
109     NEWROOT <- ROOT <- Ntip + 1
110     Nnode <- phy$Nnode
111     Nedge <- dim(phy$edge)[1]
112     if (subtree) {
113         trim.internal <- TRUE
114         tr <- reorder(phy, "pruningwise")
115         N <- .C("node_depth", as.integer(Ntip), as.integer(Nnode),
116                 as.integer(tr$edge[, 1]), as.integer(tr$edge[, 2]),
117                 as.integer(Nedge), double(Ntip + Nnode),
118                 DUP = FALSE, PACKAGE = "ape")[[6]]
119     }
120     wbl <- !is.null(phy$edge.length)
121     edge1 <- phy$edge[, 1] # local copies
122     edge2 <- phy$edge[, 2] #
123     keep <- !logical(Nedge)
124
125     ## delete the terminal edges given by `tip':
126     keep[match(tip, edge2)] <- FALSE
127
128     if (trim.internal) {
129         ints <- edge2 > Ntip
130         ## delete the internal edges that do not have anymore
131         ## descendants (ie, they are in the 2nd col of `edge' but
132         ## not in the 1st one)
133         repeat {
134             sel <- !(edge2 %in% edge1[keep]) & ints & keep
135             if (!sum(sel)) break
136             keep[sel] <- FALSE
137         }
138         if (subtree) {
139             ## keep the subtending edge(s):
140             subt <- edge1 %in% edge1[keep] & edge1 %in% edge1[!keep]
141             keep[subt] <- TRUE
142         }
143         if (root.edge && wbl) {
144             degree <- tabulate(edge1[keep])
145             if (degree[ROOT] == 1) {
146                 j <- integer(0) # will store the indices of the edges below the new root
147                 repeat {
148                     i <- which(edge1 == NEWROOT & keep)
149                     j <- c(i, j)
150                     NEWROOT <- edge2[i]
151                     degree <- tabulate(edge1[keep])
152                     if (degree[NEWROOT] > 1) break
153                 }
154                 keep[j] <- FALSE
155                 if (length(j) > root.edge) j <- 1:root.edge
156                 NewRootEdge <- sum(phy$edge.length[j])
157                 if (length(j) < root.edge && !is.null(phy$root.edge))
158                     NewRootEdge <- NewRootEdge + phy$root.edge
159                 phy$root.edge <- NewRootEdge
160             }
161         }
162     }
163
164     if (!root.edge) phy$root.edge <- NULL
165
166     ## drop the edges
167     phy$edge <- phy$edge[keep, ]
168     if (wbl) phy$edge.length <- phy$edge.length[keep]
169
170     ## find the new terminal edges (works whatever 'subtree' and 'trim.internal'):
171     TERMS <- !(phy$edge[, 2] %in% phy$edge[, 1])
172
173     ## get the old No. of the nodes and tips that become tips:
174     oldNo.ofNewTips <- phy$edge[TERMS, 2]
175
176     ## in case some tips are dropped but kept because of 'subtree = TRUE':
177     if (subtree) {
178         i <- which(tip %in% oldNo.ofNewTips)
179         if (length(i)) {
180             phy$tip.label[tip[i]] <- "[1_tip]"
181             tip <- tip[-i]
182         }
183     }
184
185     n <- length(oldNo.ofNewTips) # the new number of tips in the tree
186
187     ## the tips may not be sorted in increasing order in the
188     ## 2nd col of edge, so no need to reorder $tip.label
189     phy$edge[TERMS, 2] <- rank(phy$edge[TERMS, 2])
190     phy$tip.label <- phy$tip.label[-tip]
191
192     ## make new tip labels if necessary:
193     if (subtree || !trim.internal) {
194         ## get the numbers of the nodes that become tips:
195         node2tip <- oldNo.ofNewTips[oldNo.ofNewTips > Ntip]
196         new.tip.label <- if (subtree) {
197             paste("[", N[node2tip], "_tips]", sep = "")
198         } else {
199             if (is.null(phy$node.label)) rep("NA", length(node2tip))
200             else phy$node.label[node2tip - Ntip]
201         }
202         if (!is.null(phy$node.label))
203             phy$node.label <- phy$node.label[-(node2tip - Ntip)]
204         phy$tip.label <- c(phy$tip.label, new.tip.label)
205     }
206
207     ## update node.label if needed:
208     if (!is.null(phy$node.label))
209         phy$node.label <- phy$node.label[sort(unique(phy$edge[, 1])) - Ntip]
210
211     phy$Nnode <- dim(phy$edge)[1] - n + 1L # update phy$Nnode
212
213     ## The block below renumbers the nodes so that they conform
214     ## to the "phylo" format -- same as in root()
215     newNb <- integer(n + phy$Nnode)
216     newNb[NEWROOT] <- n + 1L
217     sndcol <- phy$edge[, 2] > n
218     ## executed from right to left, so newNb is modified before phy$edge:
219     phy$edge[sndcol, 2] <- newNb[phy$edge[sndcol, 2]] <-
220         (n + 2):(n + phy$Nnode)
221     phy$edge[, 1] <- newNb[phy$edge[, 1]]
222     storage.mode(phy$edge) <- "integer"
223     collapse.singles(phy)
224 }
225