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 data model is encoded as a so-called property graph. Nodes represent entities such as Projects, Cases, Diagnoses, Files (various kinds), and Annotations. The relationships between these entities are maintained as edges. Both nodes and edges may have Properties that supply instance details.
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] "e588f035feefee17f562b3a1bc2816c49a2b2b19"
## 
## $data_release
## [1] "Data Release 16.0 - March 26, 2019"
## 
## $status
## [1] "OK"
## 
## $tag
## [1] "1.20.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 tibbles.
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
##  $ df78950c-198d-4afc-a81e-8cfa4dea0abd:'data.frame':    1 obs. of  22 variables:
##   ..$ classification_of_tumor          : chr "not reported"
##   ..$ last_known_disease_status        : chr "not reported"
##   ..$ updated_datetime                 : chr "2018-09-10T15:06:41.433124-05:00"
##   ..$ primary_diagnosis                : chr "Adenocarcinoma, NOS"
##   ..$ submitter_id                     : chr "TCGA-L5-A8NF_diagnosis"
##   ..$ tumor_stage                      : chr "stage iva"
##   .. [list output truncated]
##  $ cccd0568-e267-44a2-8421-6471736bb45f:'data.frame':    1 obs. of  22 variables:
##   ..$ classification_of_tumor          : chr "not reported"
##   ..$ last_known_disease_status        : chr "not reported"
##   ..$ updated_datetime                 : chr "2018-09-10T15:06:41.433124-05:00"
##   ..$ primary_diagnosis                : chr "Adenocarcinoma, NOS"
##   ..$ submitter_id                     : chr "TCGA-R6-A8W8_diagnosis"
##   ..$ tumor_stage                      : chr "not reported"
##   .. [list output truncated]
##  $ ab2755d2-5bd9-4e2f-b255-e813afc8d268:'data.frame':    1 obs. of  22 variables:
##   ..$ classification_of_tumor          : chr "not reported"
##   ..$ last_known_disease_status        : chr "not reported"
##   ..$ updated_datetime                 : chr "2018-09-10T15:06:41.433124-05:00"
##   ..$ primary_diagnosis                : chr "Adenocarcinoma, NOS"
##   ..$ submitter_id                     : chr "TCGA-2H-A9GL_diagnosis"
##   ..$ tumor_stage                      : chr "stage iii"
##   .. [list output truncated]
##  $ 91e37435-16a1-4077-b342-83826480fdc9:'data.frame':    1 obs. of  22 variables:
##   ..$ classification_of_tumor          : chr "not reported"
##   ..$ last_known_disease_status        : chr "not reported"
##   ..$ updated_datetime                 : chr "2018-09-05T23:16:44.045570-05:00"
##   ..$ primary_diagnosis                : chr "Pheochromocytoma, malignant"
##   ..$ submitter_id                     : chr "TCGA-TT-A6YN_diagnosis"
##   ..$ tumor_stage                      : chr "not reported"
##   .. [list output truncated]
##  $ f49438a7-c5e5-47e7-bef4-5e6881a0732a:'data.frame':    1 obs. of  22 variables:
##   ..$ classification_of_tumor          : chr "not reported"
##   ..$ last_known_disease_status        : chr "not reported"
##   ..$ updated_datetime                 : chr "2018-09-05T23:16:44.045570-05:00"
##   ..$ primary_diagnosis                : chr "Pheochromocytoma, NOS"
##   ..$ submitter_id                     : chr "TCGA-RW-A689_diagnosis"
##   ..$ tumor_stage                      : chr "not reported"
##   .. [list output truncated]
##  $ 117eb38b-0c62-4336-a866-fa5bd013256a:'data.frame':    1 obs. of  22 variables:
##   ..$ classification_of_tumor          : chr "not reported"
##   ..$ last_known_disease_status        : chr "not reported"
##   ..$ updated_datetime                 : chr "2018-09-06T18:42:58.220092-05:00"
##   ..$ primary_diagnosis                : chr "Adenocarcinoma, NOS"
##   ..$ submitter_id                     : chr "TCGA-BR-4279_diagnosis"
##   ..$ tumor_stage                      : chr "stage ii"
##   .. [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] 45
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 "phs000466" NA NA ...
##  $ disease_type          :List of 10
##   ..$ TCGA-READ   : chr [1:2] "Cystic, Mucinous and Serous Neoplasms" "Adenomas and Adenocarcinomas"
##   ..$ TARGET-CCSK : chr "Clear Cell Sarcoma of the Kidney"
##   ..$ TCGA-MESO   : chr "Mesothelial Neoplasms"
##   ..$ TCGA-CHOL   : chr "Adenomas and Adenocarcinomas"
##   ..$ NCICCR-DLBCL: chr "Lymphoid Neoplasm Diffuse Large B-cell Lymphoma"
##   ..$ TARGET-WT   : chr "High-Risk Wilms Tumor"
##   ..$ TCGA-TGCT   : chr "Germ Cell Neoplasms"
##   ..$ TCGA-PRAD   : chr [1:3] "Cystic, Mucinous and Serous Neoplasms" "Adenomas and Adenocarcinomas" "Ductal and Lobular Neoplasms"
##   ..$ TCGA-LAML   : chr "Myeloid Leukemias"
##   ..$ TCGA-ESCA   : chr [1:3] "Cystic, Mucinous and Serous Neoplasms" "Adenomas and Adenocarcinomas" "Squamous Cell Neoplasms"
##  $ 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-READ   : chr [1:5] "Rectosigmoid junction" "Unknown" "Rectum" "Colon" ...
##   ..$ TARGET-CCSK : chr "Kidney"
##   ..$ TCGA-MESO   : chr [1:2] "Heart, mediastinum, and pleura" "Bronchus and lung"
##   ..$ TCGA-CHOL   : chr [1:3] "Other and unspecified parts of biliary tract" "Gallbladder" "Liver and intrahepatic bile ducts"
##   ..$ NCICCR-DLBCL: chr "Lymph Nodes"
##   ..$ TARGET-WT   : chr "Kidney"
##   ..$ TCGA-TGCT   : chr "Testis"
##   ..$ TCGA-PRAD   : chr "Prostate gland"
##   ..$ TCGA-LAML   : chr "Hematopoietic and reticuloendothelial systems"
##   ..$ TCGA-ESCA   : chr [1:2] "Esophagus" "Stomach"
##  $ project_id            : chr [1:10] "TCGA-READ" "TARGET-CCSK" "TCGA-MESO" "TCGA-CHOL" ...
##  $ id                    : chr [1:10] "TCGA-READ" "TARGET-CCSK" "TCGA-MESO" "TCGA-CHOL" ...
##  $ name                  : chr [1:10] "Rectum Adenocarcinoma" "Clear Cell Sarcoma of the Kidney" "Mesothelioma" "Cholangiocarcinoma" ...
##  - 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] 45
# 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] "origin"                         "pairs_on_diff_chr"             
## [27] "plate_name"                     "plate_well"                    
## [29] "platform"                       "proportion_base_mismatch"      
## [31] "proportion_coverage_10X"        "proportion_coverage_30X"       
## [33] "proportion_reads_duplicated"    "proportion_reads_mapped"       
## [35] "proportion_targets_no_coverage" "read_pair_number"              
## [37] "revision"                       "state"                         
## [39] "state_comment"                  "submitter_id"                  
## [41] "tags"                           "total_reads"                   
## [43] "type"                           "updated_datetime"
# The number of fields available for files endpoint
length(available_fields('files'))
## [1] 734
# 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" "disease_type"            
##  [9] "index_date"               "lost_to_followup"        
## [11] "portion_ids"              "primary_site"            
## [13] "sample_ids"               "slide_ids"               
## [15] "state"                    "submitter_aliquot_ids"   
## [17] "submitter_analyte_ids"    "submitter_id"            
## [19] "submitter_portion_ids"    "submitter_sample_ids"    
## [21] "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] 365463
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     31524
## 3       TCGA-LUAD     17052
## 4       TCGA-UCEC     16174
## 5       TCGA-HNSC     15289
## 6         TCGA-OV     15166
## 7       TCGA-THCA     14421
## 8       TCGA-LUSC     15324
## 9        TCGA-LGG     14728
## 10      TCGA-KIRC     15091
## 11      TCGA-PRAD     14288
## 12      TCGA-COAD     14279
## 13       TCGA-GBM     11996
## 14      TCGA-SKCM     12725
## 15      TCGA-STAD     12846
## 16      TCGA-BLCA     11711
## 17      TCGA-LIHC     10815
## 18      TCGA-CESC      8595
## 19      TCGA-KIRP      8507
## 20      TCGA-SARC      7494
## 21      TCGA-PAAD      5307
## 22      TCGA-ESCA      5271
## 23      TCGA-PCPG      5035
## 24      TCGA-READ      4925
## 25      TCGA-TGCT      4284
## 26      TCGA-LAML      4434
## 27      TCGA-THYM      3445
## 28     TARGET-NBL      2809
## 29       TCGA-ACC      2547
## 30      TCGA-KICH      2325
## 31      TCGA-MESO      2339
## 32     TARGET-AML      2459
## 33        CPTAC-3      4268
## 34       TCGA-UVM      2180
## 35       TCGA-UCS      1659
## 36      TARGET-WT      1408
## 37  TARGET-ALL-P3      2502
## 38      TCGA-CHOL      1349
## 39      TCGA-DLBC      1229
## 40   NCICCR-DLBCL       957
## 41      TARGET-OS        66
## 42      TARGET-RT       381
## 43    CTSP-DLBCL1        89
## 44    TARGET-CCSK        15
## 45 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_TOKENGDC_TOKEN_FILE.gdc_tokenAs 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-NBL      1127
## 3       TCGA-BRCA      1098
## 4      TARGET-AML       988
## 5       TARGET-WT       652
## 6        TCGA-GBM       617
## 7         TCGA-OV       608
## 8       TCGA-LUAD       585
## 9       TCGA-UCEC       560
## 10      TCGA-KIRC       537
## 11      TCGA-HNSC       528
## 12       TCGA-LGG       516
## 13      TCGA-THCA       507
## 14      TCGA-LUSC       504
## 15      TCGA-PRAD       500
## 16   NCICCR-DLBCL       489
## 17      TCGA-SKCM       470
## 18      TCGA-COAD       461
## 19      TCGA-STAD       443
## 20      TCGA-BLCA       412
## 21      TARGET-OS       381
## 22      TCGA-LIHC       377
## 23        CPTAC-3       322
## 24      TCGA-CESC       307
## 25      TCGA-KIRP       291
## 26      TCGA-SARC       261
## 27      TCGA-LAML       200
## 28      TCGA-ESCA       185
## 29      TCGA-PAAD       185
## 30      TCGA-PCPG       179
## 31      TCGA-READ       172
## 32      TCGA-TGCT       150
## 33  TARGET-ALL-P3       131
## 34      TCGA-THYM       124
## 35      TCGA-KICH       113
## 36       TCGA-ACC        92
## 37      TCGA-MESO        87
## 38       TCGA-UVM        80
## 39      TARGET-RT        75
## 40      TCGA-DLBC        58
## 41       TCGA-UCS        57
## 42      TCGA-CHOL        51
## 43    CTSP-DLBCL1        45
## 44    TARGET-CCSK        13
## 45 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] 3367
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] "2018-09-06T13:49:20.245333-05:00" "2018-11-28T10:07:18.958951-06:00" "2018-09-06T13:49:20.245333-05:00" "2018-09-06T13:49:20.245333-05:00" ...
head(ids(resp))
## [1] "7bcf6d42-eb8f-4aa2-825b-38671d4f9e3a"
## [2] "4a032bad-e726-48f2-8f39-e3acc109cc91"
## [3] "26573441-eedb-4364-966c-e7f803deef19"
## [4] "d9fd2724-7db0-4af3-ac14-217bdfa5203f"
## [5] "6cdc0d53-f813-4101-81e5-9bee68270536"
## [6] "cb9f5e50-f49d-4899-8895-9367afcc1015"
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()
## $analysis.workflow_type
##               key doc_count
## 1  HTSeq - Counts       174
## 2    HTSeq - FPKM       174
## 3 HTSeq - FPKM-UQ       174
# 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.0 (2019-04-26)
## Platform: x86_64-pc-linux-gnu (64-bit)
## Running under: Ubuntu 18.04.2 LTS
## 
## Matrix products: default
## BLAS:   /home/biocbuild/bbs-3.9-bioc/R/lib/libRblas.so
## LAPACK: /home/biocbuild/bbs-3.9-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.1.1            GenomicDataCommons_1.8.0
## [3] magrittr_1.5             knitr_1.22              
## [5] BiocStyle_2.12.0        
## 
## loaded via a namespace (and not attached):
##  [1] SummarizedExperiment_1.14.0 tidyselect_0.2.5           
##  [3] xfun_0.6                    purrr_0.3.2                
##  [5] lattice_0.20-38             colorspace_1.4-1           
##  [7] htmltools_0.3.6             stats4_3.6.0               
##  [9] yaml_2.2.0                  rlang_0.3.4                
## [11] pillar_1.3.1                withr_2.1.2                
## [13] glue_1.3.1                  BiocParallel_1.18.0        
## [15] rappdirs_0.3.1              BiocGenerics_0.30.0        
## [17] plyr_1.8.4                  matrixStats_0.54.0         
## [19] GenomeInfoDbData_1.2.1      stringr_1.4.0              
## [21] zlibbioc_1.30.0             munsell_0.5.0              
## [23] gtable_0.3.0                evaluate_0.13              
## [25] labeling_0.3                Biobase_2.44.0             
## [27] IRanges_2.18.0              GenomeInfoDb_1.20.0        
## [29] parallel_3.6.0              curl_3.3                   
## [31] Rcpp_1.0.1                  readr_1.3.1                
## [33] scales_1.0.0                BiocManager_1.30.4         
## [35] DelayedArray_0.10.0         S4Vectors_0.22.0           
## [37] jsonlite_1.6                XVector_0.24.0             
## [39] hms_0.4.2                   digest_0.6.18              
## [41] stringi_1.4.3               bookdown_0.9               
## [43] dplyr_0.8.0.1               GenomicRanges_1.36.0       
## [45] grid_3.6.0                  tools_3.6.0                
## [47] bitops_1.0-6                lazyeval_0.2.2             
## [49] RCurl_1.95-4.12             tibble_2.1.1               
## [51] crayon_1.3.4                pkgconfig_2.0.2            
## [53] Matrix_1.2-17               xml2_1.2.0                 
## [55] assertthat_0.2.1            rmarkdown_1.12             
## [57] httr_1.4.0                  R6_2.4.0                   
## [59] compiler_3.6.0
S3 object-oriented programming paradigm is used.