raxml | R Documentation |
Provides an interface to the C program RAxML (see Reference section) for maximum likelihood estimation of tree topology and/or branch lengths, rapid and conventional non-parametric bootstrapping, mapping splits onto individual topologies, and a lot more. See the RAxML manual for details, especially if you are a new user of RAxML.
raxml(
DNAbin,
m = "GTRCAT",
f,
N,
p,
b,
x,
k,
weights,
partitions,
outgroup,
backbone = NULL,
file = paste0("fromR_", Sys.Date()),
exec,
threads
)
DNAbin |
A matrix of DNA sequences of class |
m |
A vector of mode |
f |
A vector of mode |
N |
Either of mode |
p |
Integer, setting a random seed for the parsimony starting trees. |
b |
Integer, setting a random seed for bootstrapping. |
x |
Integer, setting a random seed for rapid bootstrapping. |
k |
Logical, if |
weights |
A vector of mode |
partitions |
A data frame giving the partitions of the alignment. |
outgroup |
A vector of mode |
backbone |
A |
file |
A vector of mode |
exec |
A vector of mode |
threads |
Integer, giving the number of parallel threads to use (PTHREADS only). |
There are some limitations of this wrapper compared to RAxML run directly from the command line.
Only DNA is allowed as data type.
Option f
can only take a limited number of values
(d
, a
).
RAxML needs the specification of random seeds for parsimony estimation of
starting trees and for bootstrap resampling. The corresponding argument
names in raxml
are identical to the flags used by RAxML (-p
,
-b
, and -x
). If you choose not to give any values,
raxml
will generate a (different) value for each required random
seed every time it is called. Be aware that set.seed
will work
only for p
, but not for b
or x
.
A list with a variable number of elements, depending on the analysis chosen:
"info" | RAxML log file as character string |
"bestTree" | MLE of tree |
"bipartitions" | MLE of tree annotated with bootstrap proportions |
"bootstrap" | bootstrapped trees |
RAxML is a C program and the source code is not contained in this
package. This means that in order to run this function you will need to
install RAxML yourself. See
https://cme.h-its.org/exelixis/web/software/raxml/ for the most recent
documentation and source code of RAxML. Depending on where you chose to
install RAxML, you need to adjust the exec
argument.
raxml
was last tested and running fine on Mac OS X with RAxML 8.0.29.
Please be aware that calling third-party software from within R is a
platform-specific process and I cannot guarantee that raxml
will
behave properly on any system.
(in chronological order)
Stamatakis, A., T. Ludwig and H. Meier. 2004. RAxML-III: A fast program for maximum likelihood-based inference of large phylogenetic trees. Bioinformatics 1: 1–8.
Stamatakis, A. 2006. RAxML-VI-HPC: Maximum likelihood-based phylogenetic analyses with thousands of taxa and mixed models. Bioinformatics 22: 2688–2690.
Stamatakis, A., P. Hoover, and J. Rougemont. 2008. A rapid bootstrap algorithm for the RAxML web-servers. Syst. Biol. 75: 758–771.
Pattengale, N. D., M. Alipour, O. R. P. Bininda-Emonds, B. M. E. Moret, and A. Stamatakis. 2010. How many bootstrap replicates are necessary? Journal of Computational Biology 17: 337-354.
Stamatakis, A. 2014. RAxML Version 8: A tool for phylogenetic analysis and post-analysis of large phylogenies. Bioinformatics Advance Access.
raxml.partitions
to store partitioning information in a
data frame suitable for input as partitions
argument in raxml
.
## bark beetle sequences
data(ips.cox1)
data(ips.16S)
data(ips.28S)
ips <- cbind(ips.cox1, ips.16S, ips.28S,
fill.with.gaps = TRUE)
exec <- "/Applications/RAxML-code/standard-RAxML/raxmlHPC-PTHREADS-AVX"
w <- sample(1:5, ncol(ips.cox1), replace = TRUE)
## Not run:
# Simple tree search with GTRCAT and GTRGAMMA
tr <- raxml(ips.cox1, f = "d", N = 2, p = 1234,
exec = exec) # -1743.528461
tr <- raxml(ips.cox1, m = "GTRGAMMA", f = "d", N = 2, p = 1234,
exec = exec)
# Applying weights to columns
tr <- raxml(ips.cox1, f = "d", N = 2, p = 1234,
weights = w, exec = exec) # -1743.528461
# Rapid bootstrap
tr <- raxml(ips.cox1, m = "GTRGAMMA",
f = "a", N = 10, p = 1234, x = 1234,
exec = exec)
# Rapid bootstrap with automatic halt
tr <- raxml(ips.cox1, m = "GTRGAMMA",
f = "a", N = "autoMRE", p = 1234, x = 1234,
exec = exec)
## End(Not run)
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.