]> git.donarmstrong.com Git - ape.git/commitdiff
added corPagel + gave version # 2.2
authorparadis <paradis@6e262413-ae40-0410-9e79-b911bd7a66b7>
Fri, 18 Apr 2008 10:02:18 +0000 (10:02 +0000)
committerparadis <paradis@6e262413-ae40-0410-9e79-b911bd7a66b7>
Fri, 18 Apr 2008 10:02:18 +0000 (10:02 +0000)
git-svn-id: https://svn.mpl.ird.fr/ape/dev/ape@27 6e262413-ae40-0410-9e79-b911bd7a66b7

Changes
DESCRIPTION
R/PGLS.R
man/ape-internal.Rd
man/corClasses.Rd
man/corPagel.Rd [new file with mode: 0644]

diff --git a/Changes b/Changes
index 3eaf5b7a8c52a99ccc0cc4263ec0652c0e9dd5f3..c3b05a4d59dc0e509bce577a046b954c2417b080 100644 (file)
--- a/Changes
+++ b/Changes
@@ -1,4 +1,4 @@
-               CHANGES IN APE VERSION 2.1-4
+               CHANGES IN APE VERSION 2.2
 
 
 NEW FEATURES
@@ -8,6 +8,9 @@ NEW FEATURES
       subtreeplot), and to return the graphical coordinates of tree
       (without plotting).
 
+    o The new function corPagel() implements the Pagel's "lambda"
+      correlation structure.
+
 
 BUG FIXES
 
index 510eef47806c3244795109bc17c91057a7156f64..3817b831973d2d0e5a5028dddd5bd6ae3dee558f 100644 (file)
@@ -1,6 +1,6 @@
 Package: ape
-Version: 2.1-4
-Date: 2008-04-17
+Version: 2.2
+Date: 2008-04-18
 Title: Analyses of Phylogenetics and Evolution
 Author: Emmanuel Paradis, Ben Bolker, Julien Claude, Hoa Sien Cuong,
   Richard Desper, Benoit Durand, Julien Dutheil, Olivier Gascuel,
index 2bc3b03887ee788c5a1c0950df32211005c0dc70..1b918e9a9574097322794c5ff5b17fd9b5277638 100644 (file)
--- a/R/PGLS.R
+++ b/R/PGLS.R
@@ -1,8 +1,8 @@
-## PGLS.R (2006-10-12)
+## PGLS.R (2008-04-18)
 
 ##   Phylogenetic Generalized Least Squares
 
-## Copyright 2004 Julien Dutheil, and 2006 Emmanuel Paradis
+## Copyright 2004 Julien Dutheil, and 2006-2008 Emmanuel Paradis
 
 ## This file is part of the R-package `ape'.
 ## See the file ../COPYING for licensing issues.
@@ -196,3 +196,43 @@ compute.brlen <- function(phy, method = "Grafen", power = 1, ...)
         return(phy)
     }
 }
+
+### by EP:
+
+corPagel <- function(value, phy, form = ~1, fixed = FALSE)
+{
+    if (value < 0 || value > 1)
+        stop("the value of lambda must be between 0 and 1.")
+    if (!("phylo" %in% class(phy)))
+        stop('object "phy" is not of class "phylo"')
+    attr(value, "formula") <- form
+    attr(value, "fixed") <- fixed
+    attr(value, "tree") <- phy
+    class(value) <- c("corPagel", "corPhyl", "corStruct")
+    value
+}
+
+corMatrix.corPagel <-
+    function(object, covariate = getCovariate(object), corr = TRUE, ...)
+{
+    if (!any(attr(object, "index")))
+        stop("object has not been initialized")
+    mat <- vcv.phylo(attr(object, "tree"), cor = corr)
+    index <- attr(object, "index")
+    mat <- mat[index, index]
+    tmp <- diag(mat)
+    mat <- object[1]*mat
+    diag(mat) <- tmp
+    mat
+}
+
+coef.corPagel <- function(object, unconstrained = TRUE, ...)
+{
+    if (unconstrained) {
+        if (attr(object, "fixed")) return(numeric(0))
+        else return(object[1])
+    }
+    aux <- object[1]
+    names(aux) <- "lambda"
+    aux
+}
index 242afb5774b604760dcc92b7584922076fdee5f5..441d39f9b1fb8b522b35a072211ebd860d819a57 100644 (file)
@@ -44,6 +44,8 @@
 \alias{floating.pie.asp}
 \alias{checkLabel}
 \alias{getMRCA}
