Recently the TCGA data has been moved from the DCC server to The National Cancer Institute (NCI) Genomic Data Commons (GDC) Data Portal In this version of the package, we rewrote all the functions that were acessing the old TCGA server to GDC.
The GDC, which receives, processes, harmonizes, and distributes clinical, biospecimen, and genomic data from multiple cancer research programs, has data from the following programs:
The big change is that the GDC data is harmonized against GRCh38. However, not all data has been harmonized yet. The old TCGA data can be acessed through GDC legacy Archive, in which the majority of data can be found.
More information about the project can be found in GDC FAQS
The functions TCGAquery
, TCGAdownload
, TCGAPrepare
, TCGAquery_maf
, TCGAquery_clinical
, were replaced by GDCquery
, GDCdownload
, GDCprepare
, GDCquery_maf
, GDCquery_clinical
.
And it can acess both the GDC and GDC Legacy Archive.
Note: Not all the examples in this vignette were updated.
Motivation: The Cancer Genome Atlas (TCGA) provides us with an enormous collection of data sets, not only spanning a large number of cancers but also a large number of experimental platforms. Even though the data can be accessed and downloaded from the database, the possibility to analyse these downloaded data directly in one single R package has not yet been available.
TCGAbiolinks consists of three parts or levels. Firstly, we provide different options to query and download from TCGA relevant data from all currently platforms and their subsequent pre-processing for commonly used bio-informatics (tools) packages in Bioconductor or CRAN. Secondly, the package allows to integrate different data types and it can be used for different types of analyses dealing with all platforms such as diff.expression, network inference or survival analysis, etc, and then it allows to visualize the obtained results. Thirdly we added a social level where a researcher can found a similar intereset in a bioinformatic community, and allows both to find a validation of results in literature in pubmed and also to retrieve questions and answers from site such as support.bioconductor.org, biostars.org, stackoverflow,etc.
This document describes how to search, download and analyze TCGA data using the TCGAbiolinks
package.
To install use the code below.
source("https://bioconductor.org/biocLite.R")
biocLite("TCGAbiolinks")
For a Graphical User Interface, please see TCGAbiolinksGUI. The GUI in under review and will soon be available in Bioconductor repository.
Please cite TCGAbiolinks package:
Related publications to this package:
Also, if you have used ELMER analysis please cite:
GDCquery
: Searching TCGA open-access dataGDCquery
: Searching GDC data for downloadYou can easily search GDC data using the GDCquery
function.
Using a summary of filters as used in the TCGA portal, the function works with the following arguments:
The next subsections will detail each of the search arguments. Below, we show some search examples:
#---------------------------------------------------------------
# For available entries and combinations please se table below
#---------------------------------------------------------------
# Gene expression aligned against Hg38
query <- GDCquery(project = "TARGET-AML",
data.category = "Transcriptome Profiling",
data.type = "Gene Expression Quantification",
workflow.type = "HTSeq - Counts")
# All DNA methylation data for TCGA-GBM and TCGA-GBM
query.met <- GDCquery(project = c("TCGA-GBM","TCGA-LGG"),
legacy = TRUE,
data.category = "DNA methylation",
platform = c("Illumina Human Methylation 450", "Illumina Human Methylation 27"))
# Using sample type to get only Primary solid Tumor samples and Solid Tissue Normal
query.mirna <- GDCquery(project = "TCGA-ACC",
data.category = "Transcriptome Profiling",
data.type = "miRNA Expression Quantification",
sample.type = c("Primary solid Tumor","Solid Tissue Normal"))
# Example Using legacy to accessing hg19 and filtering by barcode
query <- GDCquery(project = "TCGA-GBM",
data.category = "DNA methylation",
platform = "Illumina Human Methylation 27",
legacy = TRUE,
barcode = c("TCGA-02-0047-01A-01D-0186-05","TCGA-06-2559-01A-01D-0788-05"))
# Gene expression aligned against hg19.
query.exp.hg19 <- GDCquery(project = "TCGA-GBM",
data.category = "Gene expression",
data.type = "Gene expression quantification",
platform = "Illumina HiSeq",
file.type = "normalized_results",
experimental.strategy = "RNA-Seq",
barcode = c("TCGA-14-0736-02A-01R-2005-01", "TCGA-06-0211-02A-02R-2005-01"),
legacy = TRUE)
# Searching idat file for DNA methylation
query <- GDCquery(project = "TCGA-OV",
data.category = "Raw microarray data",
data.type = "Raw intensities",
experimental.strategy = "Methylation array",
legacy = TRUE,
file.type = ".idat",
platform = "Illumina Human Methylation 450")
The list of projects is below:
disease_type.0 | project_id | name | primary_site.0 | released | state | dbgap_accession_number | id | tumor |
---|---|---|---|---|---|---|---|---|
Cholangiocarcinoma | TCGA-CHOL | Cholangiocarcinoma | Bile Duct | True | legacy | NA | TCGA-CHOL | CHOL |
Ovarian Serous Cystadenocarcinoma | TCGA-OV | Ovarian Serous Cystadenocarcinoma | Ovary | True | legacy | NA | TCGA-OV | OV |
Liver Hepatocellular Carcinoma | TCGA-LIHC | Liver Hepatocellular Carcinoma | Liver | True | legacy | NA | TCGA-LIHC | LIHC |
Head and Neck Squamous Cell Carcinoma | TCGA-HNSC | Head and Neck Squamous Cell Carcinoma | Head and Neck | True | legacy | NA | TCGA-HNSC | HNSC |
Colon Adenocarcinoma | TCGA-COAD | Colon Adenocarcinoma | Colorectal | True | legacy | NA | TCGA-COAD | COAD |
Adrenocortical Carcinoma | TCGA-ACC | Adrenocortical Carcinoma | Adrenal Gland | True | legacy | NA | TCGA-ACC | ACC |
Mesothelioma | TCGA-MESO | Mesothelioma | Pleura | True | legacy | NA | TCGA-MESO | MESO |
Acute Myeloid Leukemia | TARGET-AML | Acute Myeloid Leukemia | Blood | True | legacy | phs000465 | TARGET-AML | AML |
Neuroblastoma | TARGET-NBL | Neuroblastoma | Nervous System | True | legacy | phs000467 | TARGET-NBL | NBL |
Brain Lower Grade Glioma | TCGA-LGG | Brain Lower Grade Glioma | Brain | True | legacy | NA | TCGA-LGG | LGG |
Pancreatic Adenocarcinoma | TCGA-PAAD | Pancreatic Adenocarcinoma | Pancreas | True | legacy | NA | TCGA-PAAD | PAAD |
High-Risk Wilms Tumor | TARGET-WT | High-Risk Wilms Tumor | Kidney | True | legacy | phs000471 | TARGET-WT | WT |
Uterine Carcinosarcoma | TCGA-UCS | Uterine Carcinosarcoma | Uterus | True | legacy | NA | TCGA-UCS | UCS |
Uveal Melanoma | TCGA-UVM | Uveal Melanoma | Eye | True | legacy | NA | TCGA-UVM | UVM |
Stomach Adenocarcinoma | TCGA-STAD | Stomach Adenocarcinoma | Stomach | True | legacy | NA | TCGA-STAD | STAD |
Lung Squamous Cell Carcinoma | TCGA-LUSC | Lung Squamous Cell Carcinoma | Lung | True | legacy | NA | TCGA-LUSC | LUSC |
Uterine Corpus Endometrial Carcinoma | TCGA-UCEC | Uterine Corpus Endometrial Carcinoma | Uterus | True | legacy | NA | TCGA-UCEC | UCEC |
Esophageal Carcinoma | TCGA-ESCA | Esophageal Carcinoma | Esophagus | True | legacy | NA | TCGA-ESCA | ESCA |
Rhabdoid Tumor | TARGET-RT | Rhabdoid Tumor | Kidney | True | legacy | phs000470 | TARGET-RT | RT |
Osteosarcoma | TARGET-OS | Osteosarcoma | Bone | True | legacy | phs000468 | TARGET-OS | OS |
Prostate Adenocarcinoma | TCGA-PRAD | Prostate Adenocarcinoma | Prostate | True | legacy | NA | TCGA-PRAD | PRAD |
Kidney Chromophobe | TCGA-KICH | Kidney Chromophobe | Kidney | True | legacy | NA | TCGA-KICH | KICH |
Testicular Germ Cell Tumors | TCGA-TGCT | Testicular Germ Cell Tumors | Testis | True | legacy | NA | TCGA-TGCT | TGCT |
Acute Myeloid Leukemia | TCGA-LAML | Acute Myeloid Leukemia | Bone Marrow | True | legacy | NA | TCGA-LAML | LAML |
Thymoma | TCGA-THYM | Thymoma | Thymus | True | legacy | NA | TCGA-THYM | THYM |
Lymphoid Neoplasm Diffuse Large B-cell Lymphoma | TCGA-DLBC | Lymphoid Neoplasm Diffuse Large B-cell Lymphoma | Lymph Nodes | True | legacy | NA | TCGA-DLBC | DLBC |
Pheochromocytoma and Paraganglioma | TCGA-PCPG | Pheochromocytoma and Paraganglioma | Adrenal Gland | True | legacy | NA | TCGA-PCPG | PCPG |
Cervical Squamous Cell Carcinoma and Endocervical Adenocarcinoma | TCGA-CESC | Cervical Squamous Cell Carcinoma and Endocervical Adenocarcinoma | Cervix | True | legacy | NA | TCGA-CESC | CESC |
Lung Adenocarcinoma | TCGA-LUAD | Lung Adenocarcinoma | Lung | True | legacy | NA | TCGA-LUAD | LUAD |
Bladder Urothelial Carcinoma | TCGA-BLCA | Bladder Urothelial Carcinoma | Bladder | True | legacy | NA | TCGA-BLCA | BLCA |
Kidney Renal Papillary Cell Carcinoma | TCGA-KIRP | Kidney Renal Papillary Cell Carcinoma | Kidney | True | legacy | NA | TCGA-KIRP | KIRP |
Clear Cell Sarcoma of the Kidney | TARGET-CCSK | Clear Cell Sarcoma of the Kidney | Kidney | True | legacy | phs000466 | TARGET-CCSK | CCSK |
Sarcoma | TCGA-SARC | Sarcoma | Soft Tissue | True | legacy | NA | TCGA-SARC | SARC |
Thyroid Carcinoma | TCGA-THCA | Thyroid Carcinoma | Thyroid | True | legacy | NA | TCGA-THCA | THCA |
Rectum Adenocarcinoma | TCGA-READ | Rectum Adenocarcinoma | Colorectal | True | legacy | NA | TCGA-READ | READ |
Skin Cutaneous Melanoma | TCGA-SKCM | Skin Cutaneous Melanoma | Skin | True | legacy | NA | TCGA-SKCM | SKCM |
Glioblastoma Multiforme | TCGA-GBM | Glioblastoma Multiforme | Brain | True | legacy | NA | TCGA-GBM | GBM |
Kidney Renal Clear Cell Carcinoma | TCGA-KIRC | Kidney Renal Clear Cell Carcinoma | Kidney | True | legacy | NA | TCGA-KIRC | KIRC |
Breast Invasive Carcinoma | TCGA-BRCA | Breast Invasive Carcinoma | Breast | True | legacy | NA | TCGA-BRCA | BRCA |
The list of sample.type is below:
tissue.code | shortLetterCode | tissue.definition |
---|---|---|
01 | TP | Primary solid Tumor |
02 | TR | Recurrent Solid Tumor |
03 | TB | Primary Blood Derived Cancer - Peripheral Blood |
04 | TRBM | Recurrent Blood Derived Cancer - Bone Marrow |
05 | TAP | Additional - New Primary |
06 | TM | Metastatic |
07 | TAM | Additional Metastatic |
08 | THOC | Human Tumor Original Cells |
09 | TBM | Primary Blood Derived Cancer - Bone Marrow |
10 | NB | Blood Derived Normal |
11 | NT | Solid Tissue Normal |
12 | NBC | Buccal Cell Normal |
13 | NEBV | EBV Immortalized Normal |
14 | NBM | Bone Marrow Normal |
20 | CELLC | Control Analyte |
40 | TRB | Recurrent Blood Derived Cancer - Peripheral Blood |
50 | CELL | Cell Lines |
60 | XP | Primary Xenograft Tissue |
61 | XCL | Cell Line Derived Xenograft Tissue |
The other fields (data.category, data.type, workflow.type, platform, file.type) can be found below. Please, not that these tables are still incomplete.
data.category | data.type | workflow.type | platform |
---|---|---|---|
Transcriptome Profiling | Gene Expression Quantification | HTSeq - Counts | |
HTSeq - FPKM-UQ | |||
HTSeq - FPKM | |||
Isoform Expression Quantification | - | ||
miRNA Expression Quantification | - | ||
Copy Number Variation | Copy Number Segment | ||
Masked Copy Number Segment | |||
Simple Nucleotide Variation | |||
Raw Sequencing Data | |||
Biospecimen | |||
Clinical | |||
DNA Methylation | Illumina Human Methylation 450 | ||
Illumina Human Methylation 27 |
data.category | data.type | platform | file.type |
---|---|---|---|
Copy number variation | - | Affymetrix SNP Array 6.0 | nocnv_hg18.seg |
- | Affymetrix SNP Array 6.0 | hg18.seg | |
- | Affymetrix SNP Array 6.0 | nocnv_hg19.seg | |
- | Affymetrix SNP Array 6.0 | hg19.seg | |
- | Illumina HiSeq | - | |
Simple nucleotide variation | |||
Raw sequencing data | |||
Biospecimen | |||
Clinical | |||
Protein expression | MDA RPPA Core | - | |
Gene expression | Gene expression quantification | Illumina HiSeq | normalized_results |
Illumina HiSeq | results | ||
HT_HG-U133A | - | ||
AgilentG4502A_07_2 | - | ||
AgilentG4502A_07_1 | - | ||
HuEx-1_0-st-v2 | FIRMA.txt | ||
gene.txt | |||
Isoform expression quantification | |||
miRNA gene quantification | hg19.mirna | ||
hg19.mirbase20 | |||
mirna | |||
Exon junction quantification | |||
Exon quantification | |||
miRNA isoform quantification | hg19.isoform | ||
isoform | |||
DNA methylation | Illumina Human Methylation 450 | - | |
Illumina Human Methylation 27 | - | ||
Illumina DNA Methylation OMA003 CPI | - | ||
Illumina DNA Methylation OMA002 CPI | - | ||
Illumina Hi Seq | |||
Raw microarray data | Raw intensities | Illumina Human Methylation 450 | idat |
Illumina Human Methylation 27 | idat | ||
Other |
GDCquery_maf
: Mutation Annotation Format (MAF) aligned against hg38In order to download the Mutation Annotation Format (MAF) aligned against hg38, we provide the user with the GDCquery_maf
function. Briefly, it will get the open acess maf files from https://gdc-docs.nci.nih.gov/Data/Release_Notes/Data_Release_Notes/. Four separate variant calling pipelines are implemented for GDC data harmonization which are described here
The arguments that will be used to filter the data are:
acc.muse.maf <- GDCquery_Maf("ACC", pipelines = "muse")
acc.varscan2.maf <- GDCquery_Maf("ACC", pipelines = "varscan2")
acc.somaticsniper.maf <- GDCquery_Maf("ACC", pipelines = "somaticsniper")
acc.mutect.maf <- GDCquery_Maf("ACC", pipelines = "mutect")
This mutation data can be visualized through a oncoprint from the ComplexHeatmap package (Gu, Zuguang and Eils, Roland and Schlesner, Matthias 2016).
mut <- GDCquery_Maf("ACC", pipelines = "mutect2")
clin <- GDCquery_clinic("TCGA-ACC","clinical")
clin <- clin[,c("bcr_patient_barcode","disease","gender","tumor_stage","race","vital_status")]
TCGAvisualize_oncoprint(mut = mut, genes = mut$Hugo_Symbol[1:20],
filename = "oncoprint.pdf",
annotation = clin,
color=c("background"="#CCCCCC","DEL"="purple","INS"="yellow","SNP"="brown"),
rows.font.size=10,
width = 10,
heatmap.legend.side = "right",
dist.col = 0,
label.font.size = 10)
In order to download the Mutation Annotation Format (MAF) aligned against hg19, the user will need to use GDCquery
with legacy = TRUE
, GDCdownload
and GDCprepare
.
query.maf.hg19 <- GDCquery(project = "TCGA-COAD",
data.category = "Simple nucleotide variation",
data.type = "Simple somatic mutation",
access = "open",
legacy = TRUE)
# Check maf availables
knitr::kable(getResults(query.maf.hg19)[,c("created_datetime","file_name")])
created_datetime | file_name |
---|---|
2016-06-13T16:30:24.763822-05:00 | gsc_COAD_pairs.aggregated.capture.tcga.uuid.automated.somatic.maf |
NA | hgsc.bcm.edu_COAD.SOLiD_DNASeq.1.somatic.maf |
query.maf.hg19 <- GDCquery(project = "TCGA-COAD",
data.category = "Simple nucleotide variation",
data.type = "Simple somatic mutation",
access = "open",
file.type = "gsc_COAD_pairs.aggregated.capture.tcga.uuid.automated.somatic.maf",
legacy = TRUE)
GDCdownload(query.maf.hg19)
coad.mutect.maf <- GDCprepare(query.maf.hg19)
GDCquery_clinic
and GDCprepare_clinic
: Working with clinical data.In GDC database the clinical data can be retrieved in two sources:
there are two main differences:
clin.query <- GDCquery(project = "TCGA-BLCA", data.category = "Clinical", barcode = "TCGA-FD-A5C0")
json <- tryCatch(GDCdownload(clin.query),
error = function(e) GDCdownload(clin.query, method = "client"))
clinical.patient <- GDCprepare_clinic(clin.query, clinical.info = "patient")
clinical.patient.followup <- GDCprepare_clinic(clin.query, clinical.info = "follow_up")
clinical.index <- GDCquery_clinic("TCGA-BLCA")
# Example of the second difference:
clinical.patient[,c("vital_status","days_to_death","days_to_last_followup")]
## vital_status days_to_death days_to_last_followup
## 1 Alive NA 209
clinical.patient.followup[,c("vital_status","days_to_death","days_to_last_followup")]
## vital_status days_to_death days_to_last_followup
## 1 Alive NA 407
## 2 Dead 550 NA
# indexed data is equivalent to follow ups information
clinical.index[clinical.index$submitter_id=="TCGA-FD-A5C0",
c("vital_status","days_to_death","days_to_last_follow_up")]
## vital_status days_to_death days_to_last_follow_up
## 197 dead 550 407
You can retrieve clinical data using the GDCquery_clinic
function. This will get only the indexed GDC clinical data. This is the same as clicking on the download clinical buttun in the data portal. This means this function is not parsing the XML files, see GDCprepare_clinic
function.
The arguments of this function are:
Examples of use:
clin <- GDCquery_clinic("TCGA-ACC", type = "clinical", save.csv = TRUE)
clin <- GDCquery_clinic("TCGA-ACC", type = "biospecimen", save.csv = TRUE)
To parse the TCGA clinical XML files, please use GDCprepare_clinic
function. It receives as argument the query object and which clinical information we should retrieve. The possible values for the argument clinical.info
are:
query <- GDCquery(project = "TCGA-COAD",
data.category = "Clinical",
barcode = c("TCGA-RU-A8FL","TCGA-AA-3972"))
GDCdownload(query)
clinical <- GDCprepare_clinic(query, clinical.info = "patient")
clinical.drug <- GDCprepare_clinic(query, clinical.info = "drug")
clinical.radiation <- GDCprepare_clinic(query, clinical.info = "radiation")
clinical.admin <- GDCprepare_clinic(query, clinical.info = "admin")
query <- GDCquery(project = "TCGA-COAD",
data.category = "Biospecimen",
barcode = c("TCGA-RU-A8FL","TCGA-AA-3972"))
GDCdownload(query)
clinical.admin <- GDCprepare_clinic(query, clinical.info = "admin")
clinical.sample <- GDCprepare_clinic(query, clinical.info = "sample")
clinical.slide <- GDCprepare_clinic(query, clinical.info = "slide")
clinical.portion <- GDCprepare_clinic(query, clinical.info = "portion")
To get all the information for TGCA samples you can use the script below:
# This code will get all clinical indexed data from TCGA
library(TCGAbiolinks)
library(data.table)
clinical <- TCGAbiolinks:::getGDCprojects()$project_id %>%
regexPipes::grep("TCGA",value=T) %>%
sort %>%
plyr::alply(1,GDCquery_clinic, .progress = "text") %>%
rbindlist
readr::write_csv(clinical,path = paste0("all_clin_indexed.csv"))
# This code will get all clinical XML data from TCGA
getclinical <- function(proj){
message(proj)
while(1){
result = tryCatch({
query <- GDCquery(project = proj, data.category = "Clinical")
GDCdownload(query)
clinical <- GDCprepare_clinic(query, clinical.info = "patient")
for(i in c("admin","radiation","follow_up","drug","new_tumor_event")){
message(i)
aux <- GDCprepare_clinic(query, clinical.info = i)
if(is.null(aux)) next
# add suffix manually if it already exists
replicated <- which(grep("bcr_patient_barcode",colnames(aux), value = T,invert = T) %in% colnames(clinical))
colnames(aux)[replicated] <- paste0(colnames(aux)[replicated],".",i)
if(!is.null(aux)) clinical <- merge(clinical,aux,by = "bcr_patient_barcode", all = TRUE)
}
readr::write_csv(clinical,path = paste0(proj,"_clinical_from_XML.csv")) # Save the clinical data into a csv file
return(clinical)
}, error = function(e) {
message(paste0("Error clinical: ", proj))
})
}
}
clinical <- TCGAbiolinks:::getGDCprojects()$project_id %>%
regexPipes::grep("TCGA",value=T) %>% sort %>%
plyr::alply(1,getclinical, .progress = "text") %>%
rbindlist(fill = TRUE) %>% setDF %>% subset(!duplicated(clinical))
readr::write_csv(clinical,path = "all_clin_XML.csv")
# result: https://drive.google.com/open?id=0B0-8N2fjttG-WWxSVE5MSGpva1U
# Obs: this table has multiple lines for each patient, as the patient might have several followups, drug treatments,
# new tumor events etc...
Also, some functions to work with clinical data are provided.
For example the function TCGAquery_SampleTypes
will filter barcodes based on a type the argument typesample. Some of the typesamples possibilities are: TP (PRIMARY SOLID TUMOR), TR (RECURRENT SOLID TUMOR), NT (Solid Tissue Normal) etc. Please, see ?TCGAquery_SampleTypes for all the possible values.
The function TCGAquery_MatchedCoupledSampleTypes
will filter the samples that have all the typesample provided as argument. For example, if TP and TR are set as typesample, the function will return the barcodes of a patient if it has both types. So, if it has a TP, but not a TR, no barcode will be returned. If it has a TP and a TR both barcodes are returned.
An example of the function is below:
bar <- c("TCGA-G9-6378-02A-11R-1789-07", "TCGA-CH-5767-04A-11R-1789-07",
"TCGA-G9-6332-60A-11R-1789-07", "TCGA-G9-6336-01A-11R-1789-07",
"TCGA-G9-6336-11A-11R-1789-07", "TCGA-G9-7336-11A-11R-1789-07",
"TCGA-G9-7336-04A-11R-1789-07", "TCGA-G9-7336-14A-11R-1789-07",
"TCGA-G9-7036-04A-11R-1789-07", "TCGA-G9-7036-02A-11R-1789-07",
"TCGA-G9-7036-11A-11R-1789-07", "TCGA-G9-7036-03A-11R-1789-07",
"TCGA-G9-7036-10A-11R-1789-07", "TCGA-BH-A1ES-10A-11R-1789-07",
"TCGA-BH-A1F0-10A-11R-1789-07", "TCGA-BH-A0BZ-02A-11R-1789-07",
"TCGA-B6-A0WY-04A-11R-1789-07", "TCGA-BH-A1FG-04A-11R-1789-08",
"TCGA-D8-A1JS-04A-11R-2089-08", "TCGA-AN-A0FN-11A-11R-8789-08",
"TCGA-AR-A2LQ-12A-11R-8799-08", "TCGA-AR-A2LH-03A-11R-1789-07",
"TCGA-BH-A1F8-04A-11R-5789-07", "TCGA-AR-A24T-04A-55R-1789-07",
"TCGA-AO-A0J5-05A-11R-1789-07", "TCGA-BH-A0B4-11A-12R-1789-07",
"TCGA-B6-A1KN-60A-13R-1789-07", "TCGA-AO-A0J5-01A-11R-1789-07",
"TCGA-AO-A0J5-01A-11R-1789-07", "TCGA-G9-6336-11A-11R-1789-07",
"TCGA-G9-6380-11A-11R-1789-07", "TCGA-G9-6380-01A-11R-1789-07",
"TCGA-G9-6340-01A-11R-1789-07", "TCGA-G9-6340-11A-11R-1789-07")
S <- TCGAquery_SampleTypes(bar,"TP")
S2 <- TCGAquery_SampleTypes(bar,"NB")
# Retrieve multiple tissue types NOT FROM THE SAME PATIENTS
SS <- TCGAquery_SampleTypes(bar,c("TP","NB"))
# Retrieve multiple tissue types FROM THE SAME PATIENTS
SSS <- TCGAquery_MatchedCoupledSampleTypes(bar,c("NT","TP"))
TCGAquery_subtype
: Working with molecular subtypes data.The Cancer Genome Atlas (TCGA) Research Network has reported integrated genome-wide studies of various diseases. We have added some of the subtypes defined by these report in our package:
These subtypes will be automatically added in the summarizedExperiment object through TCGAprepare. But you can also use the TCGAquery_subtype
function to retrive this information.
# Check with subtypes from TCGAprepare and update examples
GBM_path_subtypes <- TCGAquery_subtype(tumor = "gbm")
LGG_path_subtypes <- TCGAquery_subtype(tumor = "lgg")
A subset of the lgg subytpe is shown below:
GDCdownload
: Downloading GDC dataYou can easily download data using the GDCdownload
function. It uses GDC transfer tool to download gdc data doing a system call. For this reason some times the update will stop to show, which does not means that the download process has stopped. Once the process has finished it will give a signal to R.
The data from query will be save in a folder: project/source/data.category (where source is harmonized or legacy)
The arguments are:
GDCquery
functionGDCdownload
: Example of usequery <- GDCquery(project = "TCGA-ACC", data.category = "Copy Number Variation",
data.type = "Copy Number Segment",
barcode = c( "TCGA-OR-A5KU-01A-11D-A29H-01", "TCGA-OR-A5JK-01A-11D-A29H-01"))
GDCdownload(query)
data <- GDCprepare(query)
#--------------------------------------
# Gene expression
#--------------------------------------
# Aligned against Hg38
# mRNA pipeline: https://gdc-docs.nci.nih.gov/Data/Bioinformatics_Pipelines/Expression_mRNA_Pipeline/
query.exp.hg38 <- GDCquery(project = "TCGA-GBM",
data.category = "Transcriptome Profiling",
data.type = "Gene Expression Quantification",
workflow.type = "HTSeq - FPKM-UQ",
barcode = c("TCGA-14-0736-02A-01R-2005-01", "TCGA-06-0211-02A-02R-2005-01"))
GDCdownload(query.exp.hg38)
expdat <- GDCprepare(query = query.exp.hg38,
save = TRUE,
save.filename = "exp.rda")
# Aligned against Hg19
query.exp.hg19 <- GDCquery(project = "TCGA-GBM",
data.category = "Gene expression",
data.type = "Gene expression quantification",
platform = "Illumina HiSeq",
file.type = "normalized_results",
experimental.strategy = "RNA-Seq",
barcode = c("TCGA-14-0736-02A-01R-2005-01", "TCGA-06-0211-02A-02R-2005-01"),
legacy = TRUE)
GDCdownload(query.exp.hg19)
data <- GDCprepare(query.exp.hg19)
#--------------------------------------
# DNA methylation data
#--------------------------------------
# DNA methylation aligned to hg38
query_met.hg38 <- GDCquery(project= "TCGA-LGG",
data.category = "DNA Methylation",
platform = "Illumina Human Methylation 450",
barcode = c("TCGA-HT-8111-01A-11D-2399-05","TCGA-HT-A5R5-01A-11D-A28N-05"))
GDCdownload(query_met.hg38)
data.hg38 <- GDCprepare(query_met.hg38)
# DNA methylation aligned to hg19
query_meth.hg19 <- GDCquery(project= "TCGA-LGG",
data.category = "DNA methylation",
platform = "Illumina Human Methylation 450",
barcode = c("TCGA-HT-8111-01A-11D-2399-05","TCGA-HT-A5R5-01A-11D-A28N-05"),
legacy = TRUE)
GDCdownload(query_meth.hg19)
data.hg19 <- GDCprepare(query_meth.hg19)
# A function to download only 20 samples
legacyPipeline <- function(project, data.category, platform, n = 20){
query <- GDCquery(project = project,
data.category = data.category,
platform = platform,
legacy = TRUE)
cases <- getResults(query, rows = 1:n, cols="cases") # Get two samples from the search
query <- GDCquery(project = project,
data.category = data.category,
platform = platform,
legacy = TRUE,
barcode = cases)
GDCdownload(query,chunks.per.download = 5)
data <- GDCprepare(query)
return(data)
}
# DNA methylation
data <- legacyPipeline("TCGA-GBM","DNA methylation","Illumina Human Methylation 27")
data <- legacyPipeline("TCGA-GBM","DNA methylation","Illumina Human Methylation 450")
data <- legacyPipeline("TCGA-GBM","DNA methylation","Illumina DNA Methylation OMA003 CPI")
data <- legacyPipeline("TCGA-GBM","DNA methylation","Illumina DNA Methylation OMA002 CPI")
#-------------------------------------------------------
# Example to idat files from TCGA projects
#-------------------------------------------------------
projects <- TCGAbiolinks:::getGDCprojects()$project_id
projects <- projects[grepl('^TCGA',projects,perl=T)]
match.file.cases.all <- NULL
for(proj in projects){
print(proj)
query <- GDCquery(project = proj,
data.category = "Raw microarray data",
data.type = "Raw intensities",
experimental.strategy = "Methylation array",
legacy = TRUE,
file.type = ".idat",
platform = "Illumina Human Methylation 450")
match.file.cases <- getResults(query,cols=c("cases","file_name"))
match.file.cases$project <- proj
match.file.cases.all <- rbind(match.file.cases.all,match.file.cases)
tryCatch(GDCdownload(query, method = "api",chunks.per.download = 20),
error = function(e) GDCdownload(query, method = "client"))
}
# This will create a map between idat file name, cases (barcode) and project
readr::write_tsv(match.file.cases.all, path = "idat_filename_case.txt")
# code to move all files to local folder
for(file in dir(".",pattern = ".idat", recursive = T)){
TCGAbiolinks::move(file,basename(file))
}
GDCprepare
: Preparing the dataThis function is still under development, it is not working for all cases. See the tables below with the status. Examples of query, download, prepare can be found in this gist.
Data.category | Data.type | Workflow Type | Status |
---|---|---|---|
Transcriptome Profiling | Gene Expression Quantification | HTSeq - Counts | Data frame or SE (losing 5% of information when mapping to genomic regions) |
HTSeq - FPKM-UQ | Returning only a (losing 5% of information when mapping to genomic regions) | ||
HTSeq - FPKM | Returning only a (losing 5% of information when mapping to genomic regions) | ||
Isoform Expression Quantification | Not needed | ||
miRNA Expression Quantification | Not needed | Returning only a dataframe for the moment | |
Copy number variation | Copy Number Segment | Returning only a dataframe for the moment | |
Masked Copy Number Segment | Returning only a dataframe for the moment | ||
Simple Nucleotide Variation | |||
Raw Sequencing Data | |||
Biospecimen | |||
Clinical |
Data.category | Data.type | Platform | file.type | Status |
---|---|---|---|---|
Transcriptome Profiling | ||||
Copy number variation | - | Affymetrix SNP Array 6.0 | nocnv_hg18.seg | Working |
- | Affymetrix SNP Array 6.0 | hg18.seg | Working | |
- | Affymetrix SNP Array 6.0 | nocnv_hg19.seg | Working | |
- | Affymetrix SNP Array 6.0 | hg19.seg | Working | |
- | Illumina HiSeq | Several | Working | |
Simple Nucleotide Variation | Simple somatic mutation | |||
Raw Sequencing Data | ||||
Biospecimen | ||||
Clinical | ||||
Protein expression | MDA RPPA Core | - | Working | |
Gene expression | Gene expression quantification | Illumina HiSeq | normalized_results | Working |
Illumina HiSeq | results | Working | ||
HT_HG-U133A | - | Working | ||
AgilentG4502A_07_2 | - | Data frame only | ||
AgilentG4502A_07_1 | - | Data frame only | ||
HuEx-1_0-st-v2 | FIRMA.txt | Not Preparing | ||
gene.txt | Not Preparing | |||
Isoform expression quantification | ||||
miRNA gene quantification | ||||
Exon junction quantification | ||||
Exon quantification | ||||
miRNA isoform quantification | ||||
DNA methylation | Illumina Human Methylation 450 | Not used | Working | |
Illumina Human Methylation 27 | Not used | Working | ||
Illumina DNA Methylation OMA003 CPI | Not used | Working | ||
Illumina DNA Methylation OMA002 CPI | Not used | Working | ||
Illumina Hi Seq | Not working | |||
Raw Microarray Data | ||||
Structural Rearrangement | ||||
Other |
You can easily read the downloaded data using the GDCprepare
function. This function will prepare the data into a SummarizedExperiment (Huber, Wolfgang and Carey, Vincent J and Gentleman, Robert and Anders, Simon and Carlson, Marc and Carvalho, Benilton S and Bravo, Hector Corrada and Davis, Sean and Gatto, Laurent and Girke, Thomas and others 2015) object for downstream analysis. For the moment this function works only with data level 3.
The arguments are:
save
is TRUE
TRUE
library(TCGAbiolinks)
# Downloading and prepare
query <- GDCquery(project = "TARGET-AML",
data.category = "Transcriptome Profiling",
data.type = "Gene Expression Quantification",
workflow.type = "HTSeq - FPKM-UQ")
GDCdownload(query)
data <- GDCprepare(query)
# Downloading and prepare using legacy
query <- GDCquery(project = "TCGA-GBM",
data.category = "Protein expression",
legacy = TRUE,
barcode = c("TCGA-OX-A56R-01A-21-A44T-20","TCGA-08-0357-01A-21-1898-20"))
GDCdownload(query)
data <- GDCprepare(query, save = TRUE,
save.filename = "gbmProteinExpression.rda",
remove.files.prepared = TRUE)
# Downloading and prepare using legacy
query <- GDCquery(project = "TCGA-GBM",
data.category = "DNA methylation",
platform = "Illumina Human Methylation 27",legacy = TRUE,
barcode = c("TCGA-02-0047-01A-01D-0186-05","TCGA-06-2559-01A-01D-0788-05"))
GDCdownload(query)
data <- GDCprepare(query, add.gistic2.mut = c("PTEN","FOXJ1"))
# To view gistic and mutation information please access the samples information matrix in the summarized Experiment object
library(SummarizedExperiment)
samples.information <- colData(data)
TCGAanalyze
: Analyze data from TCGA.You can easily analyze data using following functions:
TCGAanalyze_Preprocessing
: Preprocessing of Gene Expression data (IlluminaHiSeq_RNASeqV2)You can easily search TCGA samples, download and prepare a matrix of gene expression.
# You can define a list of samples to query and download providing relative TCGA barcodes.
listSamples <- c("TCGA-E9-A1NG-11A-52R-A14M-07","TCGA-BH-A1FC-11A-32R-A13Q-07",
"TCGA-A7-A13G-11A-51R-A13Q-07","TCGA-BH-A0DK-11A-13R-A089-07",
"TCGA-E9-A1RH-11A-34R-A169-07","TCGA-BH-A0AU-01A-11R-A12P-07",
"TCGA-C8-A1HJ-01A-11R-A13Q-07","TCGA-A7-A13D-01A-13R-A12P-07",
"TCGA-A2-A0CV-01A-31R-A115-07","TCGA-AQ-A0Y5-01A-11R-A14M-07")
# Query platform Illumina HiSeq with a list of barcode
query <- GDCquery(project = "TCGA-BRCA",
data.category = "Gene expression",
data.type = "Gene expression quantification",
experimental.strategy = "RNA-Seq",
platform = "Illumina HiSeq",
file.type = "results",
barcode = listSamples,
legacy = TRUE)
# Download a list of barcodes with platform IlluminaHiSeq_RNASeqV2
GDCdownload(query)
# Prepare expression matrix with geneID in the rows and samples (barcode) in the columns
# rsem.genes.results as values
BRCARnaseqSE <- GDCprepare(query)
BRCAMatrix <- assay(BRCARnaseqSE,"raw_count") # or BRCAMatrix <- assay(BRCARnaseqSE,"raw_count")
# For gene expression if you need to see a boxplot correlation and AAIC plot to define outliers you can run
BRCARnaseq_CorOutliers <- TCGAanalyze_Preprocessing(BRCARnaseqSE)
The result is shown below:
TCGA-A7-A13G-11A-51R-A13Q-07 | TCGA-C8-A1HJ-01A-11R-A13Q-07 | |
---|---|---|
KCNJ16|3773 | 42.0 | 1064 |
SEPHS2|22928 | 2575.0 | 10767 |
LOC100268168|100268168 | 14.0 | 22 |
LOC100271722|100271722 | 170.0 | 114 |
DUSP26|78986 | 0.0 | 2 |
ITIH4|3700 | 68.4 | 31 |
IL17C|27189 | 0.0 | 1 |
METRN|79006 | 83.0 | 45 |
ACTL7B|10880 | 0.0 | 0 |
STT3B|201595 | 4611.0 | 3130 |
The result from TCGAanalyze_Preprocessing is shown below:
TCGAanalyze_DEA
&
TCGAanalyze_LevelTab
: Differential expression analysis (DEA)Perform DEA (Differential expression analysis) to identify differentially expressed genes (DEGs) using the TCGAanalyze_DEA
function.
TCGAanalyze_DEA
performs DEA using following functions from R :
This function receives as arguments:
Next, we filter the output of dataDEGs by abs(LogFC) >=1, and uses the TCGAanalyze_LevelTab
function to create a table with DEGs (differentially expressed genes), log Fold Change (FC), false discovery rate (FDR), the gene expression level for samples in Cond1type, and Cond2type, and Delta value (the difference of gene expression between the two conditions multiplied logFC).
# Downstream analysis using gene expression data
# TCGA samples from IlluminaHiSeq_RNASeqV2 with type rsem.genes.results
# save(dataBRCA, geneInfo , file = "dataGeneExpression.rda")
library(TCGAbiolinks)
# normalization of genes
dataNorm <- TCGAanalyze_Normalization(tabDF = dataBRCA, geneInfo = geneInfo)
# quantile filter of genes
dataFilt <- TCGAanalyze_Filtering(tabDF = dataNorm,
method = "quantile",
qnt.cut = 0.25)
# selection of normal samples "NT"
samplesNT <- TCGAquery_SampleTypes(barcode = colnames(dataFilt),
typesample = c("NT"))
# selection of tumor samples "TP"
samplesTP <- TCGAquery_SampleTypes(barcode = colnames(dataFilt),
typesample = c("TP"))
# Diff.expr.analysis (DEA)
dataDEGs <- TCGAanalyze_DEA(mat1 = dataFilt[,samplesNT],
mat2 = dataFilt[,samplesTP],
Cond1type = "Normal",
Cond2type = "Tumor",
fdr.cut = 0.01 ,
logFC.cut = 1,
method = "glmLRT")
# DEGs table with expression values in normal and tumor samples
dataDEGsFiltLevel <- TCGAanalyze_LevelTab(dataDEGs,"Tumor","Normal",
dataFilt[,samplesTP],dataFilt[,samplesNT])
The result is shown below:
mRNA | logFC | FDR | Tumor | Normal | Delta |
---|---|---|---|---|---|
FN1 | 2.88 | 1.296151e-19 | 347787.48 | 41234.12 | 1001017.3 |
COL1A1 | 1.77 | 1.680844e-08 | 358010.32 | 89293.72 | 633086.3 |
C4orf7 | 5.20 | 2.826474e-50 | 87821.36 | 2132.76 | 456425.4 |
COL1A2 | 1.40 | 9.480478e-06 | 273385.44 | 91241.32 | 383242.9 |
GAPDH | 1.32 | 3.290678e-05 | 179057.44 | 63663.00 | 236255.5 |
CLEC3A | 6.79 | 7.971002e-74 | 27257.16 | 259.60 | 185158.6 |
IGFBP5 | 1.24 | 1.060717e-04 | 128186.88 | 53323.12 | 158674.6 |
CPB1 | 4.27 | 3.044021e-37 | 37001.76 | 2637.72 | 157968.8 |
CARTPT | 6.72 | 1.023371e-72 | 21700.96 | 215.16 | 145872.8 |
DCD | 7.26 | 1.047988e-80 | 19941.20 | 84.80 | 144806.3 |
TCGAanalyze_EAcomplete & TCGAvisualize_EAbarplot
: Enrichment AnalysisResearchers, in order to better understand the underlying biological processes, often want to retrieve a functional profile of a set of genes that might have an important role. This can be done by performing an enrichment analysis.
We will perform an enrichment analysis on gene sets using the TCGAanalyze_EAcomplete
function. Given a set of genes that are up-regulated under certain conditions, an enrichment analysis will identify classes of genes or proteins that are over-represented using annotations for that gene set.
To view the results you can use the TCGAvisualize_EAbarplot
function as shown below.
library(TCGAbiolinks)
# Enrichment Analysis EA
# Gene Ontology (GO) and Pathway enrichment by DEGs list
Genelist <- rownames(dataDEGsFiltLevel)
system.time(ansEA <- TCGAanalyze_EAcomplete(TFname="DEA genes Normal Vs Tumor",Genelist))
# Enrichment Analysis EA (TCGAVisualize)
# Gene Ontology (GO) and Pathway enrichment barPlot
TCGAvisualize_EAbarplot(tf = rownames(ansEA$ResBP),
GOBPTab = ansEA$ResBP,
GOCCTab = ansEA$ResCC,
GOMFTab = ansEA$ResMF,
PathTab = ansEA$ResPat,
nRGTab = Genelist,
nBar = 10)
The result is shown below:
The figure shows canonical pathways significantly overrepresented (enriched) by the DEGs (differentially expressed genes). The most statistically significant canonical pathways identified in DEGs list are listed according to their p value corrected FDR (-Log) (colored bars) and the ratio of list genes found in each pathway over the total number of genes in that pathway (Ratio, red line).
TCGAanalyze_survival
: Survival Analysis: Cox Regression and dnet packageWhen analyzing survival, different problems come up than the ones discussed so far. One question is how do we deal with subjects dropping out of a study. For example, assuming that we test a new cancer drug. While some subjects die, others may believe that the new drug is not effective, and decide to drop out of the study before the study is finished. A similar problem would be faced when we investigate how long a machine lasts before it breaks down.
Using the clinical data, it is possible to create a survival plot with the function TCGAanalyze_survival
as follows:
clin.gbm <- GDCquery_clinic("TCGA-GBM", "clinical")
TCGAanalyze_survival(clin.gbm,
"gender",
main = "TCGA Set\n GBM",height = 10, width=10)
The arguments of TCGAanalyze_survival
are:
The result is shown below:
library(TCGAbiolinks)
# Survival Analysis SA
clinical_patient_Cancer <- GDCquery_clinic("TCGA-BRCA","clinical")
dataBRCAcomplete <- log2(BRCA_rnaseqv2)
tokenStop<- 1
tabSurvKMcomplete <- NULL
for( i in 1: round(nrow(dataBRCAcomplete)/100)){
message( paste( i, "of ", round(nrow(dataBRCAcomplete)/100)))
tokenStart <- tokenStop
tokenStop <-100*i
tabSurvKM<-TCGAanalyze_SurvivalKM(clinical_patient_Cancer,
dataBRCAcomplete,
Genelist = rownames(dataBRCAcomplete)[tokenStart:tokenStop],
Survresult = F,
ThreshTop=0.67,
ThreshDown=0.33)
tabSurvKMcomplete <- rbind(tabSurvKMcomplete,tabSurvKM)
}
tabSurvKMcomplete <- tabSurvKMcomplete[tabSurvKMcomplete$pvalue < 0.01,]
tabSurvKMcomplete <- tabSurvKMcomplete[order(tabSurvKMcomplete$pvalue, decreasing=F),]
tabSurvKMcompleteDEGs <- tabSurvKMcomplete[
rownames(tabSurvKMcomplete) %in% dataDEGsFiltLevel$mRNA,
]
The result is shown below:
pvalue | Cancer Deaths | Cancer Deaths with Top | Cancer Deaths with Down | |
---|---|---|---|---|
DCTPP1 | 6.204170e-08 | 66 | 46 | 20 |
APOO | 9.390193e-06 | 65 | 49 | 16 |
LOC387646 | 1.039097e-05 | 69 | 48 | 21 |
PGK1 | 1.198577e-05 | 71 | 49 | 22 |
CCNE2 | 2.100348e-05 | 65 | 48 | 17 |
Mean Tumor Top | Mean Tumor Down | Mean Normal | |
---|---|---|---|
DCTPP1 | 13.31 | 12.01 | 11.74 |
APOO | 11.40 | 10.17 | 10.01 |
LOC387646 | 7.92 | 4.64 | 5.90 |
PGK1 | 15.66 | 14.18 | 14.28 |
CCNE2 | 11.07 | 8.23 | 7.03 |
TCGAanalyze_DMR
: Differentially methylated regions AnalysisWe will search for differentially methylated CpG sites using the TCGAanalyze_DMR
function. In order to find these regions we use the beta-values (methylation values ranging from 0.0 to 1.0) to compare two groups.
First, it calculates the difference between the mean DNA methylation of each group for each probe.
Second, it test for differential expression between two groups using the wilcoxon test adjusting by the Benjamini-Hochberg method. The default arguments was set to require a minimum absolute beta-values difference of 0.2 and am adjusted p-value of < 0.01.
After these tests, we save a volcano plot (x-axis:diff mean methylation, y-axis: statistical significance) that will help the user identify the differentially methylated CpG sites, then the results are saved in a csv file (DMR_results.groupCol.group1.group2.csv) and finally the object is returned with the calculus in the rowRanges.
The arguments of TCGAanalyze_DMR are:
data <- TCGAanalyze_DMR(data, groupCol = "methylation_subtype",
group1 = "CIMP.H",
group2="CIMP.L",
p.cut = 10^-5,
diffmean.cut = 0.25,
legend = "State",
plot.filename = "coad_CIMPHvsCIMPL_metvolcano.png")
The output will be a plot such as the figure below. The green dots are the probes that are hypomethylated in group 2 compared to group 1, while the red dots are the hypermethylated probes in group 2 compared to group 1
Also, the TCGAanalyze_DMR
function will save the plot as pdf and return the same SummarizedExperiment that was given as input with the values of p-value, p-value adjusted, diffmean and the group it belongs in the graph (non significant, hypomethylated, hypermethylated) in the rowRanges. The collumns will be (where group1 and group2 are the names of the groups):
This values can be view/acessed using the rowRanges
acessesor (rowRanges(data)
).
Observation: Calling the same function again, with the same arguments will only plot the results, as it was already calculated. If you want to have them recalculated, please set overwrite
to TRUE
or remove the calculated collumns.
TCGAvisualize
: Visualize results from analysis functions with TCGA’s data.You can easily visualize results from some following functions:
TCGAvisualize_Heatmap
: Create heatmaps with cluster barsIn order to have a better view of clusters, we normally use heatmaps. TCGAvisualize_Heatmap
will plot a heatmap and add to each sample bars representing different features. This function is a wrapper to the package ComplexHeatmap package,
The arguments of this function are:
For more information please take a look on case study #2.
The result is shown below:
TCGAvisualize_Volcano
: Create volcano plotCreates a volcano plot for DNA methylation or expression
The arguments of this function are:
For more information please take a look on case study #3.
TCGAvisualize_oncoprint
The visualize the MAF files we created a function called TCGAvisualize_oncoprint
that uses the ComplexHeatmap to create an oncoprint.
mut <- GDCquery_Maf(tumor = "ACC", pipelines = "muse")
clin <- GDCquery_clinic("TCGA-ACC","clinical")
clin <- clin[,c("bcr_patient_barcode","disease","gender","tumor_stage","race","vital_status")]
TCGAvisualize_oncoprint(mut = mut, genes = mut$Hugo_Symbol[1:20],
filename = "oncoprint.pdf",
annotation = clin,
color=c("background"="#CCCCCC","DEL"="purple","INS"="yellow","SNP"="brown"),
rows.font.size=10,
heatmap.legend.side = "right",
dist.col = 0,
width = 5,
label.font.size = 10)
The columns are the samples and the rows are genes. The upper histograms shows the number of mutation of each samples, the right histogram shows the number of patients that has that mutation (which is also represented by the percentage numbers in the left side of the plot).
TCGAvisualize_PCA
: Principal Component Analysis plot for differentially expressed genesIn order to better understand our genes, we can perform a PCA to reduce the number of dimensions of our gene set. The function TCGAvisualize_PCA
will plot the PCA for different groups.
The arguments of this function are:
# normalization of genes
dataNorm <- TCGAbiolinks::TCGAanalyze_Normalization(dataBRCA, geneInfo)
# quantile filter of genes
dataFilt <- TCGAanalyze_Filtering(tabDF = dataNorm,
method = "quantile",
qnt.cut = 0.25)
# selection of normal samples "NT"
group1 <- TCGAquery_SampleTypes(colnames(dataFilt), typesample = c("NT"))
# selection of normal samples "TP"
group2 <- TCGAquery_SampleTypes(colnames(dataFilt), typesample = c("TP"))
# Principal Component Analysis plot for ntop selected DEGs
pca <- TCGAvisualize_PCA(dataFilt,dataDEGsFiltLevel, ntopgenes = 200, group1, group2)
The result is shown below:
TCGAvisualize_meanMethylation
: Mean DNA Methylation AnalysisUsing the data and calculating the mean DNA methylation per group, it is possible to create a mean DNA methylation boxplot with the function TCGAvisualize_meanMethylation
as follows:
query <- GDCquery(project = "TCGA-GBM",
data.category = "DNA methylation",
platform = "Illumina Human Methylation 27",
legacy = TRUE,
barcode = c("TCGA-02-0058-01A-01D-0186-05", "TCGA-12-1597-01B-01D-0915-05",
"TCGA-12-0829-01A-01D-0392-05", "TCGA-06-0155-01B-01D-0521-05",
"TCGA-02-0099-01A-01D-0199-05", "TCGA-19-4068-01A-01D-1228-05",
"TCGA-19-1788-01A-01D-0595-05", "TCGA-16-0848-01A-01D-0392-05"))
GDCdownload(query, method = "api")
data <- GDCprepare(query)
# "shortLetterCode" is a column in the SummarizedExperiment::colData(data) matrix
TCGAvisualize_meanMethylation(data, groupCol = "shortLetterCode",filename = NULL)
## groups Mean Median Max Min
## 1 NT 0.5014940 0.5009496 0.5359795 0.4732053
## 2 TP 0.5014747 0.4958995 0.5298975 0.4656788
## 3 TR 0.5030028 0.4986486 0.5415967 0.4673249
## NT TP TR
## NT NA 0.9987060 0.9004057
## TP 0.9987060 NA 0.9031248
## TR 0.9004057 0.9031248 NA
The arguments of TCGAvisualize_meanMethylation
are:
GDCprepare
Other examples using the parameters:
# setting y limits: lower 0, upper 1
TCGAvisualize_meanMethylation(data,groupCol = "shortLetterCode",
filename = NULL, y.limits = c(0,1))
## groups Mean Median Max Min
## 1 NT 0.5014940 0.5009496 0.5359795 0.4732053
## 2 TP 0.5014747 0.4958995 0.5298975 0.4656788
## 3 TR 0.5030028 0.4986486 0.5415967 0.4673249
## NT TP TR
## NT NA 0.9987060 0.9004057
## TP 0.9987060 NA 0.9031248
## TR 0.9004057 0.9031248 NA
# setting y limits: lower 0
TCGAvisualize_meanMethylation(data,groupCol = "shortLetterCode",
filename = NULL, y.limits = 0)
## groups Mean Median Max Min
## 1 NT 0.5014940 0.5009496 0.5359795 0.4732053
## 2 TP 0.5014747 0.4958995 0.5298975 0.4656788
## 3 TR 0.5030028 0.4986486 0.5415967 0.4673249
## NT TP TR
## NT NA 0.9987060 0.9004057
## TP 0.9987060 NA 0.9031248
## TR 0.9004057 0.9031248 NA
# Changing shapes of jitters to show subgroups
TCGAvisualize_meanMethylation(data,
groupCol = "subtype_Pan.Glioma.DNA.Methylation.Cluster",
subgroupCol ="vital_status", filename = NULL)
## groups Mean Median Max Min
## 1 group1 0.5047761 0.4999907 0.5359795 0.4869987
## 2 group2 0.4913455 0.4921337 0.5188015 0.4656788
## 3 group3 0.5098497 0.5103386 0.5415967 0.4732053
## group1 group2 group3
## group1 NA 0.2144950 0.6538003
## group2 0.2144950 NA 0.1437124
## group3 0.6538003 0.1437124 NA
# Sorting bars by descending mean methylation
TCGAvisualize_meanMethylation(data,
groupCol = "subtype_Pan.Glioma.DNA.Methylation.Cluster",
sort="mean.desc",
filename=NULL)
## groups Mean Median Max Min
## 1 group1 0.5047761 0.4999907 0.5359795 0.4869987
## 2 group2 0.4913455 0.4921337 0.5188015 0.4656788
## 3 group3 0.5098497 0.5103386 0.5415967 0.4732053
## group1 group2 group3
## group1 NA 0.2144950 0.6538003
## group2 0.2144950 NA 0.1437124
## group3 0.6538003 0.1437124 NA
# Sorting bars by asc mean methylation
TCGAvisualize_meanMethylation(data,
groupCol = "subtype_Pan.Glioma.DNA.Methylation.Cluster",
sort = "mean.asc",
filename=NULL)
## groups Mean Median Max Min
## 1 group1 0.5047761 0.4999907 0.5359795 0.4869987
## 2 group2 0.4913455 0.4921337 0.5188015 0.4656788
## 3 group3 0.5098497 0.5103386 0.5415967 0.4732053
## group1 group2 group3
## group1 NA 0.2144950 0.6538003
## group2 0.2144950 NA 0.1437124
## group3 0.6538003 0.1437124 NA
TCGAvisualize_meanMethylation(data,
groupCol = "vital_status",
sort = "mean.asc",
filename=NULL)
## groups Mean Median Max Min
## 1 ALIVE 0.5014747 0.4958995 0.5298975 0.4656788
## 2 DEAD 0.5022484 0.4993197 0.5415967 0.4673249
## ALIVE DEAD
## ALIVE NA 0.9413677
## DEAD 0.9413677 NA
TCGAvisualize_starburst
: Integration of gene expression and DNA methylation dataThe starburst plot is proposed to combine information from two volcano plots, and is applied for a study of DNA methylation and gene expression. It first introduced in 2010 (Noushmehr, H., Weisenberger, D.J., Diefes, K., Phillips, H.S., Pujara, K., Berman, B.P., Pan, F., Pelloski, C.E., Sulman, E.P., Bhat, K.P. et al. 2010). In order to reproduce this plot, we will use the TCGAvisualize_starburst
function.
The function creates Starburst plot for comparison of DNA methylation and gene expression. The log10 (FDR-corrected P value) for DNA methylation is plotted in the x axis, and for gene expression in the y axis, for each gene. The black dashed line shows the FDR-adjusted P value of 0.01.
The arguments of this function are:
GDCprepare
or Data frame from DMR_results
file. Expected colData columns: diffmean, p.value.adj and p.value Execute volcanoPlot function in order to obtain these values for the object.TCGAanalyze_DEA
function. Expected colData columns: logFC, FDRstarburst <- TCGAvisualize_starburst(coad.SummarizeExperiment,
different.experssion.analysis.data,
group1 = "CIMP.H",
group2 = "CIMP.L",
met.p.cut = 10^-5,
exp.p.cut=10^-5,
names = TRUE)
As a result, the function will a plot the figure below and return a matrix with the Gene_symbol and it status in relation to expression (up regulated/down regulated) and to methylation (Hyper/Hypo methylated). The case study 3, shows the complete pipeline for creating this figure.
CancerProject <- "TCGA-BRCA"
DataDirectory <- paste0("../GDC/",gsub("-","_",CancerProject))
FileNameData <- paste0(DataDirectory, "_","HTSeq_Counts",".rda")
query <- GDCquery(project = CancerProject,
data.category = "Transcriptome Profiling",
data.type = "Gene Expression Quantification",
workflow.type = "HTSeq - Counts")
samplesDown <- getResults(query,cols=c("cases"))
dataSmTP <- TCGAquery_SampleTypes(barcode = samplesDown,
typesample = "TP")
dataSmNT <- TCGAquery_SampleTypes(barcode = samplesDown,
typesample = "NT")
dataSmTP_short <- dataSmTP[1:10]
dataSmNT_short <- dataSmNT[1:10]
queryDown <- GDCquery(project = CancerProject,
data.category = "Transcriptome Profiling",
data.type = "Gene Expression Quantification",
workflow.type = "HTSeq - Counts",
barcode = c(dataSmTP_short, dataSmNT_short))
GDCdownload(query = queryDown,
directory = DataDirectory)
dataPrep <- GDCprepare(query = queryDown,
save = TRUE,
directory = DataDirectory,
save.filename = FileNameData)
dataPrep <- TCGAanalyze_Preprocessing(object = dataPrep,
cor.cut = 0.6,
datatype = "HTSeq - Counts")
dataNorm <- TCGAanalyze_Normalization(tabDF = dataPrep,
geneInfo = geneInfoHT,
method = "gcContent")
boxplot(dataPrep, outline = FALSE)
boxplot(dataNorm, outline = FALSE)
dataFilt <- TCGAanalyze_Filtering(tabDF = dataNorm,
method = "quantile",
qnt.cut = 0.25)
dataDEGs <- TCGAanalyze_DEA(mat1 = dataFilt[,dataSmTP_short],
mat2 = dataFilt[,dataSmNT_short],
Cond1type = "Normal",
Cond2type = "Tumor",
fdr.cut = 0.01 ,
logFC.cut = 1,
method = "glmLRT")
require(TCGAbiolinks)
CancerProject <- "TCGA-BRCA"
DataDirectory <- paste0("../GDC/",gsub("-","_",CancerProject))
FileNameData <- paste0(DataDirectory, "_","miRNA_gene_quantification",".rda")
query.miR <- GDCquery(project = CancerProject,
data.category = "Gene expression",
data.type = "miRNA gene quantification",
file.type = "hg19.mirna",
legacy = TRUE)
samplesDown.miR <- getResults(query.miR,cols=c("cases"))
dataSmTP.miR <- TCGAquery_SampleTypes(barcode = samplesDown.miR,
typesample = "TP")
dataSmNT.miR <- TCGAquery_SampleTypes(barcode = samplesDown.miR,
typesample = "NT")
dataSmTP_short.miR <- dataSmTP.miR[1:10]
dataSmNT_short.miR <- dataSmNT.miR[1:10]
queryDown.miR <- GDCquery(project = CancerProject,
data.category = "Gene expression",
data.type = "miRNA gene quantification",
file.type = "hg19.mirna",
legacy = TRUE,
barcode = c(dataSmTP_short.miR, dataSmNT_short.miR))
GDCdownload(query = queryDown.miR,
directory = DataDirectory)
dataAssy.miR <- GDCprepare(query = queryDown.miR,
save = TRUE,
save.filename = FileNameData,
summarizedExperiment = TRUE,
directory =DataDirectory )
rownames(dataAssy.miR) <- dataAssy.miR$miRNA_ID
# using read_count's data
read_countData <- colnames(dataAssy.miR)[grep("count", colnames(dataAssy.miR))]
dataAssy.miR <- dataAssy.miR[,read_countData]
colnames(dataAssy.miR) <- gsub("read_count_","", colnames(dataAssy.miR))
dataFilt <- TCGAanalyze_Filtering(tabDF = dataAssy.miR,
method = "quantile",
qnt.cut = 0.25)
dataDEGs <- TCGAanalyze_DEA(mat1 = dataFilt[,dataSmNT_short.miR],
mat2 = dataFilt[,dataSmTP_short.miR],
Cond1type = "Normal",
Cond2type = "Tumor",
fdr.cut = 0.01 ,
logFC.cut = 1,
method = "glmLRT")
This vignette shows a complete workflow of the TCGAbiolinks package. The code is divided in 4 case study:
library(SummarizedExperiment)
library(TCGAbiolinks)
query.exp <- GDCquery(project = "TCGA-BRCA",
legacy = TRUE,
data.category = "Gene expression",
data.type = "Gene expression quantification",
platform = "Illumina HiSeq",
file.type = "results",
experimental.strategy = "RNA-Seq",
sample.type = c("Primary solid Tumor","Solid Tissue Normal"))
GDCdownload(query.exp)
brca.exp <- GDCprepare(query = query.exp, save = TRUE, save.filename = "brcaExp.rda")
# get subtype information
dataSubt <- TCGAquery_subtype(tumor = "BRCA")
# get clinical data
dataClin <- GDCquery_clinic(project = "TCGA-BRCA","clinical")
# Which samples are primary solid tumor
dataSmTP <- TCGAquery_SampleTypes(getResults(query.exp,cols="cases"),"TP")
# which samples are solid tissue normal
dataSmNT <- TCGAquery_SampleTypes(getResults(query.exp,cols="cases"),"NT")
Using TCGAnalyze_DEA
, we identified 3,390 differentially expression genes (DEG) (log fold change >=1 and FDR < 1%) between 114 normal and 1097 BRCA samples. In order to understand the underlying biological process from DEGs we performed an enrichment analysis using TCGAnalyze_EA_complete
function.
dataPrep <- TCGAanalyze_Preprocessing(object = brca.exp, cor.cut = 0.6)
dataNorm <- TCGAanalyze_Normalization(tabDF = dataPrep,
geneInfo = geneInfo,
method = "gcContent")
dataFilt <- TCGAanalyze_Filtering(tabDF = dataNorm,
method = "quantile",
qnt.cut = 0.25)
dataDEGs <- TCGAanalyze_DEA(mat1 = dataFilt[,dataSmNT],
mat2 = dataFilt[,dataSmTP],
Cond1type = "Normal",
Cond2type = "Tumor",
fdr.cut = 0.01 ,
logFC.cut = 1,
method = "glmLRT")
TCGAbiolinks outputs bar chart with the number of genes for the main categories of three ontologies (GO:biological process, GO:cellular component, and GO:molecular function, respectively).
ansEA <- TCGAanalyze_EAcomplete(TFname="DEA genes Normal Vs Tumor",
RegulonList = rownames(dataDEGs))
TCGAvisualize_EAbarplot(tf = rownames(ansEA$ResBP),
GOBPTab = ansEA$ResBP,
GOCCTab = ansEA$ResCC,
GOMFTab = ansEA$ResMF,
PathTab = ansEA$ResPat,
nRGTab = rownames(dataDEGs),
nBar = 20)
The figure resulted from the code above is shown below.
The Kaplan-Meier analysis was used to compute survival univariate curves, and
log-Ratio test was computed to assess the statistical significance by using TCGAanalyze_SurvivalKM function; starting with 3,390 DEGs genes we found 555 significantly genes with p.value <0.05.
group1 <- TCGAquery_SampleTypes(colnames(dataFilt), typesample = c("NT"))
group2 <- TCGAquery_SampleTypes(colnames(dataFilt), typesample = c("TP"))
dataSurv <- TCGAanalyze_SurvivalKM(clinical_patient = dataClin,
dataGE = dataFilt,
Genelist = rownames(dataDEGs),
Survresult = FALSE,
ThreshTop = 0.67,
ThreshDown = 0.33,
p.cut = 0.05, group1, group2)
Cox-regression analysis was used to compute survival multivariate curves, and cox p-value was computed to assess the statistical significance by using TCGAnalyze_SurvivalCoxNET function. Survival multivariate analysis found 160 significantly genes according to the cox p-value FDR 5.00e-02. From DEGs that we found to correlate significantly with survival by both univariate and multivariate analyses we analyzed the following network.
The interactome network graph was generated using STRING.,org.Hs.string version 10 (Human functional protein association network). The network graph was resized by dnet package considering only multivariate survival genes, with strong interaction (threshold = 700) we obtained a subgraphsub graph of 24 nodes and 31 edges.
require(dnet) # to change
org.Hs.string <- dRDataLoader(RData = "org.Hs.string")
TabCoxNet <- TCGAvisualize_SurvivalCoxNET(dataClin,
dataFilt,
Genelist = rownames(dataSurv),
scoreConfidence = 700,
org.Hs.string = org.Hs.string,
titlePlot = "Case Study n.1 dnet")
The figure resulted from the code above is shown below.
We focused on the analysis of LGG samples. In particular, we used TCGAbiolinks to download 293 samples with molecular subtypes. Link the complete complete code. .
library(TCGAbiolinks)
library(SummarizedExperiment)
query.exp <- GDCquery(project = "TCGA-LGG",
legacy = TRUE,
data.category = "Gene expression",
data.type = "Gene expression quantification",
platform = "Illumina HiSeq",
file.type = "results",
experimental.strategy = "RNA-Seq",
sample.type = "Primary solid Tumor")
GDCdownload(query.exp)
lgg.exp <- GDCprepare(query = query.exp, save = TRUE, save.filename = "lggExp.rda")
# get subtype information
dataSubt <- TCGAquery_subtype(tumor = "LGG")
# get indexed clinical data
dataClin <- GDCquery_clinic(project = "TCGA-LGG", "Clinical")
First, we searched for possible outliers using the TCGAanalyze_Preprocessing
function, which performs an Array Array Intensity correlation AAIC. We used all samples in expression data which contain molecular subtypes, filtering out samples without molecular information, and using only IDHmut-codel (n=85), IDHmut-non-codel (n=141) and IDHwt (n=56), NA (11), to define a square symetric matrix of pearson correlation among all samples (n=293). According to this matrix we found no samples with low correlation (cor.cut = 0.6) that can be identified as possible outliers, so we continued our analysis with 70 samples.
Second, using the TCGAanalyze_Normalization
function we normalized mRNA transcripts and miRNA, using EDASeq package. This function does use Within-lane normalization procedures to adjust for GC-content effect (or other gene-level effects) on read counts: loess robust local regression, global-scaling, and full-quantile normalization (Risso, D., Schwartz, K., Sherlock, G., & Dudoit, S. 2011) and between-lane normalization procedures to adjust for distributional differences between lanes (e.g., sequencing depth): global-scaling and full-quantile normalization (James H Bullard, Elizabeth Purdom, Kasper D Hansen and Sandrine Dudoit 2010).
Third, using the TCGAanalyze_Filtering
function we applied 3 filters removing features / mRNAs with low signal across samples obtaining 4578, 4284, 1187 mRNAs respectively.
Then we applied two Hierarchical cluster analysis on 1187 mRNAs after the three filters described above, the first cluster using as method ward.D2, and the second with ConsensusClusterPlus.
After the two clustering analysis, with cut.tree =4 we obtained n=4 espression clusters (EC).
# expression data with molecular subtypes
lgg.exp <- subset(lgg.exp, select = colData(lgg.exp)$patient %in% dataSubt$patient)
dataPrep <- TCGAanalyze_Preprocessing(object = lgg.exp,cor.cut = 0.6)
dataNorm <- TCGAanalyze_Normalization(tabDF = dataPrep,
geneInfo = geneInfo,
method = "gcContent")
datFilt1 <- TCGAanalyze_Filtering(tabDF = dataNorm,method = "varFilter")
datFilt2 <- TCGAanalyze_Filtering(tabDF = datFilt1,method = "filter1")
datFilt <- TCGAanalyze_Filtering(tabDF = datFilt2,method = "filter2")
data_Hc1 <- TCGAanalyze_Clustering(tabDF = datFilt,
method = "hclust",
methodHC = "ward.D2")
data_Hc2 <- TCGAanalyze_Clustering(tabDF = datFilt,
method = "consensus",
methodHC = "ward.D2")
For the next steps, we will merge the clinical data with with the the cluster information and we will add the subtype information.
#------ Add cluster information
cluster <- data.frame("groupsHC" = data_Hc2[[4]]$consensusClass)
cluster$groupsHC <- paste0("EC",cluster$groupsHC)
cluster$patient <- substr(colData(lgg.exp)$patient,1,12)
# Add information about gropus from consensus Cluster in clinical data
dataClin <- merge(dataClin,cluster, by.x="bcr_patient_barcode", by.y="patient")
# Merge subtype and clinical data
clin_subt <- merge(dataClin,dataSubt, by.x="bcr_patient_barcode", by.y="patient")
clin_subt_all <- merge(dataClin,dataSubt,
by.x="bcr_patient_barcode", by.y="patient", all.x = TRUE)
The next steps will be to visualize the data. First, we created the survival plot.
#----------- VISUALIZE --------------------
# plotting survival for groups EC1, EC2, EC3, EC4
TCGAanalyze_survival(data = clin_subt_all,
clusterCol = "groupsHC",
main = "TCGA kaplan meier survival plot from consensus cluster",
legend = "RNA Group",
color = c("black","red","blue","green3"),
filename = "case2_surv.png")
The result is showed below:
We will also, create a heatmap of the expression.
TCGAvisualize_Heatmap(t(datFilt),
col.metadata = clin_subt[,c("bcr_patient_barcode",
"groupsHC",
"Histology",
"IDH.codel.subtype")],
col.colors = list(
groupsHC = c("EC1"="black",
"EC2"="red",
"EC3"="blue",
"EC4"="green3"),
Histology=c("astrocytoma"="navy",
"oligoastrocytoma"="green3",
"oligodendroglioma"="red"),
IDH.codel.subtype = c("IDHmut-codel"="tomato",
"IDHmut-non-codel"="navy",
"IDHwt"="gold","NA"="white")),
sortCol = "groupsHC",
type = "expression", # sets default color
scale = "row", # use z-scores for better visualization
title = "Heatmap from concensus cluster",
filename = "case2_Heatmap.pdf",
cluster_rows = TRUE)
The result is shown below:
Finally, we will take a look in the mutation genes. We will first download the maf file with TCGAquery_maf
. In this example we will investigate the gene “ATRX”.
LGGmut <- GDCquery_Maf(tumor = "LGG", pipelines = "muse")
# Selecting gene
mRNAsel <- "ATRX"
LGGselected <- LGGmut[LGGmut$Hugo_Symbol == mRNAsel,]
dataMut <- LGGselected[!duplicated(LGGselected$Tumor_Sample_Barcode),]
dataMut$Tumor_Sample_Barcode <- substr(dataMut$Tumor_Sample_Barcode,1,12)
# Adding the Expression Cluster classification found before
dataMut <- merge(dataMut, cluster, by.y="patient", by.x="Tumor_Sample_Barcode")
dataMut <- dataMut[dataMut$Variant_Classification!=0,]
In recent years, it has been described the relationship between DNA methylation and gene expression and the study of this relationship is often difficult to accomplish.
This case study will show the steps to investigate the relationship between the two types of data.
First, we downloaded ACC DNA methylation data for HumanMethylation450k platforms, and ACC RNA expression data for Illumina HiSeq platform.
TCGAbiolinks adds by default the subtypes classification already published by researchers.
We will use this classification to do our examples. So, selected the groups CIMP-low and CIMP-high to do RNA expression and DNA methylation comparison.
library(TCGAbiolinks)
library(SummarizedExperiment)
dir.create("case3")
setwd("case3")
#-----------------------------------
# STEP 1: Search, download, prepare |
#-----------------------------------
# 1.1 - DNA methylation
# ----------------------------------
query.met <- GDCquery(project = "TCGA-ACC",
legacy = TRUE,
data.category = "DNA methylation",
platform = "Illumina Human Methylation 450")
GDCdownload(query.met)
acc.met <- GDCprepare(query = query.met,
save = TRUE,
save.filename = "accDNAmet.rda",
summarizedExperiment = TRUE)
#-----------------------------------
# 1.2 - RNA expression
# ----------------------------------
query.exp <- GDCquery(project = "TCGA-ACC",
legacy = TRUE,
data.category = "Gene expression",
data.type = "Gene expression quantification",
platform = "Illumina HiSeq",
file.type = "results")
GDCdownload(query.exp)
acc.exp <- GDCprepare(query = query.exp, save = TRUE, save.filename = "accExp.rda")
For DNA methylation, we perform a DMR (different methylated region) analysis, which will give the difference of DNA methylation for the probes of the groups and their significance value. The output can be seen in a volcano plot. Note: Depending on the number of samples this function can be very slow due to the wilcoxon test, taking from hours to days.
# na.omit
acc.met <- subset(acc.met,subset = (rowSums(is.na(assay(acc.met))) == 0))
# Volcano plot
acc.met <- TCGAanalyze_DMR(acc.met, groupCol = "subtype_MethyLevel",
group1 = "CIMP-high",
group2="CIMP-low",
p.cut = 10^-5,
diffmean.cut = 0.25,
legend = "State",
plot.filename = "CIMP-highvsCIMP-low_metvolcano.png")
The figure resulted from the code above is shown below:
For the expression analysis, we do a DEA (differential expression analysis) which will give the fold change of gene expression and their significance value.
#-------------------------------------------------
# 2.3 - DEA - Expression analysis - volcano plot
# ------------------------------------------------
acc.exp.aux <- subset(acc.exp,
select = colData(acc.exp)$subtype_MethyLevel %in% c("CIMP-high","CIMP-low"))
idx <- colData(acc.exp.aux)$subtype_MethyLevel %in% c("CIMP-high")
idx2 <- colData(acc.exp.aux)$subtype_MethyLevel %in% c("CIMP-low")
dataPrep <- TCGAanalyze_Preprocessing(object = acc.exp.aux, cor.cut = 0.6)
dataNorm <- TCGAanalyze_Normalization(tabDF = dataPrep,
geneInfo = geneInfo,
method = "gcContent")
dataFilt <- TCGAanalyze_Filtering(tabDF = dataNorm,
qnt.cut = 0.25,
method='quantile')
dataDEGs <- TCGAanalyze_DEA(mat1 = dataFilt[,idx],
mat2 = dataFilt[,idx2],
Cond1type = "CIMP-high",
Cond2type = "CIMP-low",
method = "glmLRT")
TCGAVisualize_volcano(dataDEGs$logFC,dataDEGs$FDR,
filename = "Case3_volcanoexp.png",
x.cut = 3,
y.cut = 10^-5,
names = rownames(dataDEGs),
color = c("black","red","darkgreen"),
names.size = 2,
xlab = " Gene expression fold change (Log2)",
legend = "State",
title = "Volcano plot (CIMP-high vs CIMP-low)",
width = 10)
The figure resulted from the code above is shown below:
Finally, using both previous analysis we do a starburst plot to select the genes that are Candidate Biologically Significant.
Observation: over the time, the number of samples has increased and the clinical data updated. We used only the samples that had a classification in the examples.
#------------------------------------------
# 2.4 - Starburst plot
# -----------------------------------------
# If true the argument names of the geenes in circle
# (biologically significant genes, has a change in gene
# expression and DNA methylation and respects all the thresholds)
# will be shown
# these genes are returned by the function see starburst object after the function is executed
starburst <- TCGAvisualize_starburst(acc.met, dataDEGs,"CIMP-high","CIMP-low",
filename = "starburst.png",
met.p.cut = 10^-5,
exp.p.cut = 10^-5,
diffmean.cut = 0.25,
logFC.cut = 3,
names = FALSE,
height=10,
width=15,
dpi=300)
The figure resulted from the code above is shown below:
One of the biggest problems related to anlysis of data is the preparation phase, which often consists of successive steps in order to prepare it to a format acceptable by certain algorithms and software.
With the objective of assisting users in this arduous step, TCGAbiolinks will offer to prepare the data in order to obtain the correct format for different packages.
An example of package to perform DNA methylation and RNA expression analysis is ELMER (Huber, Wolfgang and Carey, Vincent J and Gentleman, Robert and Anders, Simon and Carlson, Marc and Carvalho, Benilton S and Bravo, Hector Corrada and Davis, Sean and Gatto, Laurent and Girke, Thomas and others 2015). ELMER, which is designed to combine DNA methylation and gene expression data from human tissues to infer multi-level cis-regulatory networks. ELMER uses DNA methylation to identify enhancers, and correlates enhancer state with expression of nearby genes to identify one or more transcriptional targets. Transcription factor (TF) binding site analysis of enhancers is coupled with expression analysis of all TFs to infer upstream regulators. This package can be easily applied to TCGA public available cancer data sets and custom DNA methylation and gene expression data sets.
ELMER analyses have 5 main steps:
We will present this the study KIRC by TCGAbiolinks and ELMER integration. For more information, please consult the ELMER package.
For the DNA methylation data we will search the platform HumanMethylation450 and level 3 data with TCGAquery. After, we will download the data with TCGAdownload and the data will be prepared into a data.frame object with TCGAprepare. TCGAprepare_elmer will be used to prepare the data frame to ELMER.
library(TCGAbiolinks)
library(SummarizedExperiment)
library(ELMER)
library(parallel)
dir.create("case4")
setwd("case4")
#-----------------------------------
# STEP 1: Search, download, prepare |
#-----------------------------------
# 1.1 - DNA methylation
# ----------------------------------
query.met <- GDCquery(project = "TCGA-KIRC",
legacy = TRUE,
data.category = "DNA methylation",
platform = "Illumina Human Methylation 450")
GDCdownload(query.met)
kirc.met <- GDCprepare(query = query.met,
save = TRUE,
save.filename = "kircDNAmet.rda",
summarizedExperiment = TRUE)
kirc.met <- TCGAprepare_elmer(kirc.met,
platform = "HumanMethylation450",
save = TRUE,
met.na.cut = 0.2)
For the gene expression data we will search the platform IlluminaHiSeq_RNASeqV2 and level 3 data with TCGAquery. Next, we will download the normalized_results data with TCGAdownload and prepare it into a data.frame object with TCGAprepare. TCGAprepare_elmer will be used to prepare the data frame to ELMER.
# Step 1.2 download expression data
#-----------------------------------
# 1.2 - RNA expression
# ----------------------------------
query.exp <- GDCquery(project = "TCGA-KIRC",
legacy = TRUE,
data.category = "Gene expression",
data.type = "Gene expression quantification",
platform = "Illumina HiSeq",
file.type = "normalized_results")
GDCdownload(query.exp)
kirc.exp <- GDCprepare(query = query.exp, save = TRUE, save.filename = "kirkExp.rda")
kirc.exp <- TCGAprepare_elmer(kirc.exp,
save = TRUE,
platform = "IlluminaHiSeq_RNASeqV2")
The ELMER input is a mee object that contains the methylation matrix, the expression matrix, the probe information, the gene information and a table with a summary of the data.
#-----------------------------------
# STEP 2: ELMER integration |
#-----------------------------------
# Step 2.1 prepare mee object |
# -----------------------------------
library(ELMER)
library(parallel)
geneAnnot <- txs()
geneAnnot$GENEID <- paste0("ID",geneAnnot$GENEID)
geneInfo <- promoters(geneAnnot,upstream = 0, downstream = 0)
probe <- get.feature.probe()
mee <- fetch.mee(meth = kirc.met, exp = kirc.exp, TCGA = TRUE,
probeInfo = probe, geneInfo = geneInfo)
The code below excute steps 2-5 showed above:
direction <- c("hyper","hypo")
for (j in direction){
print(j)
dir.out <- paste0("kirc/",j)
dir.create(dir.out, recursive = TRUE)
#--------------------------------------
# STEP 3: Analysis |
#--------------------------------------
# Step 3.1: Get diff methylated probes |
#--------------------------------------
Sig.probes <- get.diff.meth(mee, cores=detectCores(),
dir.out =dir.out,
diff.dir=j,
pvalue = 0.01)
#-------------------------------------------------------------
# Step 3.2: Identify significant probe-gene pairs |
#-------------------------------------------------------------
# Collect nearby 20 genes for Sig.probes
nearGenes <- GetNearGenes(TRange=getProbeInfo(mee, probe=Sig.probes$probe),
cores=detectCores(),
geneAnnot=getGeneInfo(mee))
pair <- get.pair(mee=mee,
probes=Sig.probes$probe,
nearGenes=nearGenes,
permu.dir=paste0(dir.out,"/permu"),
dir.out=dir.out,
cores=detectCores(),
label= j,
permu.size=10000,
Pe = 0.001)
Sig.probes.paired <- fetch.pair(pair=pair,
probeInfo = getProbeInfo(mee),
geneInfo = getGeneInfo(mee))
Sig.probes.paired <-read.csv(paste0(dir.out,"/getPair.",j,".pairs.significant.csv"),
stringsAsFactors=FALSE)[,1]
#-------------------------------------------------------------
# Step 3.3: Motif enrichment analysis on the selected probes |
#-------------------------------------------------------------
if(length(Sig.probes.paired) > 0 ){
#-------------------------------------------------------------
# Step 3.3: Motif enrichment analysis on the selected probes |
#-------------------------------------------------------------
enriched.motif <- get.enriched.motif(probes=Sig.probes.paired,
dir.out=dir.out, label=j,
background.probes = probe$name)
if(length(enriched.motif) > 0){
#-------------------------------------------------------------
# Step 3.4: Identifying regulatory TFs |
#-------------------------------------------------------------
print("get.TFs")
TF <- get.TFs(mee = mee,
enriched.motif = enriched.motif,
dir.out = dir.out,
cores = detectCores(), label = j)
save(TF, enriched.motif, Sig.probes.paired,
pair, nearGenes, Sig.probes,
file=paste0(dir.out,"/ELMER_results_",j,".rda"))
}
}
}
From this analysis it is possible to verify the relationship between nearby 20 gene expression vs DNA methylation at this probe. The result of this is show by the ELMER scatter plot.
scatter.plot(mee,byProbe=list(probe=c("cg00328720"),geneNum=20), category="TN", lm_line = TRUE)
The result is shown below:
Each scatter plot showing the average DNA methylation level of sites with the UA6 motif in all KIRC samples plotted against the expression of the transcription factor ZNF677 and PEG3 respectively.
scatter.plot(mee,byTF = list(TF = c("ZNF677","PEG3"),
probe = enriched.motif[["UA6"]]), category = "TN",
save = TRUE, lm_line = TRUE, dir.out = "kirc")
The result is shown below:
The schematic plot shows probe colored in blue and the location of nearby 20 genes. The genes significantly linked to the probe were shown in red.
pair <- fetch.pair(pair="./kirc/getPair.hypo.pairs.significant.csv",
probeInfo = mee@probeInfo, geneInfo = mee@geneInfo)
schematic.plot(pair=pair, byProbe="cg15862394",save=FALSE)
The result is shown below:
The plot shows the odds ratio (x axis) for the selected motifs with OR above 1.3 and lower boundary of OR above 1.3. The range shows the 95% confidence interval for each Odds Ratio.
sessionInfo()
## R version 3.4.1 (2017-06-30)
## Platform: x86_64-pc-linux-gnu (64-bit)
## Running under: Ubuntu 16.04.3 LTS
##
## Matrix products: default
## BLAS: /home/biocbuild/bbs-3.5-bioc/R/lib/libRblas.so
## LAPACK: /home/biocbuild/bbs-3.5-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] grid parallel stats4 stats graphics grDevices utils
## [8] datasets methods base
##
## other attached packages:
## [1] png_0.1-7 TCGAbiolinks_2.5.9
## [3] bindrcpp_0.2 DT_0.2
## [5] dplyr_0.7.3 SummarizedExperiment_1.6.3
## [7] DelayedArray_0.2.7 matrixStats_0.52.2
## [9] Biobase_2.36.2 GenomicRanges_1.28.5
## [11] GenomeInfoDb_1.12.2 IRanges_2.10.3
## [13] S4Vectors_0.14.4 BiocGenerics_0.22.0
##
## loaded via a namespace (and not attached):
## [1] backports_1.1.0 circlize_0.4.1
## [3] aroma.light_3.6.0 plyr_1.8.4
## [5] selectr_0.3-1 ConsensusClusterPlus_1.40.0
## [7] lazyeval_0.2.0 splines_3.4.1
## [9] BiocParallel_1.10.1 ggplot2_2.2.1
## [11] digest_0.6.12 foreach_1.4.3
## [13] htmltools_0.3.6 viridis_0.4.0
## [15] magrittr_1.5 memoise_1.1.0
## [17] cluster_2.0.6 doParallel_1.0.10
## [19] limma_3.32.6 ComplexHeatmap_1.14.0
## [21] Biostrings_2.44.2 readr_1.1.1
## [23] annotate_1.54.0 R.utils_2.5.0
## [25] colorspace_1.3-2 blob_1.1.0
## [27] rvest_0.3.2 ggrepel_0.6.5
## [29] crayon_1.3.2 RCurl_1.95-4.8
## [31] jsonlite_1.5 roxygen2_6.0.1
## [33] genefilter_1.58.1 bindr_0.1
## [35] zoo_1.8-0 survival_2.41-3
## [37] iterators_1.0.8 glue_1.1.1
## [39] survminer_0.4.0 gtable_0.2.0
## [41] zlibbioc_1.22.0 XVector_0.16.0
## [43] GetoptLong_0.1.6 kernlab_0.9-25
## [45] shape_1.4.3 prabclus_2.2-6
## [47] DEoptimR_1.0-8 scales_0.5.0
## [49] DESeq_1.28.0 mvtnorm_1.0-6
## [51] DBI_0.7 edgeR_3.18.1
## [53] ggthemes_3.4.0 Rcpp_0.12.12
## [55] viridisLite_0.2.0 xtable_1.8-2
## [57] cmprsk_2.2-7 foreign_0.8-69
## [59] bit_1.1-12 matlab_1.0.2
## [61] mclust_5.3 km.ci_0.5-2
## [63] htmlwidgets_0.9 httr_1.3.1
## [65] RColorBrewer_1.1-2 fpc_2.1-10
## [67] modeltools_0.2-21 pkgconfig_2.0.1
## [69] XML_3.98-1.9 R.methodsS3_1.7.1
## [71] flexmix_2.3-14 nnet_7.3-12
## [73] locfit_1.5-9.1 labeling_0.3
## [75] rlang_0.1.2 reshape2_1.4.2
## [77] AnnotationDbi_1.38.2 munsell_0.4.3
## [79] tools_3.4.1 downloader_0.4
## [81] RSQLite_2.0 devtools_1.13.3
## [83] broom_0.4.2 evaluate_0.10.1
## [85] stringr_1.2.0 yaml_2.1.14
## [87] knitr_1.17 bit64_0.9-7
## [89] robustbase_0.92-7 survMisc_0.5.4
## [91] purrr_0.2.3 dendextend_1.5.2
## [93] EDASeq_2.10.0 nlme_3.1-131
## [95] whisker_0.3-2 R.oo_1.21.0
## [97] xml2_1.1.1 biomaRt_2.32.1
## [99] BiocStyle_2.4.1 compiler_3.4.1
## [101] curl_2.8.1 testthat_1.0.2
## [103] tibble_1.3.4 geneplotter_1.54.0
## [105] stringi_1.1.5 highr_0.6
## [107] GenomicFeatures_1.28.4 lattice_0.20-35
## [109] trimcluster_0.1-2 Matrix_1.2-11
## [111] commonmark_1.4 psych_1.7.8
## [113] KMsurv_0.1-5 GlobalOptions_0.0.12
## [115] data.table_1.10.4 bitops_1.0-6
## [117] rtracklayer_1.36.4 R6_2.2.2
## [119] latticeExtra_0.6-28 hwriter_1.3.2
## [121] ShortRead_1.34.1 gridExtra_2.3
## [123] codetools_0.2-15 MASS_7.3-47
## [125] assertthat_0.2.0 rprojroot_1.2
## [127] rjson_0.2.15 withr_2.0.0
## [129] GenomicAlignments_1.12.2 Rsamtools_1.28.0
## [131] mnormt_1.5-5 GenomeInfoDbData_0.99.0
## [133] diptest_0.75-7 hms_0.3
## [135] tidyr_0.7.1 class_7.3-14
## [137] rmarkdown_1.6 ggpubr_0.1.5
Cancer Genome Atlas Research Network and others. 2012a. “Comprehensive Genomic Characterization of Squamous Cell Lung Cancers.” doi:10.1038/nature11404.
———. 2012b. “Comprehensive Molecular Characterization of Human Colon and Rectal Cancer.” doi:10.1038/nature11252.
———. 2012c. “Comprehensive Molecular Portraits of Human Breast Tumours.” doi:10.1038/nature11412.
———. 2013a. “Comprehensive Molecular Characterization of Clear Cell Renal Cell Carcinoma.” doi:10.1038/nature12222.
———. 2013b. “Integrated Genomic Characterization of Endometrial Carcinoma.” doi:10.1038/nature12113.
———. 2014a. “Comprehensive Molecular Characterization of Gastric Adenocarcinoma.” doi:10.1038/nature13480.
———. 2014b. “Comprehensive Molecular Profiling of Lung Adenocarcinoma.” doi:10.1038/nature13385.
———. 2014c. “Integrated Genomic Characterization of Papillary Thyroid Carcinoma.” doi:10.1016/j.cell.2014.09.050.
———. 2015a. “Comprehensive Genomic Characterization of Head and Neck Squamous Cell Carcinomas.” doi:10.1038/nature14129.
———. 2015b. “Genomic Classification of Cutaneous Melanoma.” doi:10.1016/j.cell.2015.05.044.
———. 2015c. “The Molecular Taxonomy of Primary Prostate Cancer.” doi:10.1016/j.cell.2015.10.025.
———. 2016. “Comprehensive Pan-Genomic Characterization of Adrenocortical Carcinoma.” doi:10.1016/j.ccell.2016.04.002.
Ceccarelli, Michele and Barthel, Floris P and Malta, Tathiane M and Sabedot, Thais S and Salama, Sofie R and Murray, Bradley A and Morozova, Olena and Newton, Yulia and Radenbaugh, Amie and Pagnotta, Stefano M and others. 2016. “Molecular Profiling Reveals Biologically Discrete Subsets and Pathways of Progression in Diffuse Glioma.” doi:10.1016/j.cell.2015.12.028.
Colaprico, Antonio and Silva, Tiago C. and Olsen, Catharina and Garofano, Luciano and Cava, Claudia and Garolini, Davide and Sabedot, Thais S. and Malta, Tathiane M. and Pagnotta, Stefano M. and Castiglioni, Isabella and Ceccarelli, Michele and Bontempi, Gianluca and Noushmehr, Houtan. 2016. “TCGAbiolinks: An R/Bioconductor Package for Integrative Analysis of TCGA Data.” doi:10.1093/nar/gkv1507.
Davis, Caleb F and Ricketts, Christopher J and Wang, Min and Yang, Lixing and Cherniack, Andrew D and Shen, Hui and Buhay, Christian and Kang, Hyojin and Kim, Sang Cheol and Fahey, Catherine C and others. 2014. “The Somatic Genomic Landscape of Chromophobe Renal Cell Carcinoma.” doi:10.1016/j.ccr.2014.07.014.
Fishbein, Lauren and Leshchiner, Ignaty and Walter, Vonn and Danilova, Ludmila and Robertson, A Gordon and Johnson, Amy R and Lichtenberg, Tara M and Murray, Bradley A and Ghayee, Hans K and Else, Tobias and others. 2017. “Comprehensive Molecular Characterization of Pheochromocytoma and Paraganglioma.” doi:10.1016/j.ccell.2017.01.001.
Gu, Zuguang and Eils, Roland and Schlesner, Matthias. 2016. “Complex Heatmaps Reveal Patterns and Correlations in Multidimensional Genomic Data.” doi:10.1016/j.ccell.2016.04.002.
Huber, Wolfgang and Carey, Vincent J and Gentleman, Robert and Anders, Simon and Carlson, Marc and Carvalho, Benilton S and Bravo, Hector Corrada and Davis, Sean and Gatto, Laurent and Girke, Thomas and others. 2015. “Orchestrating High-Throughput Genomic Analysis with Bioconductor.”
James H Bullard, Elizabeth Purdom, Kasper D Hansen and Sandrine Dudoit. 2010. “Evaluation of Statistical Methods for Normalization and Differential Expression in MRNA-Seq Experiments.”
Linehan, W Marston and Spellman, Paul T and Ricketts, Christopher J and Creighton, Chad J and Fei, Suzanne S and Davis, Caleb and Wheeler, David A and Murray, Bradley A and Schmidt, Laura and Vocke, Cathy D and others. 2016. “Comprehensive Molecular Characterization of Papillary Renal-Cell Carcinoma.” doi:10.1056/NEJMoa1505917.
Noushmehr, H., Weisenberger, D.J., Diefes, K., Phillips, H.S., Pujara, K., Berman, B.P., Pan, F., Pelloski, C.E., Sulman, E.P., Bhat, K.P. et al. 2010. “Identification of a CpG Island Methylator Phenotype That Defines a Distinct Subgroup of Glioma.”
Risso, D., Schwartz, K., Sherlock, G., & Dudoit, S. 2011. “GC-Content Normalization for RNA-Seq Data.”
Silva, TC and Colaprico, A and Olsen, C and D’Angelo, F and Bontempi, G and Ceccarelli, M and Noushmehr, H. 2016. “TCGA Workflow: Analyze Cancer Genomics and Epigenomics Data Using Bioconductor Packages.” doi:10.12688/f1000research.8923.1.