mpa2meco {file2meco} | R Documentation |
Transform metagenomic classification results of 'mpa' format to 'microtable' object.
Description
Transform the classification results of mpa (MetaPhlAn) format to microtable object, such as MetaPhlAn and Kraken2/Bracken results. Kraken2/Bracken results can be obtained by merge_metaphlan_tables.py from MetaPhlAn or combine_mpa.py from KrakenTools (https://ccb.jhu.edu/software/krakentools/). The algorithm of Kraken2 determines that the abundance of a taxon is not equal to the sum of abundances of taxa in its subordinate lineage. So the default tables in taxa_abund of return microtable object are extracted from the abundances of raw file. It is totally different with the return taxa_abund of cal_abund function, which sums the abundances of taxa at different taxonomic levels based on the taxonomic table and the otu_table (i.e., taxa abundance table at a specified level, e.g., 's__').
Usage
mpa2meco(
feature_table,
sample_table = NULL,
match_table = NULL,
use_level = "s__",
rel = FALSE,
sel_same = 1,
...
)
Arguments
feature_table |
'mpa' format abundance table, see the example. |
sample_table |
default NULL; sample metadata table; If provided, must be one of the several types of formats:
1) comma seperated file with the suffix csv or tab seperated file with suffix tsv/txt;
2) Excel type file with the suffix xlsx or xls; require |
match_table |
default NULL; a two column table used to replace the sample names in abundance table; Must be two columns without column names; The first column must be raw sample names same with those in feature table, the second column must be new sample names same with the rownames in sample_table; Please also see the example files. |
use_level |
default "s__"; the prefix parsed for the otu_table and tax_table; must be one of 'd__', 'k__', 'p__', 'c__', 'o__', 'f__', 'g__' and 's__'. |
rel |
default FALSE; Whether convert the original abundance to relative abundance in generated taxa_abund list. If TRUE, all the data.frame objects in taxa_abund list have relative abundance (0-1). |
sel_same |
default 1; How to select the taxonomic information in |
... |
parameter passed to microtable$new function of microeco package, such as auto_tidy parameter. |
Value
microtable object.
Examples
library(microeco)
library(file2meco)
library(magrittr)
# use Kraken2 file stored inside the package
abund_file_path <- system.file("extdata", "example_kraken2_merge.txt", package="file2meco")
mpa2meco(abund_file_path)
# add sample information table
sample_file_path <- system.file("extdata", "example_metagenome_sample_info.tsv",
package="file2meco")
# sample names are different between abund_file_path and sample_file_path;
# use a matching table to adjust them
match_file_path <- system.file("extdata", "example_metagenome_match_table.tsv", package="file2meco")
test <- mpa2meco(abund_file_path, sample_table = sample_file_path,
match_table = match_file_path, use_level = "s__", rel = TRUE)
# make the taxonomy standard for the following analysis
test$tax_table %<>% tidy_taxonomy
test$tidy_dataset()
# calculate taxa_abund with specified level instead of raw kraken results
test1 <- clone(test)
test1$cal_abund(rel = TRUE)
identical(test$taxa_abund$Kingdom, test1$taxa_abund$Kingdom)