abouheif.moran {adephylo}R Documentation

Abouheif's test based on Moran's I

Description

The test of Abouheif (1999) is designed to detect phylogenetic autocorrelation in a quantitative trait. Pavoine et al. (2008) have shown that this tests is in fact a Moran's I test using a particular phylogenetic proximity between tips (see details). The function abouheif.moran performs basically Abouheif's test for several traits at a time, but it can incorporate other phylogenetic proximities as well.

Usage

abouheif.moran(
  x,
  W = NULL,
  method = c("oriAbouheif", "patristic", "nNodes", "Abouheif", "sumDD"),
  f = function(x) {
     1/x
 },
  nrepet = 999,
  alter = c("greater", "less", "two-sided")
)

Arguments

x

a data frame with continuous variables, or a phylo4d object (i.e. containing both a tree, and tip data). In the latter case, method argument is used to determine which proximity should be used.

W

a n by n matrix (n being the number rows in x) of phylogenetic proximities, as produced by proxTips.

method

a character string (full or unambiguously abbreviated) specifying the type of proximity to be used. By default, the proximity used is that of the original Abouheif's test. See details in proxTips for information about other methods.

f

a function to turn a distance into a proximity (see proxTips).

nrepet

number of random permutations of data for the randomization test

alter

a character string specifying the alternative hypothesis, must be one of "greater" (default), "less" or "two-sided"

Details

Note that the original Abouheif's proximity (Abouheif, 1999; Pavoine et al. 2008) unifies Moran's I and Geary'c tests (Thioulouse et al. 1995).

abouheif.moran can be used in two ways:
- providing a data.frame of traits (x) and a matrix of phylogenetic proximities (W)
- providing a phylo4d object (x) and specifying the type of proximity to be used (method).

W is a squared symmetric matrix whose terms are all positive or null.

W is firstly transformed in frequency matrix A by dividing it by the total sum of data matrix :

a_{ij} = \frac{W_{ij}}{\sum_{i=1}^{n}\sum_{j=1}^{n}W_{ij}}

The neighbouring weights is defined by the matrix D = diag(d_1,d_2, \ldots) where d_i = \sum_{j=1}^{n}W_{ij}. For each vector x of the data frame x, the test is based on the Moran statistic x^{t}Ax where x is D-centred.

Value

Returns an object of class krandtest (randomization tests from ade4), containing one Monte Carlo test for each trait.

Author(s)

Original code from ade4 (gearymoran function) by Sebastien Ollier
Adapted and maintained by Thibaut Jombart <tjombart@imperial.ac.uk>.

References

Thioulouse, J., Chessel, D. and Champely, S. (1995) Multivariate analysis of spatial patterns: a unified approach to local and global structures. Environmental and Ecological Statistics, 2, 1–14.

See Also

- gearymoran from the ade4 package
- Moran.I from the ape package for the classical Moran's I test.

Examples



if(require(ade4)&& require(ape) && require(phylobase)){
## load data
data(ungulates)
tre <- read.tree(text=ungulates$tre)
x <- phylo4d(tre, ungulates$tab)

## Abouheif's tests for each trait
myTests <- abouheif.moran(x)
myTests
plot(myTests)

## a variant using another proximity
plot(abouheif.moran(x, method="nNodes") )

## Another example

data(maples)
tre <- read.tree(text=maples$tre)
dom <- maples$tab$Dom

## Abouheif's tests for each trait (equivalent to Cmean)
W1 <- proxTips(tre,method="oriAbouheif")
abouheif.moran(dom,W1)

## Equivalence with moran.idx

W2 <- proxTips(tre,method="Abouheif")
abouheif.moran(dom,W2)
moran.idx(dom,W2) 
}


[Package adephylo version 1.1-16 Index]