]> git.donarmstrong.com Git - ape.git/blob - R/plot.phylo.R
4c23c6f0032361e6c3b6181941cbfd4c83ad3f24
[ape.git] / R / plot.phylo.R
1 ## plot.phylo.R (2012-10-20)
2
3 ##   Plot Phylogenies
4
5 ## Copyright 2002-2012 Emmanuel Paradis
6
7 ## This file is part of the R-package `ape'.
8 ## See the file ../COPYING for licensing issues.
9
10 plot.phylo <-
11     function(x, type = "phylogram", use.edge.length = TRUE,
12              node.pos = NULL, show.tip.label = TRUE,
13              show.node.label = FALSE, edge.color = "black",
14              edge.width = 1, edge.lty = 1, font = 3, cex = par("cex"),
15              adj = NULL, srt = 0, no.margin = FALSE, root.edge = FALSE,
16              label.offset = 0, underscore = FALSE, x.lim = NULL,
17              y.lim = NULL, direction = "rightwards", lab4ut = "horizontal",
18              tip.color = "black", plot = TRUE, rotate.tree = 0, ...)
19 {
20     Ntip <- length(x$tip.label)
21     if (Ntip < 2) {
22         warning("found less than 2 tips in the tree")
23         return(NULL)
24     }
25     if (any(tabulate(x$edge[, 1]) == 1))
26       stop("there are single (non-splitting) nodes in your tree; you may need to use collapse.singles()")
27
28     .nodeHeight <- function(Ntip, Nnode, edge, Nedge, yy)
29         .C("node_height", as.integer(Ntip), as.integer(Nnode),
30            as.integer(edge[, 1]), as.integer(edge[, 2]),
31            as.integer(Nedge), as.double(yy),
32            DUP = FALSE, PACKAGE = "ape")[[6]]
33
34     .nodeDepth <- function(Ntip, Nnode, edge, Nedge)
35         .C("node_depth", as.integer(Ntip), as.integer(Nnode),
36            as.integer(edge[, 1]), as.integer(edge[, 2]),
37            as.integer(Nedge), double(Ntip + Nnode),
38            DUP = FALSE, PACKAGE = "ape")[[6]]
39
40     .nodeDepthEdgelength <- function(Ntip, Nnode, edge, Nedge, edge.length)
41         .C("node_depth_edgelength", as.integer(Ntip),
42            as.integer(Nnode), as.integer(edge[, 1]),
43            as.integer(edge[, 2]), as.integer(Nedge),
44            as.double(edge.length), double(Ntip + Nnode),
45            DUP = FALSE, PACKAGE = "ape")[[7]]
46
47     Nedge <- dim(x$edge)[1]
48     Nnode <- x$Nnode
49     ROOT <- Ntip + 1
50     type <- match.arg(type, c("phylogram", "cladogram", "fan",
51                               "unrooted", "radial"))
52     direction <- match.arg(direction, c("rightwards", "leftwards",
53                                         "upwards", "downwards"))
54     if (is.null(x$edge.length)) use.edge.length <- FALSE
55
56     ## the order of the last two conditions is important:
57     if (type %in% c("unrooted", "radial") || !use.edge.length ||
58         is.null(x$root.edge) || !x$root.edge) root.edge <- FALSE
59     if (type == "fan" && root.edge) {
60         warning("drawing root edge with type = 'fan' is not yet supported")
61         root.edge <- FALSE
62     }
63
64     phyloORclado <- type %in% c("phylogram", "cladogram")
65     horizontal <- direction %in% c("rightwards", "leftwards")
66     xe <- x$edge # to save
67     if (phyloORclado) {
68         ## we first compute the y-coordinates of the tips.
69         phyOrder <- attr(x, "order")
70         ## make sure the tree is in cladewise order:
71         if (is.null(phyOrder) || phyOrder != "cladewise") {
72             x <- reorder(x) # fix from Klaus Schliep (2007-06-16)
73             if (!identical(x$edge, xe)) {
74                 ## modified from Li-San Wang's fix (2007-01-23):
75                 ereorder <- match(x$edge[, 2], xe[, 2])
76                 if (length(edge.color) > 1) {
77                     edge.color <- rep(edge.color, length.out = Nedge)
78                     edge.color <- edge.color[ereorder]
79                 }
80                 if (length(edge.width) > 1) {
81                     edge.width <- rep(edge.width, length.out = Nedge)
82                     edge.width <- edge.width[ereorder]
83                 }
84                 if (length(edge.lty) > 1) {
85                     edge.lty <- rep(edge.lty, length.out = Nedge)
86                     edge.lty <- edge.lty[ereorder]
87                 }
88             }
89         }
90 ### By contrats to ape (< 2.4), the arguments edge.color, etc., are
91 ### not elongated before being passed to segments(), except if needed
92 ### to be reordered
93         yy <- numeric(Ntip + Nnode)
94         TIPS <- x$edge[x$edge[, 2] <= Ntip, 2]
95         yy[TIPS] <- 1:Ntip
96     }
97     ## 'z' is the tree in pruningwise order used in calls to .C
98     z <- reorder(x, order = "pruningwise")
99
100     if (phyloORclado) {
101         if (is.null(node.pos)) {
102             node.pos <- 1
103             if (type == "cladogram" && !use.edge.length) node.pos <- 2
104         }
105         if (node.pos == 1)
106             yy <- .nodeHeight(Ntip, Nnode, z$edge, Nedge, yy)
107         else {
108           ## node_height_clado requires the number of descendants
109           ## for each node, so we compute `xx' at the same time
110           ans <- .C("node_height_clado", as.integer(Ntip),
111                     as.integer(Nnode), as.integer(z$edge[, 1]),
112                     as.integer(z$edge[, 2]), as.integer(Nedge),
113                     double(Ntip + Nnode), as.double(yy),
114                     DUP = FALSE, PACKAGE = "ape")
115           xx <- ans[[6]] - 1
116           yy <- ans[[7]]
117         }
118         if (!use.edge.length) {
119             if (node.pos != 2) xx <- .nodeDepth(Ntip, Nnode, z$edge, Nedge) - 1
120             xx <- max(xx) - xx
121         } else  {
122             xx <- .nodeDepthEdgelength(Ntip, Nnode, z$edge, Nedge, z$edge.length)
123         }
124     } else {
125     rotate.tree <- 2 * pi * rotate.tree/360
126     switch(type, "fan" = {
127         ## if the tips are not in the same order in tip.label
128         ## and in edge[, 2], we must reorder the angles: we
129         ## use `xx' to store temporarily the angles
130         TIPS <- x$edge[which(x$edge[, 2] <= Ntip), 2]
131         xx <- seq(0, 2*pi*(1 - 1/Ntip), 2*pi/Ntip)
132         theta <- double(Ntip)
133         theta[TIPS] <- xx
134         theta <- c(theta, numeric(Nnode))
135         theta <- .nodeHeight(Ntip, Nnode, z$edge, Nedge, theta)
136         if (use.edge.length) {
137             r <- .nodeDepthEdgelength(Ntip, Nnode, z$edge, Nedge, z$edge.length)
138         } else {
139             r <- .nodeDepth(Ntip, Nnode, z$edge, Nedge)
140             r <- 1/r
141         }
142         theta <- theta + rotate.tree
143         xx <- r*cos(theta)
144         yy <- r*sin(theta)
145     }, "unrooted" = {
146         nb.sp <- .nodeDepth(Ntip, Nnode, z$edge, Nedge)
147         XY <- if (use.edge.length)
148             unrooted.xy(Ntip, Nnode, z$edge, z$edge.length, nb.sp, rotate.tree)
149         else
150             unrooted.xy(Ntip, Nnode, z$edge, rep(1, Nedge), nb.sp, rotate.tree)
151         ## rescale so that we have only positive values
152         xx <- XY$M[, 1] - min(XY$M[, 1])
153         yy <- XY$M[, 2] - min(XY$M[, 2])
154     }, "radial" = {
155         X <- .nodeDepth(Ntip, Nnode, z$edge, Nedge)
156         X[X == 1] <- 0
157         ## radius:
158         X <- 1 - X/Ntip
159         ## angle (1st compute the angles for the tips):
160         yy <- c((1:Ntip)*2*pi/Ntip, rep(0, Nnode))
161         Y <- .nodeHeight(Ntip, Nnode, z$edge, Nedge, yy)
162         xx <- X * cos(Y + rotate.tree)
163         yy <- X * sin(Y + rotate.tree)
164     })}
165     if (phyloORclado) {
166         if (!horizontal) {
167             tmp <- yy
168             yy <- xx
169             xx <- tmp - min(tmp) + 1
170         }
171         if (root.edge) {
172             if (direction == "rightwards") xx <- xx + x$root.edge
173             if (direction == "upwards") yy <- yy + x$root.edge
174         }
175     }
176     if (no.margin) par(mai = rep(0, 4))
177     if (is.null(x.lim)) {
178         if (phyloORclado) {
179             if (horizontal) {
180                 x.lim <- c(0, NA)
181                 pin1 <- par("pin")[1] # width of the device in inches
182                 strWi <- strwidth(x$tip.label, "inches") # id. for the tip labels
183                 ## 1.04 comes from that we are using a regular axis system
184                 ## with 4% on both sides of the range of x:
185                 xx.tips <- xx[1:Ntip] * 1.04
186                 ## 'alp' is the conversion coefficient from
187                 ## user coordinates to inches:
188                 alp <- try(uniroot(function(a) max(a*xx.tips + strWi) - pin1,
189                                    c(0, 1e6))$root, silent = TRUE)
190                 ## if the above fails, give 1/3 of the device for the tip labels:
191                 if (is.character(alp)) tmp <- max(xx.tips)*1.5 else {
192                     tmp <- if (show.tip.label) max(xx.tips + strWi/alp) else max(xx.tips)
193                 }
194                 x.lim[2] <- tmp
195             } else x.lim <- c(1, Ntip)
196         } else switch(type, "fan" = {
197             if (show.tip.label) {
198                 offset <- max(nchar(x$tip.label) * 0.018 * max(yy) * cex)
199                 x.lim <- c(min(xx) - offset, max(xx) + offset)
200             } else x.lim <- c(min(xx), max(xx))
201         }, "unrooted" = {
202             if (show.tip.label) {
203                 offset <- max(nchar(x$tip.label) * 0.018 * max(yy) * cex)
204                 x.lim <- c(0 - offset, max(xx) + offset)
205             } else x.lim <- c(0, max(xx))
206         }, "radial" = {
207             if (show.tip.label) {
208                 offset <- max(nchar(x$tip.label) * 0.03 * cex)
209                 x.lim <- c(-1 - offset, 1 + offset)
210             } else x.lim <- c(-1, 1)
211         })
212     } else if (length(x.lim) == 1) {
213         x.lim <- c(0, x.lim)
214         if (phyloORclado && !horizontal) x.lim[1] <- 1
215         if (type %in% c("fan", "unrooted") && show.tip.label)
216           x.lim[1] <- -max(nchar(x$tip.label) * 0.018 * max(yy) * cex)
217         if (type == "radial")
218           x.lim[1] <-
219             if (show.tip.label) -1 - max(nchar(x$tip.label) * 0.03 * cex)
220             else -1
221     }
222     ## mirror the xx:
223     if (phyloORclado && direction == "leftwards") xx <- x.lim[2] - xx
224     if (is.null(y.lim)) {
225         if (phyloORclado) {
226             if (horizontal) y.lim <- c(1, Ntip) else {
227                 y.lim <- c(0, NA)
228                 pin2 <- par("pin")[2] # height of the device in inches
229                 strWi <- strwidth(x$tip.label, "inches")
230                 ## 1.04 comes from that we are using a regular axis system
231                 ## with 4% on both sides of the range of x:
232                 yy.tips <- yy[1:Ntip] * 1.04
233                 ## 'alp' is the conversion coefficient from
234                 ## user coordinates to inches:
235                 alp <- try(uniroot(function(a) max(a*yy.tips + strWi) - pin2,
236                                    c(0, 1e6))$root, silent = TRUE)
237                 ## if the above fails, give 1/3 of the device for the tip labels:
238                 if (is.character(alp)) tmp <- max(yy.tips)*1.5 else {
239                     tmp <- if (show.tip.label) max(yy.tips + strWi/alp) else max(yy.tips)
240                 }
241                 y.lim[2] <- tmp
242             }
243         } else switch(type, "fan" = {
244             if (show.tip.label) {
245                 offset <- max(nchar(x$tip.label) * 0.018 * max(yy) * cex)
246                 y.lim <- c(min(yy) - offset, max(yy) + offset)
247             } else y.lim <- c(min(yy), max(yy))
248         }, "unrooted" = {
249             if (show.tip.label) {
250                 offset <- max(nchar(x$tip.label) * 0.018 * max(yy) * cex)
251                 y.lim <- c(0 - offset, max(yy) + offset)
252             } else y.lim <- c(0, max(yy))
253         }, "radial" = {
254             if (show.tip.label) {
255                 offset <- max(nchar(x$tip.label) * 0.03 * cex)
256                 y.lim <- c(-1 - offset, 1 + offset)
257             } else y.lim <- c(-1, 1)
258         })
259     } else if (length(y.lim) == 1) {
260         y.lim <- c(0, y.lim)
261         if (phyloORclado && horizontal) y.lim[1] <- 1
262         if (type %in% c("fan", "unrooted") && show.tip.label)
263             y.lim[1] <- -max(nchar(x$tip.label) * 0.018 * max(yy) * cex)
264         if (type == "radial")
265             y.lim[1] <- if (show.tip.label) -1 - max(nchar(x$tip.label) * 0.018 * max(yy) * cex) else -1
266     }
267     ## mirror the yy:
268     if (phyloORclado && direction == "downwards") yy <- max(yy) - yy
269     if (phyloORclado && root.edge) {
270         if (direction == "leftwards") x.lim[2] <- x.lim[2] + x$root.edge
271         if (direction == "downwards") y.lim[2] <- y.lim[2] + x$root.edge
272     }
273     asp <- if (type %in% c("fan", "radial", "unrooted")) 1 else NA # fixes by Klaus Schliep (2008-03-28 and 2010-08-12)
274     plot(0, type = "n", xlim = x.lim, ylim = y.lim, ann = FALSE, axes = FALSE, asp = asp, ...)
275
276 if (plot) {
277     if (is.null(adj))
278         adj <- if (phyloORclado && direction == "leftwards") 1 else 0
279     if (phyloORclado && show.tip.label) {
280         MAXSTRING <- max(strwidth(x$tip.label, cex = cex))
281         loy <- 0
282         if (direction == "rightwards") {
283             lox <- label.offset + MAXSTRING * 1.05 * adj
284         }
285         if (direction == "leftwards") {
286             lox <- -label.offset - MAXSTRING * 1.05 * (1 - adj)
287             ##xx <- xx + MAXSTRING
288         }
289         if (!horizontal) {
290             psr <- par("usr")
291             MAXSTRING <- MAXSTRING * 1.09 * (psr[4] - psr[3])/(psr[2] - psr[1])
292             loy <- label.offset + MAXSTRING * 1.05 * adj
293             lox <- 0
294             srt <- 90 + srt
295             if (direction == "downwards") {
296                 loy <- -loy
297                 ##yy <- yy + MAXSTRING
298                 srt <- 180 + srt
299             }
300         }
301     }
302     if (type == "phylogram") {
303         phylogram.plot(x$edge, Ntip, Nnode, xx, yy,
304                        horizontal, edge.color, edge.width, edge.lty)
305     } else {
306         if (type == "fan") {
307             ereorder <- match(z$edge[, 2], x$edge[, 2])
308             if (length(edge.color) > 1) {
309                 edge.color <- rep(edge.color, length.out = Nedge)
310                 edge.color <- edge.color[ereorder]
311             }
312             if (length(edge.width) > 1) {
313                 edge.width <- rep(edge.width, length.out = Nedge)
314                 edge.width <- edge.width[ereorder]
315             }
316             if (length(edge.lty) > 1) {
317                 edge.lty <- rep(edge.lty, length.out = Nedge)
318                 edge.lty <- edge.lty[ereorder]
319             }
320             circular.plot(z$edge, Ntip, Nnode, xx, yy, theta,
321                           r, edge.color, edge.width, edge.lty)
322         } else
323         cladogram.plot(x$edge, xx, yy, edge.color, edge.width, edge.lty)
324     }
325     if (root.edge)
326       switch(direction,
327              "rightwards" = segments(0, yy[ROOT], x$root.edge, yy[ROOT]),
328              "leftwards" = segments(xx[ROOT], yy[ROOT], xx[ROOT] + x$root.edge, yy[ROOT]),
329              "upwards" = segments(xx[ROOT], 0, xx[ROOT], x$root.edge),
330              "downwards" = segments(xx[ROOT], yy[ROOT], xx[ROOT], yy[ROOT] + x$root.edge))
331     if (show.tip.label) {
332         if (is.expression(x$tip.label)) underscore <- TRUE
333         if (!underscore) x$tip.label <- gsub("_", " ", x$tip.label)
334
335         if (phyloORclado)
336             text(xx[1:Ntip] + lox, yy[1:Ntip] + loy, x$tip.label, adj = adj,
337                  font = font, srt = srt, cex = cex, col = tip.color)
338
339         if (type == "unrooted") {
340             if (lab4ut == "horizontal") {
341                 y.adj <- x.adj <- numeric(Ntip)
342                 sel <- abs(XY$axe) > 0.75 * pi
343                 x.adj[sel] <- -strwidth(x$tip.label)[sel] * 1.05
344                 sel <- abs(XY$axe) > pi/4 & abs(XY$axe) < 0.75 * pi
345                 x.adj[sel] <- -strwidth(x$tip.label)[sel] * (2 * abs(XY$axe)[sel] / pi - 0.5)
346                 sel <- XY$axe > pi / 4 & XY$axe < 0.75 * pi
347                 y.adj[sel] <- strheight(x$tip.label)[sel] / 2
348                 sel <- XY$axe < -pi / 4 & XY$axe > -0.75 * pi
349                 y.adj[sel] <- -strheight(x$tip.label)[sel] * 0.75
350                 text(xx[1:Ntip] + x.adj * cex, yy[1:Ntip] + y.adj * cex,
351                      x$tip.label, adj = c(adj, 0), font = font,
352                      srt = srt, cex = cex, col = tip.color)
353             } else { # if lab4ut == "axial"
354                 adj <- abs(XY$axe) > pi/2
355                 srt <- 180 * XY$axe / pi
356                 srt[adj] <- srt[adj] - 180
357                 adj <- as.numeric(adj)
358                 xx.tips <- xx[1:Ntip]
359                 yy.tips <- yy[1:Ntip]
360                 if (label.offset) {
361                     xx.tips <- xx.tips + label.offset * cos(XY$axe)
362                     yy.tips <- yy.tips + label.offset * sin(XY$axe)
363                 }
364                 ## `srt' takes only a single value, so can't vectorize this:
365                 ## (and need to 'elongate' these vectors:)
366                 font <- rep(font, length.out = Ntip)
367                 tip.color <- rep(tip.color, length.out = Ntip)
368                 cex <- rep(cex, length.out = Ntip)
369                 for (i in 1:Ntip)
370                     text(xx.tips[i], yy.tips[i], cex = cex[i],
371                          x$tip.label[i], adj = adj[i], font = font[i],
372                          srt = srt[i], col = tip.color[i])
373             }
374         }
375         if (type %in% c("fan", "radial")) {
376             xx.tips <- xx[1:Ntip]
377             yy.tips <- yy[1:Ntip]
378             ## using atan2 considerably facilitates things compared to acos...
379             angle <- atan2(yy.tips, xx.tips) # in radians
380             if (label.offset) {
381                 xx.tips <- xx.tips + label.offset * cos(angle)
382                 yy.tips <- yy.tips + label.offset * sin(angle)
383             }
384             s <- xx.tips < 0
385             angle <- angle * 180/pi # switch to degrees
386             angle[s] <- angle[s] + 180
387             adj <- as.numeric(s)
388             ## `srt' takes only a single value, so can't vectorize this:
389             ## (and need to 'elongate' these vectors:)
390             font <- rep(font, length.out = Ntip)
391             tip.color <- rep(tip.color, length.out = Ntip)
392             cex <- rep(cex, length.out = Ntip)
393             for (i in 1:Ntip)
394                 text(xx.tips[i], yy.tips[i], x$tip.label[i], font = font[i],
395                      cex = cex[i], srt = angle[i], adj = adj[i],
396                      col = tip.color[i])
397         }
398     }
399     if (show.node.label)
400         text(xx[ROOT:length(xx)] + label.offset, yy[ROOT:length(yy)],
401              x$node.label, adj = adj, font = font, srt = srt, cex = cex)
402 }
403     L <- list(type = type, use.edge.length = use.edge.length,
404               node.pos = node.pos, show.tip.label = show.tip.label,
405               show.node.label = show.node.label, font = font,
406               cex = cex, adj = adj, srt = srt, no.margin = no.margin,
407               label.offset = label.offset, x.lim = x.lim, y.lim = y.lim,
408               direction = direction, tip.color = tip.color,
409               Ntip = Ntip, Nnode = Nnode)
410     assign("last_plot.phylo", c(L, list(edge = xe, xx = xx, yy = yy)),
411            envir = .PlotPhyloEnv)
412     invisible(L)
413 }
414
415 phylogram.plot <- function(edge, Ntip, Nnode, xx, yy, horizontal,
416                            edge.color, edge.width, edge.lty)
417 {
418     nodes <- (Ntip + 1):(Ntip + Nnode)
419     if (!horizontal) {
420         tmp <- yy
421         yy <- xx
422         xx <- tmp
423     }
424     ## un trait vertical a chaque noeud...
425     x0v <- xx[nodes]
426     y0v <- y1v <- numeric(Nnode)
427     ## store the index of each node in the 1st column of edge:
428     NodeInEdge1 <- vector("list", Nnode)
429     for (i in nodes) {
430         ii <- i - Ntip
431         j <- NodeInEdge1[[ii]] <- which(edge[, 1] == i)
432         tmp <- range(yy[edge[j, 2]])
433         y0v[ii] <- tmp[1]
434         y1v[ii] <- tmp[2]
435     }
436     ## ... et un trait horizontal partant de chaque tip et chaque noeud
437     ##  vers la racine
438     x0h <- xx[edge[, 1]]
439     x1h <- xx[edge[, 2]]
440     y0h <- yy[edge[, 2]]
441
442     nc <- length(edge.color)
443     nw <- length(edge.width)
444     nl <- length(edge.lty)
445
446     if (nc + nw + nl == 3) {
447         color.v <- edge.color
448         width.v <- edge.width
449         lty.v <- edge.lty
450     } else {
451         Nedge <- dim(edge)[1]
452         edge.color <- rep(edge.color, length.out = Nedge)
453         edge.width <- rep(edge.width, length.out = Nedge)
454         edge.lty <- rep(edge.lty, length.out = Nedge)
455         DF <- data.frame(edge.color, edge.width, edge.lty, stringsAsFactors = FALSE)
456         color.v <- rep("black", Nnode)
457         width.v <- rep(1, Nnode)
458         lty.v <- rep(1, Nnode)
459         for (i in 1:Nnode) {
460             br <- NodeInEdge1[[i]]
461             if (length(br) > 2) {
462                 x <- unique(DF[br, 1])
463                 if (length(x) == 1) color.v[i] <- x
464                 x <- unique(DF[br, 2])
465                 if (length(x) == 1) width.v[i] <- x
466                 x <- unique(DF[br, 3])
467                 if (length(x) == 1) lty.v[i] <- x
468             } else {
469                 A <- br[1]
470                 B <- br[2]
471                 if (any(DF[A, ] != DF[B, ])) {
472                     color.v[i] <- edge.color[B]
473                     width.v[i] <- edge.width[B]
474                     lty.v[i] <- edge.lty[B]
475                     ## add a new line:
476                     y0v <- c(y0v, y0v[i])
477                     y1v <- c(y1v, yy[i + Ntip])
478                     x0v <- c(x0v, x0v[i])
479                     color.v <- c(color.v, edge.color[A])
480                     width.v <- c(width.v, edge.width[A])
481                     lty.v <- c(lty.v, edge.lty[A])
482                     ## shorten the line:
483                     y0v[i] <- yy[i + Ntip]
484                 } else {
485                     color.v[i] <- edge.color[A]
486                     width.v[i] <- edge.width[A]
487                     lty.v[i] <- edge.lty[A]
488                 }
489             }
490         }
491     }
492
493     if (horizontal) {
494         segments(x0h, y0h, x1h, y0h, col = edge.color, lwd = edge.width, lty = edge.lty) # draws horizontal lines
495         segments(x0v, y0v, x0v, y1v, col = color.v, lwd = width.v, lty = lty.v) # draws vertical lines
496     } else {
497         segments(y0h, x0h, y0h, x1h, col = edge.color, lwd = edge.width, lty = edge.lty) # draws vertical lines
498         segments(y0v, x0v, y1v, x0v, col = color.v, lwd = width.v, lty = lty.v) # draws horizontal lines
499     }
500 }
501
502 cladogram.plot <- function(edge, xx, yy, edge.color, edge.width, edge.lty)
503     segments(xx[edge[, 1]], yy[edge[, 1]], xx[edge[, 2]], yy[edge[, 2]],
504              col = edge.color, lwd = edge.width, lty = edge.lty)
505
506 circular.plot <- function(edge, Ntip, Nnode, xx, yy, theta,
507                           r, edge.color, edge.width, edge.lty)
508 ### 'edge' must be in pruningwise order
509 {
510     r0 <- r[edge[, 1]]
511     r1 <- r[edge[, 2]]
512     theta0 <- theta[edge[, 2]]
513     costheta0 <- cos(theta0)
514     sintheta0 <- sin(theta0)
515
516     x0 <- r0 * costheta0
517     y0 <- r0 * sintheta0
518     x1 <- r1 * costheta0
519     y1 <- r1 * sintheta0
520
521     segments(x0, y0, x1, y1, col = edge.color, lwd = edge.width, lty = edge.lty)
522
523     tmp <- which(diff(edge[, 1]) != 0)
524     start <- c(1, tmp + 1)
525     Nedge <- dim(edge)[1]
526     end <- c(tmp, Nedge)
527
528     ## function dispatching the features to the arcs
529     foo <- function(edge.feat, default) {
530         if (length(edge.feat) == 1) return(rep(edge.feat, Nnode))
531         else {
532             edge.feat <- rep(edge.feat, length.out = Nedge)
533             feat.arc <- rep(default, Nnode)
534             for (k in 1:Nnode) {
535                 tmp <- edge.feat[start[k]]
536                 if (tmp == edge.feat[end[k]]) feat.arc[k] <- tmp
537             }
538         }
539         feat.arc
540     }
541     co <- foo(edge.color, "black")
542     lw <- foo(edge.width, 1)
543     ly <- foo(edge.lty, 1)
544
545     for (k in 1:Nnode) {
546         i <- start[k]
547         j <- end[k]
548         X <- rep(r[edge[i, 1]], 100)
549         Y <- seq(theta[edge[i, 2]], theta[edge[j, 2]], length.out = 100)
550         lines(X*cos(Y), X*sin(Y), col = co[k], lwd = lw[k], lty = ly[k])
551     }
552 }
553
554 unrooted.xy <- function(Ntip, Nnode, edge, edge.length, nb.sp, rotate.tree)
555 {
556     foo <- function(node, ANGLE, AXIS) {
557         ind <- which(edge[, 1] == node)
558         sons <- edge[ind, 2]
559         start <- AXIS - ANGLE/2
560         for (i in 1:length(sons)) {
561             h <- edge.length[ind[i]]
562             angle[sons[i]] <<- alpha <- ANGLE*nb.sp[sons[i]]/nb.sp[node]
563             axis[sons[i]] <<- beta <- start + alpha/2
564             start <- start + alpha
565             xx[sons[i]] <<- h*cos(beta) + xx[node]
566             yy[sons[i]] <<- h*sin(beta) + yy[node]
567         }
568         for (i in sons)
569             if (i > Ntip) foo(i, angle[i], axis[i])
570     }
571     Nedge <- dim(edge)[1]
572     yy <- xx <- numeric(Ntip + Nnode)
573     ## `angle': the angle allocated to each node wrt their nb of tips
574     ## `axis': the axis of each branch
575     axis <- angle <- numeric(Ntip + Nnode)
576     ## start with the root...
577     foo(Ntip + 1L, 2*pi, 0 + rotate.tree)
578
579     M <- cbind(xx, yy)
580     axe <- axis[1:Ntip] # the axis of the terminal branches (for export)
581     axeGTpi <- axe > pi
582     ## insures that returned angles are in [-PI, +PI]:
583     axe[axeGTpi] <- axe[axeGTpi] - 2*pi
584     list(M = M, axe = axe)
585 }
586
587 node.depth <- function(phy)
588 {
589     n <- length(phy$tip.label)
590     m <- phy$Nnode
591     N <- dim(phy$edge)[1]
592     phy <- reorder(phy, order = "pruningwise")
593     .C("node_depth", as.integer(n), as.integer(m),
594        as.integer(phy$edge[, 1]), as.integer(phy$edge[, 2]),
595        as.integer(N), double(n + m), DUP = FALSE, PACKAGE = "ape")[[6]]
596 }
597
598 node.depth.edgelength <- function(phy)
599 {
600     n <- length(phy$tip.label)
601     m <- phy$Nnode
602     N <- dim(phy$edge)[1]
603     phy <- reorder(phy, order = "pruningwise")
604     .C("node_depth_edgelength", as.integer(n), as.integer(n),
605        as.integer(phy$edge[, 1]), as.integer(phy$edge[, 2]),
606        as.integer(N), as.double(phy$edge.length), double(n + m),
607        DUP = FALSE, PACKAGE = "ape")[[7]]
608 }
609
610 node.height <- function(phy)
611 {
612     n <- length(phy$tip.label)
613     m <- phy$Nnode
614     N <- dim(phy$edge)[1]
615     phy <- reorder(phy, order = "pruningwise")
616
617     e1 <- phy$edge[, 1]
618     e2 <- phy$edge[, 2]
619
620     yy <- numeric(n + m)
621     TIPS <- e2[e2 <= n]
622     yy[TIPS] <- 1:n
623
624     .C("node_height", as.integer(n), as.integer(m),
625        as.integer(e1), as.integer(e2), as.integer(N),
626        as.double(yy), DUP = FALSE, PACKAGE = "ape")[[6]]
627 }
628
629 node.height.clado <- function(phy)
630 {
631     n <- length(phy$tip.label)
632     m <- phy$Nnode
633     N <- dim(phy$edge)[1]
634     phy <- reorder(phy, order = "pruningwise")
635
636     e1 <- phy$edge[, 1]
637     e2 <- phy$edge[, 2]
638
639     yy <- numeric(n + m)
640     TIPS <- e2[e2 <= n]
641     yy[TIPS] <- 1:n
642
643     .C("node_height_clado", as.integer(n), as.integer(m),
644        as.integer(e1), as.integer(e2), as.integer(N),
645        double(n + m), as.double(yy), DUP = FALSE,
646        PACKAGE = "ape")[[7]]
647 }
648
649 plot.multiPhylo <- function(x, layout = 1, ...)
650 {
651     if (layout > 1)
652       layout(matrix(1:layout, ceiling(sqrt(layout)), byrow = TRUE))
653     else layout(matrix(1))
654     if (!par("ask")) {
655         par(ask = TRUE)
656         on.exit(par(ask = FALSE))
657     }
658     for (i in 1:length(x)) plot(x[[i]], ...)
659 }
660
661 trex <- function(phy, title = TRUE, subbg = "lightyellow3",
662                  return.tree = FALSE, ...)
663 {
664     lastPP <- get("last_plot.phylo", envir = .PlotPhyloEnv)
665     devmain <- dev.cur() # where the main tree is plotted
666
667     restore <- function() {
668         dev.set(devmain)
669         assign("last_plot.phylo", lastPP, envir = .PlotPhyloEnv)
670     }
671
672     on.exit(restore())
673     NEW <- TRUE
674     cat("Click close to a node. Right-click to exit.\n")
675     repeat {
676         x <- identify.phylo(phy, quiet = TRUE)
677         if (is.null(x)) return(invisible(NULL)) else {
678             x <- x$nodes
679             if (is.null(x)) cat("Try again!\n") else {
680                 if (NEW) {
681                     dev.new()
682                     par(bg = subbg)
683                     devsub <- dev.cur()
684                     NEW <- FALSE
685                 } else dev.set(devsub)
686
687                 tr <- extract.clade(phy, x)
688                 plot(tr, ...)
689                 if (is.character(title)) title(title)
690                 else if (title) {
691                      tl <-
692                          if (is.null(phy$node.label))
693                          paste("From node #", x, sep = "")
694                          else paste("From", phy$node.label[x - Ntip(phy)])
695                      title(tl)
696                 }
697                 if (return.tree) return(tr)
698                 restore()
699             }
700         }
701     }
702 }
703
704 kronoviz <- function(x, layout = length(x), horiz = TRUE, ...)
705 {
706     par(mar = rep(0.5, 4), oma = rep(2, 4))
707     rts <- sapply(x, function(x) branching.times(x)[1])
708     maxrts <- max(rts)
709     lim <- cbind(rts - maxrts, rts)
710     Ntree <- length(x)
711     Ntips <- sapply(x, Ntip)
712     if (horiz) {
713         nrow <- layout
714         w <- 1
715         h <- Ntips
716     } else {
717         nrow <- 1
718         w <- Ntips
719         h <- 1
720     }
721     layout(matrix(1:layout, nrow), widths = w, heights = h)
722     if (layout > Ntree && !par("ask")) {
723         par(ask = TRUE)
724         on.exit(par(ask = FALSE))
725     }
726     if (horiz) {
727         for (i in 1:Ntree)
728             plot(x[[i]], x.lim = lim[i, ], ...)
729     } else {
730         for (i in 1:Ntree)
731             plot(x[[i]], y.lim = lim[i, ], direction = "u", ...)
732     }
733     axisPhylo(if (horiz) 1 else 4) # better if the deepest tree is last ;)
734 }