]> git.donarmstrong.com Git - ape.git/commitdiff
some updates for ape 3.0-6
authorparadis <paradis@6e262413-ae40-0410-9e79-b911bd7a66b7>
Tue, 2 Oct 2012 04:00:29 +0000 (04:00 +0000)
committerparadis <paradis@6e262413-ae40-0410-9e79-b911bd7a66b7>
Tue, 2 Oct 2012 04:00:29 +0000 (04:00 +0000)
git-svn-id: https://svn.mpl.ird.fr/ape/dev/ape@196 6e262413-ae40-0410-9e79-b911bd7a66b7

17 files changed:
NEWS
R/CDF.birth.death.R
R/DNA.R
R/dbd.R
R/evonet.R
R/ladderize.R
R/me.R
R/pic.R
R/plot.phylo.R
R/read.nexus.R
R/read.tree.R
R/rtree.R
man/ape-internal.Rd
man/dbd.Rd
man/node.depth.Rd
man/rTraitCont.Rd
src/plot_phylo.c

diff --git a/NEWS b/NEWS
index 46af03fc06ef1fcf9e8a61d37760c26251265397..729e3ed930e14e4404d261978374069c41576c6c 100644 (file)
--- a/NEWS
+++ b/NEWS
@@ -5,13 +5,18 @@ NEW FEATURES
 
     o reorder.phylo() has a new order, "postorder", and a new option
       index.only = TRUE to return only the vector of indices (the tree
-      is unmodified; see ?reorder.phylo for details).
+      is unmodified, see ?reorder.phylo for details).
+
+    o The three new functions node.depth.length, node.height, and
+      node.height.clado make some internal code available from R. See
+      ?node.depth (which was already available and documented) for
+      details.
 
 
 BUG FIXES
 
     o reorder(, "pruningwise") made R crash if the rows of the edge
-      matrix are in random order.
+      matrix are in random order: this is now fixed.
 
 
 OTHER CHANGES
@@ -22,6 +27,14 @@ OTHER CHANGES
       visible for small trees (n < 1000) but this can be more than
       1000 faster for big trees (n >= 1e4).
 
+    o The attribute "order" of the objects of class "phylo" is now
+      strongly recommended, though not mandatory. Most functions in
+      ape should return a tree with this attribute correctly set.
+
+    o dbd() is now vectorized on both arguments 'x' (number of species
+      in clade) and 't' (clade age) to make likelihood calculations
+      easier and faster.
+
 
 
                CHANGES IN APE VERSION 3.0-5
index 847c37dcd20dd60ecc4649a725ff27ccfe27fd59..7ba5f5eb175d078472c6a8c7ccc1a818c1c67d57 100644 (file)
@@ -1,9 +1,9 @@
-## CDF.birth.death.R (2010-09-27)
+## CDF.birth.death.R (2012-09-14)
 
-##    Functions to simulate and fit
+##    Functions to Simulate and Fit
 ##  Time-Dependent Birth-Death Models
 
-## Copyright 2010 Emmanuel Paradis
+## Copyright 2010-2012 Emmanuel Paradis
 
 ## This file is part of the R-package `ape'.
 ## See the file ../COPYING for licensing issues.
@@ -158,7 +158,9 @@ CDF.birth.death <-
 
     if (!case %in% c(1, 3, 6)) Pi <- Vectorize(Pi)
 
-    denom <- if (fast) integrateTrapeze(Pi, 0, Tmax) else integrate(Pi, 0, Tmax)$value
+    denom <-
+        if (fast) integrateTrapeze(Pi, 0, Tmax)
+        else integrate(Pi, 0, Tmax)$value
     n <- length(x)
     p <- numeric(n)
     if (fast) {
@@ -178,6 +180,7 @@ CDF.birth.death <-
     phy <- list(edge = edge, edge.length = edge.length,
                 tip.label = paste("t", 1:(i + 1), sep = ""), Nnode = i)
     class(phy) <- "phylo"
+    attr(phy, "order") <- "cladewise"
     phy
 }
 
diff --git a/R/DNA.R b/R/DNA.R
index 2bd752060805a657374051f7493775f5d402ceb4..b426dbd6126728c02992f994fe408a593d2a4f27 100644 (file)
--- a/R/DNA.R
+++ b/R/DNA.R
@@ -1,4 +1,4 @@
-## DNA.R (2012-06-19)
+## DNA.R (2012-09-13)
 
 ##   Manipulations and Comparisons of DNA Sequences
 
