Skip to content

Commit

Permalink
update to phylo.to.map to permit more than one observation per species
Browse files Browse the repository at this point in the history
  • Loading branch information
liamrevell committed Mar 8, 2019
1 parent 56c02b0 commit 9c8c0f0
Show file tree
Hide file tree
Showing 5 changed files with 37 additions and 24 deletions.
8 changes: 4 additions & 4 deletions DESCRIPTION
@@ -1,6 +1,6 @@
Package: phytools
Version: 0.6-70
Date: 2019-3-07
Version: 0.6-71
Date: 2019-3-08
Title: Phylogenetic Tools for Comparative Biology (and Other Things)
Author: Liam J. Revell
Maintainer: Liam J. Revell <liam.revell@umb.edu>
Expand All @@ -13,6 +13,6 @@ ZipData: no
Description: A wide range of functions for phylogenetic analysis. Functionality is concentrated in phylogenetic comparative biology, but also includes a diverse array of methods for visualizing, manipulating, reading or writing, and even inferring phylogenetic trees and data. Included among the functions in phylogenetic comparative biology are various for ancestral state reconstruction, model-fitting, simulation of phylogenies and data, and multivariate analysis. There are a broad range of plotting methods for phylogenies and comparative data which include, but are not restricted to, methods for mapping trait evolution on trees, for projecting trees into phenotypic space or a geographic map, and for visualizing correlated speciation between trees. Finally, there are a number of functions for reading, writing, analyzing, inferring, simulating, and manipulating phylogenetic trees and comparative data not covered by other packages. For instance, there are functions for randomly or non-randomly attaching species or clades to a phylogeny, for estimating supertrees or consensus phylogenies from a set, for simulating trees and phylogenetic data under a range of models, and for a wide variety of other manipulations and analyses that phylogenetic biologists might find useful in their research.
License: GPL (>= 2)
URL: http://github.com/liamrevell/phytools
Packaged: 2019-3-07 12:00:00 EST
Packaged: 2019-3-08 12:00:00 EST
Repository:
Date/Publication: 2019-3-07 12:00:00 EST
Date/Publication: 2019-3-08 12:00:00 EST
4 changes: 2 additions & 2 deletions NAMESPACE
Expand Up @@ -156,7 +156,7 @@ S3method(plot, fitpolyMk)
importFrom(animation, ani.options, ani.record, ani.replay, saveVideo)
importFrom(ape, .PlotPhyloEnv, .uncompressTipLabel, ace, all.equal.phylo, as.DNAbin, as.phylo, bind.tree, branching.times, collapse.singles)
importFrom(ape, consensus, compute.brlen, cophenetic.phylo, corBrownian, di2multi, dist.dna, dist.nodes, drop.tip, edgelabels, extract.clade)
importFrom(ape, getMRCA, is.binary.phylo)
importFrom(ape, getMRCA, is.binary)
importFrom(ape, is.monophyletic, is.rooted, is.ultrametric, ladderize, matexpo, mrca, multi2di)
importFrom(ape, nodelabels, Ntip, pic, plot.phylo, prop.part, read.tree, reorder.phylo, root, rotate, rtree, stree, tiplabels)
importFrom(ape, unroot, vcv, vcv.phylo, write.tree)
Expand All @@ -170,7 +170,7 @@ importFrom(scatterplot3d, scatterplot3d)
importFrom(stats, cophenetic, reorder, runif, setNames, optim, dexp, dnorm, rnorm, var, nlminb, biplot, pchisq, rexp, optimize, logLik, as.dist, pnorm, median, cor, aggregate, rgamma, dgamma, lm, rbinom, rgeom, pt, anova, p.adjust, screeplot, rchisq, optimHess, rmultinom, sd, dunif, dist, plogis, density, coef, cmdscale)
importFrom(methods, hasArg)
importFrom(graphics, strwidth, par, segments, locator, lines, text, strheight, symbols, plot, layout, plot.new, title, axis, points, plot.window, polygon, rect, curve, image, mtext, arrows, barplot, boxplot, contour, hist, abline, legend, grid)
importFrom(utils, flush.console, installed.packages, str, capture.output)
importFrom(utils, flush.console, head, installed.packages, str, capture.output)
importFrom(grDevices, palette, colors, dev.hold, dev.flush, rainbow, heat.colors, gray, colorRampPalette, dev.new, rgb,
pdf, dev.off, col2rgb)
importFrom(combinat, permn)
Expand Down
45 changes: 29 additions & 16 deletions R/phylo.to.map.R
@@ -1,5 +1,5 @@
## function depends on phytools (& dependencies) and maps (& dependencies)
## written by Liam J. Revell 2013, 2017
## written by Liam J. Revell 2013, 2017, 2019

