]> git.donarmstrong.com Git - ape.git/blob - R/ace.R
last changes for ape 2.4-1
[ape.git] / R / ace.R
1 ## ace.R (2009-11-12)
2
3 ##   Ancestral Character Estimation
4
5 ## Copyright 2005-2009 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 argument "x" and the tip labels of the tree
29 did not match: the former were ignored in the analysis.')
30     }
31     obj <- list()
32     if (kappa != 1) phy$edge.length <- phy$edge.length^kappa
33     if (type == "continuous") {
34         if (method == "pic") {
35             if (model != "BM")
36               stop('the "pic" method can be used only with model = "BM".')
37             ## See pic.R for some annotations.
38             phy <- reorder(phy, "pruningwise")
39             phenotype <- numeric(nb.tip + nb.node)
40             phenotype[1:nb.tip] <- if (is.null(names(x))) x else x[phy$tip.label]
41             contr <- var.con <- numeric(nb.node)
42             ans <- .C("pic", as.integer(nb.tip), as.integer(nb.node),
43                       as.integer(phy$edge[, 1]), as.integer(phy$edge[, 2]),
44                       as.double(phy$edge.length), as.double(phenotype),
45                       as.double(contr), as.double(var.con),
46                       as.integer(CI), as.integer(scaled),
47                       PACKAGE = "ape")
48             obj$ace <- ans[[6]][-(1:nb.tip)]
49             names(obj$ace) <- (nb.tip + 1):(nb.tip + nb.node)
50             if (CI) {
51                 se <- sqrt(ans[[8]])
52                 CI95 <- matrix(NA, nb.node, 2)
53                 CI95[, 1] <- obj$ace + se * qnorm(0.025)
54                 CI95[, 2] <- obj$ace - se * qnorm(0.025)
55                 obj$CI95 <- CI95
56             }
57         }
58         if (method == "ML") {
59             if (model == "BM") {
60                 tip <- phy$edge[, 2] <= nb.tip
61                 dev.BM <- function(p) {
62                     if (p[1] < 0) return(1e100) # in case sigma² is negative
63                     x1 <- p[-1][phy$edge[, 1] - nb.tip]
64                     x2 <- numeric(length(x1))
65                     x2[tip] <- x[phy$edge[tip, 2]]
66                     x2[!tip] <- p[-1][phy$edge[!tip, 2] - nb.tip]
67                     -2 * (-sum((x1 - x2)^2/phy$edge.length)/(2*p[1]) -
68                           nb.node * log(p[1]))
69                 }
70                 out <- nlm(function(p) dev.BM(p),
71                            p = c(1, rep(mean(x), nb.node)), hessian = TRUE)
72                 obj$loglik <- -out$minimum / 2
73                 obj$ace <- out$estimate[-1]
74                 names(obj$ace) <- (nb.tip + 1):(nb.tip + nb.node)
75                 se <- sqrt(diag(solve(out$hessian)))
76                 obj$sigma2 <- c(out$estimate[1], se[1])
77                 se <- se[-1]
78                 if (CI) {
79                     CI95 <- matrix(NA, nb.node, 2)
80                     CI95[, 1] <- obj$ace + se * qt(0.025, nb.node)
81                     CI95[, 2] <- obj$ace - se * qt(0.025, nb.node)
82                     obj$CI95 <- CI95
83                 }
84             }
85         }
86         if (method == "GLS") {
87             if (is.null(corStruct))
88               stop('you must give a correlation structure if method = "GLS".')
89             if (class(corStruct)[1] == "corMartins")
90               M <- corStruct[1] * dist.nodes(phy)
91             if (class(corStruct)[1] == "corGrafen")
92               phy <- compute.brlen(attr(corStruct, "tree"),
93                                    method = "Grafen",
94                                    power = exp(corStruct[1]))
95             if (class(corStruct)[1] %in% c("corBrownian", "corGrafen")) {
96                 dis <- dist.nodes(attr(corStruct, "tree"))
97                 MRCA <- mrca(attr(corStruct, "tree"), full = TRUE)
98                 M <- dis[as.character(nb.tip + 1), MRCA]
99                 dim(M) <- rep(sqrt(length(M)), 2)
100             }
101             varAY <- M[-(1:nb.tip), 1:nb.tip]
102             varA <- M[-(1:nb.tip), -(1:nb.tip)]
103             V <- corMatrix(Initialize(corStruct, data.frame(x)),
104                            corr = FALSE)
105             invV <- solve(V)
106             o <- gls(x ~ 1, data.frame(x), correlation = corStruct)
107             GM <- o$coefficients
108             obj$ace <- drop(varAY %*% invV %*% (x - GM) + GM)
109             names(obj$ace) <- (nb.tip + 1):(nb.tip + nb.node)
110             if (CI) {
111                 CI95 <- matrix(NA, nb.node, 2)
112                 se <- sqrt((varA - varAY %*% invV %*% t(varAY))[cbind(1:nb.node, 1:nb.node)])
113                 CI95[, 1] <- obj$ace + se * qnorm(0.025)
114                 CI95[, 2] <- obj$ace - se * qnorm(0.025)
115                 obj$CI95 <- CI95
116             }
117         }
118     } else { # type == "discrete"
119         if (method != "ML")
120           stop("only ML estimation is possible for discrete characters.")
121         if (!is.factor(x)) x <- factor(x)
122         nl <- nlevels(x)
123         lvls <- levels(x)
124         x <- as.integer(x)
125         if (is.character(model)) {
126             rate <- matrix(NA, nl, nl)
127             if (model == "ER") np <- rate[] <- 1
128             if (model == "ARD") {
129                 np <- nl*(nl - 1)
130                 rate[col(rate) != row(rate)] <- 1:np
131             }
132             if (model == "SYM") {
133                 np <- nl * (nl - 1)/2
134                 rate[col(rate) < row(rate)] <- 1:np
135                 rate <- t(rate)
136                 rate[col(rate) < row(rate)] <- 1:np
137             }
138         } else {
139             if (ncol(model) != nrow(model))
140               stop("the matrix given as `model' is not square")
141             if (ncol(model) != nl)
142               stop("the matrix `model' must have as many rows
143 as the number of categories in `x'")
144             rate <- model
145             np <- max(rate)
146         }
147         index.matrix <- rate
148         index.matrix[cbind(1:nl, 1:nl)] <- NA
149         rate[cbind(1:nl, 1:nl)] <- 0
150         rate[rate == 0] <- np + 1 # to avoid 0's since we will use this an numeric indexing
151
152         liks <- matrix(0, nb.tip + nb.node, nl)
153         TIPS <- 1:nb.tip
154         liks[cbind(TIPS, x)] <- 1
155         phy <- reorder(phy, "pruningwise")
156
157         Q <- matrix(0, nl, nl)
158         ## from Rich FitzJohn:
159         comp <- numeric(nb.tip + nb.node) # Storage...
160         dev <- function(p, output.liks = FALSE) {
161             Q[] <- c(p, 0)[rate]
162             diag(Q) <- -rowSums(Q)
163             for (i  in seq(from = 1, by = 2, length.out = nb.node)) {
164                 j <- i + 1
165                 anc <- phy$edge[i, 1]
166                 des1 <- phy$edge[i, 2]
167                 des2 <- phy$edge[j, 2]
168                 v.l <- matexpo(Q * phy$edge.length[i]) %*% liks[des1, ]
169                 v.r <- matexpo(Q * phy$edge.length[j]) %*% liks[des2, ]
170                 v <- v.l * v.r
171                 comp[anc] <- sum(v)
172                 liks[anc, ] <- v/comp[anc]
173             }
174             if (output.liks) return(liks[-TIPS, ])
175             -2 * sum(log(comp[-TIPS]))
176         }
177         out <- nlminb(rep(ip, length.out = np), function(p) dev(p),
178                       lower = rep(0, np), upper = rep(Inf, np))
179         obj$loglik <- -out$objective/2
180         obj$rates <- out$par
181         oldwarn <- options("warn")
182         options(warn = -1)
183         h <- nlm(function(p) dev(p), p = obj$rates, iterlim = 1,
184                  stepmax = 0, hessian = TRUE)$hessian
185         options(oldwarn)
186         if (any(h == 0))
187           warning("The likelihood gradient seems flat in at least one dimension (gradient null):\ncannot compute the standard-errors of the transition rates.\n")
188         else obj$se <- sqrt(diag(solve(h)))
189         obj$index.matrix <- index.matrix
190         if (CI) {
191             obj$lik.anc <- dev(obj$rates, TRUE)
192             colnames(obj$lik.anc) <- lvls
193         }
194     }
195     obj$call <- match.call()
196     class(obj) <- "ace"
197     obj
198 }
199
200 logLik.ace <- function(object, ...) object$loglik
201
202 deviance.ace <- function(object, ...) -2*object$loglik
203
204 AIC.ace <- function(object, ..., k = 2)
205 {
206     if (is.null(object$loglik)) return(NULL)
207     ## Trivial test of "type"; may need to be improved
208     ## if other models are included in ace(type = "c")
209     np <- if (!is.null(object$sigma2)) 1 else length(object$rates)
210     -2*object$loglik + np*k
211 }
212
213 ### by BB:
214 anova.ace <- function(object, ...)
215 {
216     X <- c(list(object), list(...))
217     df <- sapply(lapply(X, "[[", "rates"), length)
218     ll <- sapply(X, "[[", "loglik")
219     ## check if models are in correct order?
220     dev <- c(NA, 2*diff(ll))
221     ddf <- c(NA, diff(df))
222     table <- data.frame(ll, df, ddf, dev,
223                         pchisq(dev, ddf, lower.tail = FALSE))
224     dimnames(table) <- list(1:length(X), c("Log lik.", "Df",
225                                            "Df change", "Deviance",
226                                            "Pr(>|Chi|)"))
227     structure(table, heading = "Likelihood Ratio Test Table",
228               class = c("anova", "data.frame"))
229 }