@@ -295,7 +295,7 @@ base.freq <- function(x, freq = FALSE, all = FALSE)
 {
     if (is.list(x)) x <- unlist(x)
     n <- length(x)
-    BF <-.C("BaseProportion", x, n, double(17),
+    BF <-.C("BaseProportion", x, as.integer(n), double(17),
             DUP = FALSE, NAOK = TRUE, PACKAGE = "ape")[[3]]
     names(BF) <- c("a", "c", "g", "t", "r", "m", "w", "s",
                    "k", "y", "v", "h", "d", "b", "n", "-", "?")
@@ -351,8 +351,8 @@ seg.sites <- function(x)
         n <- n[1]
     }
     if (n == 1) return(integer(0))
-    ans <- .C("SegSites", x, n, s, integer(s),
-              DUP = FALSE, NAOK = TRUE, PACKAGE = "ape")
+    ans <- .C("SegSites", x, as.integer(n), as.integer(s),
+              integer(s), DUP = FALSE, NAOK = TRUE, PACKAGE = "ape")
     which(as.logical(ans[[4]]))
 }
 
@@ -397,10 +397,10 @@ dist.dna <- function(x, model = "K80", variance = FALSE, gamma = FALSE,
     var <- if (variance) double(Ndist) else 0
     if (!gamma) gamma <- alpha <- 0
     else alpha <- gamma <- 1
-    d <- .C("dist_dna", x, n, s, imod, double(Ndist), BF,
-            as.integer(pairwise.deletion), as.integer(variance),
-            var, as.integer(gamma), alpha, DUP = FALSE, NAOK = TRUE,
-            PACKAGE = "ape")
+    d <- .C("dist_dna", x, as.integer(n), as.integer(s), imod,
+            double(Ndist), BF, as.integer(pairwise.deletion),
+            as.integer(variance), var, as.integer(gamma),
+            alpha, DUP = FALSE, NAOK = TRUE, PACKAGE = "ape")
     if (variance) var <- d[[9]]
     d <- d[[5]]
     if (imod == 11) {
diff --git a/R/dbd.R b/R/dbd.R
index 3a325dc6fcee654815446b6aa57c19e59b48893a..b0beb7a4078af31a2ba1d9d77d956c11968a8d86 100644 (file)
--- a/R/dbd.R
+++ b/R/dbd.R
@@ -1,4 +1,4 @@
-## dbd.R (2012-03-19)
+## dbd.R (2012-09-14)
 
 ##   Probability Density Under Birth--Death Models
 
@@ -27,10 +27,6 @@ dbd <- function(x, lambda, mu, t, conditional = FALSE, log = FALSE)
         mu <- mu[1]
         warning("only the first value of 'mu' was considered")
     }
-    if (length(t) > 1) {
-        t <- t[1]
-        warning("only the first value of 't' was considered")
-    }
 
     if (mu == 0) return(dyule(x, lambda, t, log))
 
index 9dbaaa3d75db4e71f7cb1cdc7b5414c27baf0a4d..7f6b45f7e4ad136c4ea1507cb88a0fe2d6241d89 100644 (file)
@@ -1,8 +1,8 @@
-## evonet.R (2011-06-09)
+## evonet.R (2012-09-14)
 
 ##   Evolutionary Networks
 
-## Copyright 2011 Emmanuel Paradis
+## Copyright 2011-2012 Emmanuel Paradis
 
 ## This file is part of the R-package `ape'.
 ## See the file ../COPYING for licensing issues.
@@ -67,7 +67,7 @@ as.networx.evonet <- function(x, weight = NA, ...)
 {
     if (any(x$reticulation <= Ntip(x)))
         stop("some tips are involved in reticulations: cannot convert to \"networx\"")
-    x <- reorder(x, "p")
+    x <- reorder(x, "pruningwise")
     ned <- Nedge(x)
     nrt <- nrow(x$reticulation)
     x$edge <- rbind(x$edge, x$reticulation)
index 8d2e94a240010eee9a17eb6b763777b66744128d..596a8fdd7aca261b89e89009388d5f54c0aec8a1 100644 (file)
@@ -23,7 +23,7 @@ ladderize <- function(phy, right = TRUE)
         start <- start[o]
         neworder[newpos] <<- start
         for (i in 1:Nclade)
-          if (desc[i] > nb.tip) foo(desc[i], end[i], newpos[i] + 1)
+            if (desc[i] > nb.tip) foo(desc[i], end[i], newpos[i] + 1)
     }
     nb.tip <- length(phy$tip.label)
     nb.node <- phy$Nnode
@@ -37,6 +37,6 @@ ladderize <- function(phy, right = TRUE)
     foo(nb.tip + 1, nb.edge, 1)
     phy$edge <- phy$edge[neworder, ]
     if (!is.null(phy$edge.length))
-      phy$edge.length <- phy$edge.length[neworder]
+        phy$edge.length <- phy$edge.length[neworder]
     phy
 }
