cBioPortal: The R interface to the cBioPortal API Data Service

Description Usage Arguments Value API Metadata Patient Data Molecular Profiles Mutation Data Sample Data Gene Panels Examples

View source: R/cBioPortal.R

Description

This section of the documentation lists the functions that allow users to access the cBioPortal API. The main representation of the API can be obtained from the 'cBioPortal' function. The supporting functions listed here give access to specific parts of the API and allow the user to explore the API with individual calls. Many of the functions here are listed for documentation purposes and are recommended for advanced usage only. Users should only need to use the 'cBioPortalData' main function to obtain data.

Usage

 1
 2
 3
 4
 5
 6
 7
 8
 9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
cBioPortal(
  hostname = "www.cbioportal.org",
  protocol = "https",
  api. = "/api/api-docs"
)

getStudies(api)

clinicalData(api, studyId = NA_character_)

molecularProfiles(
  api,
  studyId = NA_character_,
  projection = c("SUMMARY", "ID", "DETAILED", "META")
)

mutationData(
  api,
  molecularProfileIds = NA_character_,
  entrezGeneIds = NULL,
  sampleIds = NULL
)

molecularData(
  api,
  molecularProfileIds = NA_character_,
  entrezGeneIds = NULL,
  sampleIds = NULL
)

searchOps(api, keyword)

geneTable(api, pageSize = 1000, pageNumber = 0, ...)

samplesInSampleLists(api, sampleListIds = NA_character_)

sampleLists(api, studyId = NA_character_)

allSamples(api, studyId = NA_character_)

genePanels(api)

getGenePanel(api, genePanelId = NA_character_)

genePanelMolecular(
  api,
  molecularProfileId = NA_character_,
  sampleListId = NULL,
  sampleIds = NULL
)

getGenePanelMolecular(api, molecularProfileIds = NA_character_, sampleIds)

getSampleInfo(
  api,
  studyId = NA_character_,
  sampleListIds = NULL,
  projection = c("SUMMARY", "ID", "DETAILED", "META")
)

getDataByGenePanel(
  api,
  studyId = NA_character_,
  genePanelId = NA_character_,
  molecularProfileIds = NULL,
  sampleListId = NULL,
  sampleIds = NULL
)

Arguments

hostname

character(1) The internet location of the service (default: 'www.cbioportal.org')

protocol

character(1) The internet protocol used to access the hostname (default: 'https')

api.

character(1) The directory location of the API protocol within the hostname (default: '/api/api-docs')

api

An API object of class 'cBioPortal' from the 'cBioPortal' function

studyId

character(1) Indicates the "studyId" as taken from 'getStudies'

projection

character(default: "SUMMARY") Specify the projection type for data retrieval for details see API documentation

molecularProfileIds

character() A vector of molecular profile IDs

entrezGeneIds

numeric() A vector indicating entrez gene IDs

sampleIds

character() Sample identifiers

keyword

character(1) Keyword or pattern for searching through available operations

pageSize

numeric(1) The number of rows in the table to return

pageNumber

numeric(1) The pagination page number

...

Additional arguments to lower level API functions

sampleListIds

character() A vector of 'sampleListId' as obtained from 'sampleLists'

genePanelId

character(1) Identifies the gene panel, as obtained from the 'genePanels' function

molecularProfileId

character(1) Indicates a molecular profile ID

sampleListId

character(1) A sample list identifier as obtained from 'sampleLists()“

Value

cBioPortal: An API object of class 'cBioPortal'

cBioPortalData: A data object of class 'MultiAssayExperiment'

API Metadata

* getStudies - Obtain a table of studies and associated metadata

* searchOps - Search through API operations with a keyword

* geneTable - Get a table of all genes by 'entrezGeneId' or 'hugoGeneSymbol'

* sampleLists - obtain all 'sampleListIds' for a particular 'studyId'

* allSamples - obtain all samples within a particular 'studyId'

* genePanels - Show all available gene panels

Patient Data

* clinicalData - Obtain clinical data for a particular study identifier ('studyId')

Molecular Profiles

* molecularProfiles - Produce a molecular profiles dataset for a given study identifier ('studyId')

* molecularData - Produce a dataset of molecular profile data based on 'molecularProfileId', 'entrezGeneIds', and 'sampleIds'

Mutation Data

* mutationData - Produce a dataset of mutation data using 'molecularProfileId', 'entrezGeneIds', and 'sampleIds'

Sample Data

* samplesInSampleLists - get all samples associated with a 'sampleListId'

* getSampleInfo - Obtain sample metadata for a particular 'studyId' or 'sampleListId'

Gene Panels

* getGenePanels - Obtain the gene panel for a particular 'genePanelId'

* genePanelMolecular - get gene panel data for a paricular 'molecularProfileId' and 'sampleListId' combination

* getGenePanelMolecular - get gene panel data for a combination of 'molecularProfileId' and 'sampleListId' vectors

* getDataByGenePanel - Download data for a gene panel and 'molecularProfileId' combination, optionally a 'sampleListId' can be provided.

Examples

 1
 2
 3
 4
 5
 6
 7
 8
 9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
cbio <- cBioPortal()

getStudies(api = cbio)

searchOps(api = cbio, keyword = "molecular")

## obtain clinical data
acc_clin <- clinicalData(api = cbio, studyId = "acc_tcga")
acc_clin

molecularProfiles(api = cbio, studyId = "acc_tcga")

genePanels(cbio)

(gp <- getGenePanel(cbio, "AmpliSeq"))

muts <- mutationData(
    api = cbio,
    molecularProfileIds = "acc_tcga_mutations",
    entrezGeneIds = 1:1000,
    sampleIds = c("TCGA-OR-A5J1-01", "TCGA-OR-A5J2-01")
)
exps <- molecularData(
    api = cbio,
    molecularProfileIds = c("acc_tcga_rna_seq_v2_mrna", "acc_tcga_rppa"),
    entrezGeneIds = 1:1000,
    sampleIds = c("TCGA-OR-A5J1-01", "TCGA-OR-A5J2-01")
)

sampleLists(api = cbio, studyId = "acc_tcga")

samplesInSampleLists(
    api = cbio,
    sampleListIds = c("acc_tcga_rppa", "acc_tcga_cnaseq")
)

genePanels(api = cbio)

getGenePanel(api = cbio, genePanelId = "IMPACT341")


getDataByGenePanel(cbio, studyId = "acc_tcga", genePanelId = "IMPACT341",
   molecularProfileId = "acc_tcga_rppa", sampleListId = "acc_tcga_rppa")

cBioPortalData documentation built on April 17, 2021, 6:07 p.m.