]> git.donarmstrong.com Git - ape.git/commitdiff
some modifs for ape 3.0-8
authorparadis <paradis@6e262413-ae40-0410-9e79-b911bd7a66b7>
Thu, 31 Jan 2013 11:58:31 +0000 (11:58 +0000)
committerparadis <paradis@6e262413-ae40-0410-9e79-b911bd7a66b7>
Thu, 31 Jan 2013 11:58:31 +0000 (11:58 +0000)
git-svn-id: https://svn.mpl.ird.fr/ape/dev/ape@206 6e262413-ae40-0410-9e79-b911bd7a66b7

DESCRIPTION
NAMESPACE
NEWS
R/ace.R
R/read.dna.R
man/ace.Rd
man/read.dna.Rd

index bed8b3da3e33c96b41a999a5c53de11fda480683..925db387bb5f600162ed55aafeed6331b60fe54c 100644 (file)
@@ -1,12 +1,12 @@
 Package: ape
-Version: 3.0-7
-Date: 2013-01-19
+Version: 3.0-8
+Date: 2013-01-31
 Title: Analyses of Phylogenetics and Evolution
 Author: Emmanuel Paradis, Ben Bolker, Julien Claude, Hoa Sien Cuong, Richard Desper, Benoit Durand, Julien Dutheil, Olivier Gascuel, Christoph Heibl, Daniel Lawson, Vincent Lefort, Pierre Legendre, Jim Lemon, Yvonnick Noel, Johan Nylander, Rainer Opgen-Rhein, Andrei-Alin Popescu, Klaus Schliep, Korbinian Strimmer, Damien de Vienne
 Maintainer: Emmanuel Paradis <Emmanuel.Paradis@ird.fr>
 Depends: R (>= 2.6.0)
-Suggests: gee
-Imports: gee, nlme, lattice
+Suggests: gee, expm
+Imports: gee, nlme, lattice, expm
 ZipData: no
 Description: ape provides functions for reading, writing, plotting, and manipulating phylogenetic trees, analyses of comparative data in a phylogenetic framework, ancestral character analyses, analyses of diversification and macroevolution, computing distances from allelic and nucleotide data, reading and writing nucleotide sequences, and several tools such as Mantel's test, minimum spanning tree, generalized skyline plots, graphical exploration of phylogenetic data (alex, trex, kronoviz), estimation of absolute evolutionary rates and clock-like trees using mean path lengths and penalized likelihood. Phylogeny estimation can be done with the NJ, BIONJ, ME, MVR, SDM, and triangle methods, and several methods handling incomplete distance matrices (NJ*, BIONJ*, MVR*, and the corresponding triangle method). Some functions call external applications (PhyML, Clustal, T-Coffee, Muscle) whose results are returned into R.
 License: GPL (>= 2)
index 6951d4220e1129fbdf6f8cf220e533cee2853e34..08f08996e6a84b74d053c6a599df7b790300b02a 100644 (file)
--- a/NAMESPACE
+++ b/NAMESPACE
@@ -5,6 +5,7 @@ exportPattern(".+")
 import(gee, nlme)
 importFrom(lattice, xyplot, panel.lines, panel.points)
 importFrom(stats, as.hclust, cophenetic, reorder)
+importFrom(expm, expm)
 
 S3method(print, phylo)
 S3method(plot, phylo)
diff --git a/NEWS b/NEWS
index 6ca4ae5ac27602474b382f7c026c3d1a533c008e..5f1d7589aafb9c97738201ad10ea1874273099a7 100644 (file)
--- a/NEWS
+++ b/NEWS
@@ -1,3 +1,18 @@
+               CHANGES IN APE VERSION 3.0-8
+
+
+NEW FEATURES
+
+    o ace() gains a new option 'use.expm' to use expm() from the
+      package of the same name in place of matexpo().
+
+
+BUG FIXES
+
+    o read.dna(, "fasta") may add '\r' in labels: this is fixed.
+
+
+
                CHANGES IN APE VERSION 3.0-7
 
 
diff --git a/R/ace.R b/R/ace.R
index 55e32c35ebc1e3c2580ce0af1be4d04650c53de4..2b9230ec7c9d43e7063ea9add8f24e8cb87a88e3 100644 (file)
--- a/R/ace.R
+++ b/R/ace.R
@@ -1,8 +1,8 @@
-## ace.R (2012-06-25)
+## ace.R (2013-01-31)
 
 ##   Ancestral Character Estimation
 
-## Copyright 2005-2012 Emmanuel Paradis and Ben Bolker
+## Copyright 2005-2013 Emmanuel Paradis and Ben Bolker
 
 ## This file is part of the R-package `ape'.
 ## See the file ../COPYING for licensing issues.