diff --git a/R/me.R b/R/me.R
index 1006c6432ec0343115d1fb0003d93fe8d59fd656..533c40edab6d7b52079943f7e9ee57e1c668d6af 100644 (file)
--- a/R/me.R
+++ b/R/me.R
@@ -1,4 +1,4 @@
-## me.R (2012-04-30)
+## me.R (2012-09-14)
 
 ##   Tree Estimation Based on Minimum Evolution Algorithm
 
@@ -20,10 +20,12 @@ fastme.bal <- function(X, nni = TRUE, spr = TRUE, tbr = TRUE)
     labels <- attr(X, "Labels")
     if (is.null(labels)) labels <- as.character(1:N)
     labels <- labels[ans[[3]]]
-    structure(list(edge =  cbind(ans[[7]], ans[[8]]),
-                   edge.length = ans[[9]],
-                   tip.label = labels, Nnode = N - 2L),
-              class = "phylo")
+    obj <- list(edge =  cbind(ans[[7]], ans[[8]]),
+                edge.length = ans[[9]],
+                tip.label = labels, Nnode = N - 2L)
+    class(obj) <- "phylo"
+    attr(obj, "order") <- "cladewise"
+    obj
 }
 
 fastme.ols <- function(X, nni = TRUE)
@@ -37,10 +39,12 @@ fastme.ols <- function(X, nni = TRUE)
     labels <- attr(X, "Labels")
     if (is.null(labels)) labels <- as.character(1:N)
     labels <- labels[ans[[3]]]
-    structure(list(edge =  cbind(ans[[5]], ans[[6]]),
-                   edge.length = ans[[7]],
-                   tip.label = labels, Nnode = N - 2L),
-              class = "phylo")
+    obj <- list(edge =  cbind(ans[[5]], ans[[6]]),
+                edge.length = ans[[7]],
+                tip.label = labels, Nnode = N - 2L)
+    class(obj) <- "phylo"
+    attr(obj, "order") <- "cladewise"
+    obj
 }
 
 bionj <- function(X)
diff --git a/R/pic.R b/R/pic.R
index 816a09b0a6b29b4de083c6f681ea2fe90e0b3b77..b1cb20b17921668f98e4375a95093d4ad3fda44b 100644 (file)
--- a/R/pic.R
+++ b/R/pic.R
@@ -1,8 +1,8 @@
-## pic.R (2011-03-01)
+## pic.R (2012-09-11)
 
 ##   Phylogenetically Independent Contrasts
 
-## Copyright 2002-2011 Emmanuel Paradis
+## Copyright 2002-2012 Emmanuel Paradis
 
 ## This file is part of the R-package `ape'.
 ## See the file ../COPYING for licensing issues.
@@ -22,7 +22,7 @@ pic <- function(x, phy, scaled = TRUE, var.contrasts = FALSE)
     if (any(is.na(x)))
       stop("missing data in 'x': you may consider removing the species with missing data from your tree with the function 'drop.tip'.")
 
-    phy <- reorder(phy, "pruningwise")
+    phy <- reorder(phy, "postorder")
     phenotype <- numeric(nb.tip + nb.node)
 
     if (is.null(names(x))) {
index 3fd0637295d84063afd81cbe2bea7d5509bbe7d4..862291de6ccd885bddeaa5a40f09a70f87a8a1ee 100644 (file)
@@ -1,4 +1,4 @@
-## plot.phylo.R (2012-03-22)
+## plot.phylo.R (2012-10-02)
 
 ##   Plot Phylogenies
 
@@ -595,6 +595,57 @@ node.depth <- function(phy)
        as.integer(N), double(n + m), DUP = FALSE, PACKAGE = "ape")[[6]]
 }
 
