Contents

Installation

Please use the devel version of the AnVIL Bioconductor package.

library(cBioPortalData)
library(AnVIL)

Introduction

This vignette is for users / developers who would like to learn more about the available in cBioPortalData and possibly hit other endpoints in the cBioPortal API implementation. The functionality demonstrated here is used internally by the package to create integrative representations of study datasets.

API representation

Obtaining the cBioPortal API representation object

(cbio <- cBioPortal())
## service: cBioPortal
## tags(); use cbioportal$<tab completion>:
## # A tibble: 61 x 3
##    tag           operation                  summary                             
##    <chr>         <chr>                      <chr>                               
##  1 Cancer Types  getAllCancerTypesUsingGET  Get all cancer types                
##  2 Cancer Types  getCancerTypeUsingGET      Get a cancer type                   
##  3 Clinical Att… fetchClinicalAttributesUs… Fetch clinical attributes           
##  4 Clinical Att… getAllClinicalAttributesI… Get all clinical attributes in the …
##  5 Clinical Att… getAllClinicalAttributesU… Get all clinical attributes         
##  6 Clinical Att… getClinicalAttributeInStu… Get specified clinical attribute    
##  7 Clinical Data fetchAllClinicalDataInStu… Fetch clinical data by patient IDs …
##  8 Clinical Data fetchClinicalDataUsingPOST Fetch clinical data by patient IDs …
##  9 Clinical Data getAllClinicalDataInStudy… Get all clinical data in a study    
## 10 Clinical Data getAllClinicalDataOfPatie… Get all clinical data of a patient …
## # … with 51 more rows
## tag values:
##   Cancer Types, Clinical Attributes, Clinical Data, Copy Number
##   Segments, Discrete Copy Number Alterations, Gene Panels, Generic
##   Assays, Genes, Molecular Data, Molecular Profiles, Mutations,
##   Patients, Sample Lists, Samples, Structural Variants, Studies,
##   Treatments
## schemas():
##   AlleleSpecificCopyNumber, AndedPatientTreatmentFilters,
##   AndedSampleTreatmentFilters, CancerStudy, CancerStudyTags
##   # ... with 55 more elements

Operations

Check available tags, operations, and descriptions as a tibble:

tags(cbio)
## # A tibble: 61 x 3
##    tag           operation                  summary                             
##    <chr>         <chr>                      <chr>                               
##  1 Cancer Types  getAllCancerTypesUsingGET  Get all cancer types                
##  2 Cancer Types  getCancerTypeUsingGET      Get a cancer type                   
##  3 Clinical Att… fetchClinicalAttributesUs… Fetch clinical attributes           
##  4 Clinical Att… getAllClinicalAttributesI… Get all clinical attributes in the …
##  5 Clinical Att… getAllClinicalAttributesU… Get all clinical attributes         
##  6 Clinical Att… getClinicalAttributeInStu… Get specified clinical attribute    
##  7 Clinical Data fetchAllClinicalDataInStu… Fetch clinical data by patient IDs …
##  8 Clinical Data fetchClinicalDataUsingPOST Fetch clinical data by patient IDs …
##  9 Clinical Data getAllClinicalDataInStudy… Get all clinical data in a study    
## 10 Clinical Data getAllClinicalDataOfPatie… Get all clinical data of a patient …
## # … with 51 more rows
head(tags(cbio)$operation)
## [1] "getAllCancerTypesUsingGET"              
## [2] "getCancerTypeUsingGET"                  
## [3] "fetchClinicalAttributesUsingPOST"       
## [4] "getAllClinicalAttributesInStudyUsingGET"
## [5] "getAllClinicalAttributesUsingGET"       
## [6] "getClinicalAttributeInStudyUsingGET"

Searching through the API

