]> git.donarmstrong.com Git - ape.git/commitdiff
fix in read.dna and change in write.dna
authorparadis <paradis@6e262413-ae40-0410-9e79-b911bd7a66b7>
Thu, 3 May 2012 08:41:44 +0000 (08:41 +0000)
committerparadis <paradis@6e262413-ae40-0410-9e79-b911bd7a66b7>
Thu, 3 May 2012 08:41:44 +0000 (08:41 +0000)
git-svn-id: https://svn.mpl.ird.fr/ape/dev/ape@188 6e262413-ae40-0410-9e79-b911bd7a66b7

DESCRIPTION
NEWS
R/DNA.R
R/read.dna.R
R/write.dna.R
man/compar.cheverud.Rd
man/dist.topo.Rd
man/howmanytrees.Rd
man/pcoa.Rd
man/read.dna.Rd
man/varcomp.Rd

index 36e0a8e0088960e213bedb240ef62c32c1d0d5e2..0f031c13dacf5c207506a9c113cd3ea00e52140a 100644 (file)
@@ -1,6 +1,6 @@
 Package: ape
 Version: 3.0-4
-Date: 2012-04-30
+Date: 2012-05-03
 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>
diff --git a/NEWS b/NEWS
index 6d74722e49d33a07139128b1f6bb68be65f6e75e..e4235a9a0cf397de73f5f19215f6401525e98f96 100644 (file)
--- a/NEWS
+++ b/NEWS
@@ -3,17 +3,25 @@
 
 BUG FIXES
 
-    o read.dna() failed to read interleaved files if the first line
-      uses tabulations instead of white spaces.
+    o read.dna() failed to read Phylip files if the first line used
+      tabulations instead of white spaces.
+
+    o read.dna() failed to read Phylip or Clustal files with less than
+      10 nucleotides failed. (See other changes in this function
+      below.)
 
 OTHER CHANGES
 