+node.depth.edgelength <- function(phy)
+{
+    n <- length(phy$tip.label)
+    m <- phy$Nnode
+    N <- dim(phy$edge)[1]
+    phy <- reorder(phy, order = "pruningwise")
+    .C("node_depth_edgelength", as.integer(n), as.integer(n),
+       as.integer(phy$edge[, 1]), as.integer(phy$edge[, 2]),
+       as.integer(N), as.double(phy$edge.length), double(n + m),
+       DUP = FALSE, PACKAGE = "ape")[[7]]
+}
+
+node.height <- function(phy)
+{
+    n <- length(phy$tip.label)
+    m <- phy$Nnode
+    N <- dim(phy$edge)[1]
+    phy <- reorder(phy, order = "pruningwise")
+
+    e1 <- phy$edge[, 1]
+    e2 <- phy$edge[, 2]
+
+    yy <- numeric(n + m)
+    TIPS <- e2[e2 <= n]
+    yy[TIPS] <- 1:n
+
+    .C("node_height", as.integer(n), as.integer(m),
+       as.integer(e1), as.integer(e2), as.integer(N),
+       as.double(yy), DUP = FALSE, PACKAGE = "ape")[[6]]
+}
+
+node.height.clado <- function(phy)
+{
+    n <- length(phy$tip.label)
+    m <- phy$Nnode
+    N <- dim(phy$edge)[1]
+    phy <- reorder(phy, order = "pruningwise")
+
+    e1 <- phy$edge[, 1]
+    e2 <- phy$edge[, 2]
+
+    yy <- numeric(n + m)
+    TIPS <- e2[e2 <= n]
+    yy[TIPS] <- 1:n
+
+    .C("node_height_clado", as.integer(n), as.integer(m),
+       as.integer(e1), as.integer(e2), as.integer(N),
+       double(n + m), as.double(yy), DUP = FALSE,
+       PACKAGE = "ape")[[7]]
+}
+
 plot.multiPhylo <- function(x, layout = 1, ...)
 {
     if (layout > 1)
index 0268e6833a448853d78822c12ef9e8bbe4cb0904..2f4177e5aed51ebab7d48969b6fda7e7db86b421 100644 (file)
@@ -1,4 +1,4 @@
-## read.nexus.R (2012-02-09)
+## read.nexus.R (2012-09-28)
 
 ##   Read Tree File in Nexus Format
 
@@ -16,6 +16,7 @@
     names(phy) <- nms
     if (all(phy$node.label == "")) phy$node.label <- NULL
     class(phy) <- "phylo"
+    attr(phy, "order") <- "cladewise"
     phy
 }
 
@@ -97,6 +98,7 @@ clado.build <- function(tp)
         if (all(obj$node.label == "NA")) NULL
         else gsub("^NA", "", obj$node.label)
     class(obj) <- "phylo"
+    attr(obj, "order") <- "cladewise"
     obj
 }
 
index cba51d10fd7a5641568b12480006bc6c7635078e..47a432c39793e6b58ceef307e4248d588b1c6249 100644 (file)
@@ -1,8 +1,8 @@
-## read.tree.R (2010-09-27)
+## read.tree.R (2012-09-14)
 
 ##   Read Tree Files in Parenthetic Format
 
-## Copyright 2002-2010 Emmanuel Paradis, Daniel Lawson and Klaus Schliep
+## Copyright 2002-2012 Emmanuel Paradis, Daniel Lawson and Klaus Schliep
 
 ## This file is part of the R-package `ape'.
 ## See the file ../COPYING for licensing issues.
@@ -95,6 +95,7 @@ tree.build <- function(tp)
     if (any(nzchar(node.label))) obj$node.label <- node.label
     if (!is.na(root.edge)) obj$root.edge <- root.edge
     class(obj) <- "phylo"
+    attr(obj, "order") <- "cladewise"
     obj
 }
 
index 0df359677a0a78678e6b7b8dfe59c8a38d68a345..eb00b60c087212a8905657ac5a6aa40b3fb5618c 100644 (file)
--- a/R/rtree.R
+++ b/R/rtree.R
@@ -1,4 +1,4 @@
-## rtree.R (2012-02-09)
+## rtree.R (2012-09-14)
 
 ##   Generates Trees
 
@@ -100,6 +100,7 @@ rtree <- function(n, rooted = TRUE, tip.label = NULL, br = runif, ...)
     }
     phy$Nnode <- n - 2L + rooted
     class(phy) <- "phylo"
+    attr(phy, "order") <- "cladewise"
     phy
 }
 
index 838a3a3b706fd7851cf7f8e41a592e0df06d6638..c923093a662be4fd910b4e55874e946afb12a773 100644 (file)
@@ -13,8 +13,6 @@
 \alias{circular.plot}
 \alias{unrooted.plot}
 \alias{unrooted.xy}
-\alias{node.height}
-\alias{node.height.clado}
 \alias{birth.step}
 \alias{death.step}
 \alias{ht.move}