searchOps(cbio, "clinical")
## [1] "getAllClinicalAttributesUsingGET"          
## [2] "fetchClinicalAttributesUsingPOST"          
## [3] "fetchClinicalDataUsingPOST"                
## [4] "getAllClinicalAttributesInStudyUsingGET"   
## [5] "getClinicalAttributeInStudyUsingGET"       
## [6] "getAllClinicalDataInStudyUsingGET"         
## [7] "fetchAllClinicalDataInStudyUsingPOST"      
## [8] "getAllClinicalDataOfPatientInStudyUsingGET"
## [9] "getAllClinicalDataOfSampleInStudyUsingGET"

Studies

Get the list of studies available:

getStudies(cbio)
## # A tibble: 309 x 13
##    name     shortName  description     publicStudy pmid  citation  groups status
##    <chr>    <chr>      <chr>           <lgl>       <chr> <chr>     <chr>   <int>
##  1 Oral Sq… Head & ne… Comprehensive … TRUE        2361… Pickerin… ""          0
##  2 Hepatoc… HCC (Inse… Whole-exome se… TRUE        2582… Schulze … "PUBL…      0
##  3 Uveal M… UM (QIMR)  Whole-genome o… TRUE        2668… Johansso… "PUBL…      0
##  4 Neurobl… NBL (AMC)  Whole genome s… TRUE        2236… Molenaar… "PUBL…      0
##  5 Nasopha… NPC (Sing… Whole exome se… TRUE        2495… Lin et a… "PUBL…      0
##  6 Neurobl… NBL (Colo… Whole-genome s… TRUE        2646… Peifer e… ""          0
##  7 Myelody… MDS (Toky… Whole exome se… TRUE        2190… Yoshida … ""          0
##  8 Insulin… Panet (Sh… Whole exome se… TRUE        2432… Cao et a… ""          0
##  9 Pleural… PLMESO (N… Whole-exome se… TRUE        2548… Guo et a… ""          0
## 10 Pilocyt… PAST (Nat… Whole-genome s… TRUE        2381… Jones et… "PUBL…      0
## # … with 299 more rows, and 5 more variables: importDate <chr>,
## #   allSampleCount <int>, studyId <chr>, cancerTypeId <chr>,
## #   referenceGenome <chr>

Clinical Data

Obtain the clinical data for a particular study:

clinicalData(cbio, "acc_tcga")
## # A tibble: 92 x 21
##    uniqueSampleKey   uniquePatientKey   sampleId  patientId studyId CANCER_TYPE 
##    <chr>             <chr>              <chr>     <chr>     <chr>   <chr>       
##  1 VENHQS1PUi1BNUox… VENHQS1PUi1BNUoxO… TCGA-OR-… TCGA-OR-… acc_tc… Adrenocorti…
##  2 VENHQS1PUi1BNUoy… VENHQS1PUi1BNUoyO… TCGA-OR-… TCGA-OR-… acc_tc… Adrenocorti…
##  3 VENHQS1PUi1BNUoz… VENHQS1PUi1BNUozO… TCGA-OR-… TCGA-OR-… acc_tc… Adrenocorti…
##  4 VENHQS1PUi1BNUo0… VENHQS1PUi1BNUo0O… TCGA-OR-… TCGA-OR-… acc_tc… Adrenocorti…
##  5 VENHQS1PUi1BNUo1… VENHQS1PUi1BNUo1O… TCGA-OR-… TCGA-OR-… acc_tc… Adrenocorti…
##  6 VENHQS1PUi1BNUo2… VENHQS1PUi1BNUo2O… TCGA-OR-… TCGA-OR-… acc_tc… Adrenocorti…
##  7 VENHQS1PUi1BNUo3… VENHQS1PUi1BNUo3O… TCGA-OR-… TCGA-OR-… acc_tc… Adrenocorti…
##  8 VENHQS1PUi1BNUo4… VENHQS1PUi1BNUo4O… TCGA-OR-… TCGA-OR-… acc_tc… Adrenocorti…
##  9 VENHQS1PUi1BNUo5… VENHQS1PUi1BNUo5O… TCGA-OR-… TCGA-OR-… acc_tc… Adrenocorti…
## 10 VENHQS1PUi1BNUpB… VENHQS1PUi1BNUpBO… TCGA-OR-… TCGA-OR-… acc_tc… Adrenocorti…
## # … with 82 more rows, and 15 more variables: CANCER_TYPE_DETAILED <chr>,
## #   DAYS_TO_COLLECTION <chr>, FRACTION_GENOME_ALTERED <chr>, IS_FFPE <chr>,
## #   MUTATION_COUNT <chr>, OCT_EMBEDDED <chr>, ONCOTREE_CODE <chr>,
## #   OTHER_SAMPLE_ID <chr>, PATHOLOGY_REPORT_FILE_NAME <chr>,
## #   PATHOLOGY_REPORT_UUID <chr>, SAMPLE_INITIAL_WEIGHT <chr>,
## #   SAMPLE_TYPE <chr>, SAMPLE_TYPE_ID <chr>, SOMATIC_STATUS <chr>,
## #   VIAL_NUMBER <chr>

