PGLS_fossil | R Documentation |
The function performs pgls for non-ultrametric trees using a
variety of evolutionary models or RRphylo
rates to change the
tree correlation structure.
PGLS_fossil(modform,data=NULL,tree=NULL,RR=NULL,GItransform=FALSE,
intercept=FALSE,model="BM",method=NULL,permutation="none",...)
modform |
the formula for the regression model. |
data |
a data.frame or list including response and predictor variables
as named in |
tree |
a phylogenetic tree to be indicated for any model except if
|
RR |
the result of |
GItransform |
logical indicating whether the PGLS approach as in Garland and Ives (2000) must be applied. |
intercept |
under |
model |
an evolutionary model as indicated in
|
method |
optional argument to be passed to |
permutation |
passed to |
... |
further argument passed to |
The function is meant to work with both univariate or multivariate
data, both low- or high-dimensional. In the first case, PGLS_fossil
uses phylolm
, returning an object of class "phylolm". In the latter,
regression coefficients are estimated by mvgls
, and statistical
significance is obtained by means of permutations within manova.gls
.
In this case, PGLS_fossil
returns a list including the output of
both analyses. In all cases, for both univariate or multivariate data, if
GItransform = TRUE
the functions returns a standard lm
output. In the latter case, the output additionally includes the result of
manova applied on the multivariate linear model.
Fitted pgls parameters and significance. The class of the output object depends on input data (see details).
Silvia Castiglione, Pasquale Raia, Carmela Serio, Gabriele Sansalone, Giorgia Girardi
Garland, Jr, T., & Ives, A. R. (2000). Using the past to predict the present: confidence intervals for regression equations in phylogenetic comparative methods. The American Naturalist, 155: 346-364. doi:10.1086/303327
RRphylo
vignette;
mvgls
; manova.gls
;phylolm
## Not run:
library(ape)
library(phytools)
cc<- 2/parallel::detectCores()
rtree(100)->tree
fastBM(tree)->resp
fastBM(tree,nsim=3)->resp.multi
fastBM(tree)->pred1
fastBM(tree)->pred2
PGLS_fossil(modform=resp~pred1+pred2,tree=tree)->pgls_noRR
PGLS_fossil(modform=resp~pred1+pred2,tree=tree,GItransform=TRUE)->GIpgls_noRR
RRphylo(tree,resp,clus=cc)->RR
PGLS_fossil(modform=resp~pred1+pred2,tree=tree,RR=RR)->pgls_RR
PGLS_fossil(modform=resp~pred1+pred2,tree=tree,RR=RR,GItransform=TRUE)->GIpgls_RR
# To derive log-likelihood and AIC for outputs of PGLS_fossil applied on univariate
# response variables the function AIC can be applied
AIC(pgls_noRR)
AIC(pgls_RR)
AIC(GIpgls_noRR)
AIC(GIpgls_RR)
PGLS_fossil(modform=resp.multi~pred1+pred2,tree=tree)->pgls2_noRR
PGLS_fossil(modform=resp.multi~pred1+pred2,tree=tree,GItransform=TRUE)->GIpgls2_noRR
# To evaluate statistical significance of multivariate models, the '$manova'
# object must be inspected
pgls2_noRR$manova
summary(GIpgls2_noRR$manova)
RRphylo(tree,resp.multi,clus=cc)->RR
PGLS_fossil(modform=resp.multi~pred1+pred2,tree=tree,RR=RR)->pgls2_RR
PGLS_fossil(modform=resp.multi~pred1+pred2,tree=tree,RR=RR,GItransform=TRUE)->GIpgls2_RR
# To evaluate statistical significance of multivariate models, the '$manova'
# object must be inspected
pgls2_noRR$manova
summary(GIpgls2_noRR$manova)
pgls2_RR$manova
summary(GIpgls2_RR$manova)
logLik(pgls2_noRR$pgls)
logLik(pgls2_RR$pgls)
## End(Not run)
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.