cBioPortalData: API Reference Guide for Devs
Marcel Ramos & Levi Waldron
March 25, 2026
Source:vignettes/cBioPortalRClient.Rmd
cBioPortalRClient.RmdIntroduction
The cBioPortal for Cancer Genomics website is a great resource for interactive exploration of study datasets. However, it does not easily allow the analyst to obtain and further analyze the data.
We’ve developed the cBioPortalData package to fill this
need to programmatically access the data resources available on the
cBioPortal.
The cBioPortalData package provides an R interface for
accessing the cBioPortal study data within the Bioconductor
ecosystem.
It downloads study data from the cBioPortal API (https://cbioportal.org/api) and uses Bioconductor infrastructure to cache and represent the data.
We use the MultiAssayExperiment
(@Ramos2017-er) package to integrate,
represent, and coordinate multiple experiments for the studies availble
in the cBioPortal. This package in conjunction with
curatedTCGAData give access to a large trove of publicly
available bioinformatic data. Please see our publication here (@Ramos2020-ya).
We demonstrate common use cases of cBioPortalData and
curatedTCGAData during Bioconductor conference workshops.
Overview
This vignette is for users / developers who would like to learn more
about the available data in cBioPortalData and to learn how
to 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.
Note. To avoid overloading the API service, the API was designed to only query a part of the study data. Therefore, the user is required to enter either a set of genes of interest or a gene panel identifier.
API representation
Obtaining the cBioPortal API representation object
(cbio <- cBioPortal())## service: cBioPortal
## host: www.cbioportal.org
## tags(); use cbioportal$<tab completion>:
## # A tibble: 65 × 3
## tag operation summary
## <chr> <chr> <chr>
## 1 Cancer Types getAllCancerTypesUsingGET Get all cance…
## 2 Cancer Types getCancerTypeUsingGET Get a cancer …
## 3 Clinical Attributes fetchClinicalAttributesUsingPOST Fetch clinica…
## 4 Clinical Attributes getAllClinicalAttributesInStudyUsingGET Get all clini…
## 5 Clinical Attributes getAllClinicalAttributesUsingGET Get all clini…
## 6 Clinical Attributes getClinicalAttributeInStudyUsingGET Get specified…
## 7 Clinical Data fetchAllClinicalDataInStudyUsingPOST Fetch clinica…
## 8 Clinical Data fetchClinicalDataUsingPOST Fetch clinica…
## 9 Clinical Data getAllClinicalDataInStudyUsingGET Get all clini…
## 10 Clinical Data getAllClinicalDataOfPatientInStudyUsingGET Get all clini…
## # ℹ 55 more rows
## tag values:
## Cancer Types, Clinical Attributes, Clinical Data, Copy Number
## Segments, Discrete Copy Number Alterations, Gene Panel Data, Gene
## Panels, Generic Assay Data, Generic Assays, Genes, Info, Molecular
## Data, Molecular Profiles, Mutations, Patients, Sample Lists, Samples,
## Server running status, Studies, Treatments
## schemas():
## AlleleSpecificCopyNumber, AlterationFilter,
## AndedPatientTreatmentFilters, AndedSampleTreatmentFilters,
## CancerStudy
## # ... with 62 more elements
Operations
Check available tags, operations, and descriptions as a
tibble:
tags(cbio)## # A tibble: 65 × 3
## tag operation summary
## <chr> <chr> <chr>
## 1 Cancer Types getAllCancerTypesUsingGET Get all cance…
## 2 Cancer Types getCancerTypeUsingGET Get a cancer …
## 3 Clinical Attributes fetchClinicalAttributesUsingPOST Fetch clinica…
## 4 Clinical Attributes getAllClinicalAttributesInStudyUsingGET Get all clini…
## 5 Clinical Attributes getAllClinicalAttributesUsingGET Get all clini…
## 6 Clinical Attributes getClinicalAttributeInStudyUsingGET Get specified…
## 7 Clinical Data fetchAllClinicalDataInStudyUsingPOST Fetch clinica…
## 8 Clinical Data fetchClinicalDataUsingPOST Fetch clinica…
## 9 Clinical Data getAllClinicalDataInStudyUsingGET Get all clini…
## 10 Clinical Data getAllClinicalDataOfPatientInStudyUsingGET Get all clini…
## # ℹ 55 more rows
## [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: 535 × 13
## name description publicStudy pmid citation groups status importDate
## <chr> <chr> <lgl> <chr> <chr> <chr> <int> <chr>
## 1 Mixed cfDNA … "IMPACT se… TRUE 3405… Tsui et… "TSUI… 0 2026-01-0…
## 2 Uterine Endo… "CPTAC Ute… TRUE NA NA "PUBL… 0 2026-01-1…
## 3 Mantle Cell … "Whole exo… TRUE 2414… Beà et… "" 0 2026-01-0…
## 4 Anaplastic T… "Whole-gen… TRUE 3841… Zeng, P… "" 0 2026-01-1…
## 5 Stomach Aden… "Whole gen… TRUE 2481… Wang et… "PUBL… 0 2026-01-1…
## 6 Non-Small Ce… "Whole exo… TRUE 2576… Rivzi e… "PUBL… 0 2026-01-0…
## 7 NCI-60 Cell … "Whole-exo… TRUE 2280… Reinhol… "" 0 2026-01-0…
## 8 Thyroid Carc… "TCGA Thyr… TRUE NA NA "PUBL… 0 2026-01-1…
## 9 Retinoblasto… "Cell-free… TRUE 3263… Kothari… "TSUI… 0 2026-01-1…
## 10 Diffuse Larg… "Targeted … TRUE 3849… Zhu, Me… "" 0 2026-01-0…
## # ℹ 525 more rows
## # ℹ 5 more variables: allSampleCount <int>, readPermission <lgl>,
## # studyId <chr>, cancerTypeId <chr>, referenceGenome <chr>
Clinical Data
Obtain the clinical data for a particular study:
clinicalData(cbio, "acc_tcga")## # A tibble: 92 × 85
## patientId AGE AJCC_PATHOLOGIC_TUMOR_STAGE ATYPICAL_MITOTIC_FIGURES
## <chr> <chr> <chr> <chr>
## 1 TCGA-OR-A5J1 58 Stage II Atypical Mitotic Figures Abse…
## 2 TCGA-OR-A5J2 44 Stage IV Atypical Mitotic Figures Pres…
## 3 TCGA-OR-A5J3 23 Stage III Atypical Mitotic Figures Abse…
## 4 TCGA-OR-A5J4 23 Stage IV Atypical Mitotic Figures Abse…
## 5 TCGA-OR-A5J5 30 Stage III Atypical Mitotic Figures Pres…
## 6 TCGA-OR-A5J6 29 Stage II Atypical Mitotic Figures Abse…
## 7 TCGA-OR-A5J7 30 Stage III Atypical Mitotic Figures Pres…
## 8 TCGA-OR-A5J8 66 Stage III Atypical Mitotic Figures Pres…
## 9 TCGA-OR-A5J9 22 Stage II Atypical Mitotic Figures Abse…
## 10 TCGA-OR-A5JA 53 Stage IV Atypical Mitotic Figures Pres…
## # ℹ 82 more rows
## # ℹ 81 more variables: CAPSULAR_INVASION <chr>, CLIN_M_STAGE <chr>,
## # CT_SCAN_PREOP_RESULTS <chr>,
## # CYTOPLASM_PRESENCE_LESS_THAN_EQUAL_25_PERCENT <chr>,
## # DAYS_TO_INITIAL_PATHOLOGIC_DIAGNOSIS <chr>, DFS_MONTHS <chr>,
## # DFS_STATUS <chr>, DIFFUSE_ARCHITECTURE <chr>, ETHNICITY <chr>,
## # FORM_COMPLETION_DATE <chr>, HISTOLOGICAL_DIAGNOSIS <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_gistic"
## [2] "acc_tcga_linear_CNA"
## [3] "acc_tcga_methylation_hm450"
## [4] "acc_tcga_mutations"
## [5] "acc_tcga_rna_seq_v2_mrna"
## [6] "acc_tcga_rna_seq_v2_mrna_median_Zscores"
## [7] "acc_tcga_rna_seq_v2_mrna_median_all_sample_Zscores"
## [8] "acc_tcga_rppa"
## [9] "acc_tcga_rppa_Zscores"
Molecular Profile Data
The data for a molecular profile can be obtained with prior knowledge
of available entrezGeneIds:
molecularData(cbio, molecularProfileIds = "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 × 8
## uniqueSampleKey uniquePatientKey entrezGeneId molecularProfileId sampleId
## <chr> <chr> <int> <chr> <chr>
## 1 VENHQS1PUi1BNUoxLTA… VENHQS1PUi1BNUo… 1 acc_tcga_rna_seq_… TCGA-OR…
## 2 VENHQS1PUi1BNUoxLTA… VENHQS1PUi1BNUo… 2 acc_tcga_rna_seq_… TCGA-OR…
## 3 VENHQS1PUi1BNUoyLTA… VENHQS1PUi1BNUo… 1 acc_tcga_rna_seq_… TCGA-OR…
## 4 VENHQS1PUi1BNUoyLTA… VENHQS1PUi1BNUo… 2 acc_tcga_rna_seq_… TCGA-OR…
## # ℹ 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 × 3
## entrezGeneId hugoGeneSymbol type
## <int> <chr> <chr>
## 1 1 A1BG protein-coding
## 2 2 A2M protein-coding
## 3 13 AADAC protein-coding
## 4 14 AAMP protein-coding
## 5 15 AANAT protein-coding
## 6 16 AARS1 protein-coding
## 7 17 AAVS1 other
## 8 18 ABAT protein-coding
## 9 19 ABCA1 protein-coding
## 10 20 ABCA2 protein-coding
## # ℹ 990 more rows
Gene Panels
genePanels(cbio)## # A tibble: 69 × 2
## description genePanelId
## <chr> <chr>
## 1 Targeted (129 genes) sequencing of cfDNA via MSK-ACCESS on Illum… ACCESS129
## 2 Targeted (27 cancer genes) sequencing of adenoid cystic carcinom… ACYC_FMI_27
## 3 ARCHER-HEME gene panel (199 genes) ARCHER-HEM…
## 4 ARCHER-SOLID Gene Panel (62 genes) ARCHER-SOL…
## 5 ARCHER-SOLID-CV4 Gene Panel (123 genes) ARCHER-SOL…
## 6 Targeted panel of 232 genes. Agilent
## 7 Targeted panel of 8 genes. AmpliSeq
## 8 Targeted panel of 509 genes. Bait_UCSF5…
## 9 DFCI-ONCOPANEL-1, Number of Genes - 304 DFCI-ONCOP…
## 10 DFCI-ONCOPANEL-2, Number of Genes - 326 DFCI-ONCOP…
## # ℹ 59 more rows
getGenePanel(cbio, "IMPACT341")## # A tibble: 341 × 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
## # ℹ 331 more rows
Molecular Gene Panels
genePanelMolecular
gprppa <- genePanelMolecular(cbio,
molecularProfileId = "acc_tcga_rppa",
sampleListId = "acc_tcga_all")
gprppa## # A tibble: 92 × 7
## uniqueSampleKey uniquePatientKey molecularProfileId sampleId patientId
## <chr> <chr> <chr> <chr> <chr>
## 1 VENHQS1PUi1BNUoxLTAxO… VENHQS1PUi1BNUo… acc_tcga_rppa TCGA-OR… TCGA-OR-…
## 2 VENHQS1PUi1BNUoyLTAxO… VENHQS1PUi1BNUo… acc_tcga_rppa TCGA-OR… TCGA-OR-…
## 3 VENHQS1PUi1BNUozLTAxO… VENHQS1PUi1BNUo… acc_tcga_rppa TCGA-OR… TCGA-OR-…
## 4 VENHQS1PUi1BNUo0LTAxO… VENHQS1PUi1BNUo… acc_tcga_rppa TCGA-OR… TCGA-OR-…
## 5 VENHQS1PUi1BNUo1LTAxO… VENHQS1PUi1BNUo… acc_tcga_rppa TCGA-OR… TCGA-OR-…
## 6 VENHQS1PUi1BNUo2LTAxO… VENHQS1PUi1BNUo… acc_tcga_rppa TCGA-OR… TCGA-OR-…
## 7 VENHQS1PUi1BNUo3LTAxO… VENHQS1PUi1BNUo… acc_tcga_rppa TCGA-OR… TCGA-OR-…
## 8 VENHQS1PUi1BNUo4LTAxO… VENHQS1PUi1BNUo… acc_tcga_rppa TCGA-OR… TCGA-OR-…
## 9 VENHQS1PUi1BNUo5LTAxO… VENHQS1PUi1BNUo… acc_tcga_rppa TCGA-OR… TCGA-OR-…
## 10 VENHQS1PUi1BNUpBLTAxO… VENHQS1PUi1BNUp… acc_tcga_rppa TCGA-OR… TCGA-OR-…
## # ℹ 82 more rows
## # ℹ 2 more variables: studyId <chr>, profiled <lgl>
getGenePanelMolecular
getGenePanelMolecular(cbio,
molecularProfileIds = c("acc_tcga_rppa", "acc_tcga_gistic"),
sampleIds = allSamples(cbio, "acc_tcga")$sampleId
)## # A tibble: 184 × 7
## uniqueSampleKey uniquePatientKey molecularProfileId sampleId patientId
## <chr> <chr> <chr> <chr> <chr>
## 1 VENHQS1PUi1BNUoxLTAxO… VENHQS1PUi1BNUo… acc_tcga_gistic TCGA-OR… TCGA-OR-…
## 2 VENHQS1PUi1BNUoyLTAxO… VENHQS1PUi1BNUo… acc_tcga_gistic TCGA-OR… TCGA-OR-…
## 3 VENHQS1PUi1BNUozLTAxO… VENHQS1PUi1BNUo… acc_tcga_gistic TCGA-OR… TCGA-OR-…
## 4 VENHQS1PUi1BNUo0LTAxO… VENHQS1PUi1BNUo… acc_tcga_gistic TCGA-OR… TCGA-OR-…
## 5 VENHQS1PUi1BNUo1LTAxO… VENHQS1PUi1BNUo… acc_tcga_gistic TCGA-OR… TCGA-OR-…
## 6 VENHQS1PUi1BNUo2LTAxO… VENHQS1PUi1BNUo… acc_tcga_gistic TCGA-OR… TCGA-OR-…
## 7 VENHQS1PUi1BNUo3LTAxO… VENHQS1PUi1BNUo… acc_tcga_gistic TCGA-OR… TCGA-OR-…
## 8 VENHQS1PUi1BNUo4LTAxO… VENHQS1PUi1BNUo… acc_tcga_gistic TCGA-OR… TCGA-OR-…
## 9 VENHQS1PUi1BNUo5LTAxO… VENHQS1PUi1BNUo… acc_tcga_gistic TCGA-OR… TCGA-OR-…
## 10 VENHQS1PUi1BNUpBLTAxO… VENHQS1PUi1BNUp… acc_tcga_gistic TCGA-OR… TCGA-OR-…
## # ℹ 174 more rows
## # ℹ 2 more variables: studyId <chr>, profiled <lgl>
getDataByGenes
getDataByGenes(cbio, "acc_tcga", genePanelId = "IMPACT341",
molecularProfileIds = "acc_tcga_rppa", sampleListId = "acc_tcga_rppa")## $acc_tcga_rppa
## # A tibble: 2,622 × 9
## uniqueSampleKey uniquePatientKey entrezGeneId molecularProfileId sampleId
## <chr> <chr> <int> <chr> <chr>
## 1 VENHQS1PUi1BNUoyLT… VENHQS1PUi1BNUo… 10000 acc_tcga_rppa TCGA-OR…
## 2 VENHQS1PUi1BNUoyLT… VENHQS1PUi1BNUo… 5604 acc_tcga_rppa TCGA-OR…
## 3 VENHQS1PUi1BNUoyLT… VENHQS1PUi1BNUo… 5290 acc_tcga_rppa TCGA-OR…
## 4 VENHQS1PUi1BNUoyLT… VENHQS1PUi1BNUo… 253260 acc_tcga_rppa TCGA-OR…
## 5 VENHQS1PUi1BNUoyLT… VENHQS1PUi1BNUo… 7157 acc_tcga_rppa TCGA-OR…
## 6 VENHQS1PUi1BNUoyLT… VENHQS1PUi1BNUo… 595 acc_tcga_rppa TCGA-OR…
## 7 VENHQS1PUi1BNUoyLT… VENHQS1PUi1BNUo… 898 acc_tcga_rppa TCGA-OR…
## 8 VENHQS1PUi1BNUoyLT… VENHQS1PUi1BNUo… 1499 acc_tcga_rppa TCGA-OR…
## 9 VENHQS1PUi1BNUoyLT… VENHQS1PUi1BNUo… 1956 acc_tcga_rppa TCGA-OR…
## 10 VENHQS1PUi1BNUoyLT… VENHQS1PUi1BNUo… 10111 acc_tcga_rppa TCGA-OR…
## # ℹ 2,612 more rows
## # ℹ 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 × 5
## category name description sampleListId studyId
## <chr> <chr> <chr> <chr> <chr>
## 1 all_cases_with_mutation_and_cna_and_mr… Comp… Samples wi… acc_tcga_3w… acc_tc…
## 2 all_cases_in_study All … All sample… acc_tcga_all acc_tc…
## 3 all_cases_with_cna_data Samp… Samples wi… acc_tcga_cna acc_tc…
## 4 all_cases_with_mutation_and_cna_data Samp… Samples wi… acc_tcga_cn… acc_tc…
## 5 all_cases_with_methylation_data Samp… Samples wi… acc_tcga_me… acc_tc…
## 6 all_cases_with_methylation_data Samp… Samples wi… acc_tcga_me… acc_tc…
## 7 all_cases_with_mrna_rnaseq_data Samp… Samples wi… acc_tcga_rn… acc_tc…
## 8 all_cases_with_rppa_data Samp… Samples pr… acc_tcga_rp… acc_tc…
## 9 all_cases_with_mutation_data Samp… Samples wi… acc_tcga_se… 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 × 6
## uniqueSampleKey uniquePatientKey sampleType sampleId patientId studyId
## <chr> <chr> <chr> <chr> <chr> <chr>
## 1 VENHQS1PUi1BNUoxLTAxO… VENHQS1PUi1BNUo… Primary S… TCGA-OR… TCGA-OR-… acc_tc…
## 2 VENHQS1PUi1BNUoyLTAxO… VENHQS1PUi1BNUo… Primary S… TCGA-OR… TCGA-OR-… acc_tc…
## 3 VENHQS1PUi1BNUozLTAxO… VENHQS1PUi1BNUo… Primary S… TCGA-OR… TCGA-OR-… acc_tc…
## 4 VENHQS1PUi1BNUo0LTAxO… VENHQS1PUi1BNUo… Primary S… TCGA-OR… TCGA-OR-… acc_tc…
## 5 VENHQS1PUi1BNUo1LTAxO… VENHQS1PUi1BNUo… Primary S… TCGA-OR… TCGA-OR-… acc_tc…
## 6 VENHQS1PUi1BNUo2LTAxO… VENHQS1PUi1BNUo… Primary S… TCGA-OR… TCGA-OR-… acc_tc…
## 7 VENHQS1PUi1BNUo3LTAxO… VENHQS1PUi1BNUo… Primary S… TCGA-OR… TCGA-OR-… acc_tc…
## 8 VENHQS1PUi1BNUo4LTAxO… VENHQS1PUi1BNUo… Primary S… TCGA-OR… TCGA-OR-… acc_tc…
## 9 VENHQS1PUi1BNUo5LTAxO… VENHQS1PUi1BNUo… Primary S… TCGA-OR… TCGA-OR-… acc_tc…
## 10 VENHQS1PUi1BNUpBLTAxO… VENHQS1PUi1BNUp… Primary S… TCGA-OR… TCGA-OR-… acc_tc…
## # ℹ 82 more rows
Info on Samples
getSampleInfo(cbio, studyId = "acc_tcga",
sampleListIds = c("acc_tcga_rppa", "acc_tcga_cna"))## # A tibble: 136 × 6
## uniqueSampleKey uniquePatientKey sampleType sampleId patientId studyId
## <chr> <chr> <chr> <chr> <chr> <chr>
## 1 VENHQS1PUi1BNUoyLTAxO… VENHQS1PUi1BNUo… Primary S… TCGA-OR… TCGA-OR-… acc_tc…
## 2 VENHQS1PUi1BNUozLTAxO… VENHQS1PUi1BNUo… Primary S… TCGA-OR… TCGA-OR-… acc_tc…
## 3 VENHQS1PUi1BNUo2LTAxO… VENHQS1PUi1BNUo… Primary S… TCGA-OR… TCGA-OR-… acc_tc…
## 4 VENHQS1PUi1BNUo3LTAxO… VENHQS1PUi1BNUo… Primary S… TCGA-OR… TCGA-OR-… acc_tc…
## 5 VENHQS1PUi1BNUo4LTAxO… VENHQS1PUi1BNUo… Primary S… TCGA-OR… TCGA-OR-… acc_tc…
## 6 VENHQS1PUi1BNUo5LTAxO… VENHQS1PUi1BNUo… Primary S… TCGA-OR… TCGA-OR-… acc_tc…
## 7 VENHQS1PUi1BNUpBLTAxO… VENHQS1PUi1BNUp… Primary S… TCGA-OR… TCGA-OR-… acc_tc…
## 8 VENHQS1PUi1BNUpQLTAxO… VENHQS1PUi1BNUp… Primary S… TCGA-OR… TCGA-OR-… acc_tc…
## 9 VENHQS1PUi1BNUpSLTAxO… VENHQS1PUi1BNUp… Primary S… TCGA-OR… TCGA-OR-… acc_tc…
## 10 VENHQS1PUi1BNUpTLTAxO… VENHQS1PUi1BNUp… Primary S… TCGA-OR… TCGA-OR-… acc_tc…
## # ℹ 126 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: 2026-03-25 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
## R Under development (unstable) (2026-03-22 r89674)
## Platform: x86_64-pc-linux-gnu
## Running under: Ubuntu 24.04.4 LTS
##
## Matrix products: default
## BLAS: /usr/lib/x86_64-linux-gnu/openblas-pthread/libblas.so.3
## LAPACK: /usr/lib/x86_64-linux-gnu/openblas-pthread/libopenblasp-r0.3.26.so; LAPACK version 3.12.0
##
## locale:
## [1] LC_CTYPE=en_US.UTF-8 LC_NUMERIC=C
## [3] LC_TIME=en_US.UTF-8 LC_COLLATE=en_US.UTF-8
## [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
##
## time zone: Etc/UTC
## tzcode source: system (glibc)
##
## attached base packages:
## [1] stats4 stats graphics grDevices utils datasets methods
## [8] base
##
## other attached packages:
## [1] cBioPortalData_2.22.3 MultiAssayExperiment_1.37.2
## [3] SummarizedExperiment_1.41.1 Biobase_2.71.0
## [5] GenomicRanges_1.63.1 Seqinfo_1.1.0
## [7] IRanges_2.45.0 S4Vectors_0.49.0
## [9] BiocGenerics_0.57.0 generics_0.1.4
## [11] MatrixGenerics_1.23.0 matrixStats_1.5.0
## [13] AnVIL_1.23.8 AnVILBase_1.5.1
## [15] dplyr_1.2.0 BiocStyle_2.39.0
##
## loaded via a namespace (and not attached):
## [1] DBI_1.3.0 bitops_1.0-9
## [3] httr2_1.2.2 formatR_1.14
## [5] rlang_1.1.7 magrittr_2.0.4
## [7] otel_0.2.0 compiler_4.6.0
## [9] RSQLite_2.4.6 GenomicFeatures_1.63.1
## [11] png_0.1-9 systemfonts_1.3.2
## [13] vctrs_0.7.2 rvest_1.0.5
## [15] stringr_1.6.0 crayon_1.5.3
## [17] pkgconfig_2.0.3 fastmap_1.2.0
## [19] dbplyr_2.5.2 XVector_0.51.0
## [21] utf8_1.2.6 Rsamtools_2.27.1
## [23] promises_1.5.0 rmarkdown_2.30
## [25] tzdb_0.5.0 UCSC.utils_1.7.1
## [27] ragg_1.5.2 purrr_1.2.1
## [29] bit_4.6.0 xfun_0.57
## [31] cachem_1.1.0 cigarillo_1.1.0
## [33] GenomeInfoDb_1.47.2 jsonlite_2.0.0
## [35] blob_1.3.0 later_1.4.8
## [37] DelayedArray_0.37.0 BiocParallel_1.45.0
## [39] parallel_4.6.0 R6_2.6.1
## [41] bslib_0.10.0 stringi_1.8.7
## [43] rtracklayer_1.71.3 jquerylib_0.1.4
## [45] Rcpp_1.1.1 bookdown_0.46
## [47] knitr_1.51 readr_2.2.0
## [49] BiocBaseUtils_1.13.0 httpuv_1.6.17
## [51] Matrix_1.7-5 tidyselect_1.2.1
## [53] abind_1.4-8 yaml_2.3.12
## [55] codetools_0.2-20 miniUI_0.1.2
## [57] curl_7.0.0 lattice_0.22-9
## [59] tibble_3.3.1 withr_3.0.2
## [61] KEGGREST_1.51.1 shiny_1.13.0
## [63] evaluate_1.0.5 desc_1.4.3
## [65] lambda.r_1.2.4 futile.logger_1.4.9
## [67] BiocFileCache_3.1.0 xml2_1.5.2
## [69] Biostrings_2.79.5 pillar_1.11.1
## [71] BiocManager_1.30.27 filelock_1.0.3
## [73] DT_0.34.0 TCGAutils_1.31.4
## [75] RCurl_1.98-1.18 hms_1.1.4
## [77] xtable_1.8-8 RTCGAToolbox_2.41.0
## [79] glue_1.8.0 tools_4.6.0
## [81] BiocIO_1.21.0 data.table_1.18.2.1
## [83] GenomicAlignments_1.47.0 rapiclient_0.1.8
## [85] XML_3.99-0.23 fs_2.0.1
## [87] grid_4.6.0 tidyr_1.3.2
## [89] AnnotationDbi_1.73.0 RaggedExperiment_1.35.0
## [91] RJSONIO_2.0.0 restfulr_0.0.16
## [93] cli_3.6.5 rappdirs_0.3.4
## [95] textshaping_1.0.5 futile.options_1.0.1
## [97] GenomicDataCommons_1.35.1 S4Arrays_1.11.1
## [99] GCPtools_1.1.0 sass_0.4.10
## [101] digest_0.6.39 SparseArray_1.11.11
## [103] rjson_0.2.23 htmlwidgets_1.6.4
## [105] memoise_2.0.1 htmltools_0.5.9
## [107] pkgdown_2.2.0 lifecycle_1.0.5
## [109] httr_1.4.8 mime_0.13
## [111] bit64_4.6.0-1