@@ -21,7 +21,8 @@
 
 ace <- function(x, phy, type = "continuous", method = "ML", CI = TRUE,
                 model = if (type == "continuous") "BM" else "ER",
-                scaled = TRUE, kappa = 1, corStruct = NULL, ip = 0.1)
+                scaled = TRUE, kappa = 1, corStruct = NULL, ip = 0.1,
+                use.expm = FALSE)
 {
     if (!inherits(phy, "phylo"))
         stop('object "phy" is not of class "phylo"')
@@ -190,6 +191,8 @@ ace <- function(x, phy, type = "continuous", method = "ML", CI = TRUE,
         liks[cbind(TIPS, x)] <- 1
         phy <- reorder(phy, "pruningwise")
 
+        E <- if (use.expm) expm::expm else ape::matexpo
+
         Q <- matrix(0, nl, nl)
         dev <- function(p, output.liks = FALSE) {
             if (any(is.nan(p)) || any(is.infinite(p))) return(1e50)
@@ -202,8 +205,8 @@ ace <- function(x, phy, type = "continuous", method = "ML", CI = TRUE,
                 anc <- phy$edge[i, 1]
                 des1 <- phy$edge[i, 2]
                 des2 <- phy$edge[j, 2]
-                v.l <- matexpo(Q * phy$edge.length[i]) %*% liks[des1, ]
-                v.r <- matexpo(Q * phy$edge.length[j]) %*% liks[des2, ]
+                v.l <- E(Q * phy$edge.length[i]) %*% liks[des1, ]
+                v.r <- E(Q * phy$edge.length[j]) %*% liks[des2, ]
                 v <- v.l * v.r
                 comp[anc] <- sum(v)
                 liks[anc, ] <- v/comp[anc]
index f898f3bc8a2df3c28dc3697d8edf18edb70b76f6..fb8b54b7c7cc70746686ce6b5b8c8272b06f3509 100644 (file)
@@ -1,8 +1,8 @@
-## read.dna.R (2013-01-04)
+## read.dna.R (2013-01-31)
 
 ##   Read DNA Sequences in a File
 
-## Copyright 2003-2012 Emmanuel Paradis
+## Copyright 2003-2013 Emmanuel Paradis
 
 ## This file is part of the R-package `ape'.
 ## See the file ../COPYING for licensing issues.
@@ -11,10 +11,8 @@ read.FASTA <- function(file)
 {
     sz <- file.info(file)$size
     x <- readBin(file, "raw", sz)
-    if (Sys.info()[1] == "Windows") {
-        icr <- which(x == as.raw(0x0d)) # CR
-        x <- x[-icr]
-    }
+    icr <- which(x == as.raw(0x0d)) # CR
+    x <- x[-icr]
     res <- .Call("rawStreamToDNAbin", x, PACKAGE = "ape")
     names(res) <- sub("^ +", "", names(res)) # to permit phylosim
     class(res) <- "DNAbin"
index 15bd38758e42de2c7ede2fb709125350ced28b62..101aebf2c4bb61aa863c1bffe9d8fcc1aa4964bc 100644 (file)
@@ -25,7 +25,8 @@
 \usage{
 ace(x, phy, type = "continuous", method = "ML", CI = TRUE,
     model = if (type == "continuous") "BM" else "ER",
-    scaled = TRUE, kappa = 1, corStruct = NULL, ip = 0.1)
+    scaled = TRUE, kappa = 1, corStruct = NULL, ip = 0.1,
+    use.expm = FALSE)
 \method{print}{ace}(x, digits = 4, ...)
 \method{logLik}{ace}(object, ...)
 \method{deviance}{ace}(object, ...)
@@ -56,6 +57,10 @@ ace(x, phy, type = "continuous", method = "ML", CI = TRUE,
     structure to be used (this also gives the assumed model).}
   \item{ip}{the initial value(s) used for the ML estimation procedure
     when \code{type == "discrete"} (possibly recycled).}
+  \item{use.expm}{a logical specifying whether to use the package
+    \pkg{expm} to compute the matrix exponential (relevant only if
+    \code{type = "d"}). The default is to use the function
+    \code{matexpo} from \pkg{ape} (see details).}
   \item{digits}{the number of digits to be printed.}
   \item{object}{an object of class \code{"ace"}.}
   \item{k}{a numeric value giving the penalty per estimated parameter;
@@ -108,6 +113,14 @@ ace(x, phy, type = "continuous", method = "ML", CI = TRUE,
   \code{"SYM"} is a symmetrical model (e.g., \code{matrix(c(0, 1, 2, 1,
     0, 3, 2, 3, 0), 3)}). If a short-cut is used, the number of states
   is determined from the data.
+
+  With discrete characters it is necessary to compute the exponential of
+  the rate matrix. By default (and the only possible choice until
+  \pkg{ape} 3.0-7) the function \code{\link{matexpo}} in \pkg{ape} is
+  used. If \code{use.expm = TRUE}, the function
+  \code{\link[expm]{expm}}, in the package of the same name, is
+  used. \code{matexpo} is faster but quite inaccurate for large and/or
+  asymmetric matrices. In case of doubt, use the latter.
 }
 \value{
   an object of class \code{"ace"} with the following elements:
index 4cf664ea68b2292797bae09f7c6934827d5fd45a..ca0a5c03f1c1fd9828a2c2a38bbd4e3c9d877ef9 100644 (file)
@@ -129,11 +129,11 @@ cat("CLUSTAL (ape) multiple sequence alignment", "",
 file = "exdna.txt", sep = "\n")
 ex.dna3 <- read.dna("exdna.txt", format = "clustal")
 ### ... and in FASTA format
-cat("> No305",
+cat(">No305",
 "NTTCGAAAAACACACCCACTACTAAAANTTATCAGTCACT",
-"> No304",
+">No304",
 "ATTCGAAAAACACACCCACTACTAAAAATTATCAACCACT",
-"> No306",
+">No306",
 "ATTCGAAAAACACACCCACTACTAAAAATTATCAATCACT",
 file = "exdna.txt", sep = "\n")
 ex.dna4 <- read.dna("exdna.txt", format = "fasta")