]> git.donarmstrong.com Git - ape.git/blob - R/nodelabels.R
fix rTraitDisc
[ape.git] / R / nodelabels.R
1 ## nodelabels.R (2010-01-30)
2
3 ##   Labelling Trees
4
5 ## Copyright 2004-2010 Emmanuel Paradis, 2006 Ben Bolker, and 2006 Jim Lemon
6
7 ## This file is part of the R-package `ape'.
8 ## See the file ../COPYING for licensing issues.
9
10 ## from JL:
11 ## floating.pie() from plotrix with two changes:
12 ## (1) aspect ratio fixed, so pies will appear circular
13 ##     (`radius' is the radius in user coordinates along the x axis);
14 ## (2) zero values allowed (but not negative).
15
16 floating.pie.asp <- function(xpos, ypos, x, edges = 200, radius = 1,
17                              col = NULL, startpos = 0, ...)
18 {
19     u <- par("usr")
20     user.asp <- diff(u[3:4])/diff(u[1:2])
21     p <- par("pin")
22     inches.asp <- p[2]/p[1]
23     asp <- user.asp/inches.asp
24     if (!is.numeric(x) || any(is.na(x) | x < 0)) {
25       ## browser()
26       stop("floating.pie: x values must be non-negative")
27     }
28     x <- c(0, cumsum(x)/sum(x))
29     dx <- diff(x)
30     nx <- length(dx)
31     if (is.null(col)) col <- rainbow(nx)
32     else if (length(col) < nx) col <- rep(col, nx)
33     bc <- 2 * pi * (x[1:nx] + dx/2) + startpos
34     for (i in 1:nx) {
35         n <- max(2, floor(edges * dx[i]))
36         t2p <- 2 * pi * seq(x[i], x[i + 1], length = n) + startpos
37         xc <- c(cos(t2p) * radius + xpos, xpos)
38         yc <- c(sin(t2p) * radius*asp + ypos, ypos)
39         polygon(xc, yc, col = col[i], ...)
40         ## t2p <- 2 * pi * mean(x[i + 0:1]) + startpos
41         ## xc <- cos(t2p) * radius
42         ## yc <- sin(t2p) * radius*asp
43         ## lines(c(1, 1.05) * xc, c(1, 1.05) * yc)
44     }
45     ## return(bc)
46 }
47
48 BOTHlabels <- function(text, sel, XX, YY, adj, frame, pch, thermo,
49                        pie, piecol, col, bg, ...)
50 {
51     if (missing(text)) text <- NULL
52     if (length(adj) == 1) adj <- c(adj, 0.5)
53     if (is.null(text) && is.null(pch) && is.null(thermo) && is.null(pie))
54       text <- as.character(sel)
55     frame <- match.arg(frame, c("rect", "circle", "none"))
56     args <- list(...)
57     CEX <- if ("cex" %in% names(args)) args$cex else par("cex")
58     if (frame != "none" && !is.null(text)) {
59         if (frame == "rect") {
60             width <- strwidth(text, units = "inches", cex = CEX)
61             height <- strheight(text, units = "inches", cex = CEX)
62             if ("srt" %in% names(args)) {
63                 args$srt <- args$srt %% 360 # just in case srt >= 360
64                 if (args$srt == 90 || args$srt == 270) {
65                     tmp <- width
66                     width <- height
67                     height <- tmp
68                 } else if (args$srt != 0)
69                   warning("only right angle rotation of frame is supported;\n         try  `frame = \"n\"' instead.\n")
70             }
71             width <- xinch(width)
72             height <- yinch(height)
73             xl <- XX - width*adj[1] - xinch(0.03)
74             xr <- xl + width + xinch(0.03)
75             yb <- YY - height*adj[2] - yinch(0.02)
76             yt <- yb + height + yinch(0.05)
77             rect(xl, yb, xr, yt, col = bg)
78         }
79         if (frame == "circle") {
80             radii <- 0.8*apply(cbind(strheight(text, units = "inches", cex = CEX),
81                                      strwidth(text, units = "inches", cex = CEX)), 1, max)
82             symbols(XX, YY, circles = radii, inches = max(radii), add = TRUE, bg = bg)
83         }
84     }
85     if (!is.null(thermo)) {
86         parusr <- par("usr")
87         width <- CEX * (parusr[2] - parusr[1]) / 40
88         height <- CEX * (parusr[4] - parusr[3]) / 15
89         if (is.vector(thermo)) thermo <- cbind(thermo, 1 - thermo)
90         thermo <- height * thermo
91         xl <- XX - width/2 + adj[1] - 0.5 # added 'adj' from Janet Young (2009-09-30)
92         xr <- xl + width
93         yb <- YY - height/2 + adj[2] - 0.5
94         if (is.null(piecol)) piecol <- rainbow(ncol(thermo))
95         ## draw the first rectangle:
96         rect(xl, yb, xr, yb + thermo[, 1], border = NA, col = piecol[1])
97         for (i in 2:ncol(thermo))
98             rect(xl, yb + rowSums(thermo[, 1:(i - 1), drop = FALSE]),
99                  xr, yb + rowSums(thermo[, 1:i]),
100                  border = NA, col = piecol[i])
101         rect(xl, yb, xr, yb + height, border = "black")
102         segments(xl, YY, xl - width/5, YY)
103         segments(xr, YY, xr + width/5, YY)
104     }
105     ## from BB:
106     if (!is.null(pie)) {
107         if (is.vector(pie)) pie <- cbind(pie, 1 - pie)
108         xrad <- CEX * diff(par("usr")[1:2]) / 50
109         xrad <- rep(xrad, length(sel))
110         XX <- XX + adj[1] - 0.5
111         YY <- YY + adj[2] - 0.5
112         for (i in 1:length(sel))
113             floating.pie.asp(XX[i], YY[i], pie[i, ], radius = xrad[i], col = piecol)
114     }
115     if (!is.null(text)) text(XX, YY, text, adj = adj, col = col, ...)
116     if (!is.null(pch)) points(XX + adj[1] - 0.5, YY + adj[2] - 0.5,
117                               pch = pch, col = col, bg = bg, ...)
118 }
119
120 nodelabels <- function(text, node, adj = c(0.5, 0.5), frame = "rect",
121                        pch = NULL, thermo = NULL, pie = NULL, piecol = NULL,
122                        col = "black", bg = "lightblue", ...)
123 {
124     lastPP <- get("last_plot.phylo", envir = .PlotPhyloEnv)
125     if (missing(node)) node <- (lastPP$Ntip + 1):length(lastPP$xx)
126     XX <- lastPP$xx[node]
127     YY <- lastPP$yy[node]
128     BOTHlabels(text, node, XX, YY, adj, frame, pch, thermo,
129                pie, piecol, col, bg, ...)
130 }
131
132 tiplabels <- function(text, tip, adj = c(0.5, 0.5), frame = "rect",
133                       pch = NULL, thermo = NULL, pie = NULL, piecol = NULL,
134                       col = "black", bg = "yellow", ...)
135 {
136     lastPP <- get("last_plot.phylo", envir = .PlotPhyloEnv)
137     if (missing(tip)) tip <- 1:lastPP$Ntip
138     XX <- lastPP$xx[tip]
139     YY <- lastPP$yy[tip]
140     BOTHlabels(text, tip, XX, YY, adj, frame, pch, thermo,
141                pie, piecol, col, bg, ...)
142 }
143
144 edgelabels <- function(text, edge, adj = c(0.5, 0.5), frame = "rect",
145                       pch = NULL, thermo = NULL, pie = NULL, piecol = NULL,
146                       col = "black", bg = "lightgreen", ...)
147 {
148     lastPP <- get("last_plot.phylo", envir = .PlotPhyloEnv)
149     if (missing(edge)) {
150         sel <- 1:dim(lastPP$edge)[1]
151         subedge <- lastPP$edge
152     } else {
153         sel <- edge
154         subedge <- lastPP$edge[sel, , drop = FALSE]
155     }
156     if (lastPP$type == "phylogram") {
157         if (lastPP$direction %in% c("rightwards", "leftwards")) {
158             XX <- (lastPP$xx[subedge[, 1]] + lastPP$xx[subedge[, 2]]) / 2
159             YY <- lastPP$yy[subedge[, 2]]
160         } else {
161             XX <- lastPP$xx[subedge[, 2]]
162             YY <- (lastPP$yy[subedge[, 1]] + lastPP$yy[subedge[, 2]]) / 2
163         }
164     } else {
165         XX <- (lastPP$xx[subedge[, 1]] + lastPP$xx[subedge[, 2]]) / 2
166         YY <- (lastPP$yy[subedge[, 1]] + lastPP$yy[subedge[, 2]]) / 2
167     }
168     BOTHlabels(text, sel, XX, YY, adj, frame, pch, thermo,
169                pie, piecol, col, bg, ...)
170 }
171
172 edges <- function(nodes0, nodes1, arrows = 0, type = "classical", ...)
173 {
174     type <- match.arg(type, c("classical", "triangle", "harpoon"))
175     lastPP <- get("last_plot.phylo", envir = .PlotPhyloEnv)
176     ## we do the recycling if necessary:
177     if (length(nodes0) != length(nodes1)) {
178         tmp <- cbind(nodes0, nodes1)
179         nodes0 <- tmp[, 1]
180         nodes1 <- tmp[, 2]
181     }
182     x0 <- lastPP$xx[nodes0]
183     y0 <- lastPP$yy[nodes0]
184     x1 <- lastPP$xx[nodes1]
185     y1 <- lastPP$yy[nodes1]
186     if (arrows)
187         if (type == "classical")
188             graphics::arrows(x0, y0, x1, y1, code = arrows, ...)
189         else
190             fancyarrows(x0, y0, x1, y1, code = arrows, type = type, ...)
191     else
192         graphics::segments(x0, y0, x1, y1, ...)
193 }
194
195 fancyarrows <-
196     function(x0, y0, x1, y1, length = 0.25, angle = 30, code = 2,
197              col = par("fg"), lty = par("lty"), lwd = par("lwd"),
198              type = "triangle", ...)
199 {
200     foo <- function(x0, y0, x1, y1) {
201         ## important to correct with these parameters cause
202         ## the coordinate system will likely not be Cartesian
203         pin <- par("pin")
204         usr <- par("usr")
205         A1 <- pin[1]/diff(usr[1:2])
206         A2 <- pin[2]/diff(usr[3:4])
207         x0 <- x0 * A1
208         y0 <- y0 * A2
209         x1 <- x1 * A1
210         y1 <- y1 * A2
211         atan2(y1 - y0, x1 - x0)
212     }
213     arrow.triangle <- function(x, y) {
214         beta <- alpha - angle/2
215         xa <- xinch(length * cos(beta)) + x
216         ya <- yinch(length * sin(beta)) + y
217         beta <- beta + angle
218         xb <- xinch(length * cos(beta)) + x
219         yb <- yinch(length * sin(beta)) + y
220         n <- length(x)
221         col <- rep(col, length.out = n)
222         for (i in 1:n)
223             polygon(c(x[i], xa[i], xb[i]), c(y[i], ya[i], yb[i]),
224                     col = col[i], border = col[i])
225         list((xa + xb)/2, (ya + yb)/2)
226     }
227     arrow.harpoon <- function(x, y) {
228         beta <- alpha - angle/2
229         xa <- xinch(length * cos(beta)) + x
230         ya <- yinch(length * sin(beta)) + y
231         beta <- alpha + angle/2
232         xb <- xinch(length * cos(beta)) + x
233         yb <- yinch(length * sin(beta)) + y
234         xc <- x/2 + (xa + xb)/4
235         yc <- y/2 + (ya + yb)/4
236         n <- length(x)
237         col <- rep(col, length.out = n)
238         for (i in 1:n)
239             polygon(c(x[i], xa[i], xc[i], xb[i]),
240                     c(y[i], ya[i], yc[i], yb[i]),
241                     col = col[i], border = col[i])
242         list(xc, yc)
243     }
244
245     type <- match.arg(type, c("triangle", "harpoon"))
246     angle <- pi*angle/180 # degree -> radian
247     alpha <- foo(x0, y0, x1, y1) # angle of segment with x-axis
248     ## alpha is in [-pi, pi]
249
250     FUN <- if (type == "triangle") arrow.triangle else arrow.harpoon
251     XY0 <- if (code == 1 || code == 3) FUN(x0, y0) else list(x0, y0)
252     if (code >= 2) {
253         alpha <- (alpha + pi) %% (2 * pi)
254         XY1 <- FUN(x1, y1)
255     } else XY1 <- list(x1, y1)
256     segments(XY0[[1]], XY0[[2]], XY1[[1]], XY1[[2]], col = col, lty = lty, lwd = lwd, ...)
257 }