A Command-Line Interface for mPTP - multi-rate Poisson Tree Processes
Source:R/mptp_tbl.R
mptp_tbl.Rd
mptp_tbl()
returns species partition hypothesis estimated by mPTP software
https://github.com/Pas-Kapli/mptp.
Usage
mptp_tbl(
infile,
exe = NULL,
outfolder = NULL,
method = c("multi", "single"),
minbrlen = 1e-04,
webserver = NULL,
delimname = "mptp"
)
Source
Kapli T., Lutteropp S., Zhang J., Kobert K., Pavlidis P., Stamatakis A., Flouri T. 2016. Multi-rate Poisson tree processes for single-locus species delimitation under maximum likelihood and Markov chain Monte Carlo. Bioinformatics 33(11):1630-1638.
Arguments
- infile
Path to tree file in Newick format. Should be dichotomous and rooted.
- exe
Path to an mPTP executable.
- outfolder
Path to output folder. Default to NULL. If not specified, a temporary location is used.
- method
Which algorithm for Maximum Likelihood point-estimate to be used. Available options are:
single Single-rate PTP model. It assumes that every species evolved with the same rate.
multi Multi-rate mPTP model. It assumes that all species have different evolutionary rates.
- minbrlen
Numeric. Branch lengths smaller or equal to the value provided are ignored from computations. Default to 0.0001. Use min_brlenfor fine tuning.
- webserver
A .txt file containing mPTP results obtained from a webserver. Default to NULL.
- delimname
Character. String to rename the delimitation method in the table. Default to 'mptp'.
Value
an object of class tbl_df
Details
mptp_tbl()
relies on system to invoke mPTP software through
a command-line interface. Hence, you must have the software available as an executable file on
your system in order to use this function properly. mptp_tbl()
saves all output files in
outfolder
and imports the results generated to Environment
.
If an outfolder
is not provided by the user, then a temporary location is used.
Alternatively, mptp_tbl()
can parse a file obtained from webserver such as
https://mptp.h-its.org/.
Author
Paschalia Kapli, Sarah Lutteropp, Jiajie Zhang, Kassian Kobert, Pavlos Pavlides, Alexandros Stamatakis, Tomáš Flouri.
Examples
# \donttest{
# get path to phylogram
path_to_file <- system.file("extdata/geophagus_raxml.nwk", package = "delimtools")
# run mPTP in single threshold mode (PTP)
ptp_df <- mptp_tbl(
infile = path_to_file,
exe = "/usr/local/bin/mptp",
method = "single",
minbrlen = 0.0001,
delimname = "ptp",
outfolder = NULL
)
#> mptp 0.2.5_linux_x86_64, 16GB RAM, 4 cores
#> https://github.com/Pas-Kapli/mptp
#>
#> Parsing tree file...
#> Loaded rooted tree...
#> Number of edges greater than minimum branch length: 159 / 272
#> Score Null Model: 681.702507
#> Best score for single coalescent rate: 747.727196
#> LRT computed p-value: 0.000000
#> Writing delimitation file /tmp/RtmpRdWzeq/geophagus_raxml.nwk.mptp.single.txt ...
#> Number of delimited species: 17
#> Creating SVG delimitation file /tmp/RtmpRdWzeq/geophagus_raxml.nwk.mptp.single.svg ...
#> Done...
#>
#>
#> ℹ Warning: there are tip-to-tip distances smaller than the specified minimum branch length (0.0001).
#> Consider using `delimtools::min_brlen()` to explore branch lengths in your tree.
#>
#>
#> ℹ mPTP files are located in '/tmp/RtmpRdWzeq'.
#>
#>
# check
ptp_df
#> # A tibble: 137 × 2
#> labels ptp
#> <chr> <int>
#> 1 MZ051794.1 1
#> 2 MZ050845.1 1
#> 3 MZ051032.1 1
#> 4 MZ051706.1 1
#> 5 MZ504576.1 2
#> 6 MZ504578.1 2
#> 7 MZ504579.1 2
#> 8 MZ504588.1 2
#> 9 MZ504587.1 2
#> 10 MZ504592.1 2
#> # ℹ 127 more rows
# run mPTP in multi threshold mode (mPTP)
mptp_df <- mptp_tbl(
infile = path_to_file,
exe = "/usr/local/bin/mptp",
method = "single",
minbrlen = 0.0001,
delimname = "mptp",
outfolder = NULL
)
#> mptp 0.2.5_linux_x86_64, 16GB RAM, 4 cores
#> https://github.com/Pas-Kapli/mptp
#>
#> Parsing tree file...
#> Loaded rooted tree...
#> Number of edges greater than minimum branch length: 159 / 272
#> Score Null Model: 681.702507
#> Best score for single coalescent rate: 747.727196
#> LRT computed p-value: 0.000000
#> Writing delimitation file /tmp/RtmpRdWzeq/geophagus_raxml.nwk.mptp.single.txt ...
#> Number of delimited species: 17
#> Creating SVG delimitation file /tmp/RtmpRdWzeq/geophagus_raxml.nwk.mptp.single.svg ...
#> Done...
#>
#>
#> ℹ Warning: there are tip-to-tip distances smaller than the specified minimum branch length (0.0001).
#> Consider using `delimtools::min_brlen()` to explore branch lengths in your tree.
#>
#>
#> ℹ mPTP files are located in '/tmp/RtmpRdWzeq'.
#>
#>
# check
mptp_df
#> # A tibble: 137 × 2
#> labels mptp
#> <chr> <int>
#> 1 MZ051794.1 1
#> 2 MZ050845.1 1
#> 3 MZ051032.1 1
#> 4 MZ051706.1 1
#> 5 MZ504576.1 2
#> 6 MZ504578.1 2
#> 7 MZ504579.1 2
#> 8 MZ504588.1 2
#> 9 MZ504587.1 2
#> 10 MZ504592.1 2
#> # ℹ 127 more rows
# }