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] "6a169d494657ba07e07a7221b24431b1e91d0942"
## 
## $data_release
## [1] "Data Release 29.0 - March 31, 2021"
## 
## $status
## [1] "OK"
## 
## $tag
## [1] "3.0.0"
## 
## $version
## [1] 1And 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)##  chr [1:50] "375436b3-66ac-4d5e-b495-18a96d812a69" ...# 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] 68The 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
##  $ id                    : chr [1:10] "GENIE-MSK" "TCGA-UCEC" "TCGA-LGG" "TCGA-SARC" ...
##  $ primary_site          :List of 10
##   ..$ GENIE-MSK : chr [1:49] "Testis" "Gallbladder" "Unknown" "Other and unspecified parts of biliary tract" ...
##   ..$ TCGA-UCEC : chr [1:2] "Corpus uteri" "Uterus, NOS"
##   ..$ TCGA-LGG  : chr "Brain"
##   ..$ TCGA-SARC : chr [1:13] "Kidney" "Other and unspecified parts of tongue" "Bones, joints and articular cartilage of limbs" "Colon" ...
##   ..$ TCGA-PAAD : chr "Pancreas"
##   ..$ TCGA-ESCA : chr [1:2] "Esophagus" "Stomach"
##   ..$ TCGA-PRAD : chr "Prostate gland"
##   ..$ GENIE-VICC: chr [1:46] "Testis" "Unknown" "Other and unspecified parts of biliary tract" "Adrenal gland" ...
##   ..$ TCGA-LAML : chr "Hematopoietic and reticuloendothelial systems"
##   ..$ TCGA-KIRC : chr "Kidney"
##  $ dbgap_accession_number: logi [1:10] NA NA NA NA NA NA ...
##  $ project_id            : chr [1:10] "GENIE-MSK" "TCGA-UCEC" "TCGA-LGG" "TCGA-SARC" ...
##  $ disease_type          :List of 10
##   ..$ GENIE-MSK : chr [1:49] "Germ Cell Neoplasms" "Granular Cell Tumors and Alveolar Soft Part Sarcomas" "Immunoproliferative Diseases" "Plasma Cell Tumors" ...
##   ..$ TCGA-UCEC : chr [1:4] "Epithelial Neoplasms, NOS" "Cystic, Mucinous and Serous Neoplasms" "Adenomas and Adenocarcinomas" "Not Reported"
##   ..$ TCGA-LGG  : chr "Gliomas"
##   ..$ TCGA-SARC : chr [1:6] "Nerve Sheath Tumors" "Myomatous Neoplasms" "Fibromatous Neoplasms" "Lipomatous Neoplasms" ...
##   ..$ TCGA-PAAD : chr [1:4] "Cystic, Mucinous and Serous Neoplasms" "Ductal and Lobular Neoplasms" "Adenomas and Adenocarcinomas" "Epithelial Neoplasms, NOS"
##   ..$ TCGA-ESCA : chr [1:3] "Cystic, Mucinous and Serous Neoplasms" "Squamous Cell Neoplasms" "Adenomas and Adenocarcinomas"
##   ..$ TCGA-PRAD : chr [1:3] "Cystic, Mucinous and Serous Neoplasms" "Ductal and Lobular Neoplasms" "Adenomas and Adenocarcinomas"
##   ..$ GENIE-VICC: chr [1:41] "Germ Cell Neoplasms" "Acinar Cell Neoplasms" "Synovial-like Neoplasms" "Plasma Cell Tumors" ...
##   ..$ TCGA-LAML : chr "Myeloid Leukemias"
##   ..$ TCGA-KIRC : chr "Adenomas and Adenocarcinomas"
##  $ name                  : chr [1:10] "AACR Project GENIE - Contributed by Memorial Sloan Kettering Cancer Center" "Uterine Corpus Endometrial Carcinoma" "Brain Lower Grade Glioma" "Sarcoma" ...
##  $ releasable            : logi [1:10] FALSE TRUE TRUE TRUE TRUE TRUE ...
##  $ state                 : chr [1:10] "open" "open" "open" "open" ...
##  $ released              : logi [1:10] TRUE TRUE TRUE TRUE TRUE TRUE ...
##  - 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] 10presults = pquery %>% results_all()
length(ids(presults))## [1] 68# includes all records
length(ids(presults)) == count(pquery)## [1] TRUEExtracting 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] "pairs_on_diff_chr"              "plate_name"                    
## [29] "plate_well"                     "platform"                      
## [31] "proportion_base_mismatch"       "proportion_coverage_10X"       
## [33] "proportion_coverage_10x"        "proportion_coverage_30X"       
## [35] "proportion_coverage_30x"        "proportion_reads_duplicated"   
## [37] "proportion_reads_mapped"        "proportion_targets_no_coverage"
## [39] "read_pair_number"               "revision"                      
## [41] "stain_type"                     "state"                         
## [43] "state_comment"                  "submitter_id"                  
## [45] "tags"                           "total_reads"                   
## [47] "tumor_ploidy"                   "tumor_purity"                  
## [49] "type"                           "updated_datetime"# The number of fields available for files endpoint
length(available_fields('files'))## [1] 981# 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] "consent_type"             "created_datetime"        
##  [7] "days_to_consent"          "days_to_lost_to_followup"
##  [9] "diagnosis_ids"            "disease_type"            
## [11] "index_date"               "lost_to_followup"        
## [13] "portion_ids"              "primary_site"            
## [15] "sample_ids"               "slide_ids"               
## [17] "state"                    "submitter_aliquot_ids"   
## [19] "submitter_analyte_ids"    "submitter_diagnosis_ids" 
## [21] "submitter_id"             "submitter_portion_ids"   
## [23] "submitter_sample_ids"     "submitter_slide_ids"     
## [25] "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$typeUsing 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] 618198To 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_10x"    
## [4] "analysis.input_files.proportion_coverage_30X"    
## [5] "analysis.input_files.proportion_coverage_30x"    
## [6] "analysis.input_files.proportion_reads_duplicated"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
##    doc_count                   key
## 1      36134                 FM-AD
## 2      33766             TCGA-BRCA
## 3      42178               CPTAC-3
## 4      36470             GENIE-MSK
## 5      18358             TCGA-LUAD
## 6      17277             TCGA-UCEC
## 7      16340             TCGA-HNSC
## 8      16344               TCGA-OV
## 9      15445             TCGA-THCA
## 10     29433         MMRF-COMMPASS
## 11     16368             TCGA-LUSC
## 12     15795              TCGA-LGG
## 13     28464            GENIE-DFCI
## 14     16255             TCGA-KIRC
## 15     15296             TCGA-PRAD
## 16     15338             TCGA-COAD
## 17     13089              TCGA-GBM
## 18     20772         TARGET-ALL-P2
## 19     13674             TCGA-SKCM
## 20     13739             TCGA-STAD
## 21     12513             TCGA-BLCA
## 22     11578             TCGA-LIHC
## 23      9201             TCGA-CESC
## 24      9137             TCGA-KIRP
## 25      8002             TCGA-SARC
## 26      7772            TARGET-AML
## 27      5671             TCGA-PAAD
## 28      5657             TCGA-ESCA
## 29      5378             TCGA-PCPG
## 30      5269             TCGA-READ
## 31      5796            TARGET-NBL
## 32      8981     BEATAML1.0-COHORT
## 33      8958               CPTAC-2
## 34      4605             TCGA-TGCT
## 35      4814             TCGA-LAML
## 36      3691             TCGA-THYM
## 37      6036             HCMI-CMDC
## 38      5941               CMI-MBC
## 39      5550         CGCI-HTMCP-CC
## 40      2736              TCGA-ACC
## 41      2457             TCGA-KICH
## 42      2677             TARGET-WT
## 43      4805          NCICCR-DLBCL
## 44      2518             TCGA-MESO
## 45      2340              TCGA-UVM
## 46      3113             TARGET-OS
## 47      3982         TARGET-ALL-P3
## 48      3857             GENIE-MDA
## 49      3833            GENIE-VICC
## 50      3320             GENIE-JHU
## 51      1765              TCGA-UCS
## 52      1426             TCGA-CHOL
## 53      2632             GENIE-UHN
## 54      1325             TCGA-DLBC
## 55      2477            CGCI-BLGSP
## 56      1049             TARGET-RT
## 57      1038            GENIE-GRCC
## 58       994            WCDT-MCRPC
## 59       934               CMI-ASC
## 60       801             GENIE-NKI
## 61       798              OHSU-CNL
## 62       703   ORGANOID-PANCREATIC
## 63       570               CMI-MPC
## 64       417           CTSP-DLBCL1
## 65       223 BEATAML1.0-CRENOLANIB
## 66       169           TARGET-CCSK
## 67       133         TARGET-ALL-P1
## 68        21        VAREPOP-APOLLOWe 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] 1137Asking 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] TRUEGenerating 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] 379The 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
##    doc_count                   key
## 1      18004                 FM-AD
## 2      16824             GENIE-MSK
## 3      14232            GENIE-DFCI
## 4       3857             GENIE-MDA
## 5       3320             GENIE-JHU
## 6       2632             GENIE-UHN
## 7       2146            TARGET-AML
## 8       2052            GENIE-VICC
## 9       1587         TARGET-ALL-P2
## 10      1132            TARGET-NBL
## 11      1098             TCGA-BRCA
## 12      1038            GENIE-GRCC
## 13       995         MMRF-COMMPASS
## 14       801             GENIE-NKI
## 15       795               CPTAC-3
## 16       652             TARGET-WT
## 17       617              TCGA-GBM
## 18       608               TCGA-OV
## 19       585             TCGA-LUAD
## 20       583     BEATAML1.0-COHORT
## 21       560             TCGA-UCEC
## 22       537             TCGA-KIRC
## 23       528             TCGA-HNSC
## 24       516              TCGA-LGG
## 25       507             TCGA-THCA
## 26       504             TCGA-LUSC
## 27       500             TCGA-PRAD
## 28       489          NCICCR-DLBCL
## 29       470             TCGA-SKCM
## 30       461             TCGA-COAD
## 31       443             TCGA-STAD
## 32       412             TCGA-BLCA
## 33       383             TARGET-OS
## 34       377             TCGA-LIHC
## 35       342               CPTAC-2
## 36       307             TCGA-CESC
## 37       291             TCGA-KIRP
## 38       261             TCGA-SARC
## 39       212         CGCI-HTMCP-CC
## 40       200               CMI-MBC
## 41       200             TCGA-LAML
## 42       191         TARGET-ALL-P3
## 43       185             TCGA-ESCA
## 44       185             TCGA-PAAD
## 45       179             TCGA-PCPG
## 46       176              OHSU-CNL
## 47       172             TCGA-READ
## 48       150             TCGA-TGCT
## 49       124             TCGA-THYM
## 50       120            CGCI-BLGSP
## 51       113             TCGA-KICH
## 52       101            WCDT-MCRPC
## 53        92              TCGA-ACC
## 54        87             TCGA-MESO
## 55        80             HCMI-CMDC
## 56        80              TCGA-UVM
## 57        70   ORGANOID-PANCREATIC
## 58        69             TARGET-RT
## 59        58             TCGA-DLBC
## 60        57              TCGA-UCS
## 61        56 BEATAML1.0-CRENOLANIB
## 62        51             TCGA-CHOL
## 63        45           CTSP-DLBCL1
## 64        36               CMI-ASC
## 65        30               CMI-MPC
## 66        24         TARGET-ALL-P1
## 67        13           TARGET-CCSK
## 68         7        VAREPOP-APOLLOlibrary(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] 6197cases() %>% 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] 162res = resp %>% results()
str(res[1],list.len=6)## List of 1
##  $ id: chr [1:162] "5cdae21d-eee5-478f-932a-0f51fcf5f031" "8c09f413-e938-4f2e-a414-84f0e7fcfe41" "d6f911b5-e895-43f8-8f86-0ac2f1bc6fae" "fa176764-a76f-44c7-b97a-cd6d21e052be" ...head(ids(resp))## [1] "5cdae21d-eee5-478f-932a-0f51fcf5f031"
## [2] "8c09f413-e938-4f2e-a414-84f0e7fcfe41"
## [3] "d6f911b5-e895-43f8-8f86-0ac2f1bc6fae"
## [4] "fa176764-a76f-44c7-b97a-cd6d21e052be"
## [5] "5d1d00c6-fcae-479e-ae1e-de76efd41d98"
## [6] "cc074b7f-d3b2-4880-902e-bf10e667b665"res = files() %>% facet('type') %>% aggregations()
res$typeggplot(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 erroropenssl to version
1.0.1 or later.
openssl, reinstall the R curl and httr packages.sessionInfo()## R version 4.1.0 (2021-05-18)
## Platform: x86_64-pc-linux-gnu (64-bit)
## Running under: Ubuntu 20.04.2 LTS
## 
## Matrix products: default
## BLAS:   /home/biocbuild/bbs-3.13-bioc/R/lib/libRblas.so
## LAPACK: /home/biocbuild/bbs-3.13-bioc/R/lib/libRlapack.so
## 
## locale:
##  [1] LC_CTYPE=en_US.UTF-8       LC_NUMERIC=C              
##  [3] LC_TIME=en_GB              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.3.3             GenomicDataCommons_1.16.0
## [3] magrittr_2.0.1            knitr_1.33               
## [5] BiocStyle_2.20.0         
## 
## loaded via a namespace (and not attached):
##  [1] Rcpp_1.0.6                  lattice_0.20-44            
##  [3] assertthat_0.2.1            digest_0.6.27              
##  [5] utf8_1.2.1                  R6_2.5.0                   
##  [7] GenomeInfoDb_1.28.0         stats4_4.1.0               
##  [9] evaluate_0.14               highr_0.9                  
## [11] httr_1.4.2                  pillar_1.6.1               
## [13] zlibbioc_1.38.0             rlang_0.4.11               
## [15] curl_4.3.1                  jquerylib_0.1.4            
## [17] magick_2.7.2                S4Vectors_0.30.0           
## [19] Matrix_1.3-3                rmarkdown_2.8              
## [21] labeling_0.4.2              readr_1.4.0                
## [23] stringr_1.4.0               RCurl_1.98-1.3             
## [25] munsell_0.5.0               DelayedArray_0.18.0        
## [27] compiler_4.1.0              xfun_0.23                  
## [29] pkgconfig_2.0.3             BiocGenerics_0.38.0        
## [31] htmltools_0.5.1.1           tidyselect_1.1.1           
## [33] SummarizedExperiment_1.22.0 tibble_3.1.2               
## [35] GenomeInfoDbData_1.2.6      bookdown_0.22              
## [37] IRanges_2.26.0              matrixStats_0.58.0         
## [39] fansi_0.4.2                 crayon_1.4.1               
## [41] dplyr_1.0.6                 withr_2.4.2                
## [43] bitops_1.0-7                rappdirs_0.3.3             
## [45] grid_4.1.0                  jsonlite_1.7.2             
## [47] gtable_0.3.0                lifecycle_1.0.0            
## [49] DBI_1.1.1                   scales_1.1.1               
## [51] stringi_1.6.2               farver_2.1.0               
## [53] XVector_0.32.0              xml2_1.3.2                 
## [55] bslib_0.2.5.1               ellipsis_0.3.2             
## [57] generics_0.1.0              vctrs_0.3.8                
## [59] tools_4.1.0                 Biobase_2.52.0             
## [61] glue_1.4.2                  purrr_0.3.4                
## [63] hms_1.1.0                   MatrixGenerics_1.4.0       
## [65] parallel_4.1.0              yaml_2.2.1                 
## [67] colorspace_2.0-1            BiocManager_1.30.15        
## [69] GenomicRanges_1.44.0        sass_0.4.0S3 object-oriented programming paradigm is used.