View source: R/shift_configuration.R
configuration_ic | R Documentation |
Computes the information criterion score for a given configuration
configuration_ic(tree, Y, shift.configuration, criterion = c("pBIC", "pBICess", "mBIC", "BIC", "AICc"), root.model = c("OUfixedRoot", "OUrandomRoot"), alpha.starting.value = NA, alpha.upper = alpha_upper_bound(tree), alpha.lower = NA, fit.OU.model = FALSE, l1ou.options = NA)
tree |
ultrametric tree of class phylo, with branch lengths, and edges in postorder. |
Y |
trait vector/matrix without missing entries. The row names of the data must be in the same order as the tip labels. |
shift.configuration |
shift positions, i.e. vector of indices of the edges where the shifts occur. |
criterion |
an information criterion (see Details). |
root.model |
an ancestral state model at the root. |
alpha.starting.value |
optional starting value for the optimization of the phylogenetic adaptation rate. |
alpha.upper |
optional upper bound for the phylogenetic adaptation rate. The default value is log(2) over the minimum length of external branches, corresponding to a half life greater or equal to the minimum external branch length. |
alpha.lower |
optional lower bound for the phylogenetic adaptation rate. |
fit.OU.model |
logical. If TRUE, it returns an object of class l1ou with all the parameters estimated. |
l1ou.options |
if provided, all the default values will be ignored. |
AICc gives the usual small-sample size modification of AIC. BIC gives the usual Bayesian information criterion, here penalizing each shift as 2 parameters. mBIC is the modified BIC proposed by Ho and Ané (2014). pBIC is the phylogenetic BIC for shifts proposed by Khabbazian et al. pBICess is a version of pBIC where the determinant term is replaced by a sum of the log of effective sample sizes (ESS), similar to the ESS proposed by Ané (2008).
Information criterion value of the given shift configuration.
Cécile Ané, 2008. "Analysis of comparative data with hierarchical autocorrelation". Annals of Applied Statistics 2(3):1078-1102.
Ho, L. S. T. and Ané, C. 2014. "Intrinsic inference difficulties for trait evolution with Ornstein-Uhlenbeck models". Methods in Ecology and Evolution. 5(11):1133-1146.
Mohammad Khabbazian, Ricardo Kriebel, Karl Rohe, and Cécile Ané (2016). "Fast and accurate detection of evolutionary shifts in Ornstein-Uhlenbeck models". Methods in Ecology and Evolution. doi:10.1111/2041-210X.12534
estimate_shift_configuration
adjust_data
data(lizard.tree, lizard.traits) lizard <- adjust_data(lizard.tree, lizard.traits[,1]) eModel <- estimate_shift_configuration(lizard$tree, lizard$Y) configuration_ic(lizard$tree, eModel$Y, eModel$shift.configuration, criterion="pBIC") ### building l1ou object out of the second best score eModel2 = configuration_ic(eModel$tree, eModel$Y, eModel$profile$configurations[[2]], fit.OU.model=TRUE, l1ou.options=eModel$l1ou.options) plot(eModel2)
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.