View source: R/reading-writing.R
reading | R Documentation |
Reading dataset(s) in the 3 tables 'dataMatrix' (or 'DM'), sampleMetadata' (or 'SM') and variableMetadata' (or 'VM') tabular format. In case of a single dataset (3 tables in the specified directory), a SummarizedExperiment instance is returned. In case of a multiple dataset (several subfolders containing 3 tables), a MultiAssayExperiment instance is created.
reading( dir.c, files.ls = NULL, subsets.vc = NA, output.c = c("exp", "set")[1], report.c = c("none", "interactive", "myfile.txt")[2] )
dir.c |
character(1): directory containing the 3 tabular files (single dataset), or containing several subdirectories with 3 tabular files (multiple datasets) |
files.ls |
list: if dir.c is set to NA, the full names of the individual files can be provided; in case of a SummarizedExperiment, the names of the list must be 'dataMatrix', 'sampleMetadata', and 'variableMetadata' with the corresponding file full names; in case of a MultiAssayExperiment, the list must consists of one such sublist per dataset |
subsets.vc |
character(): specifying a subset of the subdirectories to be included in the MultiAssayExperiment (by default, all subdirectories containing the 3 tables will be considered as datasets) |
output.c |
character(1): Either 'exp' for SummarizedExperiment (or MultiAssayExperiment), or 'set' for ExpressionSet (or MultiDataSet) output formats (the latter are deprecated) |
report.c |
character(1): File name for the printed results (call to 'sink()'); if NA (default), messages will be printed on the screen; if NULL, no verbose will be generated |
SummarizedExperiment
(one dataset) or MultiAssayExperiment
(multiple
datasets) instance containing the dataset(s)
data_dir.c <- system.file("extdata", package = "phenomis") ## 1) Single set sacurine_dir.c <- file.path(data_dir.c, "W4M00001_Sacurine-statistics") sacurine.se <- reading(sacurine_dir.c) # or sacurine.se <- reading(NA, files.ls = list(dataMatrix = file.path(sacurine_dir.c, "Galaxy1_dataMatrix.tabular"), sampleMetadata = file.path(sacurine_dir.c, "Galaxy2_sampleMetadata.tabular"), variableMetadata = file.path(sacurine_dir.c, "Galaxy3_variableMetadata.tabular"))) ## 2) Multiple sets prometis_dir.c <- file.path(data_dir.c, "prometis") prometis.mae <- reading(prometis_dir.c) metabo.mae <- reading(prometis_dir.c, subsets.vc = "metabolomics") # or prometis.mae <- reading(NA, files.ls = list(metabolomics = list(dataMatrix = file.path(prometis_dir.c, "metabolomics", "dataMatrix.tsv"), sampleMetadata = file.path(prometis_dir.c, "metabolomics", "sampleMetadata.tsv"), variableMetadata = file.path(prometis_dir.c, "metabolomics", "variableMetadata.tsv")), proteomics = list(dataMatrix = file.path(prometis_dir.c, "proteomics", "dataMatrix.tsv"), sampleMetadata = file.path(prometis_dir.c, "proteomics", "sampleMetadata.tsv"), variableMetadata = file.path(prometis_dir.c, "proteomics", "variableMetadata.tsv"))))
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.