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