phylo.to.map<-function(tree,coords,rotate=TRUE,...){
if(!inherits(tree,"phylo")) stop("tree should be an object of class \"phylo\".")
Expand All @@ -20,15 +20,19 @@ phylo.to.map<-function(tree,coords,rotate=TRUE,...){
if(hasArg(type)) type<-list(...)$type else type<-"phylogram"
if(hasArg(direction)) direction<-list(...)$direction else direction<-"downwards"
if(is.data.frame(coords)) coords<-as.matrix(coords)
if(rotate&&type=="phylogram") tree<-minRotate(tree,coords[,if(direction=="rightwards") 1 else 2])
if(rotate&&type=="phylogram"){
cc<-aggregate(coords,by=list(rownames(coords)),mean)
cc<-matrix(as.matrix(cc[,2:3]),nrow(cc),2,dimnames=list(cc[,1],colnames(cc)[2:3]))
tree<-minRotate(tree,cc[,if(direction=="rightwards") 1 else 2])
}
x<-list(tree=tree,map=map,coords=coords)
class(x)<-"phylo.to.map"
if(plot) plot.phylo.to.map(x,...)
invisible(x)
}

# function to plot object of class "phylo.to.map"
# written by Liam J. Revell 2013, 2014, 2016
## S3 method to plot object of class "phylo.to.map"
## written by Liam J. Revell 2013, 2014, 2016, 2019

plot.phylo.to.map<-function(x,type=c("phylogram","direct"),...){
type<-type[1]
Expand Down Expand Up @@ -68,6 +72,11 @@ plot.phylo.to.map<-function(x,type=c("phylogram","direct"),...){
if(length(colors)==2&&type=="phylogram"){
colors<-matrix(rep(colors,nrow(coords)),nrow(coords),2,byrow=TRUE)
rownames(colors)<-rownames(coords)
} else if(is.vector(colors)&&(length(colors)==Ntip(tree))) {
cat("MAKE IT HERE?\n")
COLS<-matrix("red",nrow(coords),2,dimnames=list(rownames(coords)))
for(i in 1:length(colors)) COLS[which(rownames(COLS)==names(colors)[i]),1:2]<-colors[i]
colors<-COLS
}
if(hasArg(direction)) direction<-list(...)$direction
else direction<-"downwards"
Expand Down Expand Up @@ -126,12 +135,14 @@ plot.phylo.to.map<-function(x,type=c("phylogram","direct"),...){
# compute start & end points of each edge
Y<-ylim[2]-nodeHeights(cw)
# plot coordinates & linking lines
coords<-coords[cw$tip.label,2:1]
for(i in 1:n) lines(c(x[i],coords[i,1]),
c(Y[which(cw$edge[,2]==i),2]-
if(from.tip) 0 else sh[i],coords[i,2]),
col=colors[cw$tip.label,][i,1],lty=lty,lwd=lwd[2])
points(coords,pch=pch,cex=cex.points[2],bg=colors[cw$tip.label,2])
coords<-coords[,2:1]
for(i in 1:nrow(coords)){
tip.i<-which(cw$tip.label==rownames(coords)[i])
lines(c(x[tip.i],coords[i,1]),c(Y[which(cw$edge[,2]==tip.i),2]-
if(from.tip) 0 else sh[tip.i],coords[i,2]),
col=colors[i,1],lty=lty,lwd=lwd[2])
}
points(coords,pch=pch,cex=cex.points[2],bg=colors[,2])
# plot vertical edges
for(i in 1:nrow(Y)) lines(rep(x[cw$edge[i,2]],2),Y[i,],
lwd=lwd[1],lend=2)
Expand Down Expand Up @@ -176,11 +187,13 @@ plot.phylo.to.map<-function(x,type=c("phylogram","direct"),...){
}
H<-nodeHeights(cw)
X<-xlim[1]+H
coords<-coords[cw$tip.label,2:1]
for(i in 1:n) lines(c(X[which(cw$edge[,2]==i),2]+
if(from.tip) 0 else sh[i],coords[i,1]),
c(y[i],coords[i,2]),col=colors[cw$tip.label,][i,1],lty=lty,lwd=lwd[2])
points(coords,pch=pch,cex=cex.points[2],bg=colors[cw$tip.label,2])
coords<-coords[,2:1]
for(i in 1:nrow(coords)){
tip.i<-which(cw$tip.label==rownames(coords)[i])
lines(c(X[which(cw$edge[,2]==tip.i),2]+if(from.tip) 0 else sh[tip.i],coords[i,1]),
c(y[tip.i],coords[i,2]),col=colors[i,1],lty=lty,lwd=lwd[2])
}
points(coords,pch=pch,cex=cex.points[2],bg=colors[,2])
for(i in 1:nrow(X)) lines(X[i,],rep(y[cw$edge[i,2]],2),lwd=lwd[1],lend=2)
for(i in 1:cw$Nnode+n) lines(X[which(cw$edge[,1]==i),1],
range(y[cw$edge[which(cw$edge[,1]==i),2]]),lwd=lwd[1],lend=2)
Expand Down Expand Up @@ -235,5 +248,5 @@ print.phylo.to.map<-function(x,...){
cat("(2) A geographic map with range:\n")
cat(paste(" ",paste(round(x$map$range[3:4],2),collapse="N, "),"N\n",sep=""))
cat(paste(" ",paste(round(x$map$range[1:2],2),collapse="W, "),"W.\n\n",sep=""))
cat(paste("(3) A table containing",nrow(x$coords),"geographic coordinates.\n\n"))
cat(paste("(3) A table containing",nrow(x$coords),"geographic coordinates (may include more than one set per species).\n\n"))
}
2 changes: 1 addition & 1 deletion man/phylo.impute.Rd
Expand Up @@ -13,7 +13,7 @@ phylo.impute(tree, X, ...)
This function performs phylogenetic imputation using Maximum Likelihood.
}
\details{
This function performs phylogenetic imputation in which the evolution of the characters in \code{X} is assumed to have occured by correlation multivariate Brownian motion. Missing values are imputed by maximizing their likelihood jointly with the parameters of the Brownian model. The function \code{\link{evol.vcv}} is used internally to compute the likelihood. Note that the \emph{Rphylopars} package (\link{https://CRAN.R-project.org/package=Rphylopars}}) also does phylogenetic imputation for multivariate trait data and it seems to be much faster.
This function performs phylogenetic imputation in which the evolution of the characters in \code{X} is assumed to have occured by correlation multivariate Brownian motion. Missing values are imputed by maximizing their likelihood jointly with the parameters of the Brownian model. The function \code{\link{evol.vcv}} is used internally to compute the likelihood. Note that the \emph{Rphylopars} package (\url{https://CRAN.R-project.org/package=Rphylopars}) also does phylogenetic imputation for multivariate trait data and it seems to be much faster.
}
\value{
An object of class \code{"phylo.impute"} consisting of a complete data frame with missing values imputed.
Expand Down
2 changes: 1 addition & 1 deletion man/phylo.to.map.Rd
Expand Up @@ -8,7 +8,7 @@ phylo.to.map(tree, coords, rotate=TRUE, ...)
}
\arguments{
\item{tree}{an object of class "phylo".}
\item{coords}{a matrix containing the latitude (in column 1) and the longitude of all tip species in the tree. The row names should be the same as \code{tree$tip.label}.}
\item{coords}{a matrix containing the latitude (in column 1) and the longitude of all tip species in the tree. The row names should be the same as \code{tree$tip.label}; however, more than one set of coordinates per species can be supplied by duplicating some row names.}
\item{rotate}{a logical value indicating whether or not to rotate nodes of the tree to better match longitudinal positions.}
\item{x}{for \code{plot.phylo.to.map}, an object of class \code{"phylo.to.map"}.}
\item{type}{a string indicating whether to map the tips of the tree onto a geographic map from a square phylogram (\code{type="phylogram"}) or to project the tree directly onto the map (\code{type="direct"}).}
Expand Down

0 comments on commit 9c8c0f0

Please sign in to comment.