Molecular Profiles

A table of molecular profiles for a particular study can be obtained by running the following:

mols <- molecularProfiles(cbio, "acc_tcga")
mols[["molecularProfileId"]]
## [1] "acc_tcga_rppa"                                     
## [2] "acc_tcga_rppa_Zscores"                             
## [3] "acc_tcga_gistic"                                   
## [4] "acc_tcga_rna_seq_v2_mrna"                          
## [5] "acc_tcga_rna_seq_v2_mrna_median_Zscores"           
## [6] "acc_tcga_linear_CNA"                               
## [7] "acc_tcga_methylation_hm450"                        
## [8] "acc_tcga_mutations"                                
## [9] "acc_tcga_rna_seq_v2_mrna_median_all_sample_Zscores"

Molecular Profile Data

The data for a molecular profile can be obtained with prior knowledge of available entrezGeneIds:

molecularData(cbio, molecularProfileId = "acc_tcga_rna_seq_v2_mrna",
    entrezGeneIds = c(1, 2),
    sampleIds = c("TCGA-OR-A5J1-01",  "TCGA-OR-A5J2-01")
)
## $acc_tcga_rna_seq_v2_mrna
## # A tibble: 4 x 8
##   uniqueSampleKey    uniquePatientKey   entrezGeneId molecularProfile… sampleId 
##   <chr>              <chr>                     <int> <chr>             <chr>    
## 1 VENHQS1PUi1BNUoxL… VENHQS1PUi1BNUoxO…            1 acc_tcga_rna_seq… TCGA-OR-…
## 2 VENHQS1PUi1BNUoxL… VENHQS1PUi1BNUoxO…            2 acc_tcga_rna_seq… TCGA-OR-…
## 3 VENHQS1PUi1BNUoyL… VENHQS1PUi1BNUoyO…            1 acc_tcga_rna_seq… TCGA-OR-…
## 4 VENHQS1PUi1BNUoyL… VENHQS1PUi1BNUoyO…            2 acc_tcga_rna_seq… TCGA-OR-…
## # … with 3 more variables: patientId <chr>, studyId <chr>, value <dbl>

Genes

All available genes

A list of all the genes provided by the API service including hugo symbols, and entrez gene IDs can be obtained by using the geneTable function:

geneTable(cbio)
## # A tibble: 1,000 x 3
##    entrezGeneId hugoGeneSymbol type          
##           <int> <chr>          <chr>         
##  1       -95835 IVNS1ABP_PT330 phosphoprotein
##  2       -95834 IVNS1ABP_PT328 phosphoprotein
##  3       -95833 IVNS1ABP_PS329 phosphoprotein
##  4       -95832 IVNS1ABP_PS277 phosphoprotein
##  5       -95831 MORC2_PS785    phosphoprotein
##  6       -95830 MORC2_PS779    phosphoprotein
##  7       -95829 MORC2_PS777    phosphoprotein
##  8       -95828 MORC2_PS743    phosphoprotein
##  9       -95827 MORC2_PS739    phosphoprotein
## 10       -95826 MORC2_PS725    phosphoprotein
## # … with 990 more rows

Gene Panels