+    o read.dna() now requires at least one space (or tab) between the
+      taxa names and the sequences (whatever the length of taxa
+      names). write.dna() now follows the same rule.
+
     o The files ape-defunct.R and ape-defunct.Rd, which have not been
       modified for almost two years, have been removed.
 
     o The C code of bionj() has been reworked: it is more stable (by
       avoiding passing character strings), slightly faster (by about
-      20%), and more accurate numerically.
+      20%), and numerically more accurate.
 
     o The C code of fastme.*() has been slightly modified and should
       be more stable by avoiding passing character strings (the
diff --git a/R/DNA.R b/R/DNA.R
index abf3a4343a33b4d456a6d1edfe549855249d55d0..70b210b14c9622212d05cbd775475e4b09d039dd 100644 (file)
--- a/R/DNA.R
+++ b/R/DNA.R
@@ -382,7 +382,7 @@ dist.dna <- function(x, model = "K80", variance = FALSE, gamma = FALSE,
     if (!pairwise.deletion) {
         keep <- .C("GlobalDeletionDNA", x, n, s,
                    rep(1L, s), PACKAGE = "ape")[[4]]
-        x <- x[,  as.logical(keep)]
+        x <- x[, as.logical(keep)]
         s <- dim(x)[2]
     }
     Ndist <- if (imod == 11) n*n else n*(n - 1)/2
index 3e2b980670d872529f78eae740215d5d1582c279..ca4da7a4daac444adbfab40c413d6ff08fa5f72d 100644 (file)
@@ -1,4 +1,4 @@
-## read.dna.R (2012-04-26)
+## read.dna.R (2012-05-03)
 
 ##   Read DNA Sequences in a File
 
@@ -8,9 +8,15 @@
 ## See the file ../COPYING for licensing issues.
 
 read.dna <- function(file, format = "interleaved", skip = 0,
-                     nlines = 0, comment.char = "#", seq.names = NULL,
+                     nlines = 0, comment.char = "#",
                      as.character = FALSE, as.matrix = NULL)
 {
+    findFirstNucleotide <- function(x) {
+        ## actually find the 1st non-blank character
+        ## just in case: pat.base <- "[-AaCcGgTtUuMmRrWwSsYyKkVvHhDdBbNn?]{10}"
+        tmp <- regexpr("[[:blank:]]+", x[1]) # consider only a single string
+        tmp[1] + attr(tmp, "match.length")
+    }
     getTaxaNames <- function(x) {
         x <- sub("^['\" ]+", "", x) # remove the leading quotes and spaces
         x <- sub("['\" ]+$", "", x) #   "     "  trailing  "     "    "
@@ -19,14 +25,14 @@ read.dna <- function(file, format = "interleaved", skip = 0,
     getNucleotide <- function(x) {
         x <- gsub(" ", "", x)
         x <- strsplit(x, NULL)
-        unlist(x)
+        tolower(unlist(x))
     }
     formats <- c("interleaved", "sequential", "fasta", "clustal")
     format <- match.arg(format, formats)
     phylip <- if (format %in% formats[1:2]) TRUE else FALSE
     X <- scan(file = file, what = "", sep = "\n", quiet = TRUE,
               skip = skip, nlines = nlines, comment.char = comment.char)
-    pat.base <- "[-AaCcGgTtUuMmRrWwSsYyKkVvHhDdBbNn?]{10}"
+
     if (phylip) {
         ## need to remove the possible leading spaces and/or tabs in the first line
         fl <- gsub("^[[:blank:]]+", "", X[1])
@@ -38,56 +44,53 @@ read.dna <- function(file, format = "interleaved", skip = 0,
         obj <- matrix("", n, s)
         X <- X[-1]
     }
-    if (format == "interleaved") {
-        start.seq <- regexpr(pat.base, X[1])[1]
-        if (is.null(seq.names))
-            seq.names <- getTaxaNames(substr(X[1:n], 1, start.seq - 1))
-        X[1:n] <- substr(X[1:n], start.seq, nchar(X[1:n]))
+    switch(format,
+    "interleaved" = {
+        start.seq <- findFirstNucleotide(X[1])
+        one2n <- 1:n
+        taxa <- getTaxaNames(substr(X[one2n], 1, start.seq - 1))
+        X[one2n] <- substr(X[one2n], start.seq, nchar(X[one2n]))
         nl <- length(X)
-        for (i in 1:n)
+        for (i in one2n)
             obj[i, ] <- getNucleotide(X[seq(i, nl, n)])
-    }
-    if (format == "sequential") {
+    },
+    "sequential" = {
         taxa <- character(n)
-        j <- 1
+        j <- 1L # line number
         for (i in 1:n) {
-            start.seq <- regexpr(pat.base, X[j])[1]
-            taxa[i] <- substr(X[j], 1, start.seq - 1)
+            start.seq <- findFirstNucleotide(X[j])
+            taxa[i] <- getTaxaNames(substr(X[j], 1, start.seq - 1))
             sequ <- getNucleotide(substr(X[j], start.seq, nchar(X[j])))
-            j <- j + 1
+            j <- j + 1L
             while (length(sequ) < s) {
                 sequ <- c(sequ, getNucleotide(X[j]))
-                j <- j + 1
+                j <- j + 1L
             }
             obj[i, ] <- sequ
         }
-        if (is.null(seq.names)) seq.names <- getTaxaNames(taxa)
-    }
-    if (format == "fasta") {
+        taxa <- getTaxaNames(taxa)
+    },
+    "fasta" = {
         start <- grep("^ {0,}>", X)
         taxa <- X[start]
+        taxa <- sub("^ {0,}>", "", taxa) # remove the hook and the spaces before
+        taxa <- getTaxaNames(taxa)
         n <- length(taxa)
         obj <- vector("list", n)
-        if (is.null(seq.names)) {
-            taxa <- sub("^ {0,}>", "", taxa) # remove the hook and the spaces before
-            seq.names <- getTaxaNames(taxa)
-        }
         start <- c(start, length(X) + 1) # this avoids the following to crash when `i = n'
         for (i in 1:n)
             obj[[i]] <- getNucleotide(X[(start[i] + 1):(start[i + 1] - 1)])
-    }
-    if (format == "clustal") {
-        X <- X[-1]
+    },
+    "clustal" = {
+        X <- X[-1] # drop the line with "Clustal bla bla..."
         ## find where the 1st sequence starts
-        start.seq <- regexpr(pat.base, X[1])[1]
+        start.seq <- findFirstNucleotide(X[1])
         ## find the lines with *********....
         nspaces <- paste("^ {", start.seq - 1, "}", sep = "", collapse = "")
         stars <- grep(nspaces, X)
         ## we now know how many sequences in the file:
         n <- stars[1] - 1
-        ## get the sequence names in the same way than "interleaved":
-        if (is.null(seq.names))
-            seq.names <- getTaxaNames(substr(X[1:n], 1, start.seq - 1))
+        taxa <- getTaxaNames(substr(X[1:n], 1, start.seq - 1))
         ## need to remove the sequence names before getting the sequences:
         X <- substr(X, start.seq, nchar(X))
         nl <- length(X)
@@ -98,13 +101,12 @@ read.dna <- function(file, format = "interleaved", skip = 0,
         obj[1, ] <- tmp
         for (i in 2:n)
             obj[i, ] <- getNucleotide(X[seq(i, nl, n + 1)])
-    }
+    })
+
     if (format != "fasta") {
-        rownames(obj) <- seq.names
-        obj <- tolower(obj)
+        rownames(obj) <- taxa
     } else {
-        names(obj) <- seq.names
-        obj <- lapply(obj, tolower)
+        names(obj) <- taxa
         LENGTHS <- unique(unlist(lapply(obj, length)))
         allSameLength <- length(LENGTHS) == 1
         if (is.logical(as.matrix)) {
@@ -115,7 +117,7 @@ read.dna <- function(file, format = "interleaved", skip = 0,
         }
         if (as.matrix) {
             obj <- matrix(unlist(obj), ncol = LENGTHS, byrow = TRUE)
-            rownames(obj) <- seq.names
+            rownames(obj) <- taxa
         }
     }
     if (!as.character) obj <- as.DNAbin(obj)
index eeedab5732219fa09d469c9495f9ed201747080e..47686db3a4b133877c0cfeb6ab85cff5c9962015 100644 (file)
@@ -1,8 +1,8 @@
-## write.dna.R (2009-05-10)
+## write.dna.R (2012-05-03)
 
 ##   Write DNA Sequences in a File
 
-## Copyright 2003-2009 Emmanuel Paradis
+## Copyright 2003-2012 Emmanuel Paradis
 
 ## This file is part of the R-package `ape'.
 ## See the file ../COPYING for licensing issues.
@@ -62,9 +62,7 @@ write.dna <- function(x, file, format = "interleaved", append = FALSE,
         }
         ## Prepare the names so that they all have the same nb of chars
         max.nc <- max(nchar(names(x)))
-        ## in case all names are 10 char long or less, sequences are
-        ## always started on the 11th column of the file
-        if (max.nc < 11) max.nc <- 9
+        ## always put a space between the sequences and the taxa names
         fmt <- paste("%-", max.nc + 1, "s", sep = "")
         names(x) <- sprintf(fmt, names(x))
     }
index d7729d021a60bbb2cffe21939123ce67a9850324..404b788135b6dbbe62ba3a000e7977a16c63e943 100644 (file)
@@ -45,8 +45,8 @@ compar.cheverud(y, W, tolerance = 1e-06, gold.tol = 1e-04)
   variables: geometric interpretations. \emph{Evolution}, \bold{55},
   2143--2160.
 
-  Harvey, P. H. and Pagel, M. D. (1991) \emph{The comparative method in
-    evolutionary biology}. Oxford University Press.
+  Harvey, P. H. and Pagel, M. D. (1991) \emph{The Comparative Method in
+    Evolutionary Biology}. Oxford University Press.
 }
 \author{Julien Dutheil \email{julien.dutheil@univ-montp2.fr}}
 \seealso{\code{\link{compar.lynch}}}
index caf1b5b6f4eaf5d1757cde044a7a843f9003aa68..d1542add7b9ac0e1233f11bbba209dd2305ff73f 100644 (file)
@@ -46,8 +46,8 @@ dist.topo(x, y, method = "PH85")
   phylogeny algorithms under equal and unequal evolutionary rates.
   \emph{Molecular Biology and Evolution}, \bold{11}, 459--468.
 
-  Nei, M. and Kumar, S. (2000) \emph{Molecular evolution and
-  phylogenetics}. Oxford: Oxford University Press.
+  Nei, M. and Kumar, S. (2000) \emph{Molecular Evolution and
+  Phylogenetics}. Oxford: Oxford University Press.
 
   Penny, D. and Hendy, M. D. (1985) The use of tree comparison
   metrics. \emph{Systemetic Zoology}, \bold{34}, 75--82.
index a576551a517e00d11567bdbd276b26e83697f058..23bc7e6e0d23b0257657c121dca2deabde9bf34b 100644 (file)
@@ -51,7 +51,7 @@ howmanytrees(n, rooted = TRUE, binary = TRUE,
   used, a named vector or matrix.
 }
 \references{
-  Felsenstein, J. (2004) \emph{Inferring phylogenies}. Sunderland:
+  Felsenstein, J. (2004) \emph{Inferring Phylogenies}. Sunderland:
   Sinauer Associates.
 }
 \author{Emmanuel Paradis}
index 8a9baac3fa3f4d9938a9a2dde3af1d1112c17b4c..c5e43eaac33e76b0e7da265a7d1abb94aaaf659f 100644 (file)
@@ -77,7 +77,7 @@ dissimilarity coefficients. \emph{Journal of Classification}, \bold{3}, 5--48.
 
 Legendre, P. and Gallagher, E. D. (2001) Ecologically meaningful transformations for ordination of species data. \emph{Oecologia}, \bold{129}, 271--280.
 
-Legendre, P. and Legendre, L. (1998) \emph{Numerical ecology, 2nd English
+Legendre, P. and Legendre, L. (1998) \emph{Numerical Ecology, 2nd English
   edition.} Amsterdam: Elsevier Science BV.
 
 Lingoes, J. C. (1971) Some boundary conditions for a monotone analysis of symmetric matrices. \emph{Psychometrika}, \bold{36}, 195--203.
index 3f350314155250032c3e5d461a1e4c13567bb337..2863ad85c32ec8c47efeca1fa852bf57a5a14315 100644 (file)
@@ -3,7 +3,7 @@
 \title{Read DNA Sequences in a File}
 \usage{
 read.dna(file, format = "interleaved", skip = 0,
-         nlines = 0, comment.char = "#", seq.names = NULL,
+         nlines = 0, comment.char = "#",
          as.character = FALSE, as.matrix = NULL)
 }
 \arguments{
@@ -19,8 +19,6 @@ read.dna(file, format = "interleaved", skip = 0,
     read untill its end).}
   \item{comment.char}{a single character, the remaining of the line
     after this character is ignored.}
-  \item{seq.names}{the names to give to each sequence; by default the
-    names read in the file are used.}
   \item{as.character}{a logical controlling whether to return the
     sequences as an object of class \code{"DNAbin"} (the default).}
   \item{as.matrix}{(used if \code{format = "fasta"}) one of the three
@@ -50,19 +48,15 @@ read.dna(file, format = "interleaved", skip = 0,
   way with blanks and line-breaks inside (with the restriction that the
   first ten nucleotides must be contiguous for the interleaved and
   sequential formats, see below). The names of the sequences are read in
-  the file unless the `seq.names' option is used. Particularities for
-  each format are detailed below.
+  the file. Particularities for each format are detailed below.
 
 \itemize{
-  \item{Interleaved:}{the function starts to read the sequences when it
-    finds 10 contiguous characters belonging to the ambiguity code of
-    the IUPAC (namely A, C, G, T, U, M, R, W, S, Y, K, V, H, D, B, and
-    N, upper- or lowercase, so you might run into trouble if you have a
-    taxa name with 10 contiguous letters among these!) All characters
-    before the sequences are taken as the taxa names after removing the
-    leading and trailing spaces (so spaces in a taxa name are
-    allowed). It is assumed that the taxa names are not repeated in the
-    subsequent blocks of nucleotides.}
+  \item{Interleaved:}{the function starts to read the sequences after it
+    finds one or more spaces (or tabulations). All characters before the
+    sequences are taken as the taxa names after removing the leading and
+    trailing spaces (so spaces in taxa names are allowed). It is assumed
+    that the taxa names are not repeated in the subsequent blocks of
+    nucleotides.}
 
   \item{Sequential:}{the same criterion than for the interleaved format
     is used to start reading the sequences and the taxa names; the
index 8eabec9c779ce46b4271a3e17412ac4e14005223..4bd30e7b8f4faffe1bb87997df918d95f7e60dad 100644 (file)
@@ -13,14 +13,14 @@ varcomp(x, scale = FALSE, cum = FALSE)
   \item{cum}{Send cumulative variance components.}
 }
 \details{
-  Variance computations is done as in Venables and Ripley's book.
+  Variance computations is done as in Venables and Ripley (2002).
 }
 \value{
   A named vector of class \code{varcomp} with estimated variance components.
 }
 \references{
-  Venables, W. N. and Ripley, B. D. (2002) \emph{Modern applied statistics
-  with S (fourth edition)}. New York: Springer-Verlag.
+  Venables, W. N. and Ripley, B. D. (2002) \emph{Modern Applied Statistics
+  with S (Fourth Edition)}. New York: Springer-Verlag.
 }
 \author{Julien Dutheil \email{julien.dutheil@univ-montp2.fr}}
 \seealso{\code{\link[nlme]{lme}}}