index ad2e6206c99a7adc1442eabc728d447542461ac0..4f5ba15b62c7c237fe6d8cc00f28e41aa8b48d14 100644 (file)
@@ -43,11 +43,12 @@ dbdTime(x, birth, death, t, conditional = FALSE,
   \code{dbdTime} is for time-varying \code{lambda} and \code{mu}
   specified as \R functions.
 
-  Only \code{dyule} is vectorized simultaneously on its three arguments
+  \code{dyule} is vectorized simultaneously on its three arguments
   \code{x}, \code{lambda}, and \code{t}, according to \R's rules of
-  recycling arguments. The two others are vectorized only on \code{x};
-  the other arguments are eventually shortened with a warning if
-  necessary.
+  recycling arguments. \code{dbd} is vectorized simultaneously \code{x}
+  and \code{t} (to make likelihood calculations easy), and
+  \code{dbdTime} is vectorized only on \code{x}; the other arguments are
+  eventually shortened with a warning if necessary.
 
   The returned value is, logically, zero for values of \code{x} out of
   range, i.e., negative or zero for \code{dyule} or if \code{conditional
index aced56a8c4952757ddc07958b0546ee4dab8fbdb..e4c7bf9bcecf6b812634ef14ad3c8d9aa5458098 100644 (file)
@@ -1,19 +1,31 @@
 \name{node.depth}
 \alias{node.depth}
-\title{Depth of Nodes and Tips}
+\alias{node.depth.length}
+\alias{node.height}
+\alias{node.height.clado}
+\title{Depth and Heights of Nodes and Tips}
 \description{
-  This function returns the depth of nodes and tips given by the number
-  of descendants (1 is returned for tips).
+  These functions return the depth or height of nodes and tips.
 }
 \usage{
 node.depth(phy)
+node.depth.length(phy)
+node.height(phy)
+node.height.clado(phy)
 }
 \arguments{
   \item{phy}{an object of class "phylo".}
 }
 \details{
-  The depth of a node is computed as the number of tips which are its
-  descendants. The value of 1 is given to the tips.
+  \code{node.depth} computes the depth of a node as the number of tips
+  which are its descendants. The value of 1 is given to the tips.
+
+  \code{node.depth.length} does the same but using branch lengths.
+
+  \code{node.height} computes the heights of nodes and tips as plotted
+  by a phylogram.
+
+  \code{node.height.clado} does the same but for a cladogram.
 }
 \value{
   A numeric vector indexed with the node numbers of the matrix `edge' of
index 9e96f9ef4b3eb548edf838403a16e200884aed12..2251c8594bb5765ecd589c8103e31a4699d94204 100644 (file)
@@ -73,7 +73,7 @@ rTraitCont(phy, model = "BM", sigma = 0.1, alpha = 1, theta = 0,
 data(bird.orders)
 rTraitCont(bird.orders) # BM with sigma = 0.1
 ### OU model with two optima:
-tr <- reorder(bird.orders, "p")
+tr <- reorder(bird.orders, "postorder")
 plot(tr)
 edgelabels()
 theta <- rep(0, Nedge(tr))
index d8666a20c83b156ba84fb2969abc124317fa5ae1..9b38d2a7f72a54b96dd5c2d859f8470e8b05449b 100644 (file)
@@ -1,6 +1,6 @@
-/* plot_phylo.c (2011-06-23) */
+/* plot_phylo.c (2012-10-01) */
 
-/* Copyright 2004-2011 Emmanuel Paradis
+/* Copyright 2004-2012 Emmanuel Paradis
 
 /* This file is part of the R-package `ape'. */
 /* See the file ../COPYING for licensing issues. */
@@ -35,7 +35,7 @@ void node_depth(int *ntip, int *nnode, int *edge1, int *edge2,
 }
 
 void node_height(int *ntip, int *nnode, int *edge1, int *edge2,
-               int *nedge, double *yy)
+                int *nedge, double *yy)
 {
     int i, n;
     double S;
@@ -77,18 +77,3 @@ void node_height_clado(int *ntip, int *nnode, int *edge1, int *edge2,
        }
     }
 }
-
-void get_single_index_integer(int *x, int *val, int *index)
-{
-       while (x[*index] != *val) (*index)++;
-       *index += 1;
-}
-
-void get_two_index_integer(int *x, int *val, int *index)
-{
-       while (x[index[0]] != *val) index[0]++;
-       index[1] = index[0] + 1;
-       while (x[index[1]] != *val) index[1]++;
-       index[0] += 1;
-       index[1] += 1;
-}