]> git.donarmstrong.com Git - ape.git/blob - R/ace.R
new option 'rotate.tree' in plot.phylo()
[ape.git] / R / ace.R
1 ## ace.R (2011-07-18)
2
3 ##   Ancestral Character Estimation
4
5 ## Copyright 2005-2011 Emmanuel Paradis and Ben Bolker
6
7 ## This file is part of the R-package `ape'.
8 ## See the file ../COPYING for licensing issues.
9
10 ace <- function(x, phy, type = "continuous", method = "ML", CI = TRUE,
11                 model = if (type == "continuous") "BM" else "ER",
12                 scaled = TRUE, kappa = 1, corStruct = NULL, ip = 0.1)
13 {
14     if (!inherits(phy, "phylo"))
15         stop('object "phy" is not of class "phylo".')
16     if (is.null(phy$edge.length))
17         stop("tree has no branch lengths")
18     type <- match.arg(type, c("continuous", "discrete"))
19     nb.tip <- length(phy$tip.label)
20     nb.node <- phy$Nnode
21     if (nb.node != nb.tip - 1)
22         stop('"phy" is not rooted AND fully dichotomous.')
23     if (length(x) != nb.tip)
24         stop("length of phenotypic and of phylogenetic data do not match.")
25     if (!is.null(names(x))) {
26         if(all(names(x) %in% phy$tip.label))
27           x <- x[phy$tip.label]
28         else warning("the names of 'x' and the tip labels of the tree do not match: the former were ignored in the analysis.")
29     }
30     obj <- list()
31     if (kappa != 1) phy$edge.length <- phy$edge.length^kappa
32     if (type == "continuous") {
33         switch(method, "REML" = {
34             minusLogLik <- function(sig2) {
35                 if (sig2 < 0) return(1e100)
36                 V <- sig2 * vcv(phy)
37                 ## next three lines borrowed from dmvnorm() in 'mvtnorm'
38                 distval <- mahalanobis(x, center = mu, cov = V)
39                 logdet <- sum(log(eigen(V, symmetric = TRUE, only.values = TRUE)$values))
40                 (nb.tip * log(2 * pi) + logdet + distval)/2
41             }
42             mu <- rep(ace(x, phy, method="pic")$ace[1], nb.tip)
43             out <- nlm(minusLogLik, 1, hessian = TRUE)
44             sigma2 <- out$estimate
45             se_sgi2 <- sqrt(1/out$hessian)
46             tip <- phy$edge[, 2] <= nb.tip
47             minus.REML.BM <- function(p) {
48                 x1 <- p[phy$edge[, 1] - nb.tip]
49                 x2 <- numeric(length(x1))
50                 x2[tip] <- x[phy$edge[tip, 2]]
51                 x2[!tip] <- p[phy$edge[!tip, 2] - nb.tip]
52                 -(-sum((x1 - x2)^2/phy$edge.length)/(2 * sigma2) -
53                   nb.node * log(sigma2))
54             }
55             out <- nlm(function(p) minus.REML.BM(p),
56                        p = rep(mu[1], nb.node), hessian = TRUE)
57             obj$resloglik <- -out$minimum
58             obj$ace <- out$estimate
59             names(obj$ace) <- (nb.tip + 1):(nb.tip + nb.node)
60             obj$sigma2 <- c(sigma2, se_sgi2)
61             if (CI) {
62                 se <- sqrt(diag(solve(out$hessian)))
63                 tmp <- se * qt(0.025, nb.node)
64                 obj$CI95 <- cbind(obj$ace + tmp, obj$ace - tmp)
65             }
66         }, "pic" = {
67             if (model != "BM")
68                 stop('the "pic" method can be used only with model = "BM".')
69             ## See pic.R for some annotations.
70             phy <- reorder(phy, "pruningwise")
71             phenotype <- numeric(nb.tip + nb.node)
72             phenotype[1:nb.tip] <- if (is.null(names(x))) x else x[phy$tip.label]
73             contr <- var.con <- numeric(nb.node)
74             ans <- .C("pic", as.integer(nb.tip), as.integer(nb.node),
75                       as.integer(phy$edge[, 1]), as.integer(phy$edge[, 2]),
76                       as.double(phy$edge.length), as.double(phenotype),
77                       as.double(contr), as.double(var.con),
78                       as.integer(CI), as.integer(scaled),
79                       PACKAGE = "ape")
80             obj$ace <- ans[[6]][-(1:nb.tip)]
81             names(obj$ace) <- (nb.tip + 1):(nb.tip + nb.node)
82             if (CI) {
83                 se <- sqrt(ans[[8]])
84                 tmp <- se * qnorm(0.025)
85                 obj$CI95 <- cbind(obj$ace + tmp, obj$ace - tmp)
86             }
87         }, "ML" = {
88             if (model == "BM") {
89                 tip <- phy$edge[, 2] <= nb.tip
90                 dev.BM <- function(p) {
91                     if (p[1] < 0) return(1e100) # in case sigma² is negative
92                     x1 <- p[-1][phy$edge[, 1] - nb.tip]
93                     x2 <- numeric(length(x1))
94                     x2[tip] <- x[phy$edge[tip, 2]]
95                     x2[!tip] <- p[-1][phy$edge[!tip, 2] - nb.tip]
96                     -2 * (-sum((x1 - x2)^2/phy$edge.length)/(2*p[1]) -
97                           nb.node * log(p[1]))
98                 }
99                 out <- nlm(function(p) dev.BM(p),
100                            p = c(1, rep(mean(x), nb.node)), hessian = TRUE)
101                 obj$loglik <- -out$minimum / 2
102                 obj$ace <- out$estimate[-1]
103                 names(obj$ace) <- (nb.tip + 1):(nb.tip + nb.node)
104                 se <- sqrt(diag(solve(out$hessian)))
105                 obj$sigma2 <- c(out$estimate[1], se[1])
106                 if (CI) {
107                     tmp <- se[-1] * qt(0.025, nb.node)
108                     obj$CI95 <- cbind(obj$ace + tmp, obj$ace - tmp)
109                 }
110             }
111         }, "GLS" = {
112             if (is.null(corStruct))
113                 stop('you must give a correlation structure if method = "GLS".')
114             if (class(corStruct)[1] == "corMartins")
115                 M <- corStruct[1] * dist.nodes(phy)
116             if (class(corStruct)[1] == "corGrafen")
117                 phy <- compute.brlen(attr(corStruct, "tree"),
118                                      method = "Grafen",
119                                      power = exp(corStruct[1]))
120             if (class(corStruct)[1] %in% c("corBrownian", "corGrafen")) {
121                 dis <- dist.nodes(attr(corStruct, "tree"))
122                 MRCA <- mrca(attr(corStruct, "tree"), full = TRUE)
123                 M <- dis[as.character(nb.tip + 1), MRCA]
124                 dim(M) <- rep(sqrt(length(M)), 2)
125             }
126             varAY <- M[-(1:nb.tip), 1:nb.tip]
127             varA <- M[-(1:nb.tip), -(1:nb.tip)]
128             V <- corMatrix(Initialize(corStruct, data.frame(x)),
129                            corr = FALSE)
130             invV <- solve(V)
131             o <- gls(x ~ 1, data.frame(x), correlation = corStruct)
132             GM <- o$coefficients
133             obj$ace <- drop(varAY %*% invV %*% (x - GM) + GM)
134             names(obj$ace) <- (nb.tip + 1):(nb.tip + nb.node)
135             if (CI) {
136                 se <- sqrt((varA - varAY %*% invV %*% t(varAY))[cbind(1:nb.node, 1:nb.node)])
137                 tmp <- se * qnorm(0.025)
138                 obj$CI95 <- cbind(obj$ace + tmp, obj$ace - tmp)
139             }
140         })
141     } else { # type == "discrete"
142         if (method != "ML")
143             stop("only ML estimation is possible for discrete characters.")
144         if (!is.factor(x)) x <- factor(x)
145         nl <- nlevels(x)
146         lvls <- levels(x)
147         x <- as.integer(x)
148         if (is.character(model)) {
149             rate <- matrix(NA, nl, nl)
150             if (model == "ER") np <- rate[] <- 1
151             if (model == "ARD") {
152                 np <- nl*(nl - 1)
153                 rate[col(rate) != row(rate)] <- 1:np
154             }
155             if (model == "SYM") {
156                 np <- nl * (nl - 1)/2
157                 sel <- col(rate) < row(rate)
158                 rate[sel] <- 1:np
159                 rate <- t(rate)
160                 rate[sel] <- 1:np
161             }
162         } else {
163             if (ncol(model) != nrow(model))
164               stop("the matrix given as 'model' is not square")
165             if (ncol(model) != nl)
166               stop("the matrix 'model' must have as many rows as the number of categories in 'x'")
167             rate <- model
168             np <- max(rate)
169         }
170         index.matrix <- rate
171         tmp <- cbind(1:nl, 1:nl)
172         index.matrix[tmp] <- NA
173         rate[tmp] <- 0
174         rate[rate == 0] <- np + 1 # to avoid 0's since we will use this as numeric indexing
175
176         liks <- matrix(0, nb.tip + nb.node, nl)
177         TIPS <- 1:nb.tip
178         liks[cbind(TIPS, x)] <- 1
179         phy <- reorder(phy, "pruningwise")
180
181         Q <- matrix(0, nl, nl)
182         dev <- function(p, output.liks = FALSE) {
183             if (any(is.nan(p)) || any(is.infinite(p))) return(1e50)
184             ## from Rich FitzJohn:
185             comp <- numeric(nb.tip + nb.node) # Storage...
186             Q[] <- c(p, 0)[rate]
187             diag(Q) <- -rowSums(Q)
188             for (i  in seq(from = 1, by = 2, length.out = nb.node)) {
189                 j <- i + 1L
190                 anc <- phy$edge[i, 1]
191                 des1 <- phy$edge[i, 2]
192                 des2 <- phy$edge[j, 2]
193                 v.l <- matexpo(Q * phy$edge.length[i]) %*% liks[des1, ]
194                 v.r <- matexpo(Q * phy$edge.length[j]) %*% liks[des2, ]
195                 v <- v.l * v.r
196                 comp[anc] <- sum(v)
197                 liks[anc, ] <- v/comp[anc]
198             }
199             if (output.liks) return(liks[-TIPS, ])
200             dev <- -2 * sum(log(comp[-TIPS]))
201             if (is.na(dev)) Inf else dev
202         }
203         out <- nlminb(rep(ip, length.out = np), function(p) dev(p),
204                       lower = rep(0, np), upper = rep(1e50, np))
205         obj$loglik <- -out$objective/2
206         obj$rates <- out$par
207         oldwarn <- options("warn")
208         options(warn = -1)
209         h <- nlm(function(p) dev(p), p = obj$rates, iterlim = 1,
210                  stepmax = 0, hessian = TRUE)$hessian
211         options(oldwarn)
212         if (any(h == 0))
213           warning("The likelihood gradient seems flat in at least one dimension (gradient null):\ncannot compute the standard-errors of the transition rates.\n")
214         else obj$se <- sqrt(diag(solve(h)))
215         obj$index.matrix <- index.matrix
216         if (CI) {
217             obj$lik.anc <- dev(obj$rates, TRUE)
218             colnames(obj$lik.anc) <- lvls
219         }
220     }
221     obj$call <- match.call()
222     class(obj) <- "ace"
223     obj
224 }
225
226 logLik.ace <- function(object, ...) object$loglik
227
228 deviance.ace <- function(object, ...) -2*object$loglik
229
230 AIC.ace <- function(object, ..., k = 2)
231 {
232     if (is.null(object$loglik)) return(NULL)
233     ## Trivial test of "type"; may need to be improved
234     ## if other models are included in ace(type = "c")
235     np <- if (!is.null(object$sigma2)) 1 else length(object$rates)
236     -2*object$loglik + np*k
237 }
238
239 ### by BB:
240 anova.ace <- function(object, ...)
241 {
242     X <- c(list(object), list(...))
243     df <- sapply(lapply(X, "[[", "rates"), length)
244     ll <- sapply(X, "[[", "loglik")
245     ## check if models are in correct order?
246     dev <- c(NA, 2*diff(ll))
247     ddf <- c(NA, diff(df))
248     table <- data.frame(ll, df, ddf, dev,
249                         pchisq(dev, ddf, lower.tail = FALSE))
250     dimnames(table) <- list(1:length(X), c("Log lik.", "Df",
251                                            "Df change", "Resid. Dev",
252                                            "Pr(>|Chi|)"))
253     structure(table, heading = "Likelihood Ratio Test Table",
254               class = c("anova", "data.frame"))
255 }
256
257 print.ace <- function(x, digits = 4, ...)
258 {
259     cat("\n    Ancestral Character Estimation\n\n")
260     cat("Call: ")
261     print(x$call)
262     cat("\n")
263     if (!is.null(x$loglik))
264         cat("    Log-likelihood:", x$loglik, "\n\n")
265     if (!is.null(x$resloglik))
266         cat("    Residual log-likelihood:", x$resloglik, "\n\n")
267     ratemat <- x$index.matrix
268     if (is.null(ratemat)) { # to be improved
269         class(x) <- NULL
270         x$resloglik <- x$loglik <- x$call <- NULL
271         print(x)
272     } else {
273         dimnames(ratemat)[1:2] <- dimnames(x$lik.anc)[2]
274         cat("Rate index matrix:\n")
275         print(ratemat, na.print = ".")
276         cat("\n")
277         npar <- length(x$rates)
278         estim <- data.frame(1:npar, round(x$rates, digits), round(x$se, digits))
279         cat("Parameter estimates:\n")
280         names(estim) <- c("rate index", "estimate", "std-err")
281         print(estim, row.names = FALSE)
282         cat("\nScaled likelihoods at the root (type '...$lik.anc' to get them for all nodes):\n")
283         print(x$lik.anc[1, ])
284     }
285 }