+\alias{plot.cophylo2}
+\alias{plot.phylo.coor}
 \title{Internal Ape Functions}
 \description{
   Internal ape functions.
index 1afb89f5571708593efae3d87c7bac3b122d8514..c7679be5ed228e968d664944b0916405e9e863da 100644 (file)
@@ -3,20 +3,26 @@
 \alias{corPhyl}
 \title{Phylogenetic Correlation Structures}
 \description{
-  Standard classes of phylogenetic correlation structures (\code{corPhyl}) available in \pkg{ape}.
-}
-\value{
-       Available standard classes:
-       \item{corBrownian}{Brownian model (Felsenstein 1985),}
-       \item{corMartins}{The covariance matrix defined in Martins and Hansen (1997),}
-       \item{corGrafen}{The covariance matrix defined in Grafen (1989).}
-       See classes documentation for reference and detailed description.
+  Classes of phylogenetic correlation structures (\code{"corPhyl"})
+  available in \pkg{ape}.
+
+  \item{corBrownian}{Brownian motion model (Felsenstein 1985)}
+  \item{corMartins}{The covariance matrix defined in Martins and Hansen
+    (1997)}
+  \item{corGrafen}{The covariance matrix defined in Grafen (1989)}
+  \item{corPagel}{The covariance matrix defined in Freckelton et al. (2002)}
+
+  See the help page of each class for references and detailed
+  description.
 }
 \seealso{
-       \code{\link[nlme]{corClasses}} and \code{\link[nlme]{gls}} in the \pkg{nlme} librarie,
-       \code{\link{corBrownian}}, \code{\link{corMartins}}, \code{\link{corGrafen}}.
+  \code{\link[nlme]{corClasses}} and \code{\link[nlme]{gls}} in the
+  \pkg{nlme} librarie, \code{\link{corBrownian}},
+  \code{\link{corMartins}}, \code{\link{corGrafen}},
+  \code{\link{corPagel}}
 }
-\author{Julien Dutheil \email{julien.dutheil@univ-montp2.fr}}
+\author{Julien Dutheil \email{julien.dutheil@univ-montp2.fr}, Emmanuel
+  Paradis}
 \examples{
 library(nlme)
 cat("((((Homo:0.21,Pongo:0.21):0.28,",
diff --git a/man/corPagel.Rd b/man/corPagel.Rd
new file mode 100644 (file)
index 0000000..7248f65
--- /dev/null
@@ -0,0 +1,50 @@
+\name{corPagel}
+\alias{corPagel}
+\alias{coef.corPagel}
+\alias{corMatrix.corPagel}
+\title{Pagel's ``lambda'' Correlation Structure}
+\usage{
+corPagel(value, phy, form = ~1, fixed = FALSE)
+\method{corMatrix}{corPagel}(object, covariate = getCovariate(object),
+                   corr = TRUE, ...)
+\method{coef}{corPagel}(object, unconstrained = TRUE, \dots)
+}
+\arguments{
+  \item{value}{the (initial) value of the parameter
+    \eqn{\gamma}{gamma}.}
+  \item{phy}{an object of class \code{"phylo"}.}
+  \item{form}{(ignored).}
+  \item{fixed}{a logical specifying whether \code{gls} should
+    estimate \eqn{\gamma}{gamma} (the default) or keep it fixed.}
+  \item{object}{an (initialized) object of class \code{"corPagel"}.}
+  \item{covariate}{(ignored).}
+  \item{corr}{a logical value specifying whether to return the
+    correlation matrix (the default) or the variance-covariance matrix.}
+  \item{unconstrained}{a logical value. If \code{TRUE} (the default),
+    the coefficients are returned in unconstrained form (the same used
+    in the optimization algorithm). If \code{FALSE} the coefficients are
+    returned in ``natural'', possibly constrained, form.}
+  \item{\dots}{further arguments passed to or from other methods.}
+}
+\description{
+  The correlation structure from the present model is derived from the
+  Brownian motion model by multiplying the off-diagonal elements (i.e.,
+  the covariances) by \eqn{\gamma}{gamma}. The variances are thus the
+  same than for a Brownian motion model.
+}
+\value{
+  an object of class \code{"corPagel"}, the coefficients from an object
+  of this class, or the correlation matrix of an initialized object of
+  this class. In most situations, only \code{corPagel()} will be called
+  by the user.
+}
+\author{Emmanuel Paradis}
+\references{
+  Freckleton, R. P., Harvey, P. H. and M. Pagel, M. (2002) Phylogenetic
+  analysis and comparative data: a test and review of evidence.
+  \emph{American Naturalist}, \bold{160}, 712--726.
+
+  Pagel, M. (1999) Inferring the historical patterns of biological
+  evolution. \emph{Nature}, \bold{401},877--884.
+}
+\keyword{models}