genePanels(cbio)
## # A tibble: 50 x 2
##    description                                                     genePanelId  
##    <chr>                                                           <chr>        
##  1 Targeted (341 cancer genes) sequencing of various tumor types … IMPACT341    
##  2 Targeted (410 cancer genes) sequencing of various tumor types … IMPACT410    
##  3 Targeted (468 cancer genes) sequencing of various tumor types … IMPACT468    
##  4 Targeted sequencing of urcc tumor via MSK-IMPACT.               IMPACT       
##  5 Targeted (300 cancer genes) sequencing of bladder urothelial c… IMPACT300    
##  6 Targeted sequencing of urcc tumor via MSK-IMPACT.               IMPACT230    
##  7 Targeted (27 cancer genes) sequencing of adenoid cystic carcin… ACYC_FMI_27  
##  8 Targeted (173 cancer genes) sequencing of breast cancers on Il… METABRIC_173 
##  9 Targeted sequencing of 504 cancer-associated genes on Illumina… DFCI_504     
## 10 DNA sequencing of 623 genes with known or potential relationsh… WUSTL-DFCI_6…
## # … with 40 more rows
getGenePanel(cbio, "IMPACT341")
## # A tibble: 341 x 2
##    entrezGeneId hugoGeneSymbol
##           <int> <chr>         
##  1           25 ABL1          
##  2        84142 ABRAXAS1      
##  3          207 AKT1          
##  4          208 AKT2          
##  5        10000 AKT3          
##  6          238 ALK           
##  7          242 ALOX12B       
##  8       139285 AMER1         
##  9          324 APC           
## 10          367 AR            
## # … with 331 more rows

Molecular Gene Panels

genePanelMolecular

gprppa <- genePanelMolecular(cbio,
    molecularProfileId = "acc_tcga_rppa",
    sampleListId = "acc_tcga_all")
gprppa
## # A tibble: 92 x 7
##    uniqueSampleKey  uniquePatientKey molecularProfil… sampleId patientId studyId
##    <chr>            <chr>            <chr>            <chr>    <chr>     <chr>  
##  1 VENHQS1PUi1BNUo… VENHQS1PUi1BNUo… acc_tcga_rppa    TCGA-OR… TCGA-OR-… acc_tc…
##  2 VENHQS1PUi1BNUo… VENHQS1PUi1BNUo… acc_tcga_rppa    TCGA-OR… TCGA-OR-… acc_tc…
##  3 VENHQS1PUi1BNUo… VENHQS1PUi1BNUo… acc_tcga_rppa    TCGA-OR… TCGA-OR-… acc_tc…
##  4 VENHQS1PUi1BNUo… VENHQS1PUi1BNUo… acc_tcga_rppa    TCGA-OR… TCGA-OR-… acc_tc…
##  5 VENHQS1PUi1BNUo… VENHQS1PUi1BNUo… acc_tcga_rppa    TCGA-OR… TCGA-OR-… acc_tc…
##  6 VENHQS1PUi1BNUo… VENHQS1PUi1BNUo… acc_tcga_rppa    TCGA-OR… TCGA-OR-… acc_tc…
##  7 VENHQS1PUi1BNUo… VENHQS1PUi1BNUo… acc_tcga_rppa    TCGA-OR… TCGA-OR-… acc_tc…
##  8 VENHQS1PUi1BNUo… VENHQS1PUi1BNUo… acc_tcga_rppa    TCGA-OR… TCGA-OR-… acc_tc…
##  9 VENHQS1PUi1BNUo… VENHQS1PUi1BNUo… acc_tcga_rppa    TCGA-OR… TCGA-OR-… acc_tc…
## 10 VENHQS1PUi1BNUp… VENHQS1PUi1BNUp… acc_tcga_rppa    TCGA-OR… TCGA-OR-… acc_tc…
## # … with 82 more rows, and 1 more variable: profiled <lgl>

getGenePanelMolecular

getGenePanelMolecular(cbio,
    molecularProfileIds = c("acc_tcga_rppa", "acc_tcga_gistic"),
    sampleIds = allSamples(cbio, "acc_tcga")$sampleId
)
## # A tibble: 184 x 7
##    uniqueSampleKey  uniquePatientKey molecularProfil… sampleId patientId studyId
##    <chr>            <chr>            <chr>            <chr>    <chr>     <chr>  
##  1 VENHQS1PUi1BNUo… VENHQS1PUi1BNUo… acc_tcga_gistic  TCGA-OR… TCGA-OR-… acc_tc…
##  2 VENHQS1PUi1BNUo… VENHQS1PUi1BNUo… acc_tcga_gistic  TCGA-OR… TCGA-OR-… acc_tc…
##  3 VENHQS1PUi1BNUo… VENHQS1PUi1BNUo… acc_tcga_gistic  TCGA-OR… TCGA-OR-… acc_tc…
##  4 VENHQS1PUi1BNUo… VENHQS1PUi1BNUo… acc_tcga_gistic  TCGA-OR… TCGA-OR-… acc_tc…
##  5 VENHQS1PUi1BNUo… VENHQS1PUi1BNUo… acc_tcga_gistic  TCGA-OR… TCGA-OR-… acc_tc…
##  6 VENHQS1PUi1BNUo… VENHQS1PUi1BNUo… acc_tcga_gistic  TCGA-OR… TCGA-OR-… acc_tc…
##  7 VENHQS1PUi1BNUo… VENHQS1PUi1BNUo… acc_tcga_gistic  TCGA-OR… TCGA-OR-… acc_tc…
##  8 VENHQS1PUi1BNUo… VENHQS1PUi1BNUo… acc_tcga_gistic  TCGA-OR… TCGA-OR-… acc_tc…
##  9 VENHQS1PUi1BNUo… VENHQS1PUi1BNUo… acc_tcga_gistic  TCGA-OR… TCGA-OR-… acc_tc…
## 10 VENHQS1PUi1BNUp… VENHQS1PUi1BNUp… acc_tcga_gistic  TCGA-OR… TCGA-OR-… acc_tc…
## # … with 174 more rows, and 1 more variable: profiled <lgl>

getDataByGenePanel

getDataByGenePanel(cbio, "acc_tcga", genePanelId = "IMPACT341",
    molecularProfileId = "acc_tcga_rppa", sampleListId = "acc_tcga_rppa")
## $acc_tcga_rppa
## # A tibble: 2,622 x 9
##    uniqueSampleKey    uniquePatientKey   entrezGeneId molecularProfil… sampleId 
##    <chr>              <chr>                     <int> <chr>            <chr>    
##  1 VENHQS1PUi1BNUoyL… VENHQS1PUi1BNUoyO…          207 acc_tcga_rppa    TCGA-OR-…
##  2 VENHQS1PUi1BNUoyL… VENHQS1PUi1BNUoyO…          208 acc_tcga_rppa    TCGA-OR-…
##  3 VENHQS1PUi1BNUoyL… VENHQS1PUi1BNUoyO…        10000 acc_tcga_rppa    TCGA-OR-…
##  4 VENHQS1PUi1BNUoyL… VENHQS1PUi1BNUoyO…          367 acc_tcga_rppa    TCGA-OR-…
##  5 VENHQS1PUi1BNUoyL… VENHQS1PUi1BNUoyO…          472 acc_tcga_rppa    TCGA-OR-…
##  6 VENHQS1PUi1BNUoyL… VENHQS1PUi1BNUoyO…         8314 acc_tcga_rppa    TCGA-OR-…
##  7 VENHQS1PUi1BNUoyL… VENHQS1PUi1BNUoyO…          596 acc_tcga_rppa    TCGA-OR-…
##  8 VENHQS1PUi1BNUoyL… VENHQS1PUi1BNUoyO…          598 acc_tcga_rppa    TCGA-OR-…
##  9 VENHQS1PUi1BNUoyL… VENHQS1PUi1BNUoyO…        10018 acc_tcga_rppa    TCGA-OR-…
## 10 VENHQS1PUi1BNUoyL… VENHQS1PUi1BNUoyO…          673 acc_tcga_rppa    TCGA-OR-…
## # … with 2,612 more rows, and 4 more variables: patientId <chr>, studyId <chr>,
## #   value <dbl>, hugoGeneSymbol <chr>

It uses the getAllGenesUsingGET function from the API.

Samples

Sample List Identifiers

To display all available sample list identifiers for a particular study ID, one can use the sampleLists function:

sampleLists(cbio, "acc_tcga")
## # A tibble: 9 x 5
##   category          name          description            sampleListId    studyId
##   <chr>             <chr>         <chr>                  <chr>           <chr>  
## 1 all_cases_with_r… Samples with… Samples protein data … acc_tcga_rppa   acc_tc…
## 2 all_cases_with_m… Samples with… Samples with mutation… acc_tcga_cnaseq acc_tc…
## 3 all_cases_in_stu… All samples   All samples (92 sampl… acc_tcga_all    acc_tc…
## 4 all_cases_with_c… Samples with… Samples with CNA data… acc_tcga_cna    acc_tc…
## 5 all_cases_with_m… Samples with… Samples with mutation… acc_tcga_seque… acc_tc…
## 6 all_cases_with_m… Samples with… Samples with methylat… acc_tcga_methy… acc_tc…
## 7 all_cases_with_m… Samples with… Samples with mRNA exp… acc_tcga_rna_s… acc_tc…
## 8 all_cases_with_m… Complete sam… Samples with mutation… acc_tcga_3way_… acc_tc…
## 9 all_cases_with_m… Samples with… Samples with methylat… acc_tcga_methy… acc_tc…

Sample Identifiers

One can obtain the barcodes / identifiers for each sample using a specific sample list identifier, in this case we want all the copy number alteration samples:

samplesInSampleLists(cbio, "acc_tcga_cna")
## CharacterList of length 1
## [["acc_tcga_cna"]] TCGA-OR-A5J1-01 TCGA-OR-A5J2-01 ... TCGA-PK-A5HC-01

This returns a CharacterList of all identifiers for each sample list identifier input:

samplesInSampleLists(cbio, c("acc_tcga_cna", "acc_tcga_cnaseq"))
## CharacterList of length 2
## [["acc_tcga_cna"]] TCGA-OR-A5J1-01 TCGA-OR-A5J2-01 ... TCGA-PK-A5HC-01
## [["acc_tcga_cnaseq"]] TCGA-OR-A5J1-01 TCGA-OR-A5J2-01 ... TCGA-PK-A5HC-01

All samples within a study ID

allSamples(cbio, "acc_tcga")
## # A tibble: 92 x 6
##    uniqueSampleKey    uniquePatientKey   sampleType  sampleId  patientId studyId
##    <chr>              <chr>              <chr>       <chr>     <chr>     <chr>  
##  1 VENHQS1PUi1BNUoxL… VENHQS1PUi1BNUoxO… Primary So… TCGA-OR-… TCGA-OR-… acc_tc…
##  2 VENHQS1PUi1BNUoyL… VENHQS1PUi1BNUoyO… Primary So… TCGA-OR-… TCGA-OR-… acc_tc…
##  3 VENHQS1PUi1BNUozL… VENHQS1PUi1BNUozO… Primary So… TCGA-OR-… TCGA-OR-… acc_tc…
##  4 VENHQS1PUi1BNUo0L… VENHQS1PUi1BNUo0O… Primary So… TCGA-OR-… TCGA-OR-… acc_tc…
##  5 VENHQS1PUi1BNUo1L… VENHQS1PUi1BNUo1O… Primary So… TCGA-OR-… TCGA-OR-… acc_tc…
##  6 VENHQS1PUi1BNUo2L… VENHQS1PUi1BNUo2O… Primary So… TCGA-OR-… TCGA-OR-… acc_tc…
##  7 VENHQS1PUi1BNUo3L… VENHQS1PUi1BNUo3O… Primary So… TCGA-OR-… TCGA-OR-… acc_tc…
##  8 VENHQS1PUi1BNUo4L… VENHQS1PUi1BNUo4O… Primary So… TCGA-OR-… TCGA-OR-… acc_tc…
##  9 VENHQS1PUi1BNUo5L… VENHQS1PUi1BNUo5O… Primary So… TCGA-OR-… TCGA-OR-… acc_tc…
## 10 VENHQS1PUi1BNUpBL… VENHQS1PUi1BNUpBO… Primary So… TCGA-OR-… TCGA-OR-… acc_tc…
## # … with 82 more rows

Info on Samples

getSampleInfo(cbio, studyId = "acc_tcga",
    sampleListIds = c("acc_tcga_rppa", "acc_tcga_gistic"))
## # A tibble: 46 x 6
##    uniqueSampleKey    uniquePatientKey   sampleType  sampleId  patientId studyId
##    <chr>              <chr>              <chr>       <chr>     <chr>     <chr>  
##  1 VENHQS1PUi1BNUoyL… VENHQS1PUi1BNUoyO… Primary So… TCGA-OR-… TCGA-OR-… acc_tc…
##  2 VENHQS1PUi1BNUozL… VENHQS1PUi1BNUozO… Primary So… TCGA-OR-… TCGA-OR-… acc_tc…
##  3 VENHQS1PUi1BNUo2L… VENHQS1PUi1BNUo2O… Primary So… TCGA-OR-… TCGA-OR-… acc_tc…
##  4 VENHQS1PUi1BNUo3L… VENHQS1PUi1BNUo3O… Primary So… TCGA-OR-… TCGA-OR-… acc_tc…
##  5 VENHQS1PUi1BNUo4L… VENHQS1PUi1BNUo4O… Primary So… TCGA-OR-… TCGA-OR-… acc_tc…
##  6 VENHQS1PUi1BNUo5L… VENHQS1PUi1BNUo5O… Primary So… TCGA-OR-… TCGA-OR-… acc_tc…
##  7 VENHQS1PUi1BNUpBL… VENHQS1PUi1BNUpBO… Primary So… TCGA-OR-… TCGA-OR-… acc_tc…
##  8 VENHQS1PUi1BNUpQL… VENHQS1PUi1BNUpQO… Primary So… TCGA-OR-… TCGA-OR-… acc_tc…
##  9 VENHQS1PUi1BNUpSL… VENHQS1PUi1BNUpSO… Primary So… TCGA-OR-… TCGA-OR-… acc_tc…
## 10 VENHQS1PUi1BNUpTL… VENHQS1PUi1BNUpTO… Primary So… TCGA-OR-… TCGA-OR-… acc_tc…
## # … with 36 more rows

Advanced Usage

The cBioPortal API representation is not limited to the functions provided in the package. Users who wish to make use of any of the endpoints provided by the API specification should use the dollar sign $ function to access the endpoints.

First the user should see the input for a particular endpoint as detailed in the API:

cbio$getGeneUsingGET
## getGeneUsingGET 
## Get a gene 
## 
## Parameters:
##   geneId (string)
##     Entrez Gene ID or Hugo Gene Symbol e.g. 1 or A1BG

Then the user can provide such input:

(resp <- cbio$getGeneUsingGET("BRCA1"))
## Response [https://www.cbioportal.org/api/genes/BRCA1]
##   Date: 2021-04-22 22:15
##   Status: 200
##   Content-Type: application/json
##   Size: 69 B

which will require the user to ‘translate’ the response using httr::content:

httr::content(resp)
## $entrezGeneId
## [1] 672
## 
## $hugoGeneSymbol
## [1] "BRCA1"
## 
## $type
## [1] "protein-coding"

Clearing the cache

For users who wish to clear the entire cBioPortalData cache, it is recommended that they use:

unlink("~/.cache/cBioPortalData/")

sessionInfo

sessionInfo()
## R version 4.0.5 (2021-03-31)
## Platform: x86_64-pc-linux-gnu (64-bit)
## Running under: Ubuntu 18.04.5 LTS
## 
## Matrix products: default
## BLAS:   /home/biocbuild/bbs-3.12-bioc/R/lib/libRblas.so
## LAPACK: /home/biocbuild/bbs-3.12-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] parallel  stats4    stats     graphics  grDevices utils     datasets 
## [8] methods   base     
## 
## other attached packages:
##  [1] cBioPortalData_2.2.11       MultiAssayExperiment_1.16.0
##  [3] SummarizedExperiment_1.20.0 Biobase_2.50.0             
##  [5] GenomicRanges_1.42.0        GenomeInfoDb_1.26.7        
##  [7] IRanges_2.24.1              S4Vectors_0.28.1           
##  [9] BiocGenerics_0.36.1         MatrixGenerics_1.2.1       
## [11] matrixStats_0.58.0          AnVIL_1.2.0                
## [13] dplyr_1.0.5                 BiocStyle_2.18.1           
## 
## loaded via a namespace (and not attached):
##  [1] bitops_1.0-6              bit64_4.0.5              
##  [3] progress_1.2.2            httr_1.4.2               
##  [5] GenomicDataCommons_1.14.0 tools_4.0.5              
##  [7] bslib_0.2.4               utf8_1.2.1               
##  [9] R6_2.5.0                  DBI_1.1.1                
## [11] withr_2.4.2               tidyselect_1.1.0         
## [13] prettyunits_1.1.1         TCGAutils_1.10.1         
## [15] bit_4.0.4                 curl_4.3                 
## [17] compiler_4.0.5            cli_2.4.0                
## [19] rvest_1.0.0               formatR_1.9              
## [21] xml2_1.3.2                DelayedArray_0.16.3      
## [23] rtracklayer_1.50.0        bookdown_0.22            
## [25] sass_0.3.1                readr_1.4.0              
## [27] askpass_1.1               rappdirs_0.3.3           
## [29] rapiclient_0.1.3          RCircos_1.2.1            
## [31] Rsamtools_2.6.0           stringr_1.4.0            
## [33] digest_0.6.27             rmarkdown_2.7            
## [35] XVector_0.30.0            pkgconfig_2.0.3          
## [37] htmltools_0.5.1.1         dbplyr_2.1.1             
## [39] fastmap_1.1.0             limma_3.46.0             
## [41] rlang_0.4.10              rstudioapi_0.13          
## [43] RSQLite_2.2.7             jquerylib_0.1.3          
## [45] generics_0.1.0            jsonlite_1.7.2           
## [47] BiocParallel_1.24.1       RCurl_1.98-1.3           
## [49] magrittr_2.0.1            GenomeInfoDbData_1.2.4   
## [51] futile.logger_1.4.3       Matrix_1.3-2             
## [53] Rcpp_1.0.6                fansi_0.4.2              
## [55] lifecycle_1.0.0           stringi_1.5.3            
## [57] yaml_2.2.1                RaggedExperiment_1.14.2  
## [59] RJSONIO_1.3-1.4           zlibbioc_1.36.0          
## [61] BiocFileCache_1.14.0      grid_4.0.5               
## [63] blob_1.2.1                crayon_1.4.1             
## [65] lattice_0.20-41           Biostrings_2.58.0        
## [67] splines_4.0.5             GenomicFeatures_1.42.3   
## [69] hms_1.0.0                 ps_1.6.0                 
## [71] knitr_1.32.9              pillar_1.6.0             
## [73] codetools_0.2-18          biomaRt_2.46.3           
## [75] futile.options_1.0.1      XML_3.99-0.6             
## [77] glue_1.4.2                evaluate_0.14            
## [79] lambda.r_1.2.4            data.table_1.14.0        
## [81] BiocManager_1.30.12       vctrs_0.3.7              
## [83] tidyr_1.1.3               openssl_1.4.3            
## [85] purrr_0.3.4               assertthat_0.2.1         
## [87] cachem_1.0.4              xfun_0.22                
## [89] survival_3.2-10           tibble_3.1.1             
## [91] RTCGAToolbox_2.20.0       GenomicAlignments_1.26.0 
## [93] AnnotationDbi_1.52.0      memoise_2.0.0            
## [95] ellipsis_0.3.1