From the Genomic Data Commons (GDC) website:
The National Cancer Institute’s (NCI’s) Genomic Data Commons (GDC) is a data sharing platform that promotes precision medicine in oncology. It is not just a database or a tool; it is an expandable knowledge network supporting the import and standardization of genomic and clinical data from cancer research programs. The GDC contains NCI-generated data from some of the largest and most comprehensive cancer genomic datasets, including The Cancer Genome Atlas (TCGA) and Therapeutically Applicable Research to Generate Effective Therapies (TARGET). For the first time, these datasets have been harmonized using a common set of bioinformatics pipelines, so that the data can be directly compared. As a growing knowledge system for cancer, the GDC also enables researchers to submit data, and harmonizes these data for import into the GDC. As more researchers add clinical and genomic data to the GDC, it will become an even more powerful tool for making discoveries about the molecular basis of cancer that may lead to better care for patients.
The data model for the GDC is complex, but it worth a quick overview and a graphical representation is included here.
The GDC API exposes these nodes and edges in a somewhat simplified set of RESTful endpoints.
This quickstart section is just meant to show basic functionality. More details of functionality are included further on in this vignette and in function-specific help.
This software is available at Bioconductor.org and can be downloaded via
BiocManager::install
.
To report bugs or problems, either
submit a new issue
or submit a bug.report(package='GenomicDataCommons')
from within R (which
will redirect you to the new issue on GitHub).
Installation can be achieved via Bioconductor’s BiocManager
package.
if (!require("BiocManager"))
install.packages("BiocManager")
BiocManager::install('GenomicDataCommons')
library(GenomicDataCommons)
The GenomicDataCommons package relies on having network
connectivity. In addition, the NCI GDC API must also be operational
and not under maintenance. Checking status
can be used to check this
connectivity and functionality.
GenomicDataCommons::status()
## $commit
## [1] "b18b2385b1e916597856067dc6437f3c20b46bca"
##
## $data_release
## [1] "Data Release 19.0 - September 17, 2019"
##
## $status
## [1] "OK"
##
## $tag
## [1] "1.22.0"
##
## $version
## [1] 1
And to check the status in code:
stopifnot(GenomicDataCommons::status()$status=="OK")
The following code builds a manifest
that can be used to guide the
download of raw data. Here, filtering finds gene expression files
quantified as raw counts using HTSeq
from ovarian cancer patients.
ge_manifest = files() %>%
filter( cases.project.project_id == 'TCGA-OV') %>%
filter( type == 'gene_expression' ) %>%
filter( analysis.workflow_type == 'HTSeq - Counts') %>%
manifest()
head(ge_manifest)
fnames = lapply(ge_manifest$id[1:20],gdcdata)
After the 379 gene expression files
specified in the query above. Using multiple processes to do the download very
significantly speeds up the transfer in many cases. On a standard 1Gb
connection, the following completes in about 30 seconds. The first time the
data are downloaded, R will ask to create a cache directory (see ?gdc_cache
for details of setting and interacting with the cache). Resulting
downloaded files will be stored in the cache directory. Future access to
the same files will be directly from the cache, alleviating multiple downloads.
fnames = lapply(ge_manifest$id[1:20],gdcdata)
If the download had included controlled-access data, the download above would
have needed to include a token
. Details are available in
the authentication section below.
Accessing clinical data is a very common task. Given a set of case_ids
,
the gdc_clinical()
function will return a list of four tibble
s.
case_ids = cases() %>% results(size=10) %>% ids()
clindat = gdc_clinical(case_ids)
names(clindat)
## [1] "demographic" "diagnoses" "exposures" "main"
head(clindat[["main"]])
head(clindat[["diagnoses"]])
The GenomicDataCommons package can access the significant
clinical, demographic, biospecimen, and annotation information
contained in the NCI GDC. The gdc_clinical()
function will often
be all that is needed, but the API and GenomicDataCommons package
make much flexibility if fine-tuning is required.
expands = c("diagnoses","annotations",
"demographic","exposures")
clinResults = cases() %>%
GenomicDataCommons::select(NULL) %>%
GenomicDataCommons::expand(expands) %>%
results(size=50)
str(clinResults[[1]],list.len=6)
## List of 50
## $ 4ca5ff6a-ace1-4349-b269-b8e3aae16785:'data.frame': 1 obs. of 28 variables:
## ..$ year_of_diagnosis : int 2007
## ..$ classification_of_tumor : chr "not reported"
## ..$ last_known_disease_status : chr "not reported"
## ..$ updated_datetime : chr "2019-08-08T16:49:47.467798-05:00"
## ..$ primary_diagnosis : chr "Clear cell adenocarcinoma, NOS"
## ..$ submitter_id : chr "TCGA-A3-3336_diagnosis"
## .. [list output truncated]
## $ 2efcefca-5e52-4b3d-a48a-9e9c73cc0141:'data.frame': 1 obs. of 28 variables:
## ..$ year_of_diagnosis : int 2005
## ..$ classification_of_tumor : chr "not reported"
## ..$ last_known_disease_status : chr "not reported"
## ..$ updated_datetime : chr "2019-08-08T16:53:29.535496-05:00"
## ..$ primary_diagnosis : chr "Clear cell adenocarcinoma, NOS"
## ..$ submitter_id : chr "TCGA-CZ-4854_diagnosis"
## .. [list output truncated]
## $ d26f3a6b-72f4-4e48-815e-57b8b211112d:'data.frame': 1 obs. of 29 variables:
## ..$ year_of_diagnosis : int 2005
## ..$ classification_of_tumor : chr "not reported"
## ..$ last_known_disease_status : chr "not reported"
## ..$ updated_datetime : chr "2019-08-08T16:50:57.630420-05:00"
## ..$ primary_diagnosis : chr "Clear cell adenocarcinoma, NOS"
## ..$ submitter_id : chr "TCGA-B0-5698_diagnosis"
## .. [list output truncated]
## $ f2801b21-5444-4cc3-a642-c60c6d82cd3d:'data.frame': 1 obs. of 28 variables:
## ..$ year_of_diagnosis : int 2006
## ..$ classification_of_tumor : chr "not reported"
## ..$ last_known_disease_status : chr "not reported"
## ..$ updated_datetime : chr "2019-08-08T16:53:29.535496-05:00"
## ..$ primary_diagnosis : chr "Clear cell adenocarcinoma, NOS"
## ..$ submitter_id : chr "TCGA-CZ-5454_diagnosis"
## .. [list output truncated]
## $ 4ee14e95-09d8-4512-b358-3e0e7cb74c9a:'data.frame': 1 obs. of 28 variables:
## ..$ year_of_diagnosis : int 2004
## ..$ classification_of_tumor : chr "not reported"
## ..$ last_known_disease_status : chr "not reported"
## ..$ updated_datetime : chr "2019-08-08T16:50:33.583003-05:00"
## ..$ primary_diagnosis : chr "Clear cell adenocarcinoma, NOS"
## ..$ submitter_id : chr "TCGA-B0-4844_diagnosis"
## .. [list output truncated]
## $ 4edff57f-4b0e-4770-beac-590da7d7232c:'data.frame': 1 obs. of 29 variables:
## ..$ year_of_diagnosis : int 2004
## ..$ classification_of_tumor : chr "not reported"
## ..$ last_known_disease_status : chr "not reported"
## ..$ updated_datetime : chr "2019-08-08T16:53:29.535496-05:00"
## ..$ primary_diagnosis : chr "Clear cell adenocarcinoma, NOS"
## ..$ submitter_id : chr "TCGA-CJ-6033_diagnosis"
## .. [list output truncated]
## [list output truncated]
# or listviewer::jsonedit(clinResults)
This package design is meant to have some similarities to the “hadleyverse” approach of dplyr. Roughly, the functionality for finding and accessing files and metadata can be divided into:
In addition, there are exhiliary functions for asking the GDC API for information about available and default fields, slicing BAM files, and downloading actual data files. Here is an overview of functionality1 See individual function and methods documentation for specific details..
projects()
cases()
files()
annotations()
filter()
facet()
select()
mapping()
available_fields()
default_fields()
grep_fields()
field_picker()
available_values()
available_expand()
results()
count()
response()
gdcdata()
transfer()
gdc_client()
aggregations()
gdc_token()
slicing()
There are two main classes of operations when working with the NCI GDC.
Both classes of operation are reviewed in detail in the following sections.
Vast amounts of metadata about cases (patients, basically), files, projects, and
so-called annotations are available via the NCI GDC API. Typically, one will
want to query metadata to either focus in on a set of files for download or
transfer or to perform so-called aggregations (pivot-tables, facets, similar
to the R table()
functionality).
Querying metadata starts with creating a “blank” query. One
will often then want to filter
the query to limit results prior
to retrieving results. The GenomicDataCommons package has
helper functions for listing fields that are available for
filtering.
In addition to fetching results, the GDC API allows faceting, or aggregating,, useful for compiling reports, generating dashboards, or building user interfaces to GDC data (see GDC web query interface for a non-R-based example).
A query of the GDC starts its life in R. Queries follow the four metadata
endpoints available at the GDC. In particular, there are four convenience
functions that each create GDCQuery
objects (actually, specific subclasses of
GDCQuery
):
projects()
cases()
files()
annotations()
pquery = projects()
The pquery
object is now an object of (S3) class, GDCQuery
(and
gdc_projects
and list
). The object contains the following elements:
projects()
function, the default fields from the GDC are used
(see default_fields()
)filter()
method and will be used to filter results on
retrieval.aggregations()
.Looking at the actual object (get used to using str()
!), note that the query
contains no results.
str(pquery)
## List of 5
## $ fields : chr [1:10] "dbgap_accession_number" "disease_type" "intended_release_date" "name" ...
## $ filters: NULL
## $ facets : NULL
## $ legacy : logi FALSE
## $ expand : NULL
## - attr(*, "class")= chr [1:3] "gdc_projects" "GDCQuery" "list"
[ GDC pagination documentation ]
With a query object available, the next step is to retrieve results from the
GDC. The GenomicDataCommons package. The most basic type of results we can get
is a simple count()
of records available that satisfy the filter criteria.
Note that we have not set any filters, so a count()
here will represent all
the project records publicly available at the GDC in the “default” archive"
pcount = count(pquery)
# or
pcount = pquery %>% count()
pcount
## [1] 53
The results()
method will fetch actual results.
presults = pquery %>% results()
These results are
returned from the GDC in JSON format and
converted into a (potentially nested) list in R. The str()
method is useful
for taking a quick glimpse of the data.
str(presults)
## List of 9
## $ dbgap_accession_number: chr [1:10] NA NA NA NA ...
## $ disease_type :List of 10
## ..$ TCGA-MESO : chr "Mesothelial Neoplasms"
## ..$ TCGA-READ : chr [1:2] "Cystic, Mucinous and Serous Neoplasms" "Adenomas and Adenocarcinomas"
## ..$ TCGA-SARC : chr [1:6] "Lipomatous Neoplasms" "Soft Tissue Tumors and Sarcomas, NOS" "Fibromatous Neoplasms" "Myomatous Neoplasms" ...
## ..$ TCGA-ACC : chr "Adenomas and Adenocarcinomas"
## ..$ TCGA-LGG : chr "Gliomas"
## ..$ TCGA-THCA : chr [1:2] "Epithelial Neoplasms, NOS" "Adenomas and Adenocarcinomas"
## ..$ TARGET-CCSK : chr "Complex Mixed and Stromal Neoplasms"
## ..$ TARGET-NBL : chr [1:2] "Neuroepitheliomatous Neoplasms" "Not Applicable"
## ..$ BEATAML1.0-CRENOLANIB: chr "Myeloid Leukemias"
## ..$ TARGET-AML : chr [1:2] "Myeloid Leukemias" "Not Applicable"
## $ releasable : logi [1:10] FALSE FALSE FALSE FALSE FALSE FALSE ...
## $ released : logi [1:10] TRUE TRUE TRUE TRUE TRUE TRUE ...
## $ state : chr [1:10] "open" "open" "open" "open" ...
## $ primary_site :List of 10
## ..$ TCGA-MESO : chr [1:2] "Heart, mediastinum, and pleura" "Bronchus and lung"
## ..$ TCGA-READ : chr [1:5] "Rectosigmoid junction" "Unknown" "Rectum" "Colon" ...
## ..$ TCGA-SARC : chr [1:13] "Corpus uteri" "Bones, joints and articular cartilage of limbs" "Other and unspecified parts of tongue" "Stomach" ...
## ..$ TCGA-ACC : chr "Adrenal gland"
## ..$ TCGA-LGG : chr "Brain"
## ..$ TCGA-THCA : chr "Thyroid gland"
## ..$ TARGET-CCSK : chr "Kidney"
## ..$ TARGET-NBL : chr [1:20] "Heart, mediastinum, and pleura" "Stomach" "Bones, joints and articular cartilage of other and unspecified sites" "Lymph nodes" ...
## ..$ BEATAML1.0-CRENOLANIB: chr "Hematopoietic and reticuloendothelial systems"
## ..$ TARGET-AML : chr [1:2] "Unknown" "Hematopoietic and reticuloendothelial systems"
## $ project_id : chr [1:10] "TCGA-MESO" "TCGA-READ" "TCGA-SARC" "TCGA-ACC" ...
## $ id : chr [1:10] "TCGA-MESO" "TCGA-READ" "TCGA-SARC" "TCGA-ACC" ...
## $ name : chr [1:10] "Mesothelioma" "Rectum Adenocarcinoma" "Sarcoma" "Adrenocortical Carcinoma" ...
## - attr(*, "row.names")= int [1:10] 1 2 3 4 5 6 7 8 9 10
## - attr(*, "class")= chr [1:3] "GDCprojectsResults" "GDCResults" "list"
A default of only 10 records are returned. We can use the size
and from
arguments to results()
to either page through results or to change the number
of results. Finally, there is a convenience method, results_all()
that will
simply fetch all the available results given a query. Note that results_all()
may take a long time and return HUGE result sets if not used carefully. Use of a
combination of count()
and results()
to get a sense of the expected data
size is probably warranted before calling results_all()
length(ids(presults))
## [1] 10
presults = pquery %>% results_all()
length(ids(presults))
## [1] 53
# includes all records
length(ids(presults)) == count(pquery)
## [1] TRUE
Extracting subsets of results or manipulating the results into a more conventional R data structure is not easily generalizable. However, the purrr, rlist, and data.tree packages are all potentially of interest for manipulating complex, nested list structures. For viewing the results in an interactive viewer, consider the listviewer package.
Central to querying and retrieving data from the GDC is the ability to specify
which fields to return, filtering by fields and values, and faceting or
aggregating. The GenomicDataCommons package includes two simple functions,
available_fields()
and default_fields()
. Each can operate on a character(1)
endpoint name (“cases”, “files”, “annotations”, or “projects”) or a GDCQuery
object.
default_fields('files')
## [1] "access" "acl"
## [3] "average_base_quality" "average_insert_size"
## [5] "average_read_length" "channel"
## [7] "chip_id" "chip_position"
## [9] "contamination" "contamination_error"
## [11] "created_datetime" "data_category"
## [13] "data_format" "data_type"
## [15] "error_type" "experimental_strategy"
## [17] "file_autocomplete" "file_id"
## [19] "file_name" "file_size"
## [21] "imaging_date" "magnification"
## [23] "md5sum" "mean_coverage"
## [25] "msi_score" "msi_status"
## [27] "origin" "pairs_on_diff_chr"
## [29] "plate_name" "plate_well"
## [31] "platform" "proportion_base_mismatch"
## [33] "proportion_coverage_10X" "proportion_coverage_30X"
## [35] "proportion_reads_duplicated" "proportion_reads_mapped"
## [37] "proportion_targets_no_coverage" "read_pair_number"
## [39] "revision" "state"
## [41] "state_comment" "submitter_id"
## [43] "tags" "total_reads"
## [45] "type" "updated_datetime"
# The number of fields available for files endpoint
length(available_fields('files'))
## [1] 850
# The first few fields available for files endpoint
head(available_fields('files'))
## [1] "access" "acl"
## [3] "analysis.analysis_id" "analysis.analysis_type"
## [5] "analysis.created_datetime" "analysis.input_files.access"
The fields to be returned by a query can be specified following a similar
paradigm to that of the dplyr package. The select()
function is a verb that
resets the fields slot of a GDCQuery
; note that this is not quite analogous to
the dplyr select()
verb that limits from already-present fields. We
completely replace the fields when using select()
on a GDCQuery
.
# Default fields here
qcases = cases()
qcases$fields
## [1] "aliquot_ids" "analyte_ids"
## [3] "case_autocomplete" "case_id"
## [5] "created_datetime" "days_to_index"
## [7] "days_to_lost_to_followup" "diagnosis_ids"
## [9] "disease_type" "index_date"
## [11] "lost_to_followup" "portion_ids"
## [13] "primary_site" "sample_ids"
## [15] "slide_ids" "state"
## [17] "submitter_aliquot_ids" "submitter_analyte_ids"
## [19] "submitter_diagnosis_ids" "submitter_id"
## [21] "submitter_portion_ids" "submitter_sample_ids"
## [23] "submitter_slide_ids" "updated_datetime"
# set up query to use ALL available fields
# Note that checking of fields is done by select()
qcases = cases() %>% GenomicDataCommons::select(available_fields('cases'))
head(qcases$fields)
## [1] "case_id" "aliquot_ids"
## [3] "analyte_ids" "annotations.annotation_id"
## [5] "annotations.case_id" "annotations.case_submitter_id"
Finding fields of interest is such a common operation that the
GenomicDataCommons includes the grep_fields()
function and the
field_picker()
widget. See the appropriate help pages for details.
The GDC API offers a feature known as aggregation or faceting. By
specifying one or more fields (of appropriate type), the GDC can
return to us a count of the number of records matching each potential
value. This is similar to the R table
method. Multiple fields can be
returned at once, but the GDC API does not have a cross-tabulation
feature; all aggregations are only on one field at a time. Results of
aggregation()
calls come back as a list of data.frames (actually,
tibbles).
# total number of files of a specific type
res = files() %>% facet(c('type','data_type')) %>% aggregations()
res$type
Using aggregations()
is an also easy way to learn the contents of individual
fields and forms the basis for faceted search pages.
[ GDC filtering
documentation ]
The GenomicDataCommons package uses a form of non-standard evaluation to specify R-like queries that are then translated into an R list. That R list is, upon calling a method that fetches results from the GDC API, translated into the appropriate JSON string. The R expression uses the formula interface as suggested by Hadley Wickham in his vignette on non-standard evaluation
It’s best to use a formula because a formula captures both the expression to evaluate and the environment where the evaluation occurs. This is important if the expression is a mixture of variables in a data frame and objects in the local environment [for example].
For the user, these details will not be too important except to note that a filter expression must begin with a “~”.
qfiles = files()
qfiles %>% count() # all files
## [1] 427407
To limit the file type, we can refer back to the section on faceting to see the possible values for the file field “type”. For example, to filter file results to only “gene_expression” files, we simply specify a filter.
qfiles = files() %>% filter( type == 'gene_expression')
# here is what the filter looks like after translation
str(get_filter(qfiles))
## List of 2
## $ op : 'scalar' chr "="
## $ content:List of 2
## ..$ field: chr "type"
## ..$ value: chr "gene_expression"
What if we want to create a filter based on the project (‘TCGA-OVCA’, for example)? Well, we have a couple of possible ways to discover available fields. The first is based on base R functionality and some intuition.
grep('pro',available_fields('files'),value=TRUE) %>%
head()
## [1] "analysis.input_files.proportion_base_mismatch"
## [2] "analysis.input_files.proportion_coverage_10X"
## [3] "analysis.input_files.proportion_coverage_30X"
## [4] "analysis.input_files.proportion_reads_duplicated"
## [5] "analysis.input_files.proportion_reads_mapped"
## [6] "analysis.input_files.proportion_targets_no_coverage"
Interestingly, the project information is “nested” inside the case. We don’t need to know that detail other than to know that we now have a few potential guesses for where our information might be in the files records. We need to know where because we need to construct the appropriate filter.
files() %>%
facet('cases.project.project_id') %>%
aggregations() %>%
head()
## $cases.project.project_id
## key doc_count
## 1 FM-AD 36134
## 2 TCGA-BRCA 31598
## 3 TCGA-LUAD 17074
## 4 TCGA-UCEC 16187
## 5 TCGA-HNSC 15298
## 6 TCGA-OV 15166
## 7 TCGA-THCA 14435
## 8 TCGA-LUSC 15340
## 9 TCGA-LGG 14741
## 10 TCGA-KIRC 15107
## 11 TCGA-PRAD 14314
## 12 TCGA-COAD 14320
## 13 TCGA-GBM 12005
## 14 TCGA-SKCM 12738
## 15 TCGA-STAD 12867
## 16 TCGA-BLCA 11721
## 17 MMRF-COMMPASS 21790
## 18 TCGA-LIHC 10832
## 19 TARGET-ALL-P2 13973
## 20 TCGA-CESC 8599
## 21 TCGA-KIRP 8541
## 22 TCGA-SARC 7494
## 23 TCGA-PAAD 5307
## 24 TCGA-ESCA 5289
## 25 TCGA-PCPG 5044
## 26 TCGA-READ 4925
## 27 CPTAC-3 9052
## 28 BEATAML1.0-COHORT 8981
## 29 TARGET-NBL 3941
## 30 TCGA-TGCT 4293
## 31 TARGET-AML 3914
## 32 TCGA-LAML 4434
## 33 TCGA-THYM 3445
## 34 TCGA-ACC 2556
## 35 TCGA-KICH 2325
## 36 NCICCR-DLBCL 4805
## 37 TCGA-MESO 2346
## 38 TARGET-WT 2115
## 39 TCGA-UVM 2180
## 40 TARGET-ALL-P3 2848
## 41 TCGA-UCS 1659
## 42 TCGA-CHOL 1354
## 43 TCGA-DLBC 1229
## 44 TARGET-OS 771
## 45 CGCI-BLGSP 1516
## 46 TARGET-RT 909
## 47 ORGANOID-PANCREATIC 703
## 48 CTSP-DLBCL1 417
## 49 HCMI-CMDC 301
## 50 BEATAML1.0-CRENOLANIB 223
## 51 TARGET-CCSK 121
## 52 TARGET-ALL-P1 109
## 53 VAREPOP-APOLLO 21
We note that cases.project.project_id
looks like it is a good fit. We also
note that TCGA-OV
is the correct project_id, not TCGA-OVCA
. Note that
unlike with dplyr and friends, the filter()
method here replaces the
filter and does not build on any previous filters.
qfiles = files() %>%
filter( cases.project.project_id == 'TCGA-OV' & type == 'gene_expression')
str(get_filter(qfiles))
## List of 2
## $ op : 'scalar' chr "and"
## $ content:List of 2
## ..$ :List of 2
## .. ..$ op : 'scalar' chr "="
## .. ..$ content:List of 2
## .. .. ..$ field: chr "cases.project.project_id"
## .. .. ..$ value: chr "TCGA-OV"
## ..$ :List of 2
## .. ..$ op : 'scalar' chr "="
## .. ..$ content:List of 2
## .. .. ..$ field: chr "type"
## .. .. ..$ value: chr "gene_expression"
qfiles %>% count()
## [1] 1137
Asking for a count()
of results given these new filter criteria gives r qfiles %>% count()
results. Filters can be chained (or nested) to
accomplish the same effect as multiple &
conditionals. The count()
below is equivalent to the &
filtering done above.
qfiles2 = files() %>%
filter( cases.project.project_id == 'TCGA-OV') %>%
filter( type == 'gene_expression')
qfiles2 %>% count()
## [1] 1137
(qfiles %>% count()) == (qfiles2 %>% count()) #TRUE
## [1] TRUE
Generating a manifest for bulk downloads is as simple as asking for the manifest from the current query.
manifest_df = qfiles %>% manifest()
head(manifest_df)
Note that we might still not be quite there. Looking at filenames, there are
suspiciously named files that might include “FPKM”, “FPKM-UQ”, or “counts”.
Another round of grep
and available_fields
, looking for “type” turned up
that the field “analysis.workflow_type” has the appropriate filter criteria.
qfiles = files() %>% filter( ~ cases.project.project_id == 'TCGA-OV' &
type == 'gene_expression' &
analysis.workflow_type == 'HTSeq - Counts')
manifest_df = qfiles %>% manifest()
nrow(manifest_df)
## [1] 379
The GDC Data Transfer Tool can be used (from R, transfer()
or from the
command-line) to orchestrate high-performance, restartable transfers of all the
files in the manifest. See the bulk downloads section for
details.
[ GDC authentication documentation ]
The GDC offers both “controlled-access” and “open” data. As of this writing, only data stored as files is “controlled-access”; that is, metadata accessible via the GDC is all “open” data and some files are “open” and some are “controlled-access”. Controlled-access data are only available after going through the process of obtaining access.
After controlled-access to one or more datasets has been granted, logging into the GDC web portal will allow you to access a GDC authentication token, which can be downloaded and then used to access available controlled-access data via the GenomicDataCommons package.
The GenomicDataCommons uses authentication tokens only for downloading
data (see transfer
and gdcdata
documentation). The package
includes a helper function, gdc_token
, that looks for the token to
be stored in one of three ways (resolved in this order):
GDC_TOKEN
GDC_TOKEN_FILE
.gdc_token
As a concrete example:
token = gdc_token()
transfer(...,token=token)
# or
transfer(...,token=get_token())
The gdcdata
function takes a character vector of one or more file
ids. A simple way of producing such a vector is to produce a
manifest
data frame and then pass in the first column, which will
contain file ids.
fnames = gdcdata(manifest_df$id[1:2],progress=FALSE)
Note that for controlled-access data, a
GDC authentication token is required. Using the
BiocParallel
package may be useful for downloading in parallel,
particularly for large numbers of smallish files.
The bulk download functionality is only efficient (as of v1.2.0 of the GDC Data Transfer Tool) for relatively large files, so use this approach only when transferring BAM files or larger VCF files, for example. Otherwise, consider using the approach shown above, perhaps in parallel.
fnames = gdcdata(manifest_df$id[3:10], access_method = 'client')
res = cases() %>% facet("project.project_id") %>% aggregations()
head(res)
## $project.project_id
## key doc_count
## 1 FM-AD 18004
## 2 TARGET-ALL-P2 1587
## 3 TARGET-NBL 1132
## 4 TCGA-BRCA 1098
## 5 TARGET-AML 1011
## 6 MMRF-COMMPASS 995
## 7 TARGET-WT 652
## 8 TCGA-GBM 617
## 9 TCGA-OV 608
## 10 TCGA-LUAD 585
## 11 BEATAML1.0-COHORT 583
## 12 TCGA-UCEC 560
## 13 TCGA-KIRC 537
## 14 TCGA-HNSC 528
## 15 TCGA-LGG 516
## 16 TCGA-THCA 507
## 17 TCGA-LUSC 504
## 18 TCGA-PRAD 500
## 19 NCICCR-DLBCL 489
## 20 TCGA-SKCM 470
## 21 TCGA-COAD 461
## 22 TCGA-STAD 443
## 23 TCGA-BLCA 412
## 24 TARGET-OS 383
## 25 TCGA-LIHC 377
## 26 CPTAC-3 322
## 27 TCGA-CESC 307
## 28 TCGA-KIRP 291
## 29 TCGA-SARC 261
## 30 TCGA-LAML 200
## 31 TARGET-ALL-P3 191
## 32 TCGA-ESCA 185
## 33 TCGA-PAAD 185
## 34 TCGA-PCPG 179
## 35 TCGA-READ 172
## 36 TCGA-TGCT 150
## 37 TCGA-THYM 124
## 38 CGCI-BLGSP 120
## 39 TCGA-KICH 113
## 40 TCGA-ACC 92
## 41 TCGA-MESO 87
## 42 TCGA-UVM 80
## 43 ORGANOID-PANCREATIC 70
## 44 TARGET-RT 69
## 45 TCGA-DLBC 58
## 46 TCGA-UCS 57
## 47 BEATAML1.0-CRENOLANIB 56
## 48 TCGA-CHOL 51
## 49 CTSP-DLBCL1 45
## 50 TARGET-ALL-P1 24
## 51 TARGET-CCSK 13
## 52 HCMI-CMDC 7
## 53 VAREPOP-APOLLO 7
library(ggplot2)
ggplot(res$project.project_id,aes(x = key, y = doc_count)) +
geom_bar(stat='identity') +
theme(axis.text.x = element_text(angle = 45, hjust = 1))
cases() %>% filter(~ project.program.name=='TARGET') %>% count()
## [1] 5062
cases() %>% filter(~ project.program.name=='TCGA') %>% count()
## [1] 11315
# The need to do the "&" here is a requirement of the
# current version of the GDC API. I have filed a feature
# request to remove this requirement.
resp = cases() %>% filter(~ project.project_id=='TCGA-BRCA' &
project.project_id=='TCGA-BRCA' ) %>%
facet('samples.sample_type') %>% aggregations()
resp$samples.sample_type
# The need to do the "&" here is a requirement of the
# current version of the GDC API. I have filed a feature
# request to remove this requirement.
resp = cases() %>% filter(~ project.project_id=='TCGA-BRCA' &
samples.sample_type=='Solid Tissue Normal') %>%
GenomicDataCommons::select(c(default_fields(cases()),'samples.sample_type')) %>%
response_all()
count(resp)
## [1] 162
res = resp %>% results()
str(res[1],list.len=6)
## List of 1
## $ updated_datetime: chr [1:162] "2019-08-06T14:16:05.792292-05:00" "2019-08-06T14:15:54.128069-05:00" "2019-08-06T14:16:05.792292-05:00" "2019-08-06T14:16:05.792292-05:00" ...
head(ids(resp))
## [1] "9ddc3e7b-8b54-4a83-8335-8053940f56c1"
## [2] "e666081b-caaf-4e8b-9446-807be71201a5"
## [3] "50619f8c-10aa-464a-a227-90a7aa6ffd43"
## [4] "cc074b7f-d3b2-4880-902e-bf10e667b665"
## [5] "5fe3c51f-65e2-4da6-8924-a6fae17c022d"
## [6] "955d4263-61f7-42e8-8a6e-772a0d6c209d"
res = files() %>% facet('type') %>% aggregations()
res$type
ggplot(res$type,aes(x = key,y = doc_count)) + geom_bar(stat='identity') +
theme(axis.text.x = element_text(angle = 45, hjust = 1))
q = files() %>%
GenomicDataCommons::select(available_fields('files')) %>%
filter(~ cases.project.project_id=='TCGA-GBM' &
data_type=='Gene Expression Quantification')
q %>% facet('analysis.workflow_type') %>% aggregations()
## list()
# so need to add another filter
file_ids = q %>% filter(~ cases.project.project_id=='TCGA-GBM' &
data_type=='Gene Expression Quantification' &
analysis.workflow_type == 'HTSeq - Counts') %>%
GenomicDataCommons::select('file_id') %>%
response_all() %>%
ids()
I need to figure out how to do slicing reproducibly in a testing environment and for vignette building.
q = files() %>%
GenomicDataCommons::select(available_fields('files')) %>%
filter(~ cases.project.project_id == 'TCGA-GBM' &
data_type == 'Aligned Reads' &
experimental_strategy == 'RNA-Seq' &
data_format == 'BAM')
file_ids = q %>% response_all() %>% ids()
bamfile = slicing(file_ids[1],regions="chr12:6534405-6538375",token=gdc_token())
library(GenomicAlignments)
aligns = readGAlignments(bamfile)
Error in curl::curl_fetch_memory(url, handle = handle) :
SSL connect error
openssl
to version
1.0.1 or later.
openssl
, reinstall the R curl
and httr
packages.sessionInfo()
## R version 3.6.1 (2019-07-05)
## Platform: x86_64-pc-linux-gnu (64-bit)
## Running under: Ubuntu 18.04.3 LTS
##
## Matrix products: default
## BLAS: /home/biocbuild/bbs-3.10-bioc/R/lib/libRblas.so
## LAPACK: /home/biocbuild/bbs-3.10-bioc/R/lib/libRlapack.so
##
## locale:
## [1] LC_CTYPE=en_US.UTF-8 LC_NUMERIC=C
## [3] LC_TIME=en_US.UTF-8 LC_COLLATE=C
## [5] LC_MONETARY=en_US.UTF-8 LC_MESSAGES=en_US.UTF-8
## [7] LC_PAPER=en_US.UTF-8 LC_NAME=C
## [9] LC_ADDRESS=C LC_TELEPHONE=C
## [11] LC_MEASUREMENT=en_US.UTF-8 LC_IDENTIFICATION=C
##
## attached base packages:
## [1] stats graphics grDevices utils datasets methods base
##
## other attached packages:
## [1] ggplot2_3.2.1 GenomicDataCommons_1.10.0
## [3] magrittr_1.5 knitr_1.25
## [5] BiocStyle_2.14.0
##
## loaded via a namespace (and not attached):
## [1] SummarizedExperiment_1.16.0 tidyselect_0.2.5
## [3] xfun_0.10 purrr_0.3.3
## [5] lattice_0.20-38 colorspace_1.4-1
## [7] vctrs_0.2.0 htmltools_0.4.0
## [9] stats4_3.6.1 yaml_2.2.0
## [11] rlang_0.4.1 pillar_1.4.2
## [13] withr_2.1.2 glue_1.3.1
## [15] BiocParallel_1.20.0 rappdirs_0.3.1
## [17] BiocGenerics_0.32.0 matrixStats_0.55.0
## [19] GenomeInfoDbData_1.2.2 stringr_1.4.0
## [21] zlibbioc_1.32.0 munsell_0.5.0
## [23] gtable_0.3.0 evaluate_0.14
## [25] labeling_0.3 Biobase_2.46.0
## [27] IRanges_2.20.0 GenomeInfoDb_1.22.0
## [29] parallel_3.6.1 curl_4.2
## [31] Rcpp_1.0.2 readr_1.3.1
## [33] scales_1.0.0 backports_1.1.5
## [35] BiocManager_1.30.9 DelayedArray_0.12.0
## [37] S4Vectors_0.24.0 jsonlite_1.6
## [39] XVector_0.26.0 hms_0.5.2
## [41] digest_0.6.22 stringi_1.4.3
## [43] bookdown_0.14 dplyr_0.8.3
## [45] GenomicRanges_1.38.0 grid_3.6.1
## [47] tools_3.6.1 bitops_1.0-6
## [49] lazyeval_0.2.2 RCurl_1.95-4.12
## [51] tibble_2.1.3 crayon_1.3.4
## [53] pkgconfig_2.0.3 zeallot_0.1.0
## [55] Matrix_1.2-17 xml2_1.2.2
## [57] assertthat_0.2.1 rmarkdown_1.16
## [59] httr_1.4.1 R6_2.4.0
## [61] compiler_3.6.1
S3
object-oriented programming paradigm is used.