This project aims to import all cBioPortal datasets as MultiAssayExperiment objects in Bioconductor. It offers some advantages over the CDGS-R package:
- The MultiAssayExperiment class explicitly links all assays to the patient clinical/pathological data
- The MultiAssayExperiment class provides a flexible API including harmonized subsetting and reshaping to convenient wide and long formats.
- It provides complete datasets, not just for subsets of genes
- It provides automatic local caching, thanks to BiocFileCache.
It is a work in progress, and due to some variation in their formats, does not yet work for all 268 (as of Dec 2019) datasets. At the time of writing, it successfully imports 93% of 200 randomly sampled datasets. Please feel free to file an issue to request prioritization of fixing any of the remaining datasets.
if (!requireNamespace("BiocManager", quietly = TRUE))
install.packages("BiocManager")
if (!requireNamespace("cBioPortalData", quietly = TRUE))
BiocManager::install("waldronlab/cBioPortalData")
library(cBioPortalData)
Flexible and granular access to cBioPortal data from
cbioportal.org/api
. This option is best used with a particular gene
panel of interest. It allows users to download sections of the data with
molecular profile and gene panel combinations within a study.
cbio <- cBioPortal()
gbm <- cBioPortalData(api = cbio, by = "hugoGeneSymbol", studyId = "gbm_tcga",
genePanelId = "IMPACT341",
molecularProfileIds = c("gbm_tcga_rppa", "gbm_tcga_mrna")
)
gbm
#> A MultiAssayExperiment object of 2 listed
#> experiments with user-defined names and respective classes.
#> Containing an ExperimentList class object of length 2:
#> [1] gbm_tcga_rppa: SummarizedExperiment with 67 rows and 244 columns
#> [2] gbm_tcga_mrna: SummarizedExperiment with 334 rows and 401 columns
#> Functionality:
#> experiments() - obtain the ExperimentList instance
#> colData() - the primary/phenotype DataFrame
#> sampleMap() - the sample coordination DataFrame
#> `$`, `[`, `[[` - extract colData columns, subset, or experiment
#> *Format() - convert into a long or wide DataFrame
#> assays() - convert ExperimentList to a SimpleList of matrices
#> exportClass() - save all data to files
This function will download a dataset from the cbioportal.org/datasets
website as a packaged tarball and serve it to users as a
MultiAssayExperiment
object. This option is good for users who are
interested in obtaining all the data for a particular study.
laml <- cBioDataPack("laml_tcga")
laml
#> A MultiAssayExperiment object of 11 listed
#> experiments with user-defined names and respective classes.
#> Containing an ExperimentList class object of length 11:
#> [1] CNA: SummarizedExperiment with 24776 rows and 191 columns
#> [2] RNA_Seq_expression_median: SummarizedExperiment with 19720 rows and 179 columns
#> [3] RNA_Seq_mRNA_median_Zscores: SummarizedExperiment with 19720 rows and 179 columns
#> [4] RNA_Seq_v2_expression_median: SummarizedExperiment with 20531 rows and 173 columns
#> [5] RNA_Seq_v2_mRNA_median_Zscores: SummarizedExperiment with 20531 rows and 173 columns
#> [6] cna_hg19.seg: RaggedExperiment with 13571 rows and 191 columns
#> [7] linear_CNA: SummarizedExperiment with 24776 rows and 191 columns
#> [8] methylation_hm27: SummarizedExperiment with 10919 rows and 194 columns
#> [9] methylation_hm450: SummarizedExperiment with 10919 rows and 194 columns
#> [10] mutations_extended: RaggedExperiment with 2584 rows and 197 columns
#> [11] mutations_mskcc: RaggedExperiment with 2584 rows and 197 columns
#> Functionality:
#> experiments() - obtain the ExperimentList instance
#> colData() - the primary/phenotype DataFrame
#> sampleMap() - the sample coordination DataFrame
#> `$`, `[`, `[[` - extract colData columns, subset, or experiment
#> *Format() - convert into a long or wide DataFrame
#> assays() - convert ExperimentList to a SimpleList of matrices
#> exportClass